Actual source code: vecs.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  1: /*
  2:    BV implemented as an array of independent Vecs

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 10:    SLEPc is free software: you can redistribute it and/or modify it under  the
 11:    terms of version 3 of the GNU Lesser General Public License as published by
 12:    the Free Software Foundation.

 14:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 15:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 16:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 17:    more details.

 19:    You  should have received a copy of the GNU Lesser General  Public  License
 20:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

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

 26: typedef struct {
 27:   Vec      *V;
 28:   PetscInt vmip;   /* Version of BVMultInPlace:
 29:        0: memory-efficient version, uses VecGetArray (default in CPU)
 30:        1: version that allocates (e-s) work vectors in every call (default in GPU) */
 31: } BV_VECS;

 35: PetscErrorCode BVMult_Vecs(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
 36: {
 38:   BV_VECS        *y = (BV_VECS*)Y->data,*x = (BV_VECS*)X->data;
 39:   PetscScalar    *q,*s=NULL;
 40:   PetscInt       i,j,ldq;

 43:   MatGetSize(Q,&ldq,NULL);
 44:   if (alpha!=1.0) {
 45:     BVAllocateWork_Private(Y,X->k-X->l);
 46:     s = Y->work;
 47:   }
 48:   MatDenseGetArray(Q,&q);
 49:   for (j=Y->l;j<Y->k;j++) {
 50:     VecScale(y->V[Y->nc+j],beta);
 51:     if (alpha!=1.0) {
 52:       for (i=X->l;i<X->k;i++) s[i-X->l] = alpha*q[i+j*ldq];
 53:     } else s = q+j*ldq+X->l;
 54:     VecMAXPY(y->V[Y->nc+j],X->k-X->l,s,x->V+X->nc+X->l);
 55:   }
 56:   MatDenseRestoreArray(Q,&q);
 57:   return(0);
 58: }

 62: PetscErrorCode BVMultVec_Vecs(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
 63: {
 65:   BV_VECS        *x = (BV_VECS*)X->data;
 66:   PetscScalar    *s=NULL;
 67:   PetscInt       i;

 70:   if (alpha!=1.0) {
 71:     BVAllocateWork_Private(X,X->k-X->l);
 72:     s = X->work;
 73:   }
 74:   VecScale(y,beta);
 75:   if (alpha!=1.0) {
 76:     for (i=0;i<X->k-X->l;i++) s[i] = alpha*q[i];
 77:   } else s = q;
 78:   VecMAXPY(y,X->k-X->l,s,x->V+X->nc+X->l);
 79:   return(0);
 80: }

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

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

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

101:   MatGetSize(Q,&ldq,NULL);
102:   MatDenseGetArray(Q,&q);
103:   /* V2 := V2*Q2 */
104:   BVMultInPlace_Vecs_Private(V,V->n,e-s,V->k,ctx->V+V->nc+s,q+s*ldq+s,PETSC_FALSE);
105:   /* V2 += V1*Q1 + V3*Q3 */
106:   for (i=s;i<e;i++) {
107:     if (s>V->l) {
108:       VecMAXPY(ctx->V[V->nc+i],s-V->l,q+i*ldq+V->l,ctx->V+V->nc+V->l);
109:     }
110:     if (V->k>e) {
111:       VecMAXPY(ctx->V[V->nc+i],V->k-e,q+i*ldq+e,ctx->V+V->nc+e);
112:     }
113:   }
114:   MatDenseRestoreArray(Q,&q);
115:   return(0);
116: }

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

123:    Version that allocates (e-s) work vectors in every call (default in GPU)
124: */
125: PetscErrorCode BVMultInPlace_Vecs_Alloc(BV V,Mat Q,PetscInt s,PetscInt e)
126: {
128:   BV_VECS        *ctx = (BV_VECS*)V->data;
129:   PetscScalar    *q;
130:   PetscInt       i,ldq;
131:   Vec            *W;

134:   MatGetSize(Q,&ldq,NULL);
135:   MatDenseGetArray(Q,&q);
136:   VecDuplicateVecs(V->t,e-s,&W);
137:   for (i=s;i<e;i++) {
138:     VecMAXPY(W[i-s],V->k-V->l,q+i*ldq+V->l,ctx->V+V->nc+V->l);
139:   }
140:   for (i=s;i<e;i++) {
141:     VecCopy(W[i-s],ctx->V[V->nc+i]);
142:   }
143:   VecDestroyVecs(e-s,&W);
144:   MatDenseRestoreArray(Q,&q);
145:   return(0);
146: }

150: /*
151:    BVMultInPlaceTranspose_Vecs - V(:,s:e-1) = V*Q'(:,s:e-1) for regular vectors.
152: */
153: PetscErrorCode BVMultInPlaceTranspose_Vecs(BV V,Mat Q,PetscInt s,PetscInt e)
154: {
156:   BV_VECS        *ctx = (BV_VECS*)V->data;
157:   PetscScalar    *q;
158:   PetscInt       i,j,ldq,n;

161:   MatGetSize(Q,&ldq,&n);
162:   MatDenseGetArray(Q,&q);
163:   /* V2 := V2*Q2' */
164:   BVMultInPlace_Vecs_Private(V,V->n,e-s,ldq,ctx->V+V->nc+s,q+s*ldq+s,PETSC_TRUE);
165:   /* V2 += V1*Q1' + V3*Q3' */
166:   for (i=s;i<e;i++) {
167:     for (j=V->l;j<s;j++) {
168:       VecAXPY(ctx->V[V->nc+i],q[i+j*ldq],ctx->V[V->nc+j]);
169:     }
170:     for (j=e;j<n;j++) {
171:       VecAXPY(ctx->V[V->nc+i],q[i+j*ldq],ctx->V[V->nc+j]);
172:     }
173:   }
174:   MatDenseRestoreArray(Q,&q);
175:   return(0);
176: }

180: PetscErrorCode BVAXPY_Vecs(BV Y,PetscScalar alpha,BV X)
181: {
183:   BV_VECS        *y = (BV_VECS*)Y->data,*x = (BV_VECS*)X->data;
184:   PetscInt       j;

187:   for (j=0;j<Y->k-Y->l;j++) {
188:     VecAXPY(y->V[Y->nc+Y->l+j],alpha,x->V[X->nc+X->l+j]);
189:   }
190:   return(0);
191: }

195: PetscErrorCode BVDot_Vecs(BV X,BV Y,Mat M)
196: {
198:   BV_VECS        *x = (BV_VECS*)X->data,*y = (BV_VECS*)Y->data;
199:   PetscScalar    *m;
200:   PetscInt       j,ldm;

203:   MatGetSize(M,&ldm,NULL);
204:   MatDenseGetArray(M,&m);
205:   for (j=X->l;j<X->k;j++) {
206:     VecMDot(x->V[X->nc+j],Y->k-Y->l,y->V+Y->nc+Y->l,m+j*ldm+Y->l);
207:   }
208:   MatDenseRestoreArray(M,&m);
209:   return(0);
210: }

214: PetscErrorCode BVDotVec_Vecs(BV X,Vec y,PetscScalar *m)
215: {
217:   BV_VECS        *x = (BV_VECS*)X->data;
218:   Vec            z = y;

221:   if (X->matrix) {
222:     BV_IPMatMult(X,y);
223:     z = X->Bx;
224:   }
225:   VecMDot(z,X->k-X->l,x->V+X->nc+X->l,m);
226:   return(0);
227: }

231: PetscErrorCode BVDotVec_Begin_Vecs(BV X,Vec y,PetscScalar *m)
232: {
234:   BV_VECS        *x = (BV_VECS*)X->data;
235:   Vec            z = y;

238:   if (X->matrix) {
239:     BV_IPMatMult(X,y);
240:     z = X->Bx;
241:   }
242:   VecMDotBegin(z,X->k-X->l,x->V+X->nc+X->l,m);
243:   return(0);
244: }

248: PetscErrorCode BVDotVec_End_Vecs(BV X,Vec y,PetscScalar *m)
249: {
251:   BV_VECS        *x = (BV_VECS*)X->data;

254:   VecMDotEnd(y,X->k-X->l,x->V+X->nc+X->l,m);
255:   return(0);
256: }

260: PetscErrorCode BVScale_Vecs(BV bv,PetscInt j,PetscScalar alpha)
261: {
263:   PetscInt       i;
264:   BV_VECS        *ctx = (BV_VECS*)bv->data;

267:   if (j<0) {
268:     for (i=bv->l;i<bv->k;i++) {
269:       VecScale(ctx->V[bv->nc+i],alpha);
270:     }
271:   } else {
272:     VecScale(ctx->V[bv->nc+j],alpha);
273:   }
274:   return(0);
275: }

279: PetscErrorCode BVNorm_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val)
280: {
282:   PetscInt       i;
283:   PetscReal      nrm;
284:   BV_VECS        *ctx = (BV_VECS*)bv->data;

287:   if (j<0) {
288:     switch (type) {
289:     case NORM_FROBENIUS:
290:       *val = 0.0;
291:       for (i=bv->l;i<bv->k;i++) {
292:         VecNorm(ctx->V[bv->nc+i],NORM_2,&nrm);
293:         *val += nrm*nrm;
294:       }
295:       *val = PetscSqrtReal(*val);
296:       break;
297:     default:
298:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Requested norm not implemented in BVVECS");
299:     }
300:   } else {
301:     VecNorm(ctx->V[bv->nc+j],type,val);
302:   }
303:   return(0);
304: }

308: PetscErrorCode BVNorm_Begin_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val)
309: {
311:   BV_VECS        *ctx = (BV_VECS*)bv->data;

314:   if (j<0) {
315:     switch (type) {
316:     default:
317:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Requested norm not implemented in BVVECS");
318:     }
319:   } else {
320:     VecNormBegin(ctx->V[bv->nc+j],type,val);
321:   }
322:   return(0);
323: }

327: PetscErrorCode BVNorm_End_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val)
328: {
330:   BV_VECS        *ctx = (BV_VECS*)bv->data;

333:   if (j<0) {
334:     switch (type) {
335:     default:
336:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Requested norm not implemented in BVVECS");
337:     }
338:   } else {
339:     VecNormEnd(ctx->V[bv->nc+j],type,val);
340:   }
341:   return(0);
342: }

346: PetscErrorCode BVMatMult_Vecs(BV V,Mat A,BV W)
347: {
349:   BV_VECS        *v = (BV_VECS*)V->data,*w = (BV_VECS*)W->data;
350:   PetscInt       j;

353:   if (V->vmm) { PetscInfo(V,"BVMatMult_Vecs: ignoring method\n"); }
354:   for (j=0;j<V->k-V->l;j++) {
355:     MatMult(A,v->V[V->nc+V->l+j],w->V[W->nc+W->l+j]);
356:   }
357:   return(0);
358: }

362: PetscErrorCode BVCopy_Vecs(BV V,BV W)
363: {
365:   BV_VECS        *v = (BV_VECS*)V->data,*w = (BV_VECS*)W->data;
366:   PetscInt       j;

369:   for (j=0;j<V->k-V->l;j++) {
370:     VecCopy(v->V[V->nc+V->l+j],w->V[W->nc+W->l+j]);
371:   }
372:   return(0);
373: }

377: PetscErrorCode BVResize_Vecs(BV bv,PetscInt m,PetscBool copy)
378: {
380:   BV_VECS        *ctx = (BV_VECS*)bv->data;
381:   Vec            *newV;
382:   PetscInt       j;
383:   char           str[50];

386:   VecDuplicateVecs(bv->t,m,&newV);
387:   PetscLogObjectParents(bv,m,newV);
388:   if (((PetscObject)bv)->name) {
389:     for (j=0;j<m;j++) {
390:       PetscSNPrintf(str,50,"%s_%D",((PetscObject)bv)->name,j);
391:       PetscObjectSetName((PetscObject)newV[j],str);
392:     }
393:   }
394:   if (copy) {
395:     for (j=0;j<PetscMin(m,bv->m);j++) {
396:       VecCopy(ctx->V[j],newV[j]);
397:     }
398:   }
399:   VecDestroyVecs(bv->m,&ctx->V);
400:   ctx->V = newV;
401:   return(0);
402: }

406: PetscErrorCode BVGetColumn_Vecs(BV bv,PetscInt j,Vec *v)
407: {
408:   BV_VECS  *ctx = (BV_VECS*)bv->data;
409:   PetscInt l;

412:   l = BVAvailableVec;
413:   bv->cv[l] = ctx->V[bv->nc+j];
414:   return(0);
415: }

419: PetscErrorCode BVGetArray_Vecs(BV bv,PetscScalar **a)
420: {
421:   PetscErrorCode    ierr;
422:   BV_VECS           *ctx = (BV_VECS*)bv->data;
423:   PetscInt          j;
424:   const PetscScalar *p;

427:   PetscMalloc((bv->nc+bv->m)*bv->n*sizeof(PetscScalar),a);
428:   for (j=0;j<bv->nc+bv->m;j++) {
429:     VecGetArrayRead(ctx->V[j],&p);
430:     PetscMemcpy(*a+j*bv->n,p,bv->n*sizeof(PetscScalar));
431:     VecRestoreArrayRead(ctx->V[j],&p);
432:   }
433:   return(0);
434: }

438: PetscErrorCode BVRestoreArray_Vecs(BV bv,PetscScalar **a)
439: {
441:   BV_VECS        *ctx = (BV_VECS*)bv->data;
442:   PetscInt       j;
443:   PetscScalar    *p;

446:   for (j=0;j<bv->nc+bv->m;j++) {
447:     VecGetArray(ctx->V[j],&p);
448:     PetscMemcpy(p,*a+j*bv->n,bv->n*sizeof(PetscScalar));
449:     VecRestoreArray(ctx->V[j],&p);
450:   }
451:   PetscFree(*a);
452:   return(0);
453: }

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

463:   PetscOptionsHead(PetscOptionsObject,"BV Vecs Options");
464:     PetscOptionsInt("-bv_vecs_vmip","Version of BVMultInPlace operation","",ctx->vmip,&ctx->vmip,NULL);
465:     if (ctx->vmip<0 || ctx->vmip>1) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Wrong version of BVMultInPlace");
466:   PetscOptionsTail();
467:   return(0);
468: }

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_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: }

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

509:   VecDestroyVecs(bv->nc+bv->m,&ctx->V);
510:   PetscFree(bv->data);
511:   return(0);
512: }

516: PETSC_EXTERN PetscErrorCode BVCreate_Vecs(BV bv)
517: {
519:   BV_VECS        *ctx;
520:   PetscInt       j;
521:   PetscBool      iscusp;
522:   char           str[50];
523:   typedef PetscErrorCode (*fmultinplace)(BV,Mat,PetscInt,PetscInt);
524:   fmultinplace   multinplace[2] = {BVMultInPlace_Vecs_ME, BVMultInPlace_Vecs_Alloc};

527:   PetscNewLog(bv,&ctx);
528:   bv->data = (void*)ctx;

530:   VecDuplicateVecs(bv->t,bv->m,&ctx->V);
531:   PetscLogObjectParents(bv,bv->m,ctx->V);
532:   if (((PetscObject)bv)->name) {
533:     for (j=0;j<bv->m;j++) {
534:       PetscSNPrintf(str,50,"%s_%D",((PetscObject)bv)->name,j);
535:       PetscObjectSetName((PetscObject)ctx->V[j],str);
536:     }
537:   }

539:   /* Default version of BVMultInPlace */
540:   PetscObjectTypeCompareAny((PetscObject)bv->t,&iscusp,VECSEQCUSP,VECMPICUSP,"");
541:   ctx->vmip = iscusp? 1: 0;

543:   /* Deferred call to setfromoptions */
544:   if (bv->defersfo) {
545:     PetscObjectOptionsBegin((PetscObject)bv);
546:     BVSetFromOptions_Vecs(PetscOptionsObject,bv);
547:     PetscOptionsEnd();
548:   }

550:   bv->ops->mult             = BVMult_Vecs;
551:   bv->ops->multvec          = BVMultVec_Vecs;
552:   bv->ops->multinplace      = multinplace[ctx->vmip];
553:   bv->ops->multinplacetrans = BVMultInPlaceTranspose_Vecs;
554:   bv->ops->axpy             = BVAXPY_Vecs;
555:   bv->ops->dot              = BVDot_Vecs;
556:   bv->ops->dotvec           = BVDotVec_Vecs;
557:   bv->ops->dotvec_begin     = BVDotVec_Begin_Vecs;
558:   bv->ops->dotvec_end       = BVDotVec_End_Vecs;
559:   bv->ops->scale            = BVScale_Vecs;
560:   bv->ops->norm             = BVNorm_Vecs;
561:   bv->ops->norm_begin       = BVNorm_Begin_Vecs;
562:   bv->ops->norm_end         = BVNorm_End_Vecs;
563:   bv->ops->matmult          = BVMatMult_Vecs;
564:   bv->ops->copy             = BVCopy_Vecs;
565:   bv->ops->resize           = BVResize_Vecs;
566:   bv->ops->getcolumn        = BVGetColumn_Vecs;
567:   bv->ops->getarray         = BVGetArray_Vecs;
568:   bv->ops->restorearray     = BVRestoreArray_Vecs;
569:   bv->ops->destroy          = BVDestroy_Vecs;
570:   bv->ops->setfromoptions   = BVSetFromOptions_Vecs;
571:   bv->ops->view             = BVView_Vecs;
572:   return(0);
573: }