Actual source code: svec.c

slepc-3.13.0 2020-03-31
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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 as a single Vec
 12: */

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

 17: PetscErrorCode BVMult_Svec(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
 18: {
 19:   PetscErrorCode    ierr;
 20:   BV_SVEC           *y = (BV_SVEC*)Y->data,*x = (BV_SVEC*)X->data;
 21:   const PetscScalar *px;
 22:   PetscScalar       *py,*q;
 23:   PetscInt          ldq;

 26:   VecGetArrayRead(x->v,&px);
 27:   VecGetArray(y->v,&py);
 28:   if (Q) {
 29:     MatGetSize(Q,&ldq,NULL);
 30:     MatDenseGetArray(Q,&q);
 31:     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);
 32:     MatDenseRestoreArray(Q,&q);
 33:   } else {
 34:     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);
 35:   }
 36:   VecRestoreArrayRead(x->v,&px);
 37:   VecRestoreArray(y->v,&py);
 38:   return(0);
 39: }

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

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

 58: PetscErrorCode BVMultInPlace_Svec(BV V,Mat Q,PetscInt s,PetscInt e)
 59: {
 61:   BV_SVEC        *ctx = (BV_SVEC*)V->data;
 62:   PetscScalar    *pv,*q;
 63:   PetscInt       ldq;

 66:   MatGetSize(Q,&ldq,NULL);
 67:   VecGetArray(ctx->v,&pv);
 68:   MatDenseGetArray(Q,&q);
 69:   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);
 70:   MatDenseRestoreArray(Q,&q);
 71:   VecRestoreArray(ctx->v,&pv);
 72:   return(0);
 73: }

 75: PetscErrorCode BVMultInPlaceTranspose_Svec(BV V,Mat Q,PetscInt s,PetscInt e)
 76: {
 78:   BV_SVEC        *ctx = (BV_SVEC*)V->data;
 79:   PetscScalar    *pv,*q;
 80:   PetscInt       ldq;

 83:   MatGetSize(Q,&ldq,NULL);
 84:   VecGetArray(ctx->v,&pv);
 85:   MatDenseGetArray(Q,&q);
 86:   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);
 87:   MatDenseRestoreArray(Q,&q);
 88:   VecRestoreArray(ctx->v,&pv);
 89:   return(0);
 90: }

 92: PetscErrorCode BVDot_Svec(BV X,BV Y,Mat M)
 93: {
 94:   PetscErrorCode    ierr;
 95:   BV_SVEC           *x = (BV_SVEC*)X->data,*y = (BV_SVEC*)Y->data;
 96:   const PetscScalar *px,*py;
 97:   PetscScalar       *m;
 98:   PetscInt          ldm;

101:   MatGetSize(M,&ldm,NULL);
102:   VecGetArrayRead(x->v,&px);
103:   VecGetArrayRead(y->v,&py);
104:   MatDenseGetArray(M,&m);
105:   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);
106:   MatDenseRestoreArray(M,&m);
107:   VecRestoreArrayRead(x->v,&px);
108:   VecRestoreArrayRead(y->v,&py);
109:   return(0);
110: }

112: PetscErrorCode BVDotVec_Svec(BV X,Vec y,PetscScalar *q)
113: {
114:   PetscErrorCode    ierr;
115:   BV_SVEC           *x = (BV_SVEC*)X->data;
116:   const PetscScalar *px,*py;
117:   PetscScalar       *qq=q;
118:   Vec               z = y;

121:   if (X->matrix) {
122:     BV_IPMatMult(X,y);
123:     z = X->Bx;
124:   }
125:   VecGetArrayRead(x->v,&px);
126:   VecGetArrayRead(z,&py);
127:   if (!q) { VecGetArray(X->buffer,&qq); }
128:   BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,qq,x->mpi);
129:   if (!q) { VecRestoreArray(X->buffer,&qq); }
130:   VecRestoreArrayRead(z,&py);
131:   VecRestoreArrayRead(x->v,&px);
132:   return(0);
133: }

135: PetscErrorCode BVDotVec_Local_Svec(BV X,Vec y,PetscScalar *m)
136: {
138:   BV_SVEC        *x = (BV_SVEC*)X->data;
139:   PetscScalar    *px,*py;
140:   Vec            z = y;

143:   if (X->matrix) {
144:     BV_IPMatMult(X,y);
145:     z = X->Bx;
146:   }
147:   VecGetArray(x->v,&px);
148:   VecGetArray(z,&py);
149:   BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,m,PETSC_FALSE);
150:   VecRestoreArray(z,&py);
151:   VecRestoreArray(x->v,&px);
152:   return(0);
153: }

155: PetscErrorCode BVScale_Svec(BV bv,PetscInt j,PetscScalar alpha)
156: {
158:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;
159:   PetscScalar    *array;

162:   VecGetArray(ctx->v,&array);
163:   if (j<0) {
164:     BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->n,array+(bv->nc+bv->l)*bv->n,alpha);
165:   } else {
166:     BVScale_BLAS_Private(bv,bv->n,array+(bv->nc+j)*bv->n,alpha);
167:   }
168:   VecRestoreArray(ctx->v,&array);
169:   return(0);
170: }

172: PetscErrorCode BVNorm_Svec(BV bv,PetscInt j,NormType type,PetscReal *val)
173: {
175:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;
176:   PetscScalar    *array;

179:   VecGetArray(ctx->v,&array);
180:   if (j<0) {
181:     BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,ctx->mpi);
182:   } else {
183:     BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,ctx->mpi);
184:   }
185:   VecRestoreArray(ctx->v,&array);
186:   return(0);
187: }

189: PetscErrorCode BVNorm_Local_Svec(BV bv,PetscInt j,NormType type,PetscReal *val)
190: {
192:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;
193:   PetscScalar    *array;

196:   VecGetArray(ctx->v,&array);
197:   if (j<0) {
198:     BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,PETSC_FALSE);
199:   } else {
200:     BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,PETSC_FALSE);
201:   }
202:   VecRestoreArray(ctx->v,&array);
203:   return(0);
204: }

206: PetscErrorCode BVMatMult_Svec(BV V,Mat A,BV W)
207: {
209:   PetscInt       j;
210:   PetscBool      flg;
211:   Mat            Vmat,Wmat;
212:   Vec            vv,ww;

215:   MatHasOperation(A,MATOP_MAT_MULT,&flg);
216:   if (V->vmm && flg) {
217:     BVGetMat(V,&Vmat);
218:     BVGetMat(W,&Wmat);
219:     MatProductCreateWithMat(A,Vmat,NULL,Wmat);
220:     MatProductSetType(Wmat,MATPRODUCT_AB);
221:     MatProductSetFromOptions(Wmat);
222:     MatProductSymbolic(Wmat);
223:     MatProductNumeric(Wmat);
224:     BVRestoreMat(V,&Vmat);
225:     BVRestoreMat(W,&Wmat);
226:   } else {
227:     for (j=0;j<V->k-V->l;j++) {
228:       BVGetColumn(V,V->l+j,&vv);
229:       BVGetColumn(W,W->l+j,&ww);
230:       MatMult(A,vv,ww);
231:       BVRestoreColumn(V,V->l+j,&vv);
232:       BVRestoreColumn(W,W->l+j,&ww);
233:     }
234:   }
235:   return(0);
236: }

238: PetscErrorCode BVCopy_Svec(BV V,BV W)
239: {
241:   BV_SVEC        *v = (BV_SVEC*)V->data,*w = (BV_SVEC*)W->data;
242:   PetscScalar    *pv,*pw,*pvc,*pwc;

245:   VecGetArray(v->v,&pv);
246:   VecGetArray(w->v,&pw);
247:   pvc = pv+(V->nc+V->l)*V->n;
248:   pwc = pw+(W->nc+W->l)*W->n;
249:   PetscArraycpy(pwc,pvc,(V->k-V->l)*V->n);
250:   VecRestoreArray(v->v,&pv);
251:   VecRestoreArray(w->v,&pw);
252:   return(0);
253: }

255: PetscErrorCode BVCopyColumn_Svec(BV V,PetscInt j,PetscInt i)
256: {
258:   BV_SVEC        *v = (BV_SVEC*)V->data;
259:   PetscScalar    *pv;

262:   VecGetArray(v->v,&pv);
263:   PetscArraycpy(pv+(V->nc+i)*V->n,pv+(V->nc+j)*V->n,V->n);
264:   VecRestoreArray(v->v,&pv);
265:   return(0);
266: }

268: PetscErrorCode BVResize_Svec(BV bv,PetscInt m,PetscBool copy)
269: {
270:   PetscErrorCode    ierr;
271:   BV_SVEC           *ctx = (BV_SVEC*)bv->data;
272:   PetscScalar       *pnew;
273:   const PetscScalar *pv;
274:   PetscInt          bs,lsplit;
275:   Vec               vnew,vpar;
276:   char              str[50];
277:   BV                parent;

280:   if (bv->issplit==2) {
281:     parent = bv->splitparent;
282:     lsplit = parent->lsplit;
283:     vpar = ((BV_SVEC*)parent->data)->v;
284:     VecGetArrayRead(vpar,&pv);
285:     VecResetArray(ctx->v);
286:     VecPlaceArray(ctx->v,pv+lsplit*bv->n);
287:     VecRestoreArrayRead(vpar,&pv);
288:   } else if (!bv->issplit) {
289:     VecGetBlockSize(bv->t,&bs);
290:     VecCreate(PetscObjectComm((PetscObject)bv->t),&vnew);
291:     VecSetType(vnew,((PetscObject)bv->t)->type_name);
292:     VecSetSizes(vnew,m*bv->n,PETSC_DECIDE);
293:     VecSetBlockSize(vnew,bs);
294:     PetscLogObjectParent((PetscObject)bv,(PetscObject)vnew);
295:     if (((PetscObject)bv)->name) {
296:       PetscSNPrintf(str,50,"%s_0",((PetscObject)bv)->name);
297:       PetscObjectSetName((PetscObject)vnew,str);
298:     }
299:     if (copy) {
300:       VecGetArrayRead(ctx->v,&pv);
301:       VecGetArray(vnew,&pnew);
302:       PetscArraycpy(pnew,pv,PetscMin(m,bv->m)*bv->n);
303:       VecRestoreArrayRead(ctx->v,&pv);
304:       VecRestoreArray(vnew,&pnew);
305:     }
306:     VecDestroy(&ctx->v);
307:     ctx->v = vnew;
308:   }
309:   return(0);
310: }

312: PetscErrorCode BVGetColumn_Svec(BV bv,PetscInt j,Vec *v)
313: {
315:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;
316:   PetscScalar    *pv;
317:   PetscInt       l;

320:   l = BVAvailableVec;
321:   VecGetArray(ctx->v,&pv);
322:   VecPlaceArray(bv->cv[l],pv+(bv->nc+j)*bv->n);
323:   return(0);
324: }

326: PetscErrorCode BVRestoreColumn_Svec(BV bv,PetscInt j,Vec *v)
327: {
329:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;
330:   PetscInt       l;

333:   l = (j==bv->ci[0])? 0: 1;
334:   VecResetArray(bv->cv[l]);
335:   VecRestoreArray(ctx->v,NULL);
336:   return(0);
337: }

339: PetscErrorCode BVGetArray_Svec(BV bv,PetscScalar **a)
340: {
342:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;

345:   VecGetArray(ctx->v,a);
346:   return(0);
347: }

349: PetscErrorCode BVRestoreArray_Svec(BV bv,PetscScalar **a)
350: {
352:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;

355:   VecRestoreArray(ctx->v,a);
356:   return(0);
357: }

359: PetscErrorCode BVGetArrayRead_Svec(BV bv,const PetscScalar **a)
360: {
362:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;

365:   VecGetArrayRead(ctx->v,a);
366:   return(0);
367: }

369: PetscErrorCode BVRestoreArrayRead_Svec(BV bv,const PetscScalar **a)
370: {
372:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;

375:   VecRestoreArrayRead(ctx->v,a);
376:   return(0);
377: }

379: PetscErrorCode BVView_Svec(BV bv,PetscViewer viewer)
380: {
381:   PetscErrorCode    ierr;
382:   BV_SVEC           *ctx = (BV_SVEC*)bv->data;
383:   PetscViewerFormat format;
384:   PetscBool         isascii;
385:   const char        *bvname,*name;

388:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
389:   if (isascii) {
390:     PetscViewerGetFormat(viewer,&format);
391:     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
392:     VecView(ctx->v,viewer);
393:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
394:       PetscObjectGetName((PetscObject)bv,&bvname);
395:       PetscObjectGetName((PetscObject)ctx->v,&name);
396:       PetscViewerASCIIPrintf(viewer,"%s=reshape(%s,%D,%D);clear %s\n",bvname,name,bv->N,bv->nc+bv->m,name);
397:       if (bv->nc) {
398:         PetscViewerASCIIPrintf(viewer,"%s=%s(:,%D:end);\n",bvname,bvname,bv->nc+1);
399:       }
400:     }
401:   } else {
402:     VecView(ctx->v,viewer);
403:   }
404:   return(0);
405: }

407: PetscErrorCode BVDestroy_Svec(BV bv)
408: {
410:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;

413:   VecDestroy(&ctx->v);
414:   VecDestroy(&bv->cv[0]);
415:   VecDestroy(&bv->cv[1]);
416:   PetscFree(bv->data);
417:   bv->cuda = PETSC_FALSE;
418:   return(0);
419: }

421: SLEPC_EXTERN PetscErrorCode BVCreate_Svec(BV bv)
422: {
423:   PetscErrorCode    ierr;
424:   BV_SVEC           *ctx;
425:   PetscInt          nloc,N,bs,tglobal=0,tlocal,lsplit;
426:   PetscBool         seq;
427:   PetscScalar       *aa,*vv;
428:   const PetscScalar *array,*ptr;
429:   char              str[50];
430:   BV                parent;
431:   Vec               vpar;
432: #if defined(PETSC_HAVE_CUDA)
433:   PetscScalar       *gpuarray,*gptr;
434: #endif

437:   PetscNewLog(bv,&ctx);
438:   bv->data = (void*)ctx;

440:   PetscObjectTypeCompareAny((PetscObject)bv->t,&bv->cuda,VECSEQCUDA,VECMPICUDA,"");
441:   PetscObjectTypeCompareAny((PetscObject)bv->t,&ctx->mpi,VECMPI,VECMPICUDA,"");

443:   PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
444:   if (!seq && !ctx->mpi && !bv->cuda) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"BVSVEC does not support the type of the provided template vector");

446:   VecGetLocalSize(bv->t,&nloc);
447:   VecGetSize(bv->t,&N);
448:   VecGetBlockSize(bv->t,&bs);
449:   tlocal  = bv->m*nloc;
450:   PetscIntMultError(bv->m,N,&tglobal);

452:   if (bv->issplit) {
453:     /* split BV: create Vec sharing the memory of the parent BV */
454:     parent = bv->splitparent;
455:     lsplit = parent->lsplit;
456:     vpar = ((BV_SVEC*)parent->data)->v;
457:     if (bv->cuda) {
458: #if defined(PETSC_HAVE_CUDA)
459:       VecCUDAGetArray(vpar,&gpuarray);
460:       gptr = (bv->issplit==1)? gpuarray: gpuarray+lsplit*nloc;
461:       VecCUDARestoreArray(vpar,&gpuarray);
462:       if (ctx->mpi) {
463:         VecCreateMPICUDAWithArray(PetscObjectComm((PetscObject)bv->t),bs,tlocal,bv->m*N,NULL,&ctx->v);
464:       } else {
465:         VecCreateSeqCUDAWithArray(PetscObjectComm((PetscObject)bv->t),bs,tlocal,NULL,&ctx->v);
466:       }
467:       VecCUDAPlaceArray(ctx->v,gptr);
468: #endif
469:     } else {
470:       VecGetArrayRead(vpar,&array);
471:       ptr = (bv->issplit==1)? array: array+lsplit*nloc;
472:       VecRestoreArrayRead(vpar,&array);
473:       if (ctx->mpi) {
474:         VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,tlocal,bv->m*N,NULL,&ctx->v);
475:       } else {
476:         VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,tlocal,NULL,&ctx->v);
477:       }
478:       VecPlaceArray(ctx->v,ptr);
479:     }
480:   } else {
481:     /* regular BV: create Vec to store the BV entries */
482:     VecCreate(PetscObjectComm((PetscObject)bv->t),&ctx->v);
483:     VecSetType(ctx->v,((PetscObject)bv->t)->type_name);
484:     VecSetSizes(ctx->v,tlocal,tglobal);
485:     VecSetBlockSize(ctx->v,bs);
486:   }
487:   PetscLogObjectParent((PetscObject)bv,(PetscObject)ctx->v);
488:   if (((PetscObject)bv)->name) {
489:     PetscSNPrintf(str,50,"%s_0",((PetscObject)bv)->name);
490:     PetscObjectSetName((PetscObject)ctx->v,str);
491:   }

493:   if (bv->Acreate) {
494:     MatDenseGetArray(bv->Acreate,&aa);
495:     VecGetArray(ctx->v,&vv);
496:     PetscArraycpy(vv,aa,tlocal);
497:     VecRestoreArray(ctx->v,&vv);
498:     MatDenseRestoreArray(bv->Acreate,&aa);
499:     MatDestroy(&bv->Acreate);
500:   }

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

505:   if (bv->cuda) {
506: #if defined(PETSC_HAVE_CUDA)
507:     bv->ops->mult             = BVMult_Svec_CUDA;
508:     bv->ops->multvec          = BVMultVec_Svec_CUDA;
509:     bv->ops->multinplace      = BVMultInPlace_Svec_CUDA;
510:     bv->ops->multinplacetrans = BVMultInPlaceTranspose_Svec_CUDA;
511:     bv->ops->dot              = BVDot_Svec_CUDA;
512:     bv->ops->dotvec           = BVDotVec_Svec_CUDA;
513:     bv->ops->dotvec_local     = BVDotVec_Local_Svec_CUDA;
514:     bv->ops->scale            = BVScale_Svec_CUDA;
515:     bv->ops->matmult          = BVMatMult_Svec_CUDA;
516:     bv->ops->copy             = BVCopy_Svec_CUDA;
517:     bv->ops->copycolumn       = BVCopyColumn_Svec_CUDA;
518:     bv->ops->resize           = BVResize_Svec_CUDA;
519:     bv->ops->getcolumn        = BVGetColumn_Svec_CUDA;
520:     bv->ops->restorecolumn    = BVRestoreColumn_Svec_CUDA;
521:     bv->ops->restoresplit     = BVRestoreSplit_Svec_CUDA;
522: #endif
523:   } else {
524:     bv->ops->mult             = BVMult_Svec;
525:     bv->ops->multvec          = BVMultVec_Svec;
526:     bv->ops->multinplace      = BVMultInPlace_Svec;
527:     bv->ops->multinplacetrans = BVMultInPlaceTranspose_Svec;
528:     bv->ops->dot              = BVDot_Svec;
529:     bv->ops->dotvec           = BVDotVec_Svec;
530:     bv->ops->dotvec_local     = BVDotVec_Local_Svec;
531:     bv->ops->scale            = BVScale_Svec;
532:     bv->ops->matmult          = BVMatMult_Svec;
533:     bv->ops->copy             = BVCopy_Svec;
534:     bv->ops->copycolumn       = BVCopyColumn_Svec;
535:     bv->ops->resize           = BVResize_Svec;
536:     bv->ops->getcolumn        = BVGetColumn_Svec;
537:     bv->ops->restorecolumn    = BVRestoreColumn_Svec;
538:   }
539:   bv->ops->norm             = BVNorm_Svec;
540:   bv->ops->norm_local       = BVNorm_Local_Svec;
541:   bv->ops->getarray         = BVGetArray_Svec;
542:   bv->ops->restorearray     = BVRestoreArray_Svec;
543:   bv->ops->getarrayread     = BVGetArrayRead_Svec;
544:   bv->ops->restorearrayread = BVRestoreArrayRead_Svec;
545:   bv->ops->destroy          = BVDestroy_Svec;
546:   if (!ctx->mpi) bv->ops->view = BVView_Svec;
547:   return(0);
548: }