Actual source code: bvmat.c
slepc-3.5.2 2014-10-10
1: /*
2: BV implemented with a dense Mat
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: Mat A;
28: PetscBool mpi;
29: } BV_MAT;
33: PetscErrorCode BVMult_Mat(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
34: {
36: BV_MAT *y = (BV_MAT*)Y->data,*x = (BV_MAT*)X->data;
37: PetscScalar *px,*py,*q;
38: PetscInt ldq;
41: MatGetSize(Q,&ldq,NULL);
42: MatDenseGetArray(x->A,&px);
43: MatDenseGetArray(y->A,&py);
44: MatDenseGetArray(Q,&q);
45: BVMult_BLAS_Private(Y,Y->n,Y->k-Y->l,X->k-X->l,ldq,alpha,px+(X->nc+X->l)*X->n,q+Y->l*ldq+X->l,beta,py+(Y->nc+Y->l)*Y->n);
46: MatDenseRestoreArray(Q,&q);
47: MatDenseRestoreArray(x->A,&px);
48: MatDenseRestoreArray(y->A,&py);
49: return(0);
50: }
54: PetscErrorCode BVMultVec_Mat(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
55: {
57: BV_MAT *x = (BV_MAT*)X->data;
58: PetscScalar *px,*py;
61: MatDenseGetArray(x->A,&px);
62: VecGetArray(y,&py);
63: BVMultVec_BLAS_Private(X,X->n,X->k-X->l,alpha,px+(X->nc+X->l)*X->n,q,beta,py);
64: MatDenseRestoreArray(x->A,&px);
65: VecRestoreArray(y,&py);
66: return(0);
67: }
71: PetscErrorCode BVMultInPlace_Mat(BV V,Mat Q,PetscInt s,PetscInt e)
72: {
74: BV_MAT *ctx = (BV_MAT*)V->data;
75: PetscScalar *pv,*q;
76: PetscInt ldq;
79: MatGetSize(Q,&ldq,NULL);
80: MatDenseGetArray(ctx->A,&pv);
81: MatDenseGetArray(Q,&q);
82: BVMultInPlace_BLAS_Private(V,V->n,V->k-V->l,ldq,s-V->l,e-V->l,pv+(V->nc+V->l)*V->n,q+V->l*ldq+V->l,PETSC_FALSE);
83: MatDenseRestoreArray(Q,&q);
84: MatDenseRestoreArray(ctx->A,&pv);
85: return(0);
86: }
90: PetscErrorCode BVMultInPlaceTranspose_Mat(BV V,Mat Q,PetscInt s,PetscInt e)
91: {
93: BV_MAT *ctx = (BV_MAT*)V->data;
94: PetscScalar *pv,*q;
95: PetscInt ldq;
98: MatGetSize(Q,&ldq,NULL);
99: MatDenseGetArray(ctx->A,&pv);
100: MatDenseGetArray(Q,&q);
101: BVMultInPlace_BLAS_Private(V,V->n,V->k-V->l,ldq,s-V->l,e-V->l,pv+(V->nc+V->l)*V->n,q+V->l*ldq+V->l,PETSC_TRUE);
102: MatDenseRestoreArray(Q,&q);
103: MatDenseRestoreArray(ctx->A,&pv);
104: return(0);
105: }
109: PetscErrorCode BVAXPY_Mat(BV Y,PetscScalar alpha,BV X)
110: {
112: BV_MAT *x = (BV_MAT*)X->data,*y = (BV_MAT*)Y->data;
113: PetscScalar *px,*py;
116: MatDenseGetArray(x->A,&px);
117: MatDenseGetArray(y->A,&py);
118: BVAXPY_BLAS_Private(Y,Y->n,Y->k-Y->l,alpha,px+(X->nc+X->l)*X->n,py+(Y->nc+Y->l)*Y->n);
119: MatDenseRestoreArray(x->A,&px);
120: MatDenseRestoreArray(y->A,&py);
121: return(0);
122: }
126: PetscErrorCode BVDot_Mat(BV X,BV Y,Mat M)
127: {
129: BV_MAT *x = (BV_MAT*)X->data,*y = (BV_MAT*)Y->data;
130: PetscScalar *px,*py,*m;
131: PetscInt ldm;
134: MatGetSize(M,&ldm,NULL);
135: MatDenseGetArray(x->A,&px);
136: MatDenseGetArray(y->A,&py);
137: MatDenseGetArray(M,&m);
138: BVDot_BLAS_Private(X,Y->k-Y->l,X->k-X->l,X->n,ldm,py+(Y->nc+Y->l)*Y->n,px+(X->nc+X->l)*X->n,m+X->l*ldm+Y->l,x->mpi);
139: MatDenseRestoreArray(M,&m);
140: MatDenseRestoreArray(x->A,&px);
141: MatDenseRestoreArray(y->A,&py);
142: return(0);
143: }
147: PetscErrorCode BVDotVec_Mat(BV X,Vec y,PetscScalar *m)
148: {
150: BV_MAT *x = (BV_MAT*)X->data;
151: PetscScalar *px,*py;
152: Vec z = y;
155: if (X->matrix) {
156: BV_IPMatMult(X,y);
157: z = X->Bx;
158: }
159: MatDenseGetArray(x->A,&px);
160: VecGetArray(z,&py);
161: BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,m,x->mpi);
162: VecRestoreArray(z,&py);
163: MatDenseRestoreArray(x->A,&px);
164: return(0);
165: }
169: PetscErrorCode BVScale_Mat(BV bv,PetscInt j,PetscScalar alpha)
170: {
172: BV_MAT *ctx = (BV_MAT*)bv->data;
173: PetscScalar *array;
176: MatDenseGetArray(ctx->A,&array);
177: if (j<0) {
178: BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->n,array+(bv->nc+bv->l)*bv->n,alpha);
179: } else {
180: BVScale_BLAS_Private(bv,bv->n,array+(bv->nc+j)*bv->n,alpha);
181: }
182: MatDenseRestoreArray(ctx->A,&array);
183: return(0);
184: }
188: PetscErrorCode BVNorm_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
189: {
191: BV_MAT *ctx = (BV_MAT*)bv->data;
192: PetscScalar *array;
195: MatDenseGetArray(ctx->A,&array);
196: if (j<0) {
197: BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,ctx->mpi);
198: } else {
199: BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,ctx->mpi);
200: }
201: MatDenseRestoreArray(ctx->A,&array);
202: return(0);
203: }
207: PetscErrorCode BVOrthogonalize_Mat(BV V,Mat R)
208: {
210: BV_MAT *ctx = (BV_MAT*)V->data;
211: PetscScalar *pv,*r=NULL;
214: if (R) { MatDenseGetArray(R,&r); }
215: MatDenseGetArray(ctx->A,&pv);
216: BVOrthogonalize_LAPACK_Private(V,V->n,V->k,pv+V->nc*V->n,r,ctx->mpi);
217: MatDenseRestoreArray(ctx->A,&pv);
218: if (R) { MatDenseRestoreArray(R,&r); }
219: return(0);
220: }
224: PetscErrorCode BVMatMult_Mat(BV V,Mat A,BV W)
225: {
227: BV_MAT *v = (BV_MAT*)V->data,*w = (BV_MAT*)W->data;
228: PetscScalar *pv,*pw;
229: PetscInt j;
232: MatDenseGetArray(v->A,&pv);
233: MatDenseGetArray(w->A,&pw);
234: for (j=0;j<V->k-V->l;j++) {
235: VecPlaceArray(V->cv[1],pv+(V->nc+V->l+j)*V->n);
236: VecPlaceArray(W->cv[1],pw+(W->nc+W->l+j)*W->n);
237: MatMult(A,V->cv[1],W->cv[1]);
238: VecResetArray(V->cv[1]);
239: VecResetArray(W->cv[1]);
240: }
241: MatDenseRestoreArray(v->A,&pv);
242: MatDenseRestoreArray(w->A,&pw);
243: return(0);
244: }
248: PetscErrorCode BVCopy_Mat(BV V,BV W)
249: {
251: BV_MAT *v = (BV_MAT*)V->data,*w = (BV_MAT*)W->data;
252: PetscScalar *pv,*pw,*pvc,*pwc;
255: MatDenseGetArray(v->A,&pv);
256: MatDenseGetArray(w->A,&pw);
257: pvc = pv+(V->nc+V->l)*V->n;
258: pwc = pw+(W->nc+W->l)*W->n;
259: PetscMemcpy(pwc,pvc,(V->k-V->l)*V->n*sizeof(PetscScalar));
260: MatDenseRestoreArray(v->A,&pv);
261: MatDenseRestoreArray(w->A,&pw);
262: return(0);
263: }
267: PetscErrorCode BVResize_Mat(BV bv,PetscInt m,PetscBool copy)
268: {
270: BV_MAT *ctx = (BV_MAT*)bv->data;
271: PetscScalar *pA,*pnew;
272: Mat A;
273: char str[50];
276: MatCreateDense(PetscObjectComm((PetscObject)bv->t),bv->n,PETSC_DECIDE,PETSC_DECIDE,m,NULL,&A);
277: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
278: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
279: PetscLogObjectParent((PetscObject)bv,(PetscObject)A);
280: if (((PetscObject)bv)->name) {
281: PetscSNPrintf(str,50,"%s_0",((PetscObject)bv)->name);
282: PetscObjectSetName((PetscObject)A,str);
283: }
284: if (copy) {
285: MatDenseGetArray(ctx->A,&pA);
286: MatDenseGetArray(A,&pnew);
287: PetscMemcpy(pnew,pA,PetscMin(m,bv->m)*bv->n*sizeof(PetscScalar));
288: MatDenseRestoreArray(ctx->A,&pA);
289: MatDenseRestoreArray(A,&pnew);
290: }
291: MatDestroy(&ctx->A);
292: ctx->A = A;
293: return(0);
294: }
298: PetscErrorCode BVGetColumn_Mat(BV bv,PetscInt j,Vec *v)
299: {
301: BV_MAT *ctx = (BV_MAT*)bv->data;
302: PetscScalar *pA;
303: PetscInt l;
306: l = BVAvailableVec;
307: MatDenseGetArray(ctx->A,&pA);
308: VecPlaceArray(bv->cv[l],pA+(bv->nc+j)*bv->n);
309: return(0);
310: }
314: PetscErrorCode BVRestoreColumn_Mat(BV bv,PetscInt j,Vec *v)
315: {
317: BV_MAT *ctx = (BV_MAT*)bv->data;
318: PetscScalar *pA;
319: PetscInt l;
322: l = (j==bv->ci[0])? 0: 1;
323: VecResetArray(bv->cv[l]);
324: MatDenseRestoreArray(ctx->A,&pA);
325: return(0);
326: }
330: PetscErrorCode BVGetArray_Mat(BV bv,PetscScalar **a)
331: {
333: BV_MAT *ctx = (BV_MAT*)bv->data;
336: MatDenseGetArray(ctx->A,a);
337: return(0);
338: }
342: PetscErrorCode BVRestoreArray_Mat(BV bv,PetscScalar **a)
343: {
345: BV_MAT *ctx = (BV_MAT*)bv->data;
348: if (a) { MatDenseRestoreArray(ctx->A,a); }
349: return(0);
350: }
354: PetscErrorCode BVView_Mat(BV bv,PetscViewer viewer)
355: {
356: PetscErrorCode ierr;
357: BV_MAT *ctx = (BV_MAT*)bv->data;
358: PetscViewerFormat format;
359: PetscBool isascii;
362: MatView(ctx->A,viewer);
363: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
364: if (isascii) {
365: PetscViewerGetFormat(viewer,&format);
366: if (format == PETSC_VIEWER_ASCII_MATLAB) {
367: PetscViewerASCIIPrintf(viewer,"%s=%s;clear %s\n",((PetscObject)bv)->name,((PetscObject)ctx->A)->name,((PetscObject)ctx->A)->name);
368: if (bv->nc) {
369: PetscViewerASCIIPrintf(viewer,"%s=%s(:,%D:end);\n",((PetscObject)bv)->name,((PetscObject)bv)->name,bv->nc+1);
370: }
371: }
372: }
373: return(0);
374: }
378: PetscErrorCode BVDestroy_Mat(BV bv)
379: {
381: BV_MAT *ctx = (BV_MAT*)bv->data;
384: MatDestroy(&ctx->A);
385: VecDestroy(&bv->cv[0]);
386: VecDestroy(&bv->cv[1]);
387: PetscFree(bv->data);
388: return(0);
389: }
393: PETSC_EXTERN PetscErrorCode BVCreate_Mat(BV bv)
394: {
396: BV_MAT *ctx;
397: PetscInt nloc,bs;
398: PetscBool seq;
399: char str[50];
402: PetscNewLog(bv,&ctx);
403: bv->data = (void*)ctx;
405: PetscObjectTypeCompare((PetscObject)bv->t,VECMPI,&ctx->mpi);
406: if (!ctx->mpi) {
407: PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
408: if (!seq) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot create a BVMAT from a non-standard template vector");
409: }
411: VecGetLocalSize(bv->t,&nloc);
412: VecGetBlockSize(bv->t,&bs);
414: MatCreateDense(PetscObjectComm((PetscObject)bv->t),nloc,PETSC_DECIDE,PETSC_DECIDE,bv->m,NULL,&ctx->A);
415: MatAssemblyBegin(ctx->A,MAT_FINAL_ASSEMBLY);
416: MatAssemblyEnd(ctx->A,MAT_FINAL_ASSEMBLY);
417: PetscLogObjectParent((PetscObject)bv,(PetscObject)ctx->A);
418: if (((PetscObject)bv)->name) {
419: PetscSNPrintf(str,50,"%s_0",((PetscObject)bv)->name);
420: PetscObjectSetName((PetscObject)ctx->A,str);
421: }
423: if (ctx->mpi) {
424: VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,NULL,&bv->cv[0]);
425: VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,NULL,&bv->cv[1]);
426: } else {
427: VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,NULL,&bv->cv[0]);
428: VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,NULL,&bv->cv[1]);
429: }
431: bv->ops->mult = BVMult_Mat;
432: bv->ops->multvec = BVMultVec_Mat;
433: bv->ops->multinplace = BVMultInPlace_Mat;
434: bv->ops->multinplacetrans = BVMultInPlaceTranspose_Mat;
435: bv->ops->axpy = BVAXPY_Mat;
436: bv->ops->dot = BVDot_Mat;
437: bv->ops->dotvec = BVDotVec_Mat;
438: bv->ops->scale = BVScale_Mat;
439: bv->ops->norm = BVNorm_Mat;
440: /*bv->ops->orthogonalize = BVOrthogonalize_Mat;*/
441: bv->ops->matmult = BVMatMult_Mat;
442: bv->ops->copy = BVCopy_Mat;
443: bv->ops->resize = BVResize_Mat;
444: bv->ops->getcolumn = BVGetColumn_Mat;
445: bv->ops->restorecolumn = BVRestoreColumn_Mat;
446: bv->ops->getarray = BVGetArray_Mat;
447: bv->ops->restorearray = BVRestoreArray_Mat;
448: bv->ops->destroy = BVDestroy_Mat;
449: if (!ctx->mpi) bv->ops->view = BVView_Mat;
450: return(0);
451: }