Actual source code: bvmat.c

slepc-3.18.3 2023-03-24
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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: {
 23:   BV_MAT            *y = (BV_MAT*)Y->data,*x = (BV_MAT*)X->data;
 24:   PetscScalar       *py;
 25:   const PetscScalar *px,*q;
 26:   PetscInt          ldq;

 28:   MatDenseGetArrayRead(x->A,&px);
 29:   MatDenseGetArray(y->A,&py);
 30:   if (Q) {
 31:     MatDenseGetLDA(Q,&ldq);
 32:     MatDenseGetArrayRead(Q,&q);
 33:     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);
 34:     MatDenseRestoreArrayRead(Q,&q);
 35:   } else 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);
 36:   MatDenseRestoreArrayRead(x->A,&px);
 37:   MatDenseRestoreArray(y->A,&py);
 38:   return 0;
 39: }

 41: PetscErrorCode BVMultVec_Mat(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
 42: {
 43:   BV_MAT            *x = (BV_MAT*)X->data;
 44:   PetscScalar       *py,*qq=q;
 45:   const PetscScalar *px;

 47:   MatDenseGetArrayRead(x->A,&px);
 48:   VecGetArray(y,&py);
 49:   if (!q) VecGetArray(X->buffer,&qq);
 50:   BVMultVec_BLAS_Private(X,X->n,X->k-X->l,alpha,px+(X->nc+X->l)*X->n,qq,beta,py);
 51:   if (!q) VecRestoreArray(X->buffer,&qq);
 52:   MatDenseRestoreArrayRead(x->A,&px);
 53:   VecRestoreArray(y,&py);
 54:   return 0;
 55: }

 57: PetscErrorCode BVMultInPlace_Mat(BV V,Mat Q,PetscInt s,PetscInt e)
 58: {
 59:   BV_MAT            *ctx = (BV_MAT*)V->data;
 60:   PetscScalar       *pv;
 61:   const PetscScalar *q;
 62:   PetscInt          ldq;

 64:   MatDenseGetLDA(Q,&ldq);
 65:   MatDenseGetArray(ctx->A,&pv);
 66:   MatDenseGetArrayRead(Q,&q);
 67:   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);
 68:   MatDenseRestoreArrayRead(Q,&q);
 69:   MatDenseRestoreArray(ctx->A,&pv);
 70:   return 0;
 71: }

 73: PetscErrorCode BVMultInPlaceHermitianTranspose_Mat(BV V,Mat Q,PetscInt s,PetscInt e)
 74: {
 75:   BV_MAT            *ctx = (BV_MAT*)V->data;
 76:   PetscScalar       *pv;
 77:   const PetscScalar *q;
 78:   PetscInt          ldq;

 80:   MatDenseGetLDA(Q,&ldq);
 81:   MatDenseGetArray(ctx->A,&pv);
 82:   MatDenseGetArrayRead(Q,&q);
 83:   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);
 84:   MatDenseRestoreArrayRead(Q,&q);
 85:   MatDenseRestoreArray(ctx->A,&pv);
 86:   return 0;
 87: }

 89: PetscErrorCode BVDot_Mat(BV X,BV Y,Mat M)
 90: {
 91:   BV_MAT            *x = (BV_MAT*)X->data,*y = (BV_MAT*)Y->data;
 92:   PetscScalar       *m;
 93:   const PetscScalar *px,*py;
 94:   PetscInt          ldm;

 96:   MatDenseGetLDA(M,&ldm);
 97:   MatDenseGetArrayRead(x->A,&px);
 98:   MatDenseGetArrayRead(y->A,&py);
 99:   MatDenseGetArray(M,&m);
100:   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);
101:   MatDenseRestoreArray(M,&m);
102:   MatDenseRestoreArrayRead(x->A,&px);
103:   MatDenseRestoreArrayRead(y->A,&py);
104:   return 0;
105: }

107: PetscErrorCode BVDotVec_Mat(BV X,Vec y,PetscScalar *q)
108: {
109:   BV_MAT            *x = (BV_MAT*)X->data;
110:   PetscScalar       *qq=q;
111:   const PetscScalar *px,*py;
112:   Vec               z = y;

114:   if (PetscUnlikely(X->matrix)) {
115:     BV_IPMatMult(X,y);
116:     z = X->Bx;
117:   }
118:   MatDenseGetArrayRead(x->A,&px);
119:   VecGetArrayRead(z,&py);
120:   if (!q) VecGetArray(X->buffer,&qq);
121:   BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,qq,x->mpi);
122:   if (!q) VecRestoreArray(X->buffer,&qq);
123:   VecRestoreArrayRead(z,&py);
124:   MatDenseRestoreArrayRead(x->A,&px);
125:   return 0;
126: }

128: PetscErrorCode BVDotVec_Local_Mat(BV X,Vec y,PetscScalar *m)
129: {
130:   BV_MAT            *x = (BV_MAT*)X->data;
131:   const PetscScalar *px,*py;
132:   Vec               z = y;

134:   if (PetscUnlikely(X->matrix)) {
135:     BV_IPMatMult(X,y);
136:     z = X->Bx;
137:   }
138:   MatDenseGetArrayRead(x->A,&px);
139:   VecGetArrayRead(z,&py);
140:   BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,m,PETSC_FALSE);
141:   VecRestoreArrayRead(z,&py);
142:   MatDenseRestoreArrayRead(x->A,&px);
143:   return 0;
144: }

146: PetscErrorCode BVScale_Mat(BV bv,PetscInt j,PetscScalar alpha)
147: {
148:   BV_MAT         *ctx = (BV_MAT*)bv->data;
149:   PetscScalar    *array;

151:   MatDenseGetArray(ctx->A,&array);
152:   if (PetscUnlikely(j<0)) BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->n,array+(bv->nc+bv->l)*bv->n,alpha);
153:   else BVScale_BLAS_Private(bv,bv->n,array+(bv->nc+j)*bv->n,alpha);
154:   MatDenseRestoreArray(ctx->A,&array);
155:   return 0;
156: }

158: PetscErrorCode BVNorm_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
159: {
160:   BV_MAT            *ctx = (BV_MAT*)bv->data;
161:   const PetscScalar *array;

163:   MatDenseGetArrayRead(ctx->A,&array);
164:   if (PetscUnlikely(j<0)) BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,ctx->mpi);
165:   else BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,ctx->mpi);
166:   MatDenseRestoreArrayRead(ctx->A,&array);
167:   return 0;
168: }

170: PetscErrorCode BVNorm_Local_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
171: {
172:   BV_MAT            *ctx = (BV_MAT*)bv->data;
173:   const PetscScalar *array;

175:   MatDenseGetArrayRead(ctx->A,&array);
176:   if (PetscUnlikely(j<0)) BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,PETSC_FALSE);
177:   else BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,PETSC_FALSE);
178:   MatDenseRestoreArrayRead(ctx->A,&array);
179:   return 0;
180: }

182: PetscErrorCode BVNormalize_Mat(BV bv,PetscScalar *eigi)
183: {
184:   BV_MAT         *ctx = (BV_MAT*)bv->data;
185:   PetscScalar    *array,*wi=NULL;

187:   MatDenseGetArray(ctx->A,&array);
188:   if (eigi) wi = eigi+bv->l;
189:   BVNormalize_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,wi,ctx->mpi);
190:   MatDenseRestoreArray(ctx->A,&array);
191:   return 0;
192: }

194: PetscErrorCode BVMatMult_Mat(BV V,Mat A,BV W)
195: {
196:   PetscInt       j;
197:   Mat            Vmat,Wmat;
198:   Vec            vv,ww;

200:   if (V->vmm) {
201:     BVGetMat(V,&Vmat);
202:     BVGetMat(W,&Wmat);
203:     MatProductCreateWithMat(A,Vmat,NULL,Wmat);
204:     MatProductSetType(Wmat,MATPRODUCT_AB);
205:     MatProductSetFromOptions(Wmat);
206:     MatProductSymbolic(Wmat);
207:     MatProductNumeric(Wmat);
208:     MatProductClear(Wmat);
209:     BVRestoreMat(V,&Vmat);
210:     BVRestoreMat(W,&Wmat);
211:   } else {
212:     for (j=0;j<V->k-V->l;j++) {
213:       BVGetColumn(V,V->l+j,&vv);
214:       BVGetColumn(W,W->l+j,&ww);
215:       MatMult(A,vv,ww);
216:       BVRestoreColumn(V,V->l+j,&vv);
217:       BVRestoreColumn(W,W->l+j,&ww);
218:     }
219:   }
220:   return 0;
221: }

223: PetscErrorCode BVCopy_Mat(BV V,BV W)
224: {
225:   BV_MAT            *v = (BV_MAT*)V->data,*w = (BV_MAT*)W->data;
226:   PetscScalar       *pw,*pwc;
227:   const PetscScalar *pv,*pvc;

229:   MatDenseGetArrayRead(v->A,&pv);
230:   MatDenseGetArray(w->A,&pw);
231:   pvc = pv+(V->nc+V->l)*V->n;
232:   pwc = pw+(W->nc+W->l)*W->n;
233:   PetscArraycpy(pwc,pvc,(V->k-V->l)*V->n);
234:   MatDenseRestoreArrayRead(v->A,&pv);
235:   MatDenseRestoreArray(w->A,&pw);
236:   return 0;
237: }

239: PetscErrorCode BVCopyColumn_Mat(BV V,PetscInt j,PetscInt i)
240: {
241:   BV_MAT         *v = (BV_MAT*)V->data;
242:   PetscScalar    *pv;

244:   MatDenseGetArray(v->A,&pv);
245:   PetscArraycpy(pv+(V->nc+i)*V->n,pv+(V->nc+j)*V->n,V->n);
246:   MatDenseRestoreArray(v->A,&pv);
247:   return 0;
248: }

250: PetscErrorCode BVResize_Mat(BV bv,PetscInt m,PetscBool copy)
251: {
252:   BV_MAT            *ctx = (BV_MAT*)bv->data;
253:   PetscScalar       *pnew;
254:   const PetscScalar *pA;
255:   Mat               A;
256:   char              str[50];

258:   MatCreateDense(PetscObjectComm((PetscObject)bv->t),bv->n,PETSC_DECIDE,PETSC_DECIDE,m,NULL,&A);
259:   if (((PetscObject)bv)->name) {
260:     PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name);
261:     PetscObjectSetName((PetscObject)A,str);
262:   }
263:   if (copy) {
264:     MatDenseGetArrayRead(ctx->A,&pA);
265:     MatDenseGetArrayWrite(A,&pnew);
266:     PetscArraycpy(pnew,pA,PetscMin(m,bv->m)*bv->n);
267:     MatDenseRestoreArrayRead(ctx->A,&pA);
268:     MatDenseRestoreArrayWrite(A,&pnew);
269:   }
270:   MatDestroy(&ctx->A);
271:   ctx->A = A;
272:   return 0;
273: }

275: PetscErrorCode BVGetColumn_Mat(BV bv,PetscInt j,Vec *v)
276: {
277:   BV_MAT         *ctx = (BV_MAT*)bv->data;
278:   PetscScalar    *pA;
279:   PetscInt       l;

281:   l = BVAvailableVec;
282:   MatDenseGetArray(ctx->A,&pA);
283:   VecPlaceArray(bv->cv[l],pA+(bv->nc+j)*bv->n);
284:   return 0;
285: }

287: PetscErrorCode BVRestoreColumn_Mat(BV bv,PetscInt j,Vec *v)
288: {
289:   BV_MAT         *ctx = (BV_MAT*)bv->data;
290:   PetscScalar    *pA;
291:   PetscInt       l;

293:   l = (j==bv->ci[0])? 0: 1;
294:   VecResetArray(bv->cv[l]);
295:   MatDenseRestoreArray(ctx->A,&pA);
296:   return 0;
297: }

299: PetscErrorCode BVGetArray_Mat(BV bv,PetscScalar **a)
300: {
301:   BV_MAT         *ctx = (BV_MAT*)bv->data;

303:   MatDenseGetArray(ctx->A,a);
304:   return 0;
305: }

307: PetscErrorCode BVRestoreArray_Mat(BV bv,PetscScalar **a)
308: {
309:   BV_MAT         *ctx = (BV_MAT*)bv->data;

311:   if (a) MatDenseRestoreArray(ctx->A,a);
312:   return 0;
313: }

315: PetscErrorCode BVGetArrayRead_Mat(BV bv,const PetscScalar **a)
316: {
317:   BV_MAT         *ctx = (BV_MAT*)bv->data;

319:   MatDenseGetArrayRead(ctx->A,a);
320:   return 0;
321: }

323: PetscErrorCode BVRestoreArrayRead_Mat(BV bv,const PetscScalar **a)
324: {
325:   BV_MAT         *ctx = (BV_MAT*)bv->data;

327:   if (a) MatDenseRestoreArrayRead(ctx->A,a);
328:   return 0;
329: }

331: PetscErrorCode BVView_Mat(BV bv,PetscViewer viewer)
332: {
333:   BV_MAT            *ctx = (BV_MAT*)bv->data;
334:   PetscViewerFormat format;
335:   PetscBool         isascii;
336:   const char        *bvname,*name;

338:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
339:   if (isascii) {
340:     PetscViewerGetFormat(viewer,&format);
341:     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return 0;
342:     MatView(ctx->A,viewer);
343:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
344:       PetscObjectGetName((PetscObject)bv,&bvname);
345:       PetscObjectGetName((PetscObject)ctx->A,&name);
346:       PetscViewerASCIIPrintf(viewer,"%s=%s;clear %s\n",bvname,name,name);
347:       if (bv->nc) PetscViewerASCIIPrintf(viewer,"%s=%s(:,%" PetscInt_FMT ":end);\n",bvname,bvname,bv->nc+1);
348:     }
349:   } else MatView(ctx->A,viewer);
350:   return 0;
351: }

353: PetscErrorCode BVDestroy_Mat(BV bv)
354: {
355:   BV_MAT         *ctx = (BV_MAT*)bv->data;

357:   MatDestroy(&ctx->A);
358:   VecDestroy(&bv->cv[0]);
359:   VecDestroy(&bv->cv[1]);
360:   PetscFree(bv->data);
361:   return 0;
362: }

364: SLEPC_EXTERN PetscErrorCode BVCreate_Mat(BV bv)
365: {
366:   BV_MAT         *ctx;
367:   PetscInt       nloc,bs,lsplit;
368:   PetscBool      seq;
369:   char           str[50];
370:   PetscScalar    *array,*ptr;
371:   BV             parent;

373:   PetscNew(&ctx);
374:   bv->data = (void*)ctx;

376:   PetscObjectTypeCompare((PetscObject)bv->t,VECMPI,&ctx->mpi);
377:   if (!ctx->mpi) {
378:     PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
380:   }

382:   VecGetLocalSize(bv->t,&nloc);
383:   VecGetBlockSize(bv->t,&bs);

385:   if (PetscUnlikely(bv->issplit)) {
386:     /* split BV: share the memory of the parent BV */
387:     parent = bv->splitparent;
388:     lsplit = parent->lsplit;
389:     MatDenseGetArray(((BV_MAT*)parent->data)->A,&array);
390:     ptr = (bv->issplit==1)? array: array+lsplit*nloc;
391:     MatDenseRestoreArray(((BV_MAT*)parent->data)->A,&array);
392:   } else {
393:     /* regular BV: allocate memory for the BV entries */
394:     ptr = NULL;
395:   }
396:   MatCreateDense(PetscObjectComm((PetscObject)bv->t),nloc,PETSC_DECIDE,PETSC_DECIDE,bv->m,ptr,&ctx->A);
397:   if (((PetscObject)bv)->name) {
398:     PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name);
399:     PetscObjectSetName((PetscObject)ctx->A,str);
400:   }

402:   if (PetscUnlikely(bv->Acreate)) {
403:     MatCopy(bv->Acreate,ctx->A,SAME_NONZERO_PATTERN);
404:     MatDestroy(&bv->Acreate);
405:   }

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

410:   bv->ops->mult             = BVMult_Mat;
411:   bv->ops->multvec          = BVMultVec_Mat;
412:   bv->ops->multinplace      = BVMultInPlace_Mat;
413:   bv->ops->multinplacetrans = BVMultInPlaceHermitianTranspose_Mat;
414:   bv->ops->dot              = BVDot_Mat;
415:   bv->ops->dotvec           = BVDotVec_Mat;
416:   bv->ops->dotvec_local     = BVDotVec_Local_Mat;
417:   bv->ops->scale            = BVScale_Mat;
418:   bv->ops->norm             = BVNorm_Mat;
419:   bv->ops->norm_local       = BVNorm_Local_Mat;
420:   bv->ops->normalize        = BVNormalize_Mat;
421:   bv->ops->matmult          = BVMatMult_Mat;
422:   bv->ops->copy             = BVCopy_Mat;
423:   bv->ops->copycolumn       = BVCopyColumn_Mat;
424:   bv->ops->resize           = BVResize_Mat;
425:   bv->ops->getcolumn        = BVGetColumn_Mat;
426:   bv->ops->restorecolumn    = BVRestoreColumn_Mat;
427:   bv->ops->getarray         = BVGetArray_Mat;
428:   bv->ops->restorearray     = BVRestoreArray_Mat;
429:   bv->ops->getarrayread     = BVGetArrayRead_Mat;
430:   bv->ops->restorearrayread = BVRestoreArrayRead_Mat;
431:   bv->ops->getmat           = BVGetMat_Default;
432:   bv->ops->restoremat       = BVRestoreMat_Default;
433:   bv->ops->destroy          = BVDestroy_Mat;
434:   if (!ctx->mpi) bv->ops->view = BVView_Mat;
435:   return 0;
436: }