Actual source code: bvmat.c

slepc-3.15.2 2021-09-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    BV implemented with a dense Mat
 12: */

 14: #include <slepc/private/bvimpl.h>

 16: typedef struct {
 17:   Mat       A;
 18:   PetscBool mpi;
 19: } BV_MAT;

 21: PetscErrorCode BVMult_Mat(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
 22: {
 24:   BV_MAT            *y = (BV_MAT*)Y->data,*x = (BV_MAT*)X->data;
 25:   PetscScalar       *py;
 26:   const PetscScalar *px,*q;
 27:   PetscInt          ldq;

 30:   MatDenseGetArrayRead(x->A,&px);
 31:   MatDenseGetArray(y->A,&py);
 32:   if (Q) {
 33:     MatGetSize(Q,&ldq,NULL);
 34:     MatDenseGetArrayRead(Q,&q);
 35:     BVMult_BLAS_Private(Y,Y->n,Y->k-Y->l,X->k-X->l,ldq,alpha,px+(X->nc+X->l)*X->n,q+Y->l*ldq+X->l,beta,py+(Y->nc+Y->l)*Y->n);
 36:     MatDenseRestoreArrayRead(Q,&q);
 37:   } else {
 38:     BVAXPY_BLAS_Private(Y,Y->n,Y->k-Y->l,alpha,px+(X->nc+X->l)*X->n,beta,py+(Y->nc+Y->l)*Y->n);
 39:   }
 40:   MatDenseRestoreArrayRead(x->A,&px);
 41:   MatDenseRestoreArray(y->A,&py);
 42:   return(0);
 43: }

 45: PetscErrorCode BVMultVec_Mat(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
 46: {
 47:   PetscErrorCode    ierr;
 48:   BV_MAT            *x = (BV_MAT*)X->data;
 49:   PetscScalar       *py,*qq=q;
 50:   const PetscScalar *px;

 53:   MatDenseGetArrayRead(x->A,&px);
 54:   VecGetArray(y,&py);
 55:   if (!q) { VecGetArray(X->buffer,&qq); }
 56:   BVMultVec_BLAS_Private(X,X->n,X->k-X->l,alpha,px+(X->nc+X->l)*X->n,qq,beta,py);
 57:   if (!q) { VecRestoreArray(X->buffer,&qq); }
 58:   MatDenseRestoreArrayRead(x->A,&px);
 59:   VecRestoreArray(y,&py);
 60:   return(0);
 61: }

 63: PetscErrorCode BVMultInPlace_Mat(BV V,Mat Q,PetscInt s,PetscInt e)
 64: {
 65:   PetscErrorCode    ierr;
 66:   BV_MAT            *ctx = (BV_MAT*)V->data;
 67:   PetscScalar       *pv;
 68:   const PetscScalar *q;
 69:   PetscInt          ldq;

 72:   MatGetSize(Q,&ldq,NULL);
 73:   MatDenseGetArray(ctx->A,&pv);
 74:   MatDenseGetArrayRead(Q,&q);
 75:   BVMultInPlace_BLAS_Private(V,V->n,V->k-V->l,ldq,s-V->l,e-V->l,pv+(V->nc+V->l)*V->n,q+V->l*ldq+V->l,PETSC_FALSE);
 76:   MatDenseRestoreArrayRead(Q,&q);
 77:   MatDenseRestoreArray(ctx->A,&pv);
 78:   return(0);
 79: }

 81: PetscErrorCode BVMultInPlaceTranspose_Mat(BV V,Mat Q,PetscInt s,PetscInt e)
 82: {
 83:   PetscErrorCode    ierr;
 84:   BV_MAT            *ctx = (BV_MAT*)V->data;
 85:   PetscScalar       *pv;
 86:   const PetscScalar *q;
 87:   PetscInt          ldq;

 90:   MatGetSize(Q,&ldq,NULL);
 91:   MatDenseGetArray(ctx->A,&pv);
 92:   MatDenseGetArrayRead(Q,&q);
 93:   BVMultInPlace_BLAS_Private(V,V->n,V->k-V->l,ldq,s-V->l,e-V->l,pv+(V->nc+V->l)*V->n,q+V->l*ldq+V->l,PETSC_TRUE);
 94:   MatDenseRestoreArrayRead(Q,&q);
 95:   MatDenseRestoreArray(ctx->A,&pv);
 96:   return(0);
 97: }

 99: PetscErrorCode BVDot_Mat(BV X,BV Y,Mat M)
100: {
102:   BV_MAT            *x = (BV_MAT*)X->data,*y = (BV_MAT*)Y->data;
103:   PetscScalar       *m;
104:   const PetscScalar *px,*py;
105:   PetscInt          ldm;

108:   MatGetSize(M,&ldm,NULL);
109:   MatDenseGetArrayRead(x->A,&px);
110:   MatDenseGetArrayRead(y->A,&py);
111:   MatDenseGetArray(M,&m);
112:   BVDot_BLAS_Private(X,Y->k-Y->l,X->k-X->l,X->n,ldm,py+(Y->nc+Y->l)*Y->n,px+(X->nc+X->l)*X->n,m+X->l*ldm+Y->l,x->mpi);
113:   MatDenseRestoreArray(M,&m);
114:   MatDenseRestoreArrayRead(x->A,&px);
115:   MatDenseRestoreArrayRead(y->A,&py);
116:   return(0);
117: }

119: PetscErrorCode BVDotVec_Mat(BV X,Vec y,PetscScalar *q)
120: {
121:   PetscErrorCode    ierr;
122:   BV_MAT            *x = (BV_MAT*)X->data;
123:   PetscScalar       *qq=q;
124:   const PetscScalar *px,*py;
125:   Vec               z = y;

128:   if (X->matrix) {
129:     BV_IPMatMult(X,y);
130:     z = X->Bx;
131:   }
132:   MatDenseGetArrayRead(x->A,&px);
133:   VecGetArrayRead(z,&py);
134:   if (!q) { VecGetArray(X->buffer,&qq); }
135:   BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,qq,x->mpi);
136:   if (!q) { VecRestoreArray(X->buffer,&qq); }
137:   VecRestoreArrayRead(z,&py);
138:   MatDenseRestoreArrayRead(x->A,&px);
139:   return(0);
140: }

142: PetscErrorCode BVDotVec_Local_Mat(BV X,Vec y,PetscScalar *m)
143: {
144:   PetscErrorCode    ierr;
145:   BV_MAT            *x = (BV_MAT*)X->data;
146:   const PetscScalar *px,*py;
147:   Vec               z = y;

150:   if (X->matrix) {
151:     BV_IPMatMult(X,y);
152:     z = X->Bx;
153:   }
154:   MatDenseGetArrayRead(x->A,&px);
155:   VecGetArrayRead(z,&py);
156:   BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,m,PETSC_FALSE);
157:   VecRestoreArrayRead(z,&py);
158:   MatDenseRestoreArrayRead(x->A,&px);
159:   return(0);
160: }

162: PetscErrorCode BVScale_Mat(BV bv,PetscInt j,PetscScalar alpha)
163: {
165:   BV_MAT         *ctx = (BV_MAT*)bv->data;
166:   PetscScalar    *array;

169:   MatDenseGetArray(ctx->A,&array);
170:   if (j<0) {
171:     BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->n,array+(bv->nc+bv->l)*bv->n,alpha);
172:   } else {
173:     BVScale_BLAS_Private(bv,bv->n,array+(bv->nc+j)*bv->n,alpha);
174:   }
175:   MatDenseRestoreArray(ctx->A,&array);
176:   return(0);
177: }

179: PetscErrorCode BVNorm_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
180: {
181:   PetscErrorCode    ierr;
182:   BV_MAT            *ctx = (BV_MAT*)bv->data;
183:   const PetscScalar *array;

186:   MatDenseGetArrayRead(ctx->A,&array);
187:   if (j<0) {
188:     BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,ctx->mpi);
189:   } else {
190:     BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,ctx->mpi);
191:   }
192:   MatDenseRestoreArrayRead(ctx->A,&array);
193:   return(0);
194: }

196: PetscErrorCode BVNorm_Local_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
197: {
198:   PetscErrorCode    ierr;
199:   BV_MAT            *ctx = (BV_MAT*)bv->data;
200:   const PetscScalar *array;

203:   MatDenseGetArrayRead(ctx->A,&array);
204:   if (j<0) {
205:     BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,PETSC_FALSE);
206:   } else {
207:     BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,PETSC_FALSE);
208:   }
209:   MatDenseRestoreArrayRead(ctx->A,&array);
210:   return(0);
211: }

213: PetscErrorCode BVNormalize_Mat(BV bv,PetscScalar *eigi)
214: {
216:   BV_MAT         *ctx = (BV_MAT*)bv->data;
217:   PetscScalar    *array,*wi=NULL;

220:   MatDenseGetArray(ctx->A,&array);
221:   if (eigi) wi = eigi+bv->l;
222:   BVNormalize_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,wi,ctx->mpi);
223:   MatDenseRestoreArray(ctx->A,&array);
224:   return(0);
225: }

227: PetscErrorCode BVMatMult_Mat(BV V,Mat A,BV W)
228: {
230:   PetscInt       j;
231:   Mat            Vmat,Wmat;
232:   Vec            vv,ww;

235:   if (V->vmm) {
236:     BVGetMat(V,&Vmat);
237:     BVGetMat(W,&Wmat);
238:     MatProductCreateWithMat(A,Vmat,NULL,Wmat);
239:     MatProductSetType(Wmat,MATPRODUCT_AB);
240:     MatProductSetFromOptions(Wmat);
241:     MatProductSymbolic(Wmat);
242:     MatProductNumeric(Wmat);
243:     MatProductClear(Wmat);
244:     BVRestoreMat(V,&Vmat);
245:     BVRestoreMat(W,&Wmat);
246:   } else {
247:     for (j=0;j<V->k-V->l;j++) {
248:       BVGetColumn(V,V->l+j,&vv);
249:       BVGetColumn(W,W->l+j,&ww);
250:       MatMult(A,vv,ww);
251:       BVRestoreColumn(V,V->l+j,&vv);
252:       BVRestoreColumn(W,W->l+j,&ww);
253:     }
254:   }
255:   return(0);
256: }

258: PetscErrorCode BVCopy_Mat(BV V,BV W)
259: {
260:   PetscErrorCode    ierr;
261:   BV_MAT            *v = (BV_MAT*)V->data,*w = (BV_MAT*)W->data;
262:   PetscScalar       *pw,*pwc;
263:   const PetscScalar *pv,*pvc;

266:   MatDenseGetArrayRead(v->A,&pv);
267:   MatDenseGetArray(w->A,&pw);
268:   pvc = pv+(V->nc+V->l)*V->n;
269:   pwc = pw+(W->nc+W->l)*W->n;
270:   PetscArraycpy(pwc,pvc,(V->k-V->l)*V->n);
271:   MatDenseRestoreArrayRead(v->A,&pv);
272:   MatDenseRestoreArray(w->A,&pw);
273:   return(0);
274: }

276: PetscErrorCode BVCopyColumn_Mat(BV V,PetscInt j,PetscInt i)
277: {
279:   BV_MAT         *v = (BV_MAT*)V->data;
280:   PetscScalar    *pv;

283:   MatDenseGetArray(v->A,&pv);
284:   PetscArraycpy(pv+(V->nc+i)*V->n,pv+(V->nc+j)*V->n,V->n);
285:   MatDenseRestoreArray(v->A,&pv);
286:   return(0);
287: }

289: PetscErrorCode BVResize_Mat(BV bv,PetscInt m,PetscBool copy)
290: {
291:   PetscErrorCode    ierr;
292:   BV_MAT            *ctx = (BV_MAT*)bv->data;
293:   PetscScalar       *pnew;
294:   const PetscScalar *pA;
295:   Mat               A;
296:   char              str[50];

299:   MatCreateDense(PetscObjectComm((PetscObject)bv->t),bv->n,PETSC_DECIDE,PETSC_DECIDE,m,NULL,&A);
300:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
301:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
302:   PetscLogObjectParent((PetscObject)bv,(PetscObject)A);
303:   if (((PetscObject)bv)->name) {
304:     PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name);
305:     PetscObjectSetName((PetscObject)A,str);
306:   }
307:   if (copy) {
308:     MatDenseGetArrayRead(ctx->A,&pA);
309:     MatDenseGetArrayWrite(A,&pnew);
310:     PetscArraycpy(pnew,pA,PetscMin(m,bv->m)*bv->n);
311:     MatDenseRestoreArrayRead(ctx->A,&pA);
312:     MatDenseRestoreArrayWrite(A,&pnew);
313:   }
314:   MatDestroy(&ctx->A);
315:   ctx->A = A;
316:   return(0);
317: }

319: PetscErrorCode BVGetColumn_Mat(BV bv,PetscInt j,Vec *v)
320: {
322:   BV_MAT         *ctx = (BV_MAT*)bv->data;
323:   PetscScalar    *pA;
324:   PetscInt       l;

327:   l = BVAvailableVec;
328:   MatDenseGetArray(ctx->A,&pA);
329:   VecPlaceArray(bv->cv[l],pA+(bv->nc+j)*bv->n);
330:   return(0);
331: }

333: PetscErrorCode BVRestoreColumn_Mat(BV bv,PetscInt j,Vec *v)
334: {
336:   BV_MAT         *ctx = (BV_MAT*)bv->data;
337:   PetscScalar    *pA;
338:   PetscInt       l;

341:   l = (j==bv->ci[0])? 0: 1;
342:   VecResetArray(bv->cv[l]);
343:   MatDenseRestoreArray(ctx->A,&pA);
344:   return(0);
345: }

347: PetscErrorCode BVGetArray_Mat(BV bv,PetscScalar **a)
348: {
350:   BV_MAT         *ctx = (BV_MAT*)bv->data;

353:   MatDenseGetArray(ctx->A,a);
354:   return(0);
355: }

357: PetscErrorCode BVRestoreArray_Mat(BV bv,PetscScalar **a)
358: {
360:   BV_MAT         *ctx = (BV_MAT*)bv->data;

363:   if (a) { MatDenseRestoreArray(ctx->A,a); }
364:   return(0);
365: }

367: PetscErrorCode BVGetArrayRead_Mat(BV bv,const PetscScalar **a)
368: {
370:   BV_MAT         *ctx = (BV_MAT*)bv->data;

373:   MatDenseGetArrayRead(ctx->A,a);
374:   return(0);
375: }

377: PetscErrorCode BVRestoreArrayRead_Mat(BV bv,const PetscScalar **a)
378: {
380:   BV_MAT         *ctx = (BV_MAT*)bv->data;

383:   if (a) { MatDenseRestoreArrayRead(ctx->A,a); }
384:   return(0);
385: }

387: PetscErrorCode BVView_Mat(BV bv,PetscViewer viewer)
388: {
389:   PetscErrorCode    ierr;
390:   BV_MAT            *ctx = (BV_MAT*)bv->data;
391:   PetscViewerFormat format;
392:   PetscBool         isascii;
393:   const char        *bvname,*name;

396:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
397:   if (isascii) {
398:     PetscViewerGetFormat(viewer,&format);
399:     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
400:     MatView(ctx->A,viewer);
401:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
402:       PetscObjectGetName((PetscObject)bv,&bvname);
403:       PetscObjectGetName((PetscObject)ctx->A,&name);
404:       PetscViewerASCIIPrintf(viewer,"%s=%s;clear %s\n",bvname,name,name);
405:       if (bv->nc) {
406:         PetscViewerASCIIPrintf(viewer,"%s=%s(:,%D:end);\n",bvname,bvname,bv->nc+1);
407:       }
408:     }
409:   } else {
410:     MatView(ctx->A,viewer);
411:   }
412:   return(0);
413: }

415: PetscErrorCode BVDestroy_Mat(BV bv)
416: {
418:   BV_MAT         *ctx = (BV_MAT*)bv->data;

421:   MatDestroy(&ctx->A);
422:   VecDestroy(&bv->cv[0]);
423:   VecDestroy(&bv->cv[1]);
424:   PetscFree(bv->data);
425:   return(0);
426: }

428: SLEPC_EXTERN PetscErrorCode BVCreate_Mat(BV bv)
429: {
431:   BV_MAT         *ctx;
432:   PetscInt       nloc,bs,lsplit;
433:   PetscBool      seq;
434:   char           str[50];
435:   PetscScalar    *array,*ptr;
436:   BV             parent;

439:   PetscNewLog(bv,&ctx);
440:   bv->data = (void*)ctx;

442:   PetscObjectTypeCompare((PetscObject)bv->t,VECMPI,&ctx->mpi);
443:   if (!ctx->mpi) {
444:     PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
445:     if (!seq) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Cannot create a BVMAT from a non-standard template vector");
446:   }

448:   VecGetLocalSize(bv->t,&nloc);
449:   VecGetBlockSize(bv->t,&bs);

451:   if (bv->issplit) {
452:     /* split BV: share the memory of the parent BV */
453:     parent = bv->splitparent;
454:     lsplit = parent->lsplit;
455:     MatDenseGetArray(((BV_MAT*)parent->data)->A,&array);
456:     ptr = (bv->issplit==1)? array: array+lsplit*nloc;
457:     MatDenseRestoreArray(((BV_MAT*)parent->data)->A,&array);
458:   } else {
459:     /* regular BV: allocate memory for the BV entries */
460:     ptr = NULL;
461:   }
462:   MatCreateDense(PetscObjectComm((PetscObject)bv->t),nloc,PETSC_DECIDE,PETSC_DECIDE,bv->m,ptr,&ctx->A);
463:   MatAssemblyBegin(ctx->A,MAT_FINAL_ASSEMBLY);
464:   MatAssemblyEnd(ctx->A,MAT_FINAL_ASSEMBLY);
465:   PetscLogObjectParent((PetscObject)bv,(PetscObject)ctx->A);
466:   if (((PetscObject)bv)->name) {
467:     PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name);
468:     PetscObjectSetName((PetscObject)ctx->A,str);
469:   }

471:   if (bv->Acreate) {
472:     MatCopy(bv->Acreate,ctx->A,SAME_NONZERO_PATTERN);
473:     MatDestroy(&bv->Acreate);
474:   }

476:   VecDuplicateEmpty(bv->t,&bv->cv[0]);
477:   VecDuplicateEmpty(bv->t,&bv->cv[1]);

479:   bv->ops->mult             = BVMult_Mat;
480:   bv->ops->multvec          = BVMultVec_Mat;
481:   bv->ops->multinplace      = BVMultInPlace_Mat;
482:   bv->ops->multinplacetrans = BVMultInPlaceTranspose_Mat;
483:   bv->ops->dot              = BVDot_Mat;
484:   bv->ops->dotvec           = BVDotVec_Mat;
485:   bv->ops->dotvec_local     = BVDotVec_Local_Mat;
486:   bv->ops->scale            = BVScale_Mat;
487:   bv->ops->norm             = BVNorm_Mat;
488:   bv->ops->norm_local       = BVNorm_Local_Mat;
489:   bv->ops->normalize        = BVNormalize_Mat;
490:   bv->ops->matmult          = BVMatMult_Mat;
491:   bv->ops->copy             = BVCopy_Mat;
492:   bv->ops->copycolumn       = BVCopyColumn_Mat;
493:   bv->ops->resize           = BVResize_Mat;
494:   bv->ops->getcolumn        = BVGetColumn_Mat;
495:   bv->ops->restorecolumn    = BVRestoreColumn_Mat;
496:   bv->ops->getarray         = BVGetArray_Mat;
497:   bv->ops->restorearray     = BVRestoreArray_Mat;
498:   bv->ops->getarrayread     = BVGetArrayRead_Mat;
499:   bv->ops->restorearrayread = BVRestoreArrayRead_Mat;
500:   bv->ops->destroy          = BVDestroy_Mat;
501:   if (!ctx->mpi) bv->ops->view = BVView_Mat;
502:   return(0);
503: }