Actual source code: vecs.c

slepc-3.14.0 2020-09-30
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 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: {
 26:   BV_VECS        *y = (BV_VECS*)Y->data,*x = (BV_VECS*)X->data;
 27:   PetscScalar    *q,*s=NULL;
 28:   PetscInt       i,j,ldq;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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