Actual source code: vecs.c

slepc-3.15.0 2021-03-31
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 as an array of independent Vecs
 12: */

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

 16: typedef struct {
 17:   Vec      *V;
 18:   PetscInt vmip;   /* Version of BVMultInPlace:
 19:        0: memory-efficient version, uses VecGetArray (default in CPU)
 20:        1: version that allocates (e-s) work vectors in every call (default in GPU) */
 21: } BV_VECS;

 23: PetscErrorCode BVMult_Vecs(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
 24: {
 25:   PetscErrorCode    ierr;
 26:   BV_VECS           *y = (BV_VECS*)Y->data,*x = (BV_VECS*)X->data;
 27:   PetscScalar       *s=NULL;
 28:   const PetscScalar *q;
 29:   PetscInt          i,j,ldq;

 32:   if (Q) {
 33:     MatGetSize(Q,&ldq,NULL);
 34:     if (alpha!=1.0) {
 35:       BVAllocateWork_Private(Y,X->k-X->l);
 36:       s = Y->work;
 37:     }
 38:     MatDenseGetArrayRead(Q,&q);
 39:     for (j=Y->l;j<Y->k;j++) {
 40:       VecScale(y->V[Y->nc+j],beta);
 41:       if (alpha!=1.0) {
 42:         for (i=X->l;i<X->k;i++) s[i-X->l] = alpha*q[i+j*ldq];
 43:       } else s = (PetscScalar*)(q+j*ldq+X->l);
 44:       VecMAXPY(y->V[Y->nc+j],X->k-X->l,s,x->V+X->nc+X->l);
 45:     }
 46:     MatDenseRestoreArrayRead(Q,&q);
 47:   } else {
 48:     for (j=0;j<Y->k-Y->l;j++) {
 49:       VecScale(y->V[Y->nc+Y->l+j],beta);
 50:       VecAXPY(y->V[Y->nc+Y->l+j],alpha,x->V[X->nc+X->l+j]);
 51:     }
 52:   }
 53:   return(0);
 54: }

 56: PetscErrorCode BVMultVec_Vecs(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
 57: {
 59:   BV_VECS        *x = (BV_VECS*)X->data;
 60:   PetscScalar    *s=NULL,*qq=q;
 61:   PetscInt       i;

 64:   if (alpha!=1.0) {
 65:     BVAllocateWork_Private(X,X->k-X->l);
 66:     s = X->work;
 67:   }
 68:   if (!q) { VecGetArray(X->buffer,&qq); }
 69:   VecScale(y,beta);
 70:   if (alpha!=1.0) {
 71:     for (i=0;i<X->k-X->l;i++) s[i] = alpha*qq[i];
 72:   } else s = qq;
 73:   VecMAXPY(y,X->k-X->l,s,x->V+X->nc+X->l);
 74:   if (!q) { VecRestoreArray(X->buffer,&qq); }
 75:   return(0);
 76: }

 78: /*
 79:    BVMultInPlace_Vecs_ME - V(:,s:e-1) = V*Q(:,s:e-1) for regular vectors.

 81:    Memory-efficient version, uses VecGetArray (default in CPU)

 83:    Writing V = [ V1 V2 V3 ] and Q(:,s:e-1) = [ Q1 Q2 Q3 ]', where V2
 84:    corresponds to the columns s:e-1, the computation is done as
 85:                   V2 := V2*Q2 + V1*Q1 + V3*Q3
 86: */
 87: PetscErrorCode BVMultInPlace_Vecs_ME(BV V,Mat Q,PetscInt s,PetscInt e)
 88: {
 89:   PetscErrorCode    ierr;
 90:   BV_VECS           *ctx = (BV_VECS*)V->data;
 91:   const PetscScalar *q;
 92:   PetscInt          i,ldq;

 95:   MatGetSize(Q,&ldq,NULL);
 96:   MatDenseGetArrayRead(Q,&q);
 97:   /* V2 := V2*Q2 */
 98:   BVMultInPlace_Vecs_Private(V,V->n,e-s,ldq,ctx->V+V->nc+s,q+s*ldq+s,PETSC_FALSE);
 99:   /* V2 += V1*Q1 + V3*Q3 */
100:   for (i=s;i<e;i++) {
101:     if (s>V->l) {
102:       VecMAXPY(ctx->V[V->nc+i],s-V->l,q+i*ldq+V->l,ctx->V+V->nc+V->l);
103:     }
104:     if (V->k>e) {
105:       VecMAXPY(ctx->V[V->nc+i],V->k-e,q+i*ldq+e,ctx->V+V->nc+e);
106:     }
107:   }
108:   MatDenseRestoreArrayRead(Q,&q);
109:   return(0);
110: }

112: /*
113:    BVMultInPlace_Vecs_Alloc - V(:,s:e-1) = V*Q(:,s:e-1) for regular vectors.

115:    Version that allocates (e-s) work vectors in every call (default in GPU)
116: */
117: PetscErrorCode BVMultInPlace_Vecs_Alloc(BV V,Mat Q,PetscInt s,PetscInt e)
118: {
119:   PetscErrorCode    ierr;
120:   BV_VECS           *ctx = (BV_VECS*)V->data;
121:   const PetscScalar *q;
122:   PetscInt          i,ldq;
123:   Vec               *W;

126:   MatGetSize(Q,&ldq,NULL);
127:   MatDenseGetArrayRead(Q,&q);
128:   VecDuplicateVecs(V->t,e-s,&W);
129:   for (i=s;i<e;i++) {
130:     VecMAXPY(W[i-s],V->k-V->l,q+i*ldq+V->l,ctx->V+V->nc+V->l);
131:   }
132:   for (i=s;i<e;i++) {
133:     VecCopy(W[i-s],ctx->V[V->nc+i]);
134:   }
135:   VecDestroyVecs(e-s,&W);
136:   MatDenseRestoreArrayRead(Q,&q);
137:   return(0);
138: }

140: /*
141:    BVMultInPlaceTranspose_Vecs - V(:,s:e-1) = V*Q'(:,s:e-1) for regular vectors.
142: */
143: PetscErrorCode BVMultInPlaceTranspose_Vecs(BV V,Mat Q,PetscInt s,PetscInt e)
144: {
145:   PetscErrorCode    ierr;
146:   BV_VECS           *ctx = (BV_VECS*)V->data;
147:   const PetscScalar *q;
148:   PetscInt          i,j,ldq,n;

151:   MatGetSize(Q,&ldq,&n);
152:   MatDenseGetArrayRead(Q,&q);
153:   /* V2 := V2*Q2' */
154:   BVMultInPlace_Vecs_Private(V,V->n,e-s,ldq,ctx->V+V->nc+s,q+s*ldq+s,PETSC_TRUE);
155:   /* V2 += V1*Q1' + V3*Q3' */
156:   for (i=s;i<e;i++) {
157:     for (j=V->l;j<s;j++) {
158:       VecAXPY(ctx->V[V->nc+i],q[i+j*ldq],ctx->V[V->nc+j]);
159:     }
160:     for (j=e;j<n;j++) {
161:       VecAXPY(ctx->V[V->nc+i],q[i+j*ldq],ctx->V[V->nc+j]);
162:     }
163:   }
164:   MatDenseRestoreArrayRead(Q,&q);
165:   return(0);
166: }

168: PetscErrorCode BVDot_Vecs(BV X,BV Y,Mat M)
169: {
171:   BV_VECS        *x = (BV_VECS*)X->data,*y = (BV_VECS*)Y->data;
172:   PetscScalar    *m;
173:   PetscInt       j,ldm;

176:   MatGetSize(M,&ldm,NULL);
177:   MatDenseGetArray(M,&m);
178:   for (j=X->l;j<X->k;j++) {
179:     VecMDot(x->V[X->nc+j],Y->k-Y->l,y->V+Y->nc+Y->l,m+j*ldm+Y->l);
180:   }
181:   MatDenseRestoreArray(M,&m);
182:   return(0);
183: }

185: PetscErrorCode BVDotVec_Vecs(BV X,Vec y,PetscScalar *q)
186: {
188:   BV_VECS        *x = (BV_VECS*)X->data;
189:   Vec            z = y;
190:   PetscScalar    *qq=q;

193:   if (X->matrix) {
194:     BV_IPMatMult(X,y);
195:     z = X->Bx;
196:   }
197:   if (!q) { VecGetArray(X->buffer,&qq); }
198:   VecMDot(z,X->k-X->l,x->V+X->nc+X->l,qq);
199:   if (!q) { VecRestoreArray(X->buffer,&qq); }
200:   return(0);
201: }

203: PetscErrorCode BVDotVec_Begin_Vecs(BV X,Vec y,PetscScalar *m)
204: {
206:   BV_VECS        *x = (BV_VECS*)X->data;
207:   Vec            z = y;

210:   if (X->matrix) {
211:     BV_IPMatMult(X,y);
212:     z = X->Bx;
213:   }
214:   VecMDotBegin(z,X->k-X->l,x->V+X->nc+X->l,m);
215:   return(0);
216: }

218: PetscErrorCode BVDotVec_End_Vecs(BV X,Vec y,PetscScalar *m)
219: {
221:   BV_VECS        *x = (BV_VECS*)X->data;

224:   VecMDotEnd(y,X->k-X->l,x->V+X->nc+X->l,m);
225:   return(0);
226: }

228: PetscErrorCode BVScale_Vecs(BV bv,PetscInt j,PetscScalar alpha)
229: {
231:   PetscInt       i;
232:   BV_VECS        *ctx = (BV_VECS*)bv->data;

235:   if (j<0) {
236:     for (i=bv->l;i<bv->k;i++) {
237:       VecScale(ctx->V[bv->nc+i],alpha);
238:     }
239:   } else {
240:     VecScale(ctx->V[bv->nc+j],alpha);
241:   }
242:   return(0);
243: }

245: PetscErrorCode BVNorm_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val)
246: {
248:   PetscInt       i;
249:   PetscReal      nrm;
250:   BV_VECS        *ctx = (BV_VECS*)bv->data;

253:   if (j<0) {
254:     if (type==NORM_FROBENIUS) {
255:       *val = 0.0;
256:       for (i=bv->l;i<bv->k;i++) {
257:         VecNorm(ctx->V[bv->nc+i],NORM_2,&nrm);
258:         *val += nrm*nrm;
259:       }
260:       *val = PetscSqrtReal(*val);
261:     } else SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not implemented in BVVECS");
262:   } else {
263:     VecNorm(ctx->V[bv->nc+j],type,val);
264:   }
265:   return(0);
266: }

268: PetscErrorCode BVNorm_Begin_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val)
269: {
271:   BV_VECS        *ctx = (BV_VECS*)bv->data;

274:   if (j<0) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not implemented in BVVECS");
275:   else {
276:     VecNormBegin(ctx->V[bv->nc+j],type,val);
277:   }
278:   return(0);
279: }

281: PetscErrorCode BVNorm_End_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val)
282: {
284:   BV_VECS        *ctx = (BV_VECS*)bv->data;

287:   if (j<0) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not implemented in BVVECS");
288:   else {
289:     VecNormEnd(ctx->V[bv->nc+j],type,val);
290:   }
291:   return(0);
292: }

294: PetscErrorCode BVMatMult_Vecs(BV V,Mat A,BV W)
295: {
297:   BV_VECS        *v = (BV_VECS*)V->data,*w = (BV_VECS*)W->data;
298:   PetscInt       j;
299:   Mat            Vmat,Wmat;

302:   if (V->vmm) {
303:     BVGetMat(V,&Vmat);
304:     BVGetMat(W,&Wmat);
305:     MatProductCreateWithMat(A,Vmat,NULL,Wmat);
306:     MatProductSetType(Wmat,MATPRODUCT_AB);
307:     MatProductSetFromOptions(Wmat);
308:     MatProductSymbolic(Wmat);
309:     MatProductNumeric(Wmat);
310:     MatProductClear(Wmat);
311:     BVRestoreMat(V,&Vmat);
312:     BVRestoreMat(W,&Wmat);
313:   } else {
314:     for (j=0;j<V->k-V->l;j++) {
315:       MatMult(A,v->V[V->nc+V->l+j],w->V[W->nc+W->l+j]);
316:     }
317:   }
318:   return(0);
319: }

321: PetscErrorCode BVCopy_Vecs(BV V,BV W)
322: {
324:   BV_VECS        *v = (BV_VECS*)V->data,*w = (BV_VECS*)W->data;
325:   PetscInt       j;

328:   for (j=0;j<V->k-V->l;j++) {
329:     VecCopy(v->V[V->nc+V->l+j],w->V[W->nc+W->l+j]);
330:   }
331:   return(0);
332: }

334: PetscErrorCode BVCopyColumn_Vecs(BV V,PetscInt j,PetscInt i)
335: {
337:   BV_VECS        *v = (BV_VECS*)V->data;

340:   VecCopy(v->V[V->nc+j],v->V[V->nc+i]);
341:   return(0);
342: }

344: PetscErrorCode BVResize_Vecs(BV bv,PetscInt m,PetscBool copy)
345: {
347:   BV_VECS        *ctx = (BV_VECS*)bv->data;
348:   Vec            *newV;
349:   PetscInt       j;
350:   char           str[50];

353:   VecDuplicateVecs(bv->t,m,&newV);
354:   PetscLogObjectParents(bv,m,newV);
355:   if (((PetscObject)bv)->name) {
356:     for (j=0;j<m;j++) {
357:       PetscSNPrintf(str,sizeof(str),"%s_%D",((PetscObject)bv)->name,j);
358:       PetscObjectSetName((PetscObject)newV[j],str);
359:     }
360:   }
361:   if (copy) {
362:     for (j=0;j<PetscMin(m,bv->m);j++) {
363:       VecCopy(ctx->V[j],newV[j]);
364:     }
365:   }
366:   VecDestroyVecs(bv->m,&ctx->V);
367:   ctx->V = newV;
368:   return(0);
369: }

371: PetscErrorCode BVGetColumn_Vecs(BV bv,PetscInt j,Vec *v)
372: {
373:   BV_VECS  *ctx = (BV_VECS*)bv->data;
374:   PetscInt l;

377:   l = BVAvailableVec;
378:   bv->cv[l] = ctx->V[bv->nc+j];
379:   return(0);
380: }

382: PetscErrorCode BVGetArray_Vecs(BV bv,PetscScalar **a)
383: {
384:   PetscErrorCode    ierr;
385:   BV_VECS           *ctx = (BV_VECS*)bv->data;
386:   PetscInt          j;
387:   const PetscScalar *p;

390:   PetscMalloc1((bv->nc+bv->m)*bv->n,a);
391:   for (j=0;j<bv->nc+bv->m;j++) {
392:     VecGetArrayRead(ctx->V[j],&p);
393:     PetscArraycpy(*a+j*bv->n,p,bv->n);
394:     VecRestoreArrayRead(ctx->V[j],&p);
395:   }
396:   return(0);
397: }

399: PetscErrorCode BVRestoreArray_Vecs(BV bv,PetscScalar **a)
400: {
402:   BV_VECS        *ctx = (BV_VECS*)bv->data;
403:   PetscInt       j;
404:   PetscScalar    *p;

407:   for (j=0;j<bv->nc+bv->m;j++) {
408:     VecGetArray(ctx->V[j],&p);
409:     PetscArraycpy(p,*a+j*bv->n,bv->n);
410:     VecRestoreArray(ctx->V[j],&p);
411:   }
412:   PetscFree(*a);
413:   return(0);
414: }

416: PetscErrorCode BVGetArrayRead_Vecs(BV bv,const PetscScalar **a)
417: {
418:   PetscErrorCode    ierr;
419:   BV_VECS           *ctx = (BV_VECS*)bv->data;
420:   PetscInt          j;
421:   const PetscScalar *p;

424:   PetscMalloc1((bv->nc+bv->m)*bv->n,(PetscScalar**)a);
425:   for (j=0;j<bv->nc+bv->m;j++) {
426:     VecGetArrayRead(ctx->V[j],&p);
427:     PetscArraycpy((PetscScalar*)*a+j*bv->n,p,bv->n);
428:     VecRestoreArrayRead(ctx->V[j],&p);
429:   }
430:   return(0);
431: }

433: PetscErrorCode BVRestoreArrayRead_Vecs(BV bv,const PetscScalar **a)
434: {

438:   PetscFree(*a);
439:   return(0);
440: }

442: /*
443:    Sets the value of vmip flag and resets ops->multinplace accordingly
444:  */
445: PETSC_STATIC_INLINE PetscErrorCode BVVecsSetVmip(BV bv,PetscInt vmip)
446: {
447:   typedef PetscErrorCode (*fmultinplace)(BV,Mat,PetscInt,PetscInt);
448:   fmultinplace multinplace[2] = {BVMultInPlace_Vecs_ME, BVMultInPlace_Vecs_Alloc};
449:   BV_VECS      *ctx = (BV_VECS*)bv->data;

452:   ctx->vmip            = vmip;
453:   bv->ops->multinplace = multinplace[vmip];
454:   return(0);
455: }

457: PetscErrorCode BVSetFromOptions_Vecs(PetscOptionItems *PetscOptionsObject,BV bv)
458: {
460:   BV_VECS        *ctx = (BV_VECS*)bv->data;

463:   PetscOptionsHead(PetscOptionsObject,"BV Vecs Options");

465:     PetscOptionsRangeInt("-bv_vecs_vmip","Version of BVMultInPlace operation","",ctx->vmip,&ctx->vmip,NULL,0,1);
466:     BVVecsSetVmip(bv,ctx->vmip);

468:   PetscOptionsTail();
469:   return(0);
470: }

472: PetscErrorCode BVView_Vecs(BV bv,PetscViewer viewer)
473: {
474:   PetscErrorCode    ierr;
475:   BV_VECS           *ctx = (BV_VECS*)bv->data;
476:   PetscInt          j;
477:   PetscViewerFormat format;
478:   PetscBool         isascii,ismatlab=PETSC_FALSE;
479:   const char        *bvname,*name;

482:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
483:   if (isascii) {
484:     PetscViewerGetFormat(viewer,&format);
485:     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
486:     if (format == PETSC_VIEWER_ASCII_MATLAB) ismatlab = PETSC_TRUE;
487:   }
488:   if (ismatlab) {
489:     PetscObjectGetName((PetscObject)bv,&bvname);
490:     PetscViewerASCIIPrintf(viewer,"%s=[];\n",bvname);
491:   }
492:   for (j=bv->nc;j<bv->nc+bv->m;j++) {
493:     VecView(ctx->V[j],viewer);
494:     if (ismatlab) {
495:       PetscObjectGetName((PetscObject)ctx->V[j],&name);
496:       PetscViewerASCIIPrintf(viewer,"%s=[%s,%s];clear %s\n",bvname,bvname,name,name);
497:     }
498:   }
499:   return(0);
500: }

502: PetscErrorCode BVDestroy_Vecs(BV bv)
503: {
505:   BV_VECS        *ctx = (BV_VECS*)bv->data;

508:   if (!bv->issplit) { VecDestroyVecs(bv->nc+bv->m,&ctx->V); }
509:   PetscFree(bv->data);
510:   return(0);
511: }

513: PetscErrorCode BVDuplicate_Vecs(BV V,BV W)
514: {
516:   BV_VECS        *ctx = (BV_VECS*)V->data;

519:   BVVecsSetVmip(W,ctx->vmip);
520:   return(0);
521: }

523: SLEPC_EXTERN PetscErrorCode BVCreate_Vecs(BV bv)
524: {
526:   BV_VECS        *ctx;
527:   PetscInt       j,lsplit;
528:   PetscBool      isgpu;
529:   char           str[50];
530:   BV             parent;
531:   Vec            *Vpar;

534:   PetscNewLog(bv,&ctx);
535:   bv->data = (void*)ctx;

537:   if (bv->issplit) {
538:     /* split BV: share the Vecs of the parent BV */
539:     parent = bv->splitparent;
540:     lsplit = parent->lsplit;
541:     Vpar   = ((BV_VECS*)parent->data)->V;
542:     ctx->V = (bv->issplit==1)? Vpar: Vpar+lsplit;
543:   } else {
544:     /* regular BV: create array of Vecs to store the BV columns */
545:     VecDuplicateVecs(bv->t,bv->m,&ctx->V);
546:     PetscLogObjectParents(bv,bv->m,ctx->V);
547:     if (((PetscObject)bv)->name) {
548:       for (j=0;j<bv->m;j++) {
549:         PetscSNPrintf(str,sizeof(str),"%s_%D",((PetscObject)bv)->name,j);
550:         PetscObjectSetName((PetscObject)ctx->V[j],str);
551:       }
552:     }
553:   }

555:   if (bv->Acreate) {
556:     for (j=0;j<bv->m;j++) {
557:       MatGetColumnVector(bv->Acreate,ctx->V[j],j);
558:     }
559:     MatDestroy(&bv->Acreate);
560:   }

562:   /* Default version of BVMultInPlace */
563:   PetscObjectTypeCompareAny((PetscObject)bv->t,&isgpu,VECSEQCUDA,VECMPICUDA,"");
564:   ctx->vmip = isgpu? 1: 0;

566:   /* Default BVMatMult method */
567:   bv->vmm = BV_MATMULT_VECS;

569:   /* Deferred call to setfromoptions */
570:   if (bv->defersfo) {
571:     PetscObjectOptionsBegin((PetscObject)bv);
572:     BVSetFromOptions_Vecs(PetscOptionsObject,bv);
573:     PetscOptionsEnd();
574:   }
575:   BVVecsSetVmip(bv,ctx->vmip);

577:   bv->ops->mult             = BVMult_Vecs;
578:   bv->ops->multvec          = BVMultVec_Vecs;
579:   bv->ops->multinplacetrans = BVMultInPlaceTranspose_Vecs;
580:   bv->ops->dot              = BVDot_Vecs;
581:   bv->ops->dotvec           = BVDotVec_Vecs;
582:   bv->ops->dotvec_begin     = BVDotVec_Begin_Vecs;
583:   bv->ops->dotvec_end       = BVDotVec_End_Vecs;
584:   bv->ops->scale            = BVScale_Vecs;
585:   bv->ops->norm             = BVNorm_Vecs;
586:   bv->ops->norm_begin       = BVNorm_Begin_Vecs;
587:   bv->ops->norm_end         = BVNorm_End_Vecs;
588:   bv->ops->matmult          = BVMatMult_Vecs;
589:   bv->ops->copy             = BVCopy_Vecs;
590:   bv->ops->copycolumn       = BVCopyColumn_Vecs;
591:   bv->ops->resize           = BVResize_Vecs;
592:   bv->ops->getcolumn        = BVGetColumn_Vecs;
593:   bv->ops->getarray         = BVGetArray_Vecs;
594:   bv->ops->restorearray     = BVRestoreArray_Vecs;
595:   bv->ops->getarrayread     = BVGetArrayRead_Vecs;
596:   bv->ops->restorearrayread = BVRestoreArrayRead_Vecs;
597:   bv->ops->destroy          = BVDestroy_Vecs;
598:   bv->ops->duplicate        = BVDuplicate_Vecs;
599:   bv->ops->setfromoptions   = BVSetFromOptions_Vecs;
600:   bv->ops->view             = BVView_Vecs;
601:   return(0);
602: }