Actual source code: contig.c
slepc-3.6.1 2015-09-03
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-2015, 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: MatGetSize(Q,&ldq,NULL);
43: MatDenseGetArray(Q,&q);
44: 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);
45: MatDenseRestoreArray(Q,&q);
46: return(0);
47: }
51: PetscErrorCode BVMultVec_Contiguous(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
52: {
54: BV_CONTIGUOUS *x = (BV_CONTIGUOUS*)X->data;
55: PetscScalar *py;
58: VecGetArray(y,&py);
59: BVMultVec_BLAS_Private(X,X->n,X->k-X->l,alpha,x->array+(X->nc+X->l)*X->n,q,beta,py);
60: VecRestoreArray(y,&py);
61: return(0);
62: }
66: PetscErrorCode BVMultInPlace_Contiguous(BV V,Mat Q,PetscInt s,PetscInt e)
67: {
69: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)V->data;
70: PetscScalar *q;
71: PetscInt ldq;
74: MatGetSize(Q,&ldq,NULL);
75: MatDenseGetArray(Q,&q);
76: 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);
77: MatDenseRestoreArray(Q,&q);
78: return(0);
79: }
83: PetscErrorCode BVMultInPlaceTranspose_Contiguous(BV V,Mat Q,PetscInt s,PetscInt e)
84: {
86: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)V->data;
87: PetscScalar *q;
88: PetscInt ldq;
91: MatGetSize(Q,&ldq,NULL);
92: MatDenseGetArray(Q,&q);
93: 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);
94: MatDenseRestoreArray(Q,&q);
95: return(0);
96: }
100: PetscErrorCode BVAXPY_Contiguous(BV Y,PetscScalar alpha,BV X)
101: {
103: BV_CONTIGUOUS *x = (BV_CONTIGUOUS*)X->data,*y = (BV_CONTIGUOUS*)Y->data;
106: BVAXPY_BLAS_Private(Y,Y->n,Y->k-Y->l,alpha,x->array+(X->nc+X->l)*X->n,y->array+(Y->nc+Y->l)*Y->n);
107: return(0);
108: }
112: PetscErrorCode BVDot_Contiguous(BV X,BV Y,Mat M)
113: {
115: BV_CONTIGUOUS *x = (BV_CONTIGUOUS*)X->data,*y = (BV_CONTIGUOUS*)Y->data;
116: PetscScalar *m;
117: PetscInt ldm;
120: MatGetSize(M,&ldm,NULL);
121: MatDenseGetArray(M,&m);
122: 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);
123: MatDenseRestoreArray(M,&m);
124: return(0);
125: }
129: PetscErrorCode BVDotVec_Contiguous(BV X,Vec y,PetscScalar *m)
130: {
131: PetscErrorCode ierr;
132: BV_CONTIGUOUS *x = (BV_CONTIGUOUS*)X->data;
133: const PetscScalar *py;
134: Vec z = y;
137: if (X->matrix) {
138: BV_IPMatMult(X,y);
139: z = X->Bx;
140: }
141: VecGetArrayRead(z,&py);
142: BVDotVec_BLAS_Private(X,X->n,X->k-X->l,x->array+(X->nc+X->l)*X->n,py,m,x->mpi);
143: VecRestoreArrayRead(z,&py);
144: return(0);
145: }
149: PetscErrorCode BVDotVec_Local_Contiguous(BV X,Vec y,PetscScalar *m)
150: {
152: BV_CONTIGUOUS *x = (BV_CONTIGUOUS*)X->data;
153: PetscScalar *py;
154: Vec z = y;
157: if (X->matrix) {
158: BV_IPMatMult(X,y);
159: z = X->Bx;
160: }
161: VecGetArray(z,&py);
162: BVDotVec_BLAS_Private(X,X->n,X->k-X->l,x->array+(X->nc+X->l)*X->n,py,m,PETSC_FALSE);
163: VecRestoreArray(z,&py);
164: return(0);
165: }
169: PetscErrorCode BVScale_Contiguous(BV bv,PetscInt j,PetscScalar alpha)
170: {
172: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
175: if (j<0) {
176: BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->n,ctx->array+(bv->nc+bv->l)*bv->n,alpha);
177: } else {
178: BVScale_BLAS_Private(bv,bv->n,ctx->array+(bv->nc+j)*bv->n,alpha);
179: }
180: return(0);
181: }
185: PetscErrorCode BVNorm_Contiguous(BV bv,PetscInt j,NormType type,PetscReal *val)
186: {
188: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
191: if (j<0) {
192: BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,ctx->array+(bv->nc+bv->l)*bv->n,type,val,ctx->mpi);
193: } else {
194: BVNorm_LAPACK_Private(bv,bv->n,1,ctx->array+(bv->nc+j)*bv->n,type,val,ctx->mpi);
195: }
196: return(0);
197: }
201: PetscErrorCode BVNorm_Local_Contiguous(BV bv,PetscInt j,NormType type,PetscReal *val)
202: {
204: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
207: if (j<0) {
208: BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,ctx->array+(bv->nc+bv->l)*bv->n,type,val,PETSC_FALSE);
209: } else {
210: BVNorm_LAPACK_Private(bv,bv->n,1,ctx->array+(bv->nc+j)*bv->n,type,val,PETSC_FALSE);
211: }
212: return(0);
213: }
217: PetscErrorCode BVOrthogonalize_Contiguous(BV V,Mat R)
218: {
220: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)V->data;
221: PetscScalar *r=NULL;
224: if (R) { MatDenseGetArray(R,&r); }
225: BVOrthogonalize_LAPACK_Private(V,V->n,V->k,ctx->array+V->nc*V->n,r,ctx->mpi);
226: if (R) { MatDenseRestoreArray(R,&r); }
227: return(0);
228: }
232: PetscErrorCode BVMatMult_Contiguous(BV V,Mat A,BV W)
233: {
235: BV_CONTIGUOUS *v = (BV_CONTIGUOUS*)V->data,*w = (BV_CONTIGUOUS*)W->data;
236: PetscScalar *pb,*pc;
237: PetscInt j,m;
238: PetscBool flg;
241: MatHasOperation(A,MATOP_MAT_MULT,&flg);
242: if (V->vmm && flg) {
243: m = V->k-V->l;
244: if (V->vmm==BV_MATMULT_MAT_SAVE) {
245: BV_AllocateMatMult(V,A,m);
246: MatDenseGetArray(V->B,&pb);
247: PetscMemcpy(pb,v->array+(V->nc+V->l)*V->n,m*V->n*sizeof(PetscScalar));
248: MatDenseRestoreArray(V->B,&pb);
249: } else { /* BV_MATMULT_MAT */
250: MatCreateDense(PetscObjectComm((PetscObject)V),V->n,PETSC_DECIDE,V->N,m,v->array+(V->nc+V->l)*V->n,&V->B);
251: }
252: if (!V->C) {
253: MatMatMultSymbolic(A,V->B,PETSC_DEFAULT,&V->C);
254: }
255: MatMatMultNumeric(A,V->B,V->C);
256: MatDenseGetArray(V->C,&pc);
257: PetscMemcpy(w->array+(W->nc+W->l)*W->n,pc,m*V->n*sizeof(PetscScalar));
258: MatDenseRestoreArray(V->C,&pc);
259: if (V->vmm==BV_MATMULT_MAT) {
260: MatDestroy(&V->B);
261: MatDestroy(&V->C);
262: }
263: } else {
264: for (j=0;j<V->k-V->l;j++) {
265: MatMult(A,v->V[V->nc+V->l+j],w->V[W->nc+W->l+j]);
266: }
267: }
268: return(0);
269: }
273: PetscErrorCode BVCopy_Contiguous(BV V,BV W)
274: {
276: BV_CONTIGUOUS *v = (BV_CONTIGUOUS*)V->data,*w = (BV_CONTIGUOUS*)W->data;
277: PetscScalar *pvc,*pwc;
280: pvc = v->array+(V->nc+V->l)*V->n;
281: pwc = w->array+(W->nc+W->l)*W->n;
282: PetscMemcpy(pwc,pvc,(V->k-V->l)*V->n*sizeof(PetscScalar));
283: return(0);
284: }
288: PetscErrorCode BVResize_Contiguous(BV bv,PetscInt m,PetscBool copy)
289: {
291: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
292: PetscInt j,bs;
293: PetscScalar *newarray;
294: Vec *newV;
295: char str[50];
298: VecGetBlockSize(bv->t,&bs);
299: PetscMalloc1(m*bv->n,&newarray);
300: PetscMemzero(newarray,m*bv->n*sizeof(PetscScalar));
301: PetscMalloc1(m,&newV);
302: for (j=0;j<m;j++) {
303: if (ctx->mpi) {
304: VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,bv->n,PETSC_DECIDE,newarray+j*bv->n,newV+j);
305: } else {
306: VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,bv->n,newarray+j*bv->n,newV+j);
307: }
308: }
309: PetscLogObjectParents(bv,m,newV);
310: if (((PetscObject)bv)->name) {
311: for (j=0;j<m;j++) {
312: PetscSNPrintf(str,50,"%s_%D",((PetscObject)bv)->name,j);
313: PetscObjectSetName((PetscObject)newV[j],str);
314: }
315: }
316: if (copy) {
317: PetscMemcpy(newarray,ctx->array,PetscMin(m,bv->m)*bv->n*sizeof(PetscScalar));
318: }
319: VecDestroyVecs(bv->m,&ctx->V);
320: ctx->V = newV;
321: PetscFree(ctx->array);
322: ctx->array = newarray;
323: return(0);
324: }
328: PetscErrorCode BVGetColumn_Contiguous(BV bv,PetscInt j,Vec *v)
329: {
330: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
331: PetscInt l;
334: l = BVAvailableVec;
335: bv->cv[l] = ctx->V[bv->nc+j];
336: return(0);
337: }
341: PetscErrorCode BVGetArray_Contiguous(BV bv,PetscScalar **a)
342: {
343: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
346: *a = ctx->array;
347: return(0);
348: }
352: PetscErrorCode BVDestroy_Contiguous(BV bv)
353: {
355: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
358: VecDestroyVecs(bv->nc+bv->m,&ctx->V);
359: PetscFree(ctx->array);
360: PetscFree(bv->data);
361: return(0);
362: }
366: PETSC_EXTERN PetscErrorCode BVCreate_Contiguous(BV bv)
367: {
369: BV_CONTIGUOUS *ctx;
370: PetscInt j,nloc,bs;
371: PetscBool seq;
372: char str[50];
375: PetscNewLog(bv,&ctx);
376: bv->data = (void*)ctx;
378: PetscObjectTypeCompare((PetscObject)bv->t,VECMPI,&ctx->mpi);
379: if (!ctx->mpi) {
380: PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
381: if (!seq) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot create a contiguous BV from a non-standard template vector");
382: }
384: VecGetLocalSize(bv->t,&nloc);
385: VecGetBlockSize(bv->t,&bs);
386: PetscMalloc1(bv->m*nloc,&ctx->array);
387: PetscMemzero(ctx->array,bv->m*nloc*sizeof(PetscScalar));
388: PetscMalloc1(bv->m,&ctx->V);
389: for (j=0;j<bv->m;j++) {
390: if (ctx->mpi) {
391: VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,ctx->array+j*nloc,ctx->V+j);
392: } else {
393: VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,ctx->array+j*nloc,ctx->V+j);
394: }
395: }
396: PetscLogObjectParents(bv,bv->m,ctx->V);
397: if (((PetscObject)bv)->name) {
398: for (j=0;j<bv->m;j++) {
399: PetscSNPrintf(str,50,"%s_%D",((PetscObject)bv)->name,j);
400: PetscObjectSetName((PetscObject)ctx->V[j],str);
401: }
402: }
404: bv->ops->mult = BVMult_Contiguous;
405: bv->ops->multvec = BVMultVec_Contiguous;
406: bv->ops->multinplace = BVMultInPlace_Contiguous;
407: bv->ops->multinplacetrans = BVMultInPlaceTranspose_Contiguous;
408: bv->ops->axpy = BVAXPY_Contiguous;
409: bv->ops->dot = BVDot_Contiguous;
410: bv->ops->dotvec = BVDotVec_Contiguous;
411: bv->ops->dotvec_local = BVDotVec_Local_Contiguous;
412: bv->ops->scale = BVScale_Contiguous;
413: bv->ops->norm = BVNorm_Contiguous;
414: bv->ops->norm_local = BVNorm_Local_Contiguous;
415: /*bv->ops->orthogonalize = BVOrthogonalize_Contiguous;*/
416: bv->ops->matmult = BVMatMult_Contiguous;
417: bv->ops->copy = BVCopy_Contiguous;
418: bv->ops->resize = BVResize_Contiguous;
419: bv->ops->getcolumn = BVGetColumn_Contiguous;
420: bv->ops->getarray = BVGetArray_Contiguous;
421: bv->ops->destroy = BVDestroy_Contiguous;
422: return(0);
423: }