Actual source code: bvmat.c
slepc-3.6.1 2015-09-03
1: /*
2: BV implemented with a dense Mat
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: 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: {
149: PetscErrorCode ierr;
150: BV_MAT *x = (BV_MAT*)X->data;
151: PetscScalar *px;
152: const PetscScalar *py;
153: Vec z = y;
156: if (X->matrix) {
157: BV_IPMatMult(X,y);
158: z = X->Bx;
159: }
160: MatDenseGetArray(x->A,&px);
161: VecGetArrayRead(z,&py);
162: BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,m,x->mpi);
163: VecRestoreArrayRead(z,&py);
164: MatDenseRestoreArray(x->A,&px);
165: return(0);
166: }
170: PetscErrorCode BVDotVec_Local_Mat(BV X,Vec y,PetscScalar *m)
171: {
173: BV_MAT *x = (BV_MAT*)X->data;
174: PetscScalar *px,*py;
175: Vec z = y;
178: if (X->matrix) {
179: BV_IPMatMult(X,y);
180: z = X->Bx;
181: }
182: MatDenseGetArray(x->A,&px);
183: VecGetArray(z,&py);
184: BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,m,PETSC_FALSE);
185: VecRestoreArray(z,&py);
186: MatDenseRestoreArray(x->A,&px);
187: return(0);
188: }
192: PetscErrorCode BVScale_Mat(BV bv,PetscInt j,PetscScalar alpha)
193: {
195: BV_MAT *ctx = (BV_MAT*)bv->data;
196: PetscScalar *array;
199: MatDenseGetArray(ctx->A,&array);
200: if (j<0) {
201: BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->n,array+(bv->nc+bv->l)*bv->n,alpha);
202: } else {
203: BVScale_BLAS_Private(bv,bv->n,array+(bv->nc+j)*bv->n,alpha);
204: }
205: MatDenseRestoreArray(ctx->A,&array);
206: return(0);
207: }
211: PetscErrorCode BVNorm_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
212: {
214: BV_MAT *ctx = (BV_MAT*)bv->data;
215: PetscScalar *array;
218: MatDenseGetArray(ctx->A,&array);
219: if (j<0) {
220: BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,ctx->mpi);
221: } else {
222: BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,ctx->mpi);
223: }
224: MatDenseRestoreArray(ctx->A,&array);
225: return(0);
226: }
230: PetscErrorCode BVNorm_Local_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
231: {
233: BV_MAT *ctx = (BV_MAT*)bv->data;
234: PetscScalar *array;
237: MatDenseGetArray(ctx->A,&array);
238: if (j<0) {
239: BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,PETSC_FALSE);
240: } else {
241: BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,PETSC_FALSE);
242: }
243: MatDenseRestoreArray(ctx->A,&array);
244: return(0);
245: }
249: PetscErrorCode BVOrthogonalize_Mat(BV V,Mat R)
250: {
252: BV_MAT *ctx = (BV_MAT*)V->data;
253: PetscScalar *pv,*r=NULL;
256: if (R) { MatDenseGetArray(R,&r); }
257: MatDenseGetArray(ctx->A,&pv);
258: BVOrthogonalize_LAPACK_Private(V,V->n,V->k,pv+V->nc*V->n,r,ctx->mpi);
259: MatDenseRestoreArray(ctx->A,&pv);
260: if (R) { MatDenseRestoreArray(R,&r); }
261: return(0);
262: }
266: PetscErrorCode BVMatMult_Mat(BV V,Mat A,BV W)
267: {
269: BV_MAT *v = (BV_MAT*)V->data,*w = (BV_MAT*)W->data;
270: PetscScalar *pv,*pw,*pb,*pc;
271: PetscInt j,m;
272: PetscBool flg;
275: MatDenseGetArray(v->A,&pv);
276: MatDenseGetArray(w->A,&pw);
277: MatHasOperation(A,MATOP_MAT_MULT,&flg);
278: if (V->vmm && flg) {
279: m = V->k-V->l;
280: if (V->vmm==BV_MATMULT_MAT_SAVE) {
281: BV_AllocateMatMult(V,A,m);
282: MatDenseGetArray(V->B,&pb);
283: PetscMemcpy(pb,pv+(V->nc+V->l)*V->n,m*V->n*sizeof(PetscScalar));
284: MatDenseRestoreArray(V->B,&pb);
285: } else { /* BV_MATMULT_MAT */
286: MatCreateDense(PetscObjectComm((PetscObject)V),V->n,PETSC_DECIDE,V->N,m,pv+(V->nc+V->l)*V->n,&V->B);
287: }
288: if (!V->C) {
289: MatMatMultSymbolic(A,V->B,PETSC_DEFAULT,&V->C);
290: }
291: MatMatMultNumeric(A,V->B,V->C);
292: MatDenseGetArray(V->C,&pc);
293: PetscMemcpy(pw+(W->nc+W->l)*W->n,pc,m*V->n*sizeof(PetscScalar));
294: MatDenseRestoreArray(V->C,&pc);
295: if (V->vmm==BV_MATMULT_MAT) {
296: MatDestroy(&V->B);
297: MatDestroy(&V->C);
298: }
299: } else {
300: for (j=0;j<V->k-V->l;j++) {
301: VecPlaceArray(V->cv[1],pv+(V->nc+V->l+j)*V->n);
302: VecPlaceArray(W->cv[1],pw+(W->nc+W->l+j)*W->n);
303: MatMult(A,V->cv[1],W->cv[1]);
304: VecResetArray(V->cv[1]);
305: VecResetArray(W->cv[1]);
306: }
307: }
308: MatDenseRestoreArray(v->A,&pv);
309: MatDenseRestoreArray(w->A,&pw);
310: return(0);
311: }
315: PetscErrorCode BVCopy_Mat(BV V,BV W)
316: {
318: BV_MAT *v = (BV_MAT*)V->data,*w = (BV_MAT*)W->data;
319: PetscScalar *pv,*pw,*pvc,*pwc;
322: MatDenseGetArray(v->A,&pv);
323: MatDenseGetArray(w->A,&pw);
324: pvc = pv+(V->nc+V->l)*V->n;
325: pwc = pw+(W->nc+W->l)*W->n;
326: PetscMemcpy(pwc,pvc,(V->k-V->l)*V->n*sizeof(PetscScalar));
327: MatDenseRestoreArray(v->A,&pv);
328: MatDenseRestoreArray(w->A,&pw);
329: return(0);
330: }
334: PetscErrorCode BVResize_Mat(BV bv,PetscInt m,PetscBool copy)
335: {
337: BV_MAT *ctx = (BV_MAT*)bv->data;
338: PetscScalar *pA,*pnew;
339: Mat A;
340: char str[50];
343: MatCreateDense(PetscObjectComm((PetscObject)bv->t),bv->n,PETSC_DECIDE,PETSC_DECIDE,m,NULL,&A);
344: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
345: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
346: PetscLogObjectParent((PetscObject)bv,(PetscObject)A);
347: if (((PetscObject)bv)->name) {
348: PetscSNPrintf(str,50,"%s_0",((PetscObject)bv)->name);
349: PetscObjectSetName((PetscObject)A,str);
350: }
351: if (copy) {
352: MatDenseGetArray(ctx->A,&pA);
353: MatDenseGetArray(A,&pnew);
354: PetscMemcpy(pnew,pA,PetscMin(m,bv->m)*bv->n*sizeof(PetscScalar));
355: MatDenseRestoreArray(ctx->A,&pA);
356: MatDenseRestoreArray(A,&pnew);
357: }
358: MatDestroy(&ctx->A);
359: ctx->A = A;
360: return(0);
361: }
365: PetscErrorCode BVGetColumn_Mat(BV bv,PetscInt j,Vec *v)
366: {
368: BV_MAT *ctx = (BV_MAT*)bv->data;
369: PetscScalar *pA;
370: PetscInt l;
373: l = BVAvailableVec;
374: MatDenseGetArray(ctx->A,&pA);
375: VecPlaceArray(bv->cv[l],pA+(bv->nc+j)*bv->n);
376: return(0);
377: }
381: PetscErrorCode BVRestoreColumn_Mat(BV bv,PetscInt j,Vec *v)
382: {
384: BV_MAT *ctx = (BV_MAT*)bv->data;
385: PetscScalar *pA;
386: PetscInt l;
389: l = (j==bv->ci[0])? 0: 1;
390: VecResetArray(bv->cv[l]);
391: MatDenseRestoreArray(ctx->A,&pA);
392: return(0);
393: }
397: PetscErrorCode BVGetArray_Mat(BV bv,PetscScalar **a)
398: {
400: BV_MAT *ctx = (BV_MAT*)bv->data;
403: MatDenseGetArray(ctx->A,a);
404: return(0);
405: }
409: PetscErrorCode BVRestoreArray_Mat(BV bv,PetscScalar **a)
410: {
412: BV_MAT *ctx = (BV_MAT*)bv->data;
415: if (a) { MatDenseRestoreArray(ctx->A,a); }
416: return(0);
417: }
421: PetscErrorCode BVView_Mat(BV bv,PetscViewer viewer)
422: {
423: PetscErrorCode ierr;
424: BV_MAT *ctx = (BV_MAT*)bv->data;
425: PetscViewerFormat format;
426: PetscBool isascii;
427: const char *bvname,*name;
430: MatView(ctx->A,viewer);
431: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
432: if (isascii) {
433: PetscViewerGetFormat(viewer,&format);
434: if (format == PETSC_VIEWER_ASCII_MATLAB) {
435: PetscObjectGetName((PetscObject)bv,&bvname);
436: PetscObjectGetName((PetscObject)ctx->A,&name);
437: PetscViewerASCIIPrintf(viewer,"%s=%s;clear %s\n",bvname,name,name);
438: if (bv->nc) {
439: PetscViewerASCIIPrintf(viewer,"%s=%s(:,%D:end);\n",bvname,bvname,bv->nc+1);
440: }
441: }
442: }
443: return(0);
444: }
448: PetscErrorCode BVDestroy_Mat(BV bv)
449: {
451: BV_MAT *ctx = (BV_MAT*)bv->data;
454: MatDestroy(&ctx->A);
455: VecDestroy(&bv->cv[0]);
456: VecDestroy(&bv->cv[1]);
457: PetscFree(bv->data);
458: return(0);
459: }
463: PETSC_EXTERN PetscErrorCode BVCreate_Mat(BV bv)
464: {
466: BV_MAT *ctx;
467: PetscInt nloc,bs;
468: PetscBool seq;
469: char str[50];
472: PetscNewLog(bv,&ctx);
473: bv->data = (void*)ctx;
475: PetscObjectTypeCompare((PetscObject)bv->t,VECMPI,&ctx->mpi);
476: if (!ctx->mpi) {
477: PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
478: if (!seq) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot create a BVMAT from a non-standard template vector");
479: }
481: VecGetLocalSize(bv->t,&nloc);
482: VecGetBlockSize(bv->t,&bs);
484: MatCreateDense(PetscObjectComm((PetscObject)bv->t),nloc,PETSC_DECIDE,PETSC_DECIDE,bv->m,NULL,&ctx->A);
485: MatAssemblyBegin(ctx->A,MAT_FINAL_ASSEMBLY);
486: MatAssemblyEnd(ctx->A,MAT_FINAL_ASSEMBLY);
487: PetscLogObjectParent((PetscObject)bv,(PetscObject)ctx->A);
488: if (((PetscObject)bv)->name) {
489: PetscSNPrintf(str,50,"%s_0",((PetscObject)bv)->name);
490: PetscObjectSetName((PetscObject)ctx->A,str);
491: }
493: if (ctx->mpi) {
494: VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,NULL,&bv->cv[0]);
495: VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,NULL,&bv->cv[1]);
496: } else {
497: VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,NULL,&bv->cv[0]);
498: VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,NULL,&bv->cv[1]);
499: }
501: bv->ops->mult = BVMult_Mat;
502: bv->ops->multvec = BVMultVec_Mat;
503: bv->ops->multinplace = BVMultInPlace_Mat;
504: bv->ops->multinplacetrans = BVMultInPlaceTranspose_Mat;
505: bv->ops->axpy = BVAXPY_Mat;
506: bv->ops->dot = BVDot_Mat;
507: bv->ops->dotvec = BVDotVec_Mat;
508: bv->ops->dotvec_local = BVDotVec_Local_Mat;
509: bv->ops->scale = BVScale_Mat;
510: bv->ops->norm = BVNorm_Mat;
511: bv->ops->norm_local = BVNorm_Local_Mat;
512: /*bv->ops->orthogonalize = BVOrthogonalize_Mat;*/
513: bv->ops->matmult = BVMatMult_Mat;
514: bv->ops->copy = BVCopy_Mat;
515: bv->ops->resize = BVResize_Mat;
516: bv->ops->getcolumn = BVGetColumn_Mat;
517: bv->ops->restorecolumn = BVRestoreColumn_Mat;
518: bv->ops->getarray = BVGetArray_Mat;
519: bv->ops->restorearray = BVRestoreArray_Mat;
520: bv->ops->destroy = BVDestroy_Mat;
521: if (!ctx->mpi) bv->ops->view = BVView_Mat;
522: return(0);
523: }