Actual source code: contig.c

slepc-3.7.2 2016-07-19
Report Typos and Errors
  1: /*
  2:    BV implemented as an array of Vecs sharing a contiguous array for elements

  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:   Vec         *V;
 28:   PetscScalar *array;
 29:   PetscBool   mpi;
 30: } BV_CONTIGUOUS;

 34: PetscErrorCode BVMult_Contiguous(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
 35: {
 37:   BV_CONTIGUOUS  *y = (BV_CONTIGUOUS*)Y->data,*x = (BV_CONTIGUOUS*)X->data;
 38:   PetscScalar    *q;
 39:   PetscInt       ldq;

 42:   if (Q) {
 43:     MatGetSize(Q,&ldq,NULL);
 44:     MatDenseGetArray(Q,&q);
 45:     BVMult_BLAS_Private(Y,Y->n,Y->k-Y->l,X->k-X->l,ldq,alpha,x->array+(X->nc+X->l)*X->n,q+Y->l*ldq+X->l,beta,y->array+(Y->nc+Y->l)*Y->n);
 46:     MatDenseRestoreArray(Q,&q);
 47:   } else {
 48:     BVAXPY_BLAS_Private(Y,Y->n,Y->k-Y->l,alpha,x->array+(X->nc+X->l)*X->n,beta,y->array+(Y->nc+Y->l)*Y->n);
 49:   }
 50:   return(0);
 51: }

 55: PetscErrorCode BVMultVec_Contiguous(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
 56: {
 58:   BV_CONTIGUOUS  *x = (BV_CONTIGUOUS*)X->data;
 59:   PetscScalar    *py;

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

 70: PetscErrorCode BVMultInPlace_Contiguous(BV V,Mat Q,PetscInt s,PetscInt e)
 71: {
 73:   BV_CONTIGUOUS  *ctx = (BV_CONTIGUOUS*)V->data;
 74:   PetscScalar    *q;
 75:   PetscInt       ldq;

 78:   MatGetSize(Q,&ldq,NULL);
 79:   MatDenseGetArray(Q,&q);
 80:   BVMultInPlace_BLAS_Private(V,V->n,V->k-V->l,ldq,s-V->l,e-V->l,ctx->array+(V->nc+V->l)*V->n,q+V->l*ldq+V->l,PETSC_FALSE);
 81:   MatDenseRestoreArray(Q,&q);
 82:   return(0);
 83: }

 87: PetscErrorCode BVMultInPlaceTranspose_Contiguous(BV V,Mat Q,PetscInt s,PetscInt e)
 88: {
 90:   BV_CONTIGUOUS  *ctx = (BV_CONTIGUOUS*)V->data;
 91:   PetscScalar    *q;
 92:   PetscInt       ldq;

 95:   MatGetSize(Q,&ldq,NULL);
 96:   MatDenseGetArray(Q,&q);
 97:   BVMultInPlace_BLAS_Private(V,V->n,V->k-V->l,ldq,s-V->l,e-V->l,ctx->array+(V->nc+V->l)*V->n,q+V->l*ldq+V->l,PETSC_TRUE);
 98:   MatDenseRestoreArray(Q,&q);
 99:   return(0);
100: }

104: PetscErrorCode BVDot_Contiguous(BV X,BV Y,Mat M)
105: {
107:   BV_CONTIGUOUS  *x = (BV_CONTIGUOUS*)X->data,*y = (BV_CONTIGUOUS*)Y->data;
108:   PetscScalar    *m;
109:   PetscInt       ldm;

112:   MatGetSize(M,&ldm,NULL);
113:   MatDenseGetArray(M,&m);
114:   BVDot_BLAS_Private(X,Y->k-Y->l,X->k-X->l,X->n,ldm,y->array+(Y->nc+Y->l)*Y->n,x->array+(X->nc+X->l)*X->n,m+X->l*ldm+Y->l,x->mpi);
115:   MatDenseRestoreArray(M,&m);
116:   return(0);
117: }

121: PetscErrorCode BVDotVec_Contiguous(BV X,Vec y,PetscScalar *m)
122: {
123:   PetscErrorCode    ierr;
124:   BV_CONTIGUOUS     *x = (BV_CONTIGUOUS*)X->data;
125:   const PetscScalar *py;
126:   Vec               z = y;

129:   if (X->matrix) {
130:     BV_IPMatMult(X,y);
131:     z = X->Bx;
132:   }
133:   VecGetArrayRead(z,&py);
134:   BVDotVec_BLAS_Private(X,X->n,X->k-X->l,x->array+(X->nc+X->l)*X->n,py,m,x->mpi);
135:   VecRestoreArrayRead(z,&py);
136:   return(0);
137: }

141: PetscErrorCode BVDotVec_Local_Contiguous(BV X,Vec y,PetscScalar *m)
142: {
144:   BV_CONTIGUOUS  *x = (BV_CONTIGUOUS*)X->data;
145:   PetscScalar    *py;
146:   Vec            z = y;

149:   if (X->matrix) {
150:     BV_IPMatMult(X,y);
151:     z = X->Bx;
152:   }
153:   VecGetArray(z,&py);
154:   BVDotVec_BLAS_Private(X,X->n,X->k-X->l,x->array+(X->nc+X->l)*X->n,py,m,PETSC_FALSE);
155:   VecRestoreArray(z,&py);
156:   return(0);
157: }

161: PetscErrorCode BVScale_Contiguous(BV bv,PetscInt j,PetscScalar alpha)
162: {
164:   BV_CONTIGUOUS  *ctx = (BV_CONTIGUOUS*)bv->data;

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

177: PetscErrorCode BVNorm_Contiguous(BV bv,PetscInt j,NormType type,PetscReal *val)
178: {
180:   BV_CONTIGUOUS  *ctx = (BV_CONTIGUOUS*)bv->data;

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

193: PetscErrorCode BVNorm_Local_Contiguous(BV bv,PetscInt j,NormType type,PetscReal *val)
194: {
196:   BV_CONTIGUOUS  *ctx = (BV_CONTIGUOUS*)bv->data;

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

209: PetscErrorCode BVOrthogonalize_Contiguous(BV V,Mat R)
210: {
212:   BV_CONTIGUOUS  *ctx = (BV_CONTIGUOUS*)V->data;
213:   PetscScalar    *r=NULL;

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

224: PetscErrorCode BVMatMult_Contiguous(BV V,Mat A,BV W)
225: {
227:   BV_CONTIGUOUS  *v = (BV_CONTIGUOUS*)V->data,*w = (BV_CONTIGUOUS*)W->data;
228:   PetscScalar    *pb,*pc;
229:   PetscInt       j,m;
230:   PetscBool      flg;

233:   MatHasOperation(A,MATOP_MAT_MULT,&flg);
234:   if (V->vmm && flg) {
235:     m = V->k-V->l;
236:     if (V->vmm==BV_MATMULT_MAT_SAVE) {
237:       BV_AllocateMatMult(V,A,m);
238:       MatDenseGetArray(V->B,&pb);
239:       PetscMemcpy(pb,v->array+(V->nc+V->l)*V->n,m*V->n*sizeof(PetscScalar));
240:       MatDenseRestoreArray(V->B,&pb);
241:     } else {  /* BV_MATMULT_MAT */
242:       MatCreateDense(PetscObjectComm((PetscObject)V),V->n,PETSC_DECIDE,V->N,m,v->array+(V->nc+V->l)*V->n,&V->B);
243:     }
244:     if (!V->C) {
245:       MatMatMultSymbolic(A,V->B,PETSC_DEFAULT,&V->C);
246:     }
247:     MatMatMultNumeric(A,V->B,V->C);
248:     MatDenseGetArray(V->C,&pc);
249:     PetscMemcpy(w->array+(W->nc+W->l)*W->n,pc,m*V->n*sizeof(PetscScalar));
250:     MatDenseRestoreArray(V->C,&pc);
251:     if (V->vmm==BV_MATMULT_MAT) {
252:       MatDestroy(&V->B);
253:       MatDestroy(&V->C);
254:     }
255:   } else {
256:     for (j=0;j<V->k-V->l;j++) {
257:       MatMult(A,v->V[V->nc+V->l+j],w->V[W->nc+W->l+j]);
258:     }
259:   }
260:   return(0);
261: }

265: PetscErrorCode BVCopy_Contiguous(BV V,BV W)
266: {
268:   BV_CONTIGUOUS  *v = (BV_CONTIGUOUS*)V->data,*w = (BV_CONTIGUOUS*)W->data;
269:   PetscScalar    *pvc,*pwc;

272:   pvc = v->array+(V->nc+V->l)*V->n;
273:   pwc = w->array+(W->nc+W->l)*W->n;
274:   PetscMemcpy(pwc,pvc,(V->k-V->l)*V->n*sizeof(PetscScalar));
275:   return(0);
276: }

280: PetscErrorCode BVResize_Contiguous(BV bv,PetscInt m,PetscBool copy)
281: {
283:   BV_CONTIGUOUS  *ctx = (BV_CONTIGUOUS*)bv->data;
284:   PetscInt       j,bs;
285:   PetscScalar    *newarray;
286:   Vec            *newV;
287:   char           str[50];

290:   VecGetBlockSize(bv->t,&bs);
291:   PetscMalloc1(m*bv->n,&newarray);
292:   PetscMemzero(newarray,m*bv->n*sizeof(PetscScalar));
293:   PetscMalloc1(m,&newV);
294:   for (j=0;j<m;j++) {
295:     if (ctx->mpi) {
296:       VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,bv->n,PETSC_DECIDE,newarray+j*bv->n,newV+j);
297:     } else {
298:       VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,bv->n,newarray+j*bv->n,newV+j);
299:     }
300:   }
301:   PetscLogObjectParents(bv,m,newV);
302:   if (((PetscObject)bv)->name) {
303:     for (j=0;j<m;j++) {
304:       PetscSNPrintf(str,50,"%s_%d",((PetscObject)bv)->name,(int)j);
305:       PetscObjectSetName((PetscObject)newV[j],str);
306:     }
307:   }
308:   if (copy) {
309:     PetscMemcpy(newarray,ctx->array,PetscMin(m,bv->m)*bv->n*sizeof(PetscScalar));
310:   }
311:   VecDestroyVecs(bv->m,&ctx->V);
312:   ctx->V = newV;
313:   PetscFree(ctx->array);
314:   ctx->array = newarray;
315:   return(0);
316: }

320: PetscErrorCode BVGetColumn_Contiguous(BV bv,PetscInt j,Vec *v)
321: {
322:   BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
323:   PetscInt      l;

326:   l = BVAvailableVec;
327:   bv->cv[l] = ctx->V[bv->nc+j];
328:   return(0);
329: }

333: PetscErrorCode BVGetArray_Contiguous(BV bv,PetscScalar **a)
334: {
335:   BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;

338:   *a = ctx->array;
339:   return(0);
340: }

344: PetscErrorCode BVGetArrayRead_Contiguous(BV bv,const PetscScalar **a)
345: {
346:   BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;

349:   *a = ctx->array;
350:   return(0);
351: }

355: PetscErrorCode BVDestroy_Contiguous(BV bv)
356: {
358:   BV_CONTIGUOUS  *ctx = (BV_CONTIGUOUS*)bv->data;

361:   VecDestroyVecs(bv->nc+bv->m,&ctx->V);
362:   PetscFree(ctx->array);
363:   PetscFree(bv->data);
364:   return(0);
365: }

369: PETSC_EXTERN PetscErrorCode BVCreate_Contiguous(BV bv)
370: {
372:   BV_CONTIGUOUS  *ctx;
373:   PetscInt       j,nloc,bs;
374:   PetscBool      seq;
375:   char           str[50];

378:   PetscNewLog(bv,&ctx);
379:   bv->data = (void*)ctx;

381:   PetscObjectTypeCompare((PetscObject)bv->t,VECMPI,&ctx->mpi);
382:   if (!ctx->mpi) {
383:     PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
384:     if (!seq) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot create a contiguous BV from a non-standard template vector");
385:   }

387:   VecGetLocalSize(bv->t,&nloc);
388:   VecGetBlockSize(bv->t,&bs);
389:   PetscMalloc1(bv->m*nloc,&ctx->array);
390:   PetscMemzero(ctx->array,bv->m*nloc*sizeof(PetscScalar));
391:   PetscMalloc1(bv->m,&ctx->V);
392:   for (j=0;j<bv->m;j++) {
393:     if (ctx->mpi) {
394:       VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,ctx->array+j*nloc,ctx->V+j);
395:     } else {
396:       VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,ctx->array+j*nloc,ctx->V+j);
397:     }
398:   }
399:   PetscLogObjectParents(bv,bv->m,ctx->V);
400:   if (((PetscObject)bv)->name) {
401:     for (j=0;j<bv->m;j++) {
402:       PetscSNPrintf(str,50,"%s_%d",((PetscObject)bv)->name,(int)j);
403:       PetscObjectSetName((PetscObject)ctx->V[j],str);
404:     }
405:   }

407:   bv->ops->mult             = BVMult_Contiguous;
408:   bv->ops->multvec          = BVMultVec_Contiguous;
409:   bv->ops->multinplace      = BVMultInPlace_Contiguous;
410:   bv->ops->multinplacetrans = BVMultInPlaceTranspose_Contiguous;
411:   bv->ops->dot              = BVDot_Contiguous;
412:   bv->ops->dotvec           = BVDotVec_Contiguous;
413:   bv->ops->dotvec_local     = BVDotVec_Local_Contiguous;
414:   bv->ops->scale            = BVScale_Contiguous;
415:   bv->ops->norm             = BVNorm_Contiguous;
416:   bv->ops->norm_local       = BVNorm_Local_Contiguous;
417:   /*bv->ops->orthogonalize    = BVOrthogonalize_Contiguous;*/
418:   bv->ops->matmult          = BVMatMult_Contiguous;
419:   bv->ops->copy             = BVCopy_Contiguous;
420:   bv->ops->resize           = BVResize_Contiguous;
421:   bv->ops->getcolumn        = BVGetColumn_Contiguous;
422:   bv->ops->getarray         = BVGetArray_Contiguous;
423:   bv->ops->getarrayread     = BVGetArrayRead_Contiguous;
424:   bv->ops->destroy          = BVDestroy_Contiguous;
425:   return(0);
426: }