Actual source code: bvmat.c

slepc-3.7.1 2016-05-27
Report Typos and Errors
  1: /*
  2:    BV implemented with a dense Mat

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2016, 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:   Mat       A;
 28:   PetscBool mpi;
 29: } BV_MAT;

 33: PetscErrorCode BVMult_Mat(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
 34: {
 36:   BV_MAT         *y = (BV_MAT*)Y->data,*x = (BV_MAT*)X->data;
 37:   PetscScalar    *px,*py,*q;
 38:   PetscInt       ldq;

 41:   MatDenseGetArray(x->A,&px);
 42:   MatDenseGetArray(y->A,&py);
 43:   if (Q) {
 44:     MatGetSize(Q,&ldq,NULL);
 45:     MatDenseGetArray(Q,&q);
 46:     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);
 47:     MatDenseRestoreArray(Q,&q);
 48:   } else {
 49:     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);
 50:   }
 51:   MatDenseRestoreArray(x->A,&px);
 52:   MatDenseRestoreArray(y->A,&py);
 53:   return(0);
 54: }

 58: PetscErrorCode BVMultVec_Mat(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
 59: {
 61:   BV_MAT         *x = (BV_MAT*)X->data;
 62:   PetscScalar    *px,*py;

 65:   MatDenseGetArray(x->A,&px);
 66:   VecGetArray(y,&py);
 67:   BVMultVec_BLAS_Private(X,X->n,X->k-X->l,alpha,px+(X->nc+X->l)*X->n,q,beta,py);
 68:   MatDenseRestoreArray(x->A,&px);
 69:   VecRestoreArray(y,&py);
 70:   return(0);
 71: }

 75: PetscErrorCode BVMultInPlace_Mat(BV V,Mat Q,PetscInt s,PetscInt e)
 76: {
 78:   BV_MAT         *ctx = (BV_MAT*)V->data;
 79:   PetscScalar    *pv,*q;
 80:   PetscInt       ldq;

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

 94: PetscErrorCode BVMultInPlaceTranspose_Mat(BV V,Mat Q,PetscInt s,PetscInt e)
 95: {
 97:   BV_MAT         *ctx = (BV_MAT*)V->data;
 98:   PetscScalar    *pv,*q;
 99:   PetscInt       ldq;

102:   MatGetSize(Q,&ldq,NULL);
103:   MatDenseGetArray(ctx->A,&pv);
104:   MatDenseGetArray(Q,&q);
105:   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);
106:   MatDenseRestoreArray(Q,&q);
107:   MatDenseRestoreArray(ctx->A,&pv);
108:   return(0);
109: }

113: PetscErrorCode BVDot_Mat(BV X,BV Y,Mat M)
114: {
116:   BV_MAT         *x = (BV_MAT*)X->data,*y = (BV_MAT*)Y->data;
117:   PetscScalar    *px,*py,*m;
118:   PetscInt       ldm;

121:   MatGetSize(M,&ldm,NULL);
122:   MatDenseGetArray(x->A,&px);
123:   MatDenseGetArray(y->A,&py);
124:   MatDenseGetArray(M,&m);
125:   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);
126:   MatDenseRestoreArray(M,&m);
127:   MatDenseRestoreArray(x->A,&px);
128:   MatDenseRestoreArray(y->A,&py);
129:   return(0);
130: }

134: PetscErrorCode BVDotVec_Mat(BV X,Vec y,PetscScalar *m)
135: {
136:   PetscErrorCode    ierr;
137:   BV_MAT            *x = (BV_MAT*)X->data;
138:   PetscScalar       *px;
139:   const PetscScalar *py;
140:   Vec               z = y;

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

157: PetscErrorCode BVDotVec_Local_Mat(BV X,Vec y,PetscScalar *m)
158: {
160:   BV_MAT         *x = (BV_MAT*)X->data;
161:   PetscScalar    *px,*py;
162:   Vec            z = y;

165:   if (X->matrix) {
166:     BV_IPMatMult(X,y);
167:     z = X->Bx;
168:   }
169:   MatDenseGetArray(x->A,&px);
170:   VecGetArray(z,&py);
171:   BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,m,PETSC_FALSE);
172:   VecRestoreArray(z,&py);
173:   MatDenseRestoreArray(x->A,&px);
174:   return(0);
175: }

179: PetscErrorCode BVScale_Mat(BV bv,PetscInt j,PetscScalar alpha)
180: {
182:   BV_MAT         *ctx = (BV_MAT*)bv->data;
183:   PetscScalar    *array;

186:   MatDenseGetArray(ctx->A,&array);
187:   if (j<0) {
188:     BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->n,array+(bv->nc+bv->l)*bv->n,alpha);
189:   } else {
190:     BVScale_BLAS_Private(bv,bv->n,array+(bv->nc+j)*bv->n,alpha);
191:   }
192:   MatDenseRestoreArray(ctx->A,&array);
193:   return(0);
194: }

198: PetscErrorCode BVNorm_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
199: {
201:   BV_MAT         *ctx = (BV_MAT*)bv->data;
202:   PetscScalar    *array;

205:   MatDenseGetArray(ctx->A,&array);
206:   if (j<0) {
207:     BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,ctx->mpi);
208:   } else {
209:     BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,ctx->mpi);
210:   }
211:   MatDenseRestoreArray(ctx->A,&array);
212:   return(0);
213: }

217: PetscErrorCode BVNorm_Local_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
218: {
220:   BV_MAT         *ctx = (BV_MAT*)bv->data;
221:   PetscScalar    *array;

224:   MatDenseGetArray(ctx->A,&array);
225:   if (j<0) {
226:     BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,PETSC_FALSE);
227:   } else {
228:     BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,PETSC_FALSE);
229:   }
230:   MatDenseRestoreArray(ctx->A,&array);
231:   return(0);
232: }

236: PetscErrorCode BVOrthogonalize_Mat(BV V,Mat R)
237: {
239:   BV_MAT         *ctx = (BV_MAT*)V->data;
240:   PetscScalar    *pv,*r=NULL;

243:   if (R) { MatDenseGetArray(R,&r); }
244:   MatDenseGetArray(ctx->A,&pv);
245:   BVOrthogonalize_LAPACK_Private(V,V->n,V->k,pv+V->nc*V->n,r,ctx->mpi);
246:   MatDenseRestoreArray(ctx->A,&pv);
247:   if (R) { MatDenseRestoreArray(R,&r); }
248:   return(0);
249: }

253: PetscErrorCode BVMatMult_Mat(BV V,Mat A,BV W)
254: {
256:   BV_MAT         *v = (BV_MAT*)V->data,*w = (BV_MAT*)W->data;
257:   PetscScalar    *pv,*pw,*pb,*pc;
258:   PetscInt       j,m;
259:   PetscBool      flg;

262:   MatDenseGetArray(v->A,&pv);
263:   MatDenseGetArray(w->A,&pw);
264:   MatHasOperation(A,MATOP_MAT_MULT,&flg);
265:   if (V->vmm && flg) {
266:     m = V->k-V->l;
267:     if (V->vmm==BV_MATMULT_MAT_SAVE) {
268:       BV_AllocateMatMult(V,A,m);
269:       MatDenseGetArray(V->B,&pb);
270:       PetscMemcpy(pb,pv+(V->nc+V->l)*V->n,m*V->n*sizeof(PetscScalar));
271:       MatDenseRestoreArray(V->B,&pb);
272:     } else {  /* BV_MATMULT_MAT */
273:       MatCreateDense(PetscObjectComm((PetscObject)V),V->n,PETSC_DECIDE,V->N,m,pv+(V->nc+V->l)*V->n,&V->B);
274:     }
275:     if (!V->C) {
276:       MatMatMultSymbolic(A,V->B,PETSC_DEFAULT,&V->C);
277:     }
278:     MatMatMultNumeric(A,V->B,V->C);
279:     MatDenseGetArray(V->C,&pc);
280:     PetscMemcpy(pw+(W->nc+W->l)*W->n,pc,m*V->n*sizeof(PetscScalar));
281:     MatDenseRestoreArray(V->C,&pc);
282:     if (V->vmm==BV_MATMULT_MAT) {
283:       MatDestroy(&V->B);
284:       MatDestroy(&V->C);
285:     }
286:   } else {
287:     for (j=0;j<V->k-V->l;j++) {
288:       VecPlaceArray(V->cv[1],pv+(V->nc+V->l+j)*V->n);
289:       VecPlaceArray(W->cv[1],pw+(W->nc+W->l+j)*W->n);
290:       MatMult(A,V->cv[1],W->cv[1]);
291:       VecResetArray(V->cv[1]);
292:       VecResetArray(W->cv[1]);
293:     }
294:   }
295:   MatDenseRestoreArray(v->A,&pv);
296:   MatDenseRestoreArray(w->A,&pw);
297:   return(0);
298: }

302: PetscErrorCode BVCopy_Mat(BV V,BV W)
303: {
305:   BV_MAT         *v = (BV_MAT*)V->data,*w = (BV_MAT*)W->data;
306:   PetscScalar    *pv,*pw,*pvc,*pwc;

309:   MatDenseGetArray(v->A,&pv);
310:   MatDenseGetArray(w->A,&pw);
311:   pvc = pv+(V->nc+V->l)*V->n;
312:   pwc = pw+(W->nc+W->l)*W->n;
313:   PetscMemcpy(pwc,pvc,(V->k-V->l)*V->n*sizeof(PetscScalar));
314:   MatDenseRestoreArray(v->A,&pv);
315:   MatDenseRestoreArray(w->A,&pw);
316:   return(0);
317: }

321: PetscErrorCode BVResize_Mat(BV bv,PetscInt m,PetscBool copy)
322: {
324:   BV_MAT         *ctx = (BV_MAT*)bv->data;
325:   PetscScalar    *pA,*pnew;
326:   Mat            A;
327:   char           str[50];

330:   MatCreateDense(PetscObjectComm((PetscObject)bv->t),bv->n,PETSC_DECIDE,PETSC_DECIDE,m,NULL,&A);
331:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
332:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
333:   PetscLogObjectParent((PetscObject)bv,(PetscObject)A);
334:   if (((PetscObject)bv)->name) {
335:     PetscSNPrintf(str,50,"%s_0",((PetscObject)bv)->name);
336:     PetscObjectSetName((PetscObject)A,str);
337:   }
338:   if (copy) {
339:     MatDenseGetArray(ctx->A,&pA);
340:     MatDenseGetArray(A,&pnew);
341:     PetscMemcpy(pnew,pA,PetscMin(m,bv->m)*bv->n*sizeof(PetscScalar));
342:     MatDenseRestoreArray(ctx->A,&pA);
343:     MatDenseRestoreArray(A,&pnew);
344:   }
345:   MatDestroy(&ctx->A);
346:   ctx->A = A;
347:   return(0);
348: }

352: PetscErrorCode BVGetColumn_Mat(BV bv,PetscInt j,Vec *v)
353: {
355:   BV_MAT         *ctx = (BV_MAT*)bv->data;
356:   PetscScalar    *pA;
357:   PetscInt       l;

360:   l = BVAvailableVec;
361:   MatDenseGetArray(ctx->A,&pA);
362:   VecPlaceArray(bv->cv[l],pA+(bv->nc+j)*bv->n);
363:   return(0);
364: }

368: PetscErrorCode BVRestoreColumn_Mat(BV bv,PetscInt j,Vec *v)
369: {
371:   BV_MAT         *ctx = (BV_MAT*)bv->data;
372:   PetscScalar    *pA;
373:   PetscInt       l;

376:   l = (j==bv->ci[0])? 0: 1;
377:   VecResetArray(bv->cv[l]);
378:   MatDenseRestoreArray(ctx->A,&pA);
379:   return(0);
380: }

384: PetscErrorCode BVGetArray_Mat(BV bv,PetscScalar **a)
385: {
387:   BV_MAT         *ctx = (BV_MAT*)bv->data;

390:   MatDenseGetArray(ctx->A,a);
391:   return(0);
392: }

396: PetscErrorCode BVRestoreArray_Mat(BV bv,PetscScalar **a)
397: {
399:   BV_MAT         *ctx = (BV_MAT*)bv->data;

402:   if (a) { MatDenseRestoreArray(ctx->A,a); }
403:   return(0);
404: }

408: PetscErrorCode BVGetArrayRead_Mat(BV bv,const PetscScalar **a)
409: {
411:   BV_MAT         *ctx = (BV_MAT*)bv->data;

414:   MatDenseGetArray(ctx->A,(PetscScalar**)a);
415:   return(0);
416: }

420: PetscErrorCode BVRestoreArrayRead_Mat(BV bv,const PetscScalar **a)
421: {
423:   BV_MAT         *ctx = (BV_MAT*)bv->data;

426:   if (a) { MatDenseRestoreArray(ctx->A,(PetscScalar**)a); }
427:   return(0);
428: }

432: PetscErrorCode BVView_Mat(BV bv,PetscViewer viewer)
433: {
434:   PetscErrorCode    ierr;
435:   BV_MAT            *ctx = (BV_MAT*)bv->data;
436:   PetscViewerFormat format;
437:   PetscBool         isascii;
438:   const char        *bvname,*name;

441:   MatView(ctx->A,viewer);
442:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
443:   if (isascii) {
444:     PetscViewerGetFormat(viewer,&format);
445:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
446:       PetscObjectGetName((PetscObject)bv,&bvname);
447:       PetscObjectGetName((PetscObject)ctx->A,&name);
448:       PetscViewerASCIIPrintf(viewer,"%s=%s;clear %s\n",bvname,name,name);
449:       if (bv->nc) {
450:         PetscViewerASCIIPrintf(viewer,"%s=%s(:,%D:end);\n",bvname,bvname,bv->nc+1);
451:       }
452:     }
453:   }
454:   return(0);
455: }

459: PetscErrorCode BVDestroy_Mat(BV bv)
460: {
462:   BV_MAT         *ctx = (BV_MAT*)bv->data;

465:   MatDestroy(&ctx->A);
466:   VecDestroy(&bv->cv[0]);
467:   VecDestroy(&bv->cv[1]);
468:   PetscFree(bv->data);
469:   return(0);
470: }

474: PETSC_EXTERN PetscErrorCode BVCreate_Mat(BV bv)
475: {
477:   BV_MAT         *ctx;
478:   PetscInt       nloc,bs;
479:   PetscBool      seq;
480:   char           str[50];

483:   PetscNewLog(bv,&ctx);
484:   bv->data = (void*)ctx;

486:   PetscObjectTypeCompare((PetscObject)bv->t,VECMPI,&ctx->mpi);
487:   if (!ctx->mpi) {
488:     PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
489:     if (!seq) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot create a BVMAT from a non-standard template vector");
490:   }

492:   VecGetLocalSize(bv->t,&nloc);
493:   VecGetBlockSize(bv->t,&bs);

495:   MatCreateDense(PetscObjectComm((PetscObject)bv->t),nloc,PETSC_DECIDE,PETSC_DECIDE,bv->m,NULL,&ctx->A);
496:   MatAssemblyBegin(ctx->A,MAT_FINAL_ASSEMBLY);
497:   MatAssemblyEnd(ctx->A,MAT_FINAL_ASSEMBLY);
498:   PetscLogObjectParent((PetscObject)bv,(PetscObject)ctx->A);
499:   if (((PetscObject)bv)->name) {
500:     PetscSNPrintf(str,50,"%s_0",((PetscObject)bv)->name);
501:     PetscObjectSetName((PetscObject)ctx->A,str);
502:   }

504:   if (ctx->mpi) {
505:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,NULL,&bv->cv[0]);
506:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,NULL,&bv->cv[1]);
507:   } else {
508:     VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,NULL,&bv->cv[0]);
509:     VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,NULL,&bv->cv[1]);
510:   }

512:   bv->ops->mult             = BVMult_Mat;
513:   bv->ops->multvec          = BVMultVec_Mat;
514:   bv->ops->multinplace      = BVMultInPlace_Mat;
515:   bv->ops->multinplacetrans = BVMultInPlaceTranspose_Mat;
516:   bv->ops->dot              = BVDot_Mat;
517:   bv->ops->dotvec           = BVDotVec_Mat;
518:   bv->ops->dotvec_local     = BVDotVec_Local_Mat;
519:   bv->ops->scale            = BVScale_Mat;
520:   bv->ops->norm             = BVNorm_Mat;
521:   bv->ops->norm_local       = BVNorm_Local_Mat;
522:   /*bv->ops->orthogonalize    = BVOrthogonalize_Mat;*/
523:   bv->ops->matmult          = BVMatMult_Mat;
524:   bv->ops->copy             = BVCopy_Mat;
525:   bv->ops->resize           = BVResize_Mat;
526:   bv->ops->getcolumn        = BVGetColumn_Mat;
527:   bv->ops->restorecolumn    = BVRestoreColumn_Mat;
528:   bv->ops->getarray         = BVGetArray_Mat;
529:   bv->ops->restorearray     = BVRestoreArray_Mat;
530:   bv->ops->getarrayread     = BVGetArrayRead_Mat;
531:   bv->ops->restorearrayread = BVRestoreArrayRead_Mat;
532:   bv->ops->destroy          = BVDestroy_Mat;
533:   if (!ctx->mpi) bv->ops->view = BVView_Mat;
534:   return(0);
535: }