Actual source code: contig.c
slepc-3.13.3 2020-06-14
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: MatProductClear(Wmat);
201: BVRestoreMat(V,&Vmat);
202: BVRestoreMat(W,&Wmat);
203: } else {
204: for (j=0;j<V->k-V->l;j++) {
205: MatMult(A,v->V[V->nc+V->l+j],w->V[W->nc+W->l+j]);
206: }
207: }
208: return(0);
209: }
211: PetscErrorCode BVCopy_Contiguous(BV V,BV W)
212: {
214: BV_CONTIGUOUS *v = (BV_CONTIGUOUS*)V->data,*w = (BV_CONTIGUOUS*)W->data;
215: PetscScalar *pvc,*pwc;
218: pvc = v->array+(V->nc+V->l)*V->n;
219: pwc = w->array+(W->nc+W->l)*W->n;
220: PetscArraycpy(pwc,pvc,(V->k-V->l)*V->n);
221: return(0);
222: }
224: PetscErrorCode BVCopyColumn_Contiguous(BV V,PetscInt j,PetscInt i)
225: {
227: BV_CONTIGUOUS *v = (BV_CONTIGUOUS*)V->data;
230: PetscArraycpy(v->array+(V->nc+i)*V->n,v->array+(V->nc+j)*V->n,V->n);
231: return(0);
232: }
234: PetscErrorCode BVResize_Contiguous(BV bv,PetscInt m,PetscBool copy)
235: {
237: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
238: PetscInt j,bs,lsplit;
239: PetscScalar *newarray;
240: Vec *newV;
241: char str[50];
242: BV parent;
245: if (bv->issplit==2) {
246: parent = bv->splitparent;
247: lsplit = parent->lsplit;
248: ctx->V = ((BV_CONTIGUOUS*)parent->data)->V+lsplit;
249: ctx->array = ((BV_CONTIGUOUS*)parent->data)->array+lsplit*bv->n;
250: } else if (!bv->issplit) {
251: VecGetBlockSize(bv->t,&bs);
252: PetscMalloc1(m*bv->n,&newarray);
253: PetscArrayzero(newarray,m*bv->n);
254: PetscMalloc1(m,&newV);
255: for (j=0;j<m;j++) {
256: if (ctx->mpi) {
257: VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,bv->n,PETSC_DECIDE,newarray+j*bv->n,newV+j);
258: } else {
259: VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,bv->n,newarray+j*bv->n,newV+j);
260: }
261: }
262: PetscLogObjectParents(bv,m,newV);
263: if (((PetscObject)bv)->name) {
264: for (j=0;j<m;j++) {
265: PetscSNPrintf(str,50,"%s_%d",((PetscObject)bv)->name,(int)j);
266: PetscObjectSetName((PetscObject)newV[j],str);
267: }
268: }
269: if (copy) {
270: PetscArraycpy(newarray,ctx->array,PetscMin(m,bv->m)*bv->n);
271: }
272: VecDestroyVecs(bv->m,&ctx->V);
273: ctx->V = newV;
274: PetscFree(ctx->array);
275: ctx->array = newarray;
276: }
277: return(0);
278: }
280: PetscErrorCode BVGetColumn_Contiguous(BV bv,PetscInt j,Vec *v)
281: {
282: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
283: PetscInt l;
286: l = BVAvailableVec;
287: bv->cv[l] = ctx->V[bv->nc+j];
288: return(0);
289: }
291: PetscErrorCode BVGetArray_Contiguous(BV bv,PetscScalar **a)
292: {
293: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
296: *a = ctx->array;
297: return(0);
298: }
300: PetscErrorCode BVGetArrayRead_Contiguous(BV bv,const PetscScalar **a)
301: {
302: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
305: *a = ctx->array;
306: return(0);
307: }
309: PetscErrorCode BVDestroy_Contiguous(BV bv)
310: {
312: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
315: if (!bv->issplit) {
316: VecDestroyVecs(bv->nc+bv->m,&ctx->V);
317: PetscFree(ctx->array);
318: }
319: PetscFree(bv->data);
320: return(0);
321: }
323: SLEPC_EXTERN PetscErrorCode BVCreate_Contiguous(BV bv)
324: {
326: BV_CONTIGUOUS *ctx;
327: PetscInt j,nloc,bs,lsplit;
328: PetscBool seq;
329: PetscScalar *aa;
330: char str[50];
331: PetscScalar *array;
332: BV parent;
333: Vec *Vpar;
336: PetscNewLog(bv,&ctx);
337: bv->data = (void*)ctx;
339: PetscObjectTypeCompare((PetscObject)bv->t,VECMPI,&ctx->mpi);
340: if (!ctx->mpi) {
341: PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
342: if (!seq) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Cannot create a contiguous BV from a non-standard template vector");
343: }
345: VecGetLocalSize(bv->t,&nloc);
346: VecGetBlockSize(bv->t,&bs);
348: if (bv->issplit) {
349: /* split BV: share memory and Vecs of the parent BV */
350: parent = bv->splitparent;
351: lsplit = parent->lsplit;
352: Vpar = ((BV_CONTIGUOUS*)parent->data)->V;
353: ctx->V = (bv->issplit==1)? Vpar: Vpar+lsplit;
354: array = ((BV_CONTIGUOUS*)parent->data)->array;
355: ctx->array = (bv->issplit==1)? array: array+lsplit*nloc;
356: } else {
357: /* regular BV: allocate memory and Vecs for the BV entries */
358: PetscCalloc1(bv->m*nloc,&ctx->array);
359: PetscMalloc1(bv->m,&ctx->V);
360: for (j=0;j<bv->m;j++) {
361: if (ctx->mpi) {
362: VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,ctx->array+j*nloc,ctx->V+j);
363: } else {
364: VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,ctx->array+j*nloc,ctx->V+j);
365: }
366: }
367: PetscLogObjectParents(bv,bv->m,ctx->V);
368: }
369: if (((PetscObject)bv)->name) {
370: for (j=0;j<bv->m;j++) {
371: PetscSNPrintf(str,50,"%s_%d",((PetscObject)bv)->name,(int)j);
372: PetscObjectSetName((PetscObject)ctx->V[j],str);
373: }
374: }
376: if (bv->Acreate) {
377: MatDenseGetArray(bv->Acreate,&aa);
378: PetscArraycpy(ctx->array,aa,bv->m*nloc);
379: MatDenseRestoreArray(bv->Acreate,&aa);
380: MatDestroy(&bv->Acreate);
381: }
383: bv->ops->mult = BVMult_Contiguous;
384: bv->ops->multvec = BVMultVec_Contiguous;
385: bv->ops->multinplace = BVMultInPlace_Contiguous;
386: bv->ops->multinplacetrans = BVMultInPlaceTranspose_Contiguous;
387: bv->ops->dot = BVDot_Contiguous;
388: bv->ops->dotvec = BVDotVec_Contiguous;
389: bv->ops->dotvec_local = BVDotVec_Local_Contiguous;
390: bv->ops->scale = BVScale_Contiguous;
391: bv->ops->norm = BVNorm_Contiguous;
392: bv->ops->norm_local = BVNorm_Local_Contiguous;
393: bv->ops->matmult = BVMatMult_Contiguous;
394: bv->ops->copy = BVCopy_Contiguous;
395: bv->ops->copycolumn = BVCopyColumn_Contiguous;
396: bv->ops->resize = BVResize_Contiguous;
397: bv->ops->getcolumn = BVGetColumn_Contiguous;
398: bv->ops->getarray = BVGetArray_Contiguous;
399: bv->ops->getarrayread = BVGetArrayRead_Contiguous;
400: bv->ops->destroy = BVDestroy_Contiguous;
401: return(0);
402: }