Actual source code: svec.c

slepc-3.8.2 2017-12-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, 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 "./svecimpl.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:   BV_SVEC        *v = (BV_SVEC*)V->data,*w = (BV_SVEC*)W->data;
210:   PetscScalar    *pv,*pw,*pb,*pc;
211:   PetscInt       j,m;
212:   PetscBool      flg;

215:   VecGetArray(v->v,&pv);
216:   VecGetArray(w->v,&pw);
217:   MatHasOperation(A,MATOP_MAT_MULT,&flg);
218:   if (V->vmm && flg) {
219:     m = V->k-V->l;
220:     if (V->vmm==BV_MATMULT_MAT_SAVE) {
221:       BV_AllocateMatMult(V,A,m);
222:       MatDenseGetArray(V->B,&pb);
223:       PetscMemcpy(pb,pv+(V->nc+V->l)*V->n,m*V->n*sizeof(PetscScalar));
224:       MatDenseRestoreArray(V->B,&pb);
225:     } else {  /* BV_MATMULT_MAT */
226:       MatCreateDense(PetscObjectComm((PetscObject)V),V->n,PETSC_DECIDE,V->N,m,pv+(V->nc+V->l)*V->n,&V->B);
227:     }
228:     if (!V->C) {
229:       MatMatMultSymbolic(A,V->B,PETSC_DEFAULT,&V->C);
230:     }
231:     MatMatMultNumeric(A,V->B,V->C);
232:     MatDenseGetArray(V->C,&pc);
233:     PetscMemcpy(pw+(W->nc+W->l)*W->n,pc,m*V->n*sizeof(PetscScalar));
234:     MatDenseRestoreArray(V->C,&pc);
235:     if (V->vmm==BV_MATMULT_MAT) {
236:       MatDestroy(&V->B);
237:       MatDestroy(&V->C);
238:     }
239:   } else {
240:     for (j=0;j<V->k-V->l;j++) {
241:       VecPlaceArray(V->cv[1],pv+(V->nc+V->l+j)*V->n);
242:       VecPlaceArray(W->cv[1],pw+(W->nc+W->l+j)*W->n);
243:       MatMult(A,V->cv[1],W->cv[1]);
244:       VecResetArray(V->cv[1]);
245:       VecResetArray(W->cv[1]);
246:     }
247:   }
248:   VecRestoreArray(v->v,&pv);
249:   VecRestoreArray(w->v,&pw);
250:   return(0);
251: }

253: PetscErrorCode BVCopy_Svec(BV V,BV W)
254: {
256:   BV_SVEC        *v = (BV_SVEC*)V->data,*w = (BV_SVEC*)W->data;
257:   PetscScalar    *pv,*pw,*pvc,*pwc;

260:   VecGetArray(v->v,&pv);
261:   VecGetArray(w->v,&pw);
262:   pvc = pv+(V->nc+V->l)*V->n;
263:   pwc = pw+(W->nc+W->l)*W->n;
264:   PetscMemcpy(pwc,pvc,(V->k-V->l)*V->n*sizeof(PetscScalar));
265:   VecRestoreArray(v->v,&pv);
266:   VecRestoreArray(w->v,&pw);
267:   return(0);
268: }

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

280:   VecGetBlockSize(bv->t,&bs);
281:   VecCreate(PetscObjectComm((PetscObject)bv->t),&vnew);
282:   VecSetType(vnew,((PetscObject)bv->t)->type_name);
283:   VecSetSizes(vnew,m*bv->n,PETSC_DECIDE);
284:   VecSetBlockSize(vnew,bs);
285:   PetscLogObjectParent((PetscObject)bv,(PetscObject)vnew);
286:   if (((PetscObject)bv)->name) {
287:     PetscSNPrintf(str,50,"%s_0",((PetscObject)bv)->name);
288:     PetscObjectSetName((PetscObject)vnew,str);
289:   }
290:   if (copy) {
291:     VecGetArray(ctx->v,&pv);
292:     VecGetArray(vnew,&pnew);
293:     PetscMemcpy(pnew,pv,PetscMin(m,bv->m)*bv->n*sizeof(PetscScalar));
294:     VecRestoreArray(ctx->v,&pv);
295:     VecRestoreArray(vnew,&pnew);
296:   }
297:   VecDestroy(&ctx->v);
298:   ctx->v = vnew;
299:   return(0);
300: }

302: PetscErrorCode BVGetColumn_Svec(BV bv,PetscInt j,Vec *v)
303: {
305:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;
306:   PetscScalar    *pv;
307:   PetscInt       l;

310:   l = BVAvailableVec;
311:   VecGetArray(ctx->v,&pv);
312:   VecPlaceArray(bv->cv[l],pv+(bv->nc+j)*bv->n);
313:   return(0);
314: }

316: PetscErrorCode BVRestoreColumn_Svec(BV bv,PetscInt j,Vec *v)
317: {
319:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;
320:   PetscInt       l;

323:   l = (j==bv->ci[0])? 0: 1;
324:   VecResetArray(bv->cv[l]);
325:   VecRestoreArray(ctx->v,NULL);
326:   return(0);
327: }

329: PetscErrorCode BVGetArray_Svec(BV bv,PetscScalar **a)
330: {
332:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;

335:   VecGetArray(ctx->v,a);
336:   return(0);
337: }

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

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

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

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

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

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

369: PetscErrorCode BVView_Svec(BV bv,PetscViewer viewer)
370: {
371:   PetscErrorCode    ierr;
372:   BV_SVEC           *ctx = (BV_SVEC*)bv->data;
373:   PetscViewerFormat format;
374:   PetscBool         isascii;
375:   const char        *bvname,*name;

378:   VecView(ctx->v,viewer);
379:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
380:   if (isascii) {
381:     PetscViewerGetFormat(viewer,&format);
382:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
383:       PetscObjectGetName((PetscObject)bv,&bvname);
384:       PetscObjectGetName((PetscObject)ctx->v,&name);
385:       PetscViewerASCIIPrintf(viewer,"%s=reshape(%s,%D,%D);clear %s\n",bvname,name,bv->N,bv->nc+bv->m,name);
386:       if (bv->nc) {
387:         PetscViewerASCIIPrintf(viewer,"%s=%s(:,%D:end);\n",bvname,bvname,bv->nc+1);
388:       }
389:     }
390:   }
391:   return(0);
392: }

394: PetscErrorCode BVDestroy_Svec(BV bv)
395: {
397:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;

400:   VecDestroy(&ctx->v);
401:   VecDestroy(&bv->cv[0]);
402:   VecDestroy(&bv->cv[1]);
403:   PetscFree(bv->data);
404:   bv->cuda = PETSC_FALSE;
405:   return(0);
406: }

408: PETSC_EXTERN PetscErrorCode BVCreate_Svec(BV bv)
409: {
411:   BV_SVEC        *ctx;
412:   PetscInt       nloc,N,bs,tlocal;
413:   PetscBool      seq;
414:   PetscScalar    *aa,*vv;
415:   char           str[50];

418:   PetscNewLog(bv,&ctx);
419:   bv->data = (void*)ctx;

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

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

427:   VecGetLocalSize(bv->t,&nloc);
428:   VecGetSize(bv->t,&N);
429:   VecGetBlockSize(bv->t,&bs);
430:   tlocal = bv->m*nloc;
431:   if (tlocal<0) SETERRQ2(PetscObjectComm((PetscObject)bv),1,"The product %D times %D overflows the size of PetscInt; consider reducing the number of columns, or use BVVECS instead",bv->m,nloc);

433:   VecCreate(PetscObjectComm((PetscObject)bv->t),&ctx->v);
434:   VecSetType(ctx->v,((PetscObject)bv->t)->type_name);
435:   VecSetSizes(ctx->v,tlocal,bv->m*N);
436:   VecSetBlockSize(ctx->v,bs);
437:   PetscLogObjectParent((PetscObject)bv,(PetscObject)ctx->v);
438:   if (((PetscObject)bv)->name) {
439:     PetscSNPrintf(str,50,"%s_0",((PetscObject)bv)->name);
440:     PetscObjectSetName((PetscObject)ctx->v,str);
441:   }

443:   if (bv->Acreate) {
444:     MatDenseGetArray(bv->Acreate,&aa);
445:     VecGetArray(ctx->v,&vv);
446:     PetscMemcpy(vv,aa,tlocal*sizeof(PetscScalar));
447:     VecRestoreArray(ctx->v,&vv);
448:     MatDenseRestoreArray(bv->Acreate,&aa);
449:     MatDestroy(&bv->Acreate);
450:   }

452:   VecDuplicateEmpty(bv->t,&bv->cv[0]);
453:   VecDuplicateEmpty(bv->t,&bv->cv[1]);
454:   VecSetType(bv->cv[0],((PetscObject)bv->t)->type_name);
455:   VecSetType(bv->cv[1],((PetscObject)bv->t)->type_name);
456:   if (bv->cuda) {
457: #if defined(PETSC_HAVE_VECCUDA)
458:     bv->ops->mult             = BVMult_Svec_CUDA;
459:     bv->ops->multvec          = BVMultVec_Svec_CUDA;
460:     bv->ops->multinplace      = BVMultInPlace_Svec_CUDA;
461:     bv->ops->multinplacetrans = BVMultInPlaceTranspose_Svec_CUDA;
462:     bv->ops->dot              = BVDot_Svec_CUDA;
463:     bv->ops->dotvec           = BVDotVec_Svec_CUDA;
464:     bv->ops->dotvec_local     = BVDotVec_Local_Svec_CUDA;
465:     bv->ops->scale            = BVScale_Svec_CUDA;
466:     bv->ops->matmult          = BVMatMult_Svec_CUDA;
467:     bv->ops->copy             = BVCopy_Svec_CUDA;
468:     bv->ops->resize           = BVResize_Svec_CUDA;
469:     bv->ops->getcolumn        = BVGetColumn_Svec_CUDA;
470:     bv->ops->restorecolumn    = BVRestoreColumn_Svec_CUDA;
471: #endif
472:   } else {
473:     bv->ops->mult             = BVMult_Svec;
474:     bv->ops->multvec          = BVMultVec_Svec;
475:     bv->ops->multinplace      = BVMultInPlace_Svec;
476:     bv->ops->multinplacetrans = BVMultInPlaceTranspose_Svec;
477:     bv->ops->dot              = BVDot_Svec;
478:     bv->ops->dotvec           = BVDotVec_Svec;
479:     bv->ops->dotvec_local     = BVDotVec_Local_Svec;
480:     bv->ops->scale            = BVScale_Svec;
481:     bv->ops->matmult          = BVMatMult_Svec;
482:     bv->ops->copy             = BVCopy_Svec;
483:     bv->ops->resize           = BVResize_Svec;
484:     bv->ops->getcolumn        = BVGetColumn_Svec;
485:     bv->ops->restorecolumn    = BVRestoreColumn_Svec;
486:   }
487:   bv->ops->norm             = BVNorm_Svec;
488:   bv->ops->norm_local       = BVNorm_Local_Svec;
489:   bv->ops->getarray         = BVGetArray_Svec;
490:   bv->ops->restorearray     = BVRestoreArray_Svec;
491:   bv->ops->getarrayread     = BVGetArrayRead_Svec;
492:   bv->ops->restorearrayread = BVRestoreArrayRead_Svec;
493:   bv->ops->destroy          = BVDestroy_Svec;
494:   if (!ctx->mpi) bv->ops->view = BVView_Svec;
495:   return(0);
496: }