Actual source code: bvmat.c

slepc-3.8.0 2017-10-20
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 with a dense Mat
 12: */

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

 16: typedef struct {
 17:   Mat       A;
 18:   PetscBool mpi;
 19: } BV_MAT;

 21: PetscErrorCode BVMult_Mat(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
 22: {
 24:   BV_MAT         *y = (BV_MAT*)Y->data,*x = (BV_MAT*)X->data;
 25:   PetscScalar    *px,*py,*q;
 26:   PetscInt       ldq;

 29:   MatDenseGetArray(x->A,&px);
 30:   MatDenseGetArray(y->A,&py);
 31:   if (Q) {
 32:     MatGetSize(Q,&ldq,NULL);
 33:     MatDenseGetArray(Q,&q);
 34:     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);
 35:     MatDenseRestoreArray(Q,&q);
 36:   } else {
 37:     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);
 38:   }
 39:   MatDenseRestoreArray(x->A,&px);
 40:   MatDenseRestoreArray(y->A,&py);
 41:   return(0);
 42: }

 44: PetscErrorCode BVMultVec_Mat(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
 45: {
 47:   BV_MAT         *x = (BV_MAT*)X->data;
 48:   PetscScalar    *px,*py,*qq=q;

 51:   MatDenseGetArray(x->A,&px);
 52:   VecGetArray(y,&py);
 53:   if (!q) { VecGetArray(X->buffer,&qq); }
 54:   BVMultVec_BLAS_Private(X,X->n,X->k-X->l,alpha,px+(X->nc+X->l)*X->n,qq,beta,py);
 55:   if (!q) { VecRestoreArray(X->buffer,&qq); }
 56:   MatDenseRestoreArray(x->A,&px);
 57:   VecRestoreArray(y,&py);
 58:   return(0);
 59: }

 61: PetscErrorCode BVMultInPlace_Mat(BV V,Mat Q,PetscInt s,PetscInt e)
 62: {
 64:   BV_MAT         *ctx = (BV_MAT*)V->data;
 65:   PetscScalar    *pv,*q;
 66:   PetscInt       ldq;

 69:   MatGetSize(Q,&ldq,NULL);
 70:   MatDenseGetArray(ctx->A,&pv);
 71:   MatDenseGetArray(Q,&q);
 72:   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);
 73:   MatDenseRestoreArray(Q,&q);
 74:   MatDenseRestoreArray(ctx->A,&pv);
 75:   return(0);
 76: }

 78: PetscErrorCode BVMultInPlaceTranspose_Mat(BV V,Mat Q,PetscInt s,PetscInt e)
 79: {
 81:   BV_MAT         *ctx = (BV_MAT*)V->data;
 82:   PetscScalar    *pv,*q;
 83:   PetscInt       ldq;

 86:   MatGetSize(Q,&ldq,NULL);
 87:   MatDenseGetArray(ctx->A,&pv);
 88:   MatDenseGetArray(Q,&q);
 89:   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);
 90:   MatDenseRestoreArray(Q,&q);
 91:   MatDenseRestoreArray(ctx->A,&pv);
 92:   return(0);
 93: }

 95: PetscErrorCode BVDot_Mat(BV X,BV Y,Mat M)
 96: {
 98:   BV_MAT         *x = (BV_MAT*)X->data,*y = (BV_MAT*)Y->data;
 99:   PetscScalar    *px,*py,*m;
100:   PetscInt       ldm;

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

114: PetscErrorCode BVDotVec_Mat(BV X,Vec y,PetscScalar *q)
115: {
116:   PetscErrorCode    ierr;
117:   BV_MAT            *x = (BV_MAT*)X->data;
118:   PetscScalar       *px,*qq=q;
119:   const PetscScalar *py;
120:   Vec               z = y;

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

137: PetscErrorCode BVDotVec_Local_Mat(BV X,Vec y,PetscScalar *m)
138: {
140:   BV_MAT         *x = (BV_MAT*)X->data;
141:   PetscScalar    *px,*py;
142:   Vec            z = y;

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

157: PetscErrorCode BVScale_Mat(BV bv,PetscInt j,PetscScalar alpha)
158: {
160:   BV_MAT         *ctx = (BV_MAT*)bv->data;
161:   PetscScalar    *array;

164:   MatDenseGetArray(ctx->A,&array);
165:   if (j<0) {
166:     BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->n,array+(bv->nc+bv->l)*bv->n,alpha);
167:   } else {
168:     BVScale_BLAS_Private(bv,bv->n,array+(bv->nc+j)*bv->n,alpha);
169:   }
170:   MatDenseRestoreArray(ctx->A,&array);
171:   return(0);
172: }

174: PetscErrorCode BVNorm_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
175: {
177:   BV_MAT         *ctx = (BV_MAT*)bv->data;
178:   PetscScalar    *array;

181:   MatDenseGetArray(ctx->A,&array);
182:   if (j<0) {
183:     BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,ctx->mpi);
184:   } else {
185:     BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,ctx->mpi);
186:   }
187:   MatDenseRestoreArray(ctx->A,&array);
188:   return(0);
189: }

191: PetscErrorCode BVNorm_Local_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
192: {
194:   BV_MAT         *ctx = (BV_MAT*)bv->data;
195:   PetscScalar    *array;

198:   MatDenseGetArray(ctx->A,&array);
199:   if (j<0) {
200:     BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,PETSC_FALSE);
201:   } else {
202:     BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,PETSC_FALSE);
203:   }
204:   MatDenseRestoreArray(ctx->A,&array);
205:   return(0);
206: }

208: PetscErrorCode BVMatMult_Mat(BV V,Mat A,BV W)
209: {
211:   BV_MAT         *v = (BV_MAT*)V->data,*w = (BV_MAT*)W->data;
212:   PetscScalar    *pv,*pw,*pb,*pc;
213:   PetscInt       j,m;
214:   PetscBool      flg;

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

255: PetscErrorCode BVCopy_Mat(BV V,BV W)
256: {
258:   BV_MAT         *v = (BV_MAT*)V->data,*w = (BV_MAT*)W->data;
259:   PetscScalar    *pv,*pw,*pvc,*pwc;

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

272: PetscErrorCode BVResize_Mat(BV bv,PetscInt m,PetscBool copy)
273: {
275:   BV_MAT         *ctx = (BV_MAT*)bv->data;
276:   PetscScalar    *pA,*pnew;
277:   Mat            A;
278:   char           str[50];

281:   MatCreateDense(PetscObjectComm((PetscObject)bv->t),bv->n,PETSC_DECIDE,PETSC_DECIDE,m,NULL,&A);
282:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
283:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
284:   PetscLogObjectParent((PetscObject)bv,(PetscObject)A);
285:   if (((PetscObject)bv)->name) {
286:     PetscSNPrintf(str,50,"%s_0",((PetscObject)bv)->name);
287:     PetscObjectSetName((PetscObject)A,str);
288:   }
289:   if (copy) {
290:     MatDenseGetArray(ctx->A,&pA);
291:     MatDenseGetArray(A,&pnew);
292:     PetscMemcpy(pnew,pA,PetscMin(m,bv->m)*bv->n*sizeof(PetscScalar));
293:     MatDenseRestoreArray(ctx->A,&pA);
294:     MatDenseRestoreArray(A,&pnew);
295:   }
296:   MatDestroy(&ctx->A);
297:   ctx->A = A;
298:   return(0);
299: }

301: PetscErrorCode BVGetColumn_Mat(BV bv,PetscInt j,Vec *v)
302: {
304:   BV_MAT         *ctx = (BV_MAT*)bv->data;
305:   PetscScalar    *pA;
306:   PetscInt       l;

309:   l = BVAvailableVec;
310:   MatDenseGetArray(ctx->A,&pA);
311:   VecPlaceArray(bv->cv[l],pA+(bv->nc+j)*bv->n);
312:   return(0);
313: }

315: PetscErrorCode BVRestoreColumn_Mat(BV bv,PetscInt j,Vec *v)
316: {
318:   BV_MAT         *ctx = (BV_MAT*)bv->data;
319:   PetscScalar    *pA;
320:   PetscInt       l;

323:   l = (j==bv->ci[0])? 0: 1;
324:   VecResetArray(bv->cv[l]);
325:   MatDenseRestoreArray(ctx->A,&pA);
326:   return(0);
327: }

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

335:   MatDenseGetArray(ctx->A,a);
336:   return(0);
337: }

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

345:   if (a) { MatDenseRestoreArray(ctx->A,a); }
346:   return(0);
347: }

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

355:   MatDenseGetArray(ctx->A,(PetscScalar**)a);
356:   return(0);
357: }

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

365:   if (a) { MatDenseRestoreArray(ctx->A,(PetscScalar**)a); }
366:   return(0);
367: }

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

378:   MatView(ctx->A,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->A,&name);
385:       PetscViewerASCIIPrintf(viewer,"%s=%s;clear %s\n",bvname,name,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_Mat(BV bv)
395: {
397:   BV_MAT         *ctx = (BV_MAT*)bv->data;

400:   MatDestroy(&ctx->A);
401:   VecDestroy(&bv->cv[0]);
402:   VecDestroy(&bv->cv[1]);
403:   PetscFree(bv->data);
404:   return(0);
405: }

407: PETSC_EXTERN PetscErrorCode BVCreate_Mat(BV bv)
408: {
410:   BV_MAT         *ctx;
411:   PetscInt       nloc,bs;
412:   PetscBool      seq;
413:   char           str[50];

416:   PetscNewLog(bv,&ctx);
417:   bv->data = (void*)ctx;

419:   PetscObjectTypeCompare((PetscObject)bv->t,VECMPI,&ctx->mpi);
420:   if (!ctx->mpi) {
421:     PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
422:     if (!seq) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Cannot create a BVMAT from a non-standard template vector");
423:   }

425:   VecGetLocalSize(bv->t,&nloc);
426:   VecGetBlockSize(bv->t,&bs);

428:   MatCreateDense(PetscObjectComm((PetscObject)bv->t),nloc,PETSC_DECIDE,PETSC_DECIDE,bv->m,NULL,&ctx->A);
429:   MatAssemblyBegin(ctx->A,MAT_FINAL_ASSEMBLY);
430:   MatAssemblyEnd(ctx->A,MAT_FINAL_ASSEMBLY);
431:   PetscLogObjectParent((PetscObject)bv,(PetscObject)ctx->A);
432:   if (((PetscObject)bv)->name) {
433:     PetscSNPrintf(str,50,"%s_0",((PetscObject)bv)->name);
434:     PetscObjectSetName((PetscObject)ctx->A,str);
435:   }

437:   if (bv->Acreate) {
438:     MatCopy(bv->Acreate,ctx->A,SAME_NONZERO_PATTERN);
439:     MatDestroy(&bv->Acreate);
440:   }

442:   if (ctx->mpi) {
443:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,NULL,&bv->cv[0]);
444:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,NULL,&bv->cv[1]);
445:   } else {
446:     VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,NULL,&bv->cv[0]);
447:     VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,NULL,&bv->cv[1]);
448:   }

450:   bv->ops->mult             = BVMult_Mat;
451:   bv->ops->multvec          = BVMultVec_Mat;
452:   bv->ops->multinplace      = BVMultInPlace_Mat;
453:   bv->ops->multinplacetrans = BVMultInPlaceTranspose_Mat;
454:   bv->ops->dot              = BVDot_Mat;
455:   bv->ops->dotvec           = BVDotVec_Mat;
456:   bv->ops->dotvec_local     = BVDotVec_Local_Mat;
457:   bv->ops->scale            = BVScale_Mat;
458:   bv->ops->norm             = BVNorm_Mat;
459:   bv->ops->norm_local       = BVNorm_Local_Mat;
460:   bv->ops->matmult          = BVMatMult_Mat;
461:   bv->ops->copy             = BVCopy_Mat;
462:   bv->ops->resize           = BVResize_Mat;
463:   bv->ops->getcolumn        = BVGetColumn_Mat;
464:   bv->ops->restorecolumn    = BVRestoreColumn_Mat;
465:   bv->ops->getarray         = BVGetArray_Mat;
466:   bv->ops->restorearray     = BVRestoreArray_Mat;
467:   bv->ops->getarrayread     = BVGetArrayRead_Mat;
468:   bv->ops->restorearrayread = BVRestoreArrayRead_Mat;
469:   bv->ops->destroy          = BVDestroy_Mat;
470:   if (!ctx->mpi) bv->ops->view = BVView_Mat;
471:   return(0);
472: }