Actual source code: svec.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    BV implemented as a single Vec

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

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

 41:   MatGetSize(Q,&ldq,NULL);
 42:   VecGetArray(x->v,&px);
 43:   VecGetArray(y->v,&py);
 44:   MatDenseGetArray(Q,&q);
 45:   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);
 46:   MatDenseRestoreArray(Q,&q);
 47:   VecRestoreArray(x->v,&px);
 48:   VecRestoreArray(y->v,&py);
 49:   return(0);
 50: }

 54: PetscErrorCode BVMultVec_Svec(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
 55: {
 57:   BV_SVEC        *x = (BV_SVEC*)X->data;
 58:   PetscScalar    *px,*py;

 61:   VecGetArray(x->v,&px);
 62:   VecGetArray(y,&py);
 63:   BVMultVec_BLAS_Private(X,X->n,X->k-X->l,alpha,px+(X->nc+X->l)*X->n,q,beta,py);
 64:   VecRestoreArray(x->v,&px);
 65:   VecRestoreArray(y,&py);
 66:   return(0);
 67: }

 71: PetscErrorCode BVMultInPlace_Svec(BV V,Mat Q,PetscInt s,PetscInt e)
 72: {
 74:   BV_SVEC        *ctx = (BV_SVEC*)V->data;
 75:   PetscScalar    *pv,*q;
 76:   PetscInt       ldq;

 79:   MatGetSize(Q,&ldq,NULL);
 80:   VecGetArray(ctx->v,&pv);
 81:   MatDenseGetArray(Q,&q);
 82:   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);
 83:   MatDenseRestoreArray(Q,&q);
 84:   VecRestoreArray(ctx->v,&pv);
 85:   return(0);
 86: }

 90: PetscErrorCode BVMultInPlaceTranspose_Svec(BV V,Mat Q,PetscInt s,PetscInt e)
 91: {
 93:   BV_SVEC        *ctx = (BV_SVEC*)V->data;
 94:   PetscScalar    *pv,*q;
 95:   PetscInt       ldq;

 98:   MatGetSize(Q,&ldq,NULL);
 99:   VecGetArray(ctx->v,&pv);
100:   MatDenseGetArray(Q,&q);
101:   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);
102:   MatDenseRestoreArray(Q,&q);
103:   VecRestoreArray(ctx->v,&pv);
104:   return(0);
105: }

109: PetscErrorCode BVAXPY_Svec(BV Y,PetscScalar alpha,BV X)
110: {
112:   BV_SVEC        *x = (BV_SVEC*)X->data,*y = (BV_SVEC*)Y->data;
113:   PetscScalar    *px,*py;

116:   VecGetArray(x->v,&px);
117:   VecGetArray(y->v,&py);
118:   BVAXPY_BLAS_Private(Y,Y->n,Y->k-Y->l,alpha,px+(X->nc+X->l)*X->n,py+(Y->nc+Y->l)*Y->n);
119:   VecRestoreArray(x->v,&px);
120:   VecRestoreArray(y->v,&py);
121:   return(0);
122: }

126: PetscErrorCode BVDot_Svec(BV X,BV Y,Mat M)
127: {
129:   BV_SVEC        *x = (BV_SVEC*)X->data,*y = (BV_SVEC*)Y->data;
130:   PetscScalar    *px,*py,*m;
131:   PetscInt       ldm;

134:   MatGetSize(M,&ldm,NULL);
135:   VecGetArray(x->v,&px);
136:   VecGetArray(y->v,&py);
137:   MatDenseGetArray(M,&m);
138:   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);
139:   MatDenseRestoreArray(M,&m);
140:   VecRestoreArray(x->v,&px);
141:   VecRestoreArray(y->v,&py);
142:   return(0);
143: }

147: PetscErrorCode BVDotVec_Svec(BV X,Vec y,PetscScalar *m)
148: {
150:   BV_SVEC        *x = (BV_SVEC*)X->data;
151:   PetscScalar    *px,*py;
152:   Vec            z = y;

155:   if (X->matrix) {
156:     BV_IPMatMult(X,y);
157:     z = X->Bx;
158:   }
159:   VecGetArray(x->v,&px);
160:   VecGetArray(z,&py);
161:   BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,m,x->mpi);
162:   VecRestoreArray(z,&py);
163:   VecRestoreArray(x->v,&px);
164:   return(0);
165: }

169: PetscErrorCode BVScale_Svec(BV bv,PetscInt j,PetscScalar alpha)
170: {
172:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;
173:   PetscScalar    *array;

176:   VecGetArray(ctx->v,&array);
177:   if (j<0) {
178:     BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->n,array+(bv->nc+bv->l)*bv->n,alpha);
179:   } else {
180:     BVScale_BLAS_Private(bv,bv->n,array+(bv->nc+j)*bv->n,alpha);
181:   }
182:   VecRestoreArray(ctx->v,&array);
183:   return(0);
184: }

188: PetscErrorCode BVNorm_Svec(BV bv,PetscInt j,NormType type,PetscReal *val)
189: {
191:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;
192:   PetscScalar    *array;

195:   VecGetArray(ctx->v,&array);
196:   if (j<0) {
197:     BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,ctx->mpi);
198:   } else {
199:     BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,ctx->mpi);
200:   }
201:   VecRestoreArray(ctx->v,&array);
202:   return(0);
203: }

207: PetscErrorCode BVOrthogonalize_Svec(BV V,Mat R)
208: {
210:   BV_SVEC        *ctx = (BV_SVEC*)V->data;
211:   PetscScalar    *pv,*r=NULL;

214:   if (R) { MatDenseGetArray(R,&r); }
215:   VecGetArray(ctx->v,&pv);
216:   BVOrthogonalize_LAPACK_Private(V,V->n,V->k,pv+V->nc*V->n,r,ctx->mpi);
217:   VecRestoreArray(ctx->v,&pv);
218:   if (R) { MatDenseRestoreArray(R,&r); }
219:   return(0);
220: }

224: PetscErrorCode BVMatMult_Svec(BV V,Mat A,BV W)
225: {
227:   BV_SVEC        *v = (BV_SVEC*)V->data,*w = (BV_SVEC*)W->data;
228:   PetscScalar    *pv,*pw;
229:   PetscInt       j;

232:   VecGetArray(v->v,&pv);
233:   VecGetArray(w->v,&pw);
234:   for (j=0;j<V->k-V->l;j++) {
235:     VecPlaceArray(V->cv[1],pv+(V->nc+V->l+j)*V->n);
236:     VecPlaceArray(W->cv[1],pw+(W->nc+W->l+j)*W->n);
237:     MatMult(A,V->cv[1],W->cv[1]);
238:     VecResetArray(V->cv[1]);
239:     VecResetArray(W->cv[1]);
240:   }
241:   VecRestoreArray(v->v,&pv);
242:   VecRestoreArray(w->v,&pw);
243:   return(0);
244: }

248: PetscErrorCode BVCopy_Svec(BV V,BV W)
249: {
251:   BV_SVEC        *v = (BV_SVEC*)V->data,*w = (BV_SVEC*)W->data;
252:   PetscScalar    *pv,*pw,*pvc,*pwc;

255:   VecGetArray(v->v,&pv);
256:   VecGetArray(w->v,&pw);
257:   pvc = pv+(V->nc+V->l)*V->n;
258:   pwc = pw+(W->nc+W->l)*W->n;
259:   PetscMemcpy(pwc,pvc,(V->k-V->l)*V->n*sizeof(PetscScalar));
260:   VecRestoreArray(v->v,&pv);
261:   VecRestoreArray(w->v,&pw);
262:   return(0);
263: }

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

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

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

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

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

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

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

338:   VecGetArray(ctx->v,a);
339:   return(0);
340: }

344: PetscErrorCode BVRestoreArray_Svec(BV bv,PetscScalar **a)
345: {
347:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;

350:   VecRestoreArray(ctx->v,a);
351:   return(0);
352: }

356: PetscErrorCode BVView_Svec(BV bv,PetscViewer viewer)
357: {
358:   PetscErrorCode    ierr;
359:   BV_SVEC           *ctx = (BV_SVEC*)bv->data;
360:   PetscViewerFormat format;
361:   PetscBool         isascii;

364:   VecView(ctx->v,viewer);
365:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
366:   if (isascii) {
367:     PetscViewerGetFormat(viewer,&format);
368:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
369:       PetscViewerASCIIPrintf(viewer,"%s=reshape(%s,%D,%D);clear %s\n",((PetscObject)bv)->name,((PetscObject)ctx->v)->name,bv->N,bv->nc+bv->m,((PetscObject)ctx->v)->name);
370:       if (bv->nc) {
371:         PetscViewerASCIIPrintf(viewer,"%s=%s(:,%D:end);\n",((PetscObject)bv)->name,((PetscObject)bv)->name,bv->nc+1);
372:       }
373:     }
374:   }
375:   return(0);
376: }

380: PetscErrorCode BVDestroy_Svec(BV bv)
381: {
383:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;

386:   VecDestroy(&ctx->v);
387:   VecDestroy(&bv->cv[0]);
388:   VecDestroy(&bv->cv[1]);
389:   PetscFree(bv->data);
390:   return(0);
391: }

395: PETSC_EXTERN PetscErrorCode BVCreate_Svec(BV bv)
396: {
398:   BV_SVEC        *ctx;
399:   PetscInt       nloc,bs;
400:   PetscBool      seq;
401:   char           str[50];

404:   PetscNewLog(bv,&ctx);
405:   bv->data = (void*)ctx;

407:   PetscObjectTypeCompare((PetscObject)bv->t,VECMPI,&ctx->mpi);
408:   if (!ctx->mpi) {
409:     PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
410:     if (!seq) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot create a BVSVEC from a non-standard template vector");
411:   }

413:   VecGetLocalSize(bv->t,&nloc);
414:   VecGetBlockSize(bv->t,&bs);

416:   VecCreate(PetscObjectComm((PetscObject)bv->t),&ctx->v);
417:   VecSetType(ctx->v,((PetscObject)bv->t)->type_name);
418:   VecSetSizes(ctx->v,bv->m*nloc,PETSC_DECIDE);
419:   VecSetBlockSize(ctx->v,bs);
420:   PetscLogObjectParent((PetscObject)bv,(PetscObject)ctx->v);
421:   if (((PetscObject)bv)->name) {
422:     PetscSNPrintf(str,50,"%s_0",((PetscObject)bv)->name);
423:     PetscObjectSetName((PetscObject)ctx->v,str);
424:   }

426:   if (ctx->mpi) {
427:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,NULL,&bv->cv[0]);
428:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,NULL,&bv->cv[1]);
429:   } else {
430:     VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,NULL,&bv->cv[0]);
431:     VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,NULL,&bv->cv[1]);
432:   }

434:   bv->ops->mult             = BVMult_Svec;
435:   bv->ops->multvec          = BVMultVec_Svec;
436:   bv->ops->multinplace      = BVMultInPlace_Svec;
437:   bv->ops->multinplacetrans = BVMultInPlaceTranspose_Svec;
438:   bv->ops->axpy             = BVAXPY_Svec;
439:   bv->ops->dot              = BVDot_Svec;
440:   bv->ops->dotvec           = BVDotVec_Svec;
441:   bv->ops->scale            = BVScale_Svec;
442:   bv->ops->norm             = BVNorm_Svec;
443:   /*bv->ops->orthogonalize    = BVOrthogonalize_Svec;*/
444:   bv->ops->matmult          = BVMatMult_Svec;
445:   bv->ops->copy             = BVCopy_Svec;
446:   bv->ops->resize           = BVResize_Svec;
447:   bv->ops->getcolumn        = BVGetColumn_Svec;
448:   bv->ops->restorecolumn    = BVRestoreColumn_Svec;
449:   bv->ops->getarray         = BVGetArray_Svec;
450:   bv->ops->restorearray     = BVRestoreArray_Svec;
451:   bv->ops->destroy          = BVDestroy_Svec;
452:   if (!ctx->mpi) bv->ops->view = BVView_Svec;
453:   return(0);
454: }