Actual source code: contig.c
slepc-3.5.2 2014-10-10
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-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: 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: {
132: BV_CONTIGUOUS *x = (BV_CONTIGUOUS*)X->data;
133: PetscScalar *py;
134: Vec z = y;
137: if (X->matrix) {
138: BV_IPMatMult(X,y);
139: z = X->Bx;
140: }
141: VecGetArray(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: VecRestoreArray(z,&py);
144: return(0);
145: }
149: PetscErrorCode BVScale_Contiguous(BV bv,PetscInt j,PetscScalar alpha)
150: {
152: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
155: if (j<0) {
156: BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->n,ctx->array+(bv->nc+bv->l)*bv->n,alpha);
157: } else {
158: BVScale_BLAS_Private(bv,bv->n,ctx->array+(bv->nc+j)*bv->n,alpha);
159: }
160: return(0);
161: }
165: PetscErrorCode BVNorm_Contiguous(BV bv,PetscInt j,NormType type,PetscReal *val)
166: {
168: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
171: if (j<0) {
172: BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,ctx->array+(bv->nc+bv->l)*bv->n,type,val,ctx->mpi);
173: } else {
174: BVNorm_LAPACK_Private(bv,bv->n,1,ctx->array+(bv->nc+j)*bv->n,type,val,ctx->mpi);
175: }
176: return(0);
177: }
181: PetscErrorCode BVOrthogonalize_Contiguous(BV V,Mat R)
182: {
184: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)V->data;
185: PetscScalar *r=NULL;
188: if (R) { MatDenseGetArray(R,&r); }
189: BVOrthogonalize_LAPACK_Private(V,V->n,V->k,ctx->array+V->nc*V->n,r,ctx->mpi);
190: if (R) { MatDenseRestoreArray(R,&r); }
191: return(0);
192: }
196: PetscErrorCode BVMatMult_Contiguous(BV V,Mat A,BV W)
197: {
199: BV_CONTIGUOUS *v = (BV_CONTIGUOUS*)V->data,*w = (BV_CONTIGUOUS*)W->data;
200: PetscInt j;
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: return(0);
207: }
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: PetscMemcpy(pwc,pvc,(V->k-V->l)*V->n*sizeof(PetscScalar));
221: return(0);
222: }
226: PetscErrorCode BVResize_Contiguous(BV bv,PetscInt m,PetscBool copy)
227: {
229: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
230: PetscInt j,bs;
231: PetscScalar *newarray;
232: Vec *newV;
233: char str[50];
236: VecGetBlockSize(bv->t,&bs);
237: PetscMalloc1(m*bv->n,&newarray);
238: PetscMemzero(newarray,m*bv->n*sizeof(PetscScalar));
239: PetscMalloc1(m,&newV);
240: for (j=0;j<m;j++) {
241: if (ctx->mpi) {
242: VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,bv->n,PETSC_DECIDE,newarray+j*bv->n,newV+j);
243: } else {
244: VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,bv->n,newarray+j*bv->n,newV+j);
245: }
246: }
247: PetscLogObjectParents(bv,m,newV);
248: if (((PetscObject)bv)->name) {
249: for (j=0;j<m;j++) {
250: PetscSNPrintf(str,50,"%s_%D",((PetscObject)bv)->name,j);
251: PetscObjectSetName((PetscObject)newV[j],str);
252: }
253: }
254: if (copy) {
255: PetscMemcpy(newarray,ctx->array,PetscMin(m,bv->m)*bv->n*sizeof(PetscScalar));
256: }
257: VecDestroyVecs(bv->m,&ctx->V);
258: ctx->V = newV;
259: PetscFree(ctx->array);
260: ctx->array = newarray;
261: return(0);
262: }
266: PetscErrorCode BVGetColumn_Contiguous(BV bv,PetscInt j,Vec *v)
267: {
268: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
269: PetscInt l;
272: l = BVAvailableVec;
273: bv->cv[l] = ctx->V[bv->nc+j];
274: return(0);
275: }
279: PetscErrorCode BVGetArray_Contiguous(BV bv,PetscScalar **a)
280: {
281: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
284: *a = ctx->array;
285: return(0);
286: }
290: PetscErrorCode BVDestroy_Contiguous(BV bv)
291: {
293: BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
296: VecDestroyVecs(bv->nc+bv->m,&ctx->V);
297: PetscFree(ctx->array);
298: PetscFree(bv->data);
299: return(0);
300: }
304: PETSC_EXTERN PetscErrorCode BVCreate_Contiguous(BV bv)
305: {
307: BV_CONTIGUOUS *ctx;
308: PetscInt j,nloc,bs;
309: PetscBool seq;
310: char str[50];
313: PetscNewLog(bv,&ctx);
314: bv->data = (void*)ctx;
316: PetscObjectTypeCompare((PetscObject)bv->t,VECMPI,&ctx->mpi);
317: if (!ctx->mpi) {
318: PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
319: if (!seq) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot create a contiguous BV from a non-standard template vector");
320: }
322: VecGetLocalSize(bv->t,&nloc);
323: VecGetBlockSize(bv->t,&bs);
324: PetscMalloc1(bv->m*nloc,&ctx->array);
325: PetscMemzero(ctx->array,bv->m*nloc*sizeof(PetscScalar));
326: PetscMalloc1(bv->m,&ctx->V);
327: for (j=0;j<bv->m;j++) {
328: if (ctx->mpi) {
329: VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,ctx->array+j*nloc,ctx->V+j);
330: } else {
331: VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,ctx->array+j*nloc,ctx->V+j);
332: }
333: }
334: PetscLogObjectParents(bv,bv->m,ctx->V);
335: if (((PetscObject)bv)->name) {
336: for (j=0;j<bv->m;j++) {
337: PetscSNPrintf(str,50,"%s_%D",((PetscObject)bv)->name,j);
338: PetscObjectSetName((PetscObject)ctx->V[j],str);
339: }
340: }
342: bv->ops->mult = BVMult_Contiguous;
343: bv->ops->multvec = BVMultVec_Contiguous;
344: bv->ops->multinplace = BVMultInPlace_Contiguous;
345: bv->ops->multinplacetrans = BVMultInPlaceTranspose_Contiguous;
346: bv->ops->axpy = BVAXPY_Contiguous;
347: bv->ops->dot = BVDot_Contiguous;
348: bv->ops->dotvec = BVDotVec_Contiguous;
349: bv->ops->scale = BVScale_Contiguous;
350: bv->ops->norm = BVNorm_Contiguous;
351: /*bv->ops->orthogonalize = BVOrthogonalize_Contiguous;*/
352: bv->ops->matmult = BVMatMult_Contiguous;
353: bv->ops->copy = BVCopy_Contiguous;
354: bv->ops->resize = BVResize_Contiguous;
355: bv->ops->getcolumn = BVGetColumn_Contiguous;
356: bv->ops->getarray = BVGetArray_Contiguous;
357: bv->ops->destroy = BVDestroy_Contiguous;
358: return(0);
359: }