Actual source code: contig.c

slepc-3.13.1 2020-04-12
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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 an array of Vecs sharing a contiguous array for elements
 12: */

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

 16: typedef struct {
 17:   Vec         *V;
 18:   PetscScalar *array;
 19:   PetscBool   mpi;
 20: } BV_CONTIGUOUS;

 22: PetscErrorCode BVMult_Contiguous(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
 23: {
 25:   BV_CONTIGUOUS  *y = (BV_CONTIGUOUS*)Y->data,*x = (BV_CONTIGUOUS*)X->data;
 26:   PetscScalar    *q;
 27:   PetscInt       ldq;

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

 41: PetscErrorCode BVMultVec_Contiguous(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
 42: {
 44:   BV_CONTIGUOUS  *x = (BV_CONTIGUOUS*)X->data;
 45:   PetscScalar    *py,*qq=q;

 48:   VecGetArray(y,&py);
 49:   if (!q) { VecGetArray(X->buffer,&qq); }
 50:   BVMultVec_BLAS_Private(X,X->n,X->k-X->l,alpha,x->array+(X->nc+X->l)*X->n,qq,beta,py);
 51:   if (!q) { VecRestoreArray(X->buffer,&qq); }
 52:   VecRestoreArray(y,&py);
 53:   return(0);
 54: }

 56: PetscErrorCode BVMultInPlace_Contiguous(BV V,Mat Q,PetscInt s,PetscInt e)
 57: {
 59:   BV_CONTIGUOUS  *ctx = (BV_CONTIGUOUS*)V->data;
 60:   PetscScalar    *q;
 61:   PetscInt       ldq;

 64:   MatGetSize(Q,&ldq,NULL);
 65:   MatDenseGetArray(Q,&q);
 66:   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);
 67:   MatDenseRestoreArray(Q,&q);
 68:   return(0);
 69: }

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

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

 86: PetscErrorCode BVDot_Contiguous(BV X,BV Y,Mat M)
 87: {
 89:   BV_CONTIGUOUS  *x = (BV_CONTIGUOUS*)X->data,*y = (BV_CONTIGUOUS*)Y->data;
 90:   PetscScalar    *m;
 91:   PetscInt       ldm;

 94:   MatGetSize(M,&ldm,NULL);
 95:   MatDenseGetArray(M,&m);
 96:   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);
 97:   MatDenseRestoreArray(M,&m);
 98:   return(0);
 99: }

101: PetscErrorCode BVDotVec_Contiguous(BV X,Vec y,PetscScalar *q)
102: {
103:   PetscErrorCode    ierr;
104:   BV_CONTIGUOUS     *x = (BV_CONTIGUOUS*)X->data;
105:   const PetscScalar *py;
106:   PetscScalar       *qq=q;
107:   Vec               z = y;

110:   if (X->matrix) {
111:     BV_IPMatMult(X,y);
112:     z = X->Bx;
113:   }
114:   VecGetArrayRead(z,&py);
115:   if (!q) { VecGetArray(X->buffer,&qq); }
116:   BVDotVec_BLAS_Private(X,X->n,X->k-X->l,x->array+(X->nc+X->l)*X->n,py,qq,x->mpi);
117:   if (!q) { VecRestoreArray(X->buffer,&qq); }
118:   VecRestoreArrayRead(z,&py);
119:   return(0);
120: }

122: PetscErrorCode BVDotVec_Local_Contiguous(BV X,Vec y,PetscScalar *m)
123: {
125:   BV_CONTIGUOUS  *x = (BV_CONTIGUOUS*)X->data;
126:   PetscScalar    *py;
127:   Vec            z = y;

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

140: PetscErrorCode BVScale_Contiguous(BV bv,PetscInt j,PetscScalar alpha)
141: {
143:   BV_CONTIGUOUS  *ctx = (BV_CONTIGUOUS*)bv->data;

146:   if (j<0) {
147:     BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->n,ctx->array+(bv->nc+bv->l)*bv->n,alpha);
148:   } else {
149:     BVScale_BLAS_Private(bv,bv->n,ctx->array+(bv->nc+j)*bv->n,alpha);
150:   }
151:   return(0);
152: }

154: PetscErrorCode BVNorm_Contiguous(BV bv,PetscInt j,NormType type,PetscReal *val)
155: {
157:   BV_CONTIGUOUS  *ctx = (BV_CONTIGUOUS*)bv->data;

160:   if (j<0) {
161:     BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,ctx->array+(bv->nc+bv->l)*bv->n,type,val,ctx->mpi);
162:   } else {
163:     BVNorm_LAPACK_Private(bv,bv->n,1,ctx->array+(bv->nc+j)*bv->n,type,val,ctx->mpi);
164:   }
165:   return(0);
166: }

168: PetscErrorCode BVNorm_Local_Contiguous(BV bv,PetscInt j,NormType type,PetscReal *val)
169: {
171:   BV_CONTIGUOUS  *ctx = (BV_CONTIGUOUS*)bv->data;

174:   if (j<0) {
175:     BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,ctx->array+(bv->nc+bv->l)*bv->n,type,val,PETSC_FALSE);
176:   } else {
177:     BVNorm_LAPACK_Private(bv,bv->n,1,ctx->array+(bv->nc+j)*bv->n,type,val,PETSC_FALSE);
178:   }
179:   return(0);
180: }

182: PetscErrorCode BVMatMult_Contiguous(BV V,Mat A,BV W)
183: {
185:   BV_CONTIGUOUS  *v = (BV_CONTIGUOUS*)V->data,*w = (BV_CONTIGUOUS*)W->data;
186:   PetscInt       j;
187:   PetscBool      flg;
188:   Mat            Vmat,Wmat;

191:   MatHasOperation(A,MATOP_MAT_MULT,&flg);
192:   if (V->vmm && flg) {
193:     BVGetMat(V,&Vmat);
194:     BVGetMat(W,&Wmat);
195:     MatProductCreateWithMat(A,Vmat,NULL,Wmat);
196:     MatProductSetType(Wmat,MATPRODUCT_AB);
197:     MatProductSetFromOptions(Wmat);
198:     MatProductSymbolic(Wmat);
199:     MatProductNumeric(Wmat);
200:     BVRestoreMat(V,&Vmat);
201:     BVRestoreMat(W,&Wmat);
202:   } else {
203:     for (j=0;j<V->k-V->l;j++) {
204:       MatMult(A,v->V[V->nc+V->l+j],w->V[W->nc+W->l+j]);
205:     }
206:   }
207:   return(0);
208: }

210: PetscErrorCode BVCopy_Contiguous(BV V,BV W)
211: {
213:   BV_CONTIGUOUS  *v = (BV_CONTIGUOUS*)V->data,*w = (BV_CONTIGUOUS*)W->data;
214:   PetscScalar    *pvc,*pwc;

217:   pvc = v->array+(V->nc+V->l)*V->n;
218:   pwc = w->array+(W->nc+W->l)*W->n;
219:   PetscArraycpy(pwc,pvc,(V->k-V->l)*V->n);
220:   return(0);
221: }

223: PetscErrorCode BVCopyColumn_Contiguous(BV V,PetscInt j,PetscInt i)
224: {
226:   BV_CONTIGUOUS  *v = (BV_CONTIGUOUS*)V->data;

229:   PetscArraycpy(v->array+(V->nc+i)*V->n,v->array+(V->nc+j)*V->n,V->n);
230:   return(0);
231: }

233: PetscErrorCode BVResize_Contiguous(BV bv,PetscInt m,PetscBool copy)
234: {
236:   BV_CONTIGUOUS  *ctx = (BV_CONTIGUOUS*)bv->data;
237:   PetscInt       j,bs,lsplit;
238:   PetscScalar    *newarray;
239:   Vec            *newV;
240:   char           str[50];
241:   BV             parent;

244:   if (bv->issplit==2) {
245:     parent = bv->splitparent;
246:     lsplit = parent->lsplit;
247:     ctx->V = ((BV_CONTIGUOUS*)parent->data)->V+lsplit;
248:     ctx->array = ((BV_CONTIGUOUS*)parent->data)->array+lsplit*bv->n;
249:   } else if (!bv->issplit) {
250:     VecGetBlockSize(bv->t,&bs);
251:     PetscMalloc1(m*bv->n,&newarray);
252:     PetscArrayzero(newarray,m*bv->n);
253:     PetscMalloc1(m,&newV);
254:     for (j=0;j<m;j++) {
255:       if (ctx->mpi) {
256:         VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,bv->n,PETSC_DECIDE,newarray+j*bv->n,newV+j);
257:       } else {
258:         VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,bv->n,newarray+j*bv->n,newV+j);
259:       }
260:     }
261:     PetscLogObjectParents(bv,m,newV);
262:     if (((PetscObject)bv)->name) {
263:       for (j=0;j<m;j++) {
264:         PetscSNPrintf(str,50,"%s_%d",((PetscObject)bv)->name,(int)j);
265:         PetscObjectSetName((PetscObject)newV[j],str);
266:       }
267:     }
268:     if (copy) {
269:       PetscArraycpy(newarray,ctx->array,PetscMin(m,bv->m)*bv->n);
270:     }
271:     VecDestroyVecs(bv->m,&ctx->V);
272:     ctx->V = newV;
273:     PetscFree(ctx->array);
274:     ctx->array = newarray;
275:   }
276:   return(0);
277: }

279: PetscErrorCode BVGetColumn_Contiguous(BV bv,PetscInt j,Vec *v)
280: {
281:   BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
282:   PetscInt      l;

285:   l = BVAvailableVec;
286:   bv->cv[l] = ctx->V[bv->nc+j];
287:   return(0);
288: }

290: PetscErrorCode BVGetArray_Contiguous(BV bv,PetscScalar **a)
291: {
292:   BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;

295:   *a = ctx->array;
296:   return(0);
297: }

299: PetscErrorCode BVGetArrayRead_Contiguous(BV bv,const PetscScalar **a)
300: {
301:   BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;

304:   *a = ctx->array;
305:   return(0);
306: }

308: PetscErrorCode BVDestroy_Contiguous(BV bv)
309: {
311:   BV_CONTIGUOUS  *ctx = (BV_CONTIGUOUS*)bv->data;

314:   if (!bv->issplit) {
315:     VecDestroyVecs(bv->nc+bv->m,&ctx->V);
316:     PetscFree(ctx->array);
317:   }
318:   PetscFree(bv->data);
319:   return(0);
320: }

322: SLEPC_EXTERN PetscErrorCode BVCreate_Contiguous(BV bv)
323: {
325:   BV_CONTIGUOUS  *ctx;
326:   PetscInt       j,nloc,bs,lsplit;
327:   PetscBool      seq;
328:   PetscScalar    *aa;
329:   char           str[50];
330:   PetscScalar    *array;
331:   BV             parent;
332:   Vec            *Vpar;

335:   PetscNewLog(bv,&ctx);
336:   bv->data = (void*)ctx;

338:   PetscObjectTypeCompare((PetscObject)bv->t,VECMPI,&ctx->mpi);
339:   if (!ctx->mpi) {
340:     PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
341:     if (!seq) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Cannot create a contiguous BV from a non-standard template vector");
342:   }

344:   VecGetLocalSize(bv->t,&nloc);
345:   VecGetBlockSize(bv->t,&bs);

347:   if (bv->issplit) {
348:     /* split BV: share memory and Vecs of the parent BV */
349:     parent = bv->splitparent;
350:     lsplit = parent->lsplit;
351:     Vpar   = ((BV_CONTIGUOUS*)parent->data)->V;
352:     ctx->V = (bv->issplit==1)? Vpar: Vpar+lsplit;
353:     array  = ((BV_CONTIGUOUS*)parent->data)->array;
354:     ctx->array = (bv->issplit==1)? array: array+lsplit*nloc;
355:   } else {
356:     /* regular BV: allocate memory and Vecs for the BV entries */
357:     PetscCalloc1(bv->m*nloc,&ctx->array);
358:     PetscMalloc1(bv->m,&ctx->V);
359:     for (j=0;j<bv->m;j++) {
360:       if (ctx->mpi) {
361:         VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,ctx->array+j*nloc,ctx->V+j);
362:       } else {
363:         VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,ctx->array+j*nloc,ctx->V+j);
364:       }
365:     }
366:     PetscLogObjectParents(bv,bv->m,ctx->V);
367:   }
368:   if (((PetscObject)bv)->name) {
369:     for (j=0;j<bv->m;j++) {
370:       PetscSNPrintf(str,50,"%s_%d",((PetscObject)bv)->name,(int)j);
371:       PetscObjectSetName((PetscObject)ctx->V[j],str);
372:     }
373:   }

375:   if (bv->Acreate) {
376:     MatDenseGetArray(bv->Acreate,&aa);
377:     PetscArraycpy(ctx->array,aa,bv->m*nloc);
378:     MatDenseRestoreArray(bv->Acreate,&aa);
379:     MatDestroy(&bv->Acreate);
380:   }

382:   bv->ops->mult             = BVMult_Contiguous;
383:   bv->ops->multvec          = BVMultVec_Contiguous;
384:   bv->ops->multinplace      = BVMultInPlace_Contiguous;
385:   bv->ops->multinplacetrans = BVMultInPlaceTranspose_Contiguous;
386:   bv->ops->dot              = BVDot_Contiguous;
387:   bv->ops->dotvec           = BVDotVec_Contiguous;
388:   bv->ops->dotvec_local     = BVDotVec_Local_Contiguous;
389:   bv->ops->scale            = BVScale_Contiguous;
390:   bv->ops->norm             = BVNorm_Contiguous;
391:   bv->ops->norm_local       = BVNorm_Local_Contiguous;
392:   bv->ops->matmult          = BVMatMult_Contiguous;
393:   bv->ops->copy             = BVCopy_Contiguous;
394:   bv->ops->copycolumn       = BVCopyColumn_Contiguous;
395:   bv->ops->resize           = BVResize_Contiguous;
396:   bv->ops->getcolumn        = BVGetColumn_Contiguous;
397:   bv->ops->getarray         = BVGetArray_Contiguous;
398:   bv->ops->getarrayread     = BVGetArrayRead_Contiguous;
399:   bv->ops->destroy          = BVDestroy_Contiguous;
400:   return(0);
401: }