Actual source code: svec.c
slepc-3.7.0 2016-05-16
1: /*
2: BV implemented as a single Vec
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2016, 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: PetscBool mpi;
29: } BV_SVEC;
33: PetscErrorCode BVMult_Svec(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
34: {
35: PetscErrorCode ierr;
36: BV_SVEC *y = (BV_SVEC*)Y->data,*x = (BV_SVEC*)X->data;
37: const PetscScalar *px;
38: PetscScalar *py,*q;
39: PetscInt ldq;
42: VecGetArrayRead(x->v,&px);
43: VecGetArray(y->v,&py);
44: if (Q) {
45: MatGetSize(Q,&ldq,NULL);
46: MatDenseGetArray(Q,&q);
47: 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);
48: MatDenseRestoreArray(Q,&q);
49: } else {
50: BVAXPY_BLAS_Private(Y,Y->n,Y->k-Y->l,alpha,px+(X->nc+X->l)*X->n,beta,py+(Y->nc+Y->l)*Y->n);
51: }
52: VecRestoreArrayRead(x->v,&px);
53: VecRestoreArray(y->v,&py);
54: return(0);
55: }
59: PetscErrorCode BVMultVec_Svec(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
60: {
62: BV_SVEC *x = (BV_SVEC*)X->data;
63: PetscScalar *px,*py;
66: VecGetArray(x->v,&px);
67: VecGetArray(y,&py);
68: BVMultVec_BLAS_Private(X,X->n,X->k-X->l,alpha,px+(X->nc+X->l)*X->n,q,beta,py);
69: VecRestoreArray(x->v,&px);
70: VecRestoreArray(y,&py);
71: return(0);
72: }
76: PetscErrorCode BVMultInPlace_Svec(BV V,Mat Q,PetscInt s,PetscInt e)
77: {
79: BV_SVEC *ctx = (BV_SVEC*)V->data;
80: PetscScalar *pv,*q;
81: PetscInt ldq;
84: MatGetSize(Q,&ldq,NULL);
85: VecGetArray(ctx->v,&pv);
86: MatDenseGetArray(Q,&q);
87: 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);
88: MatDenseRestoreArray(Q,&q);
89: VecRestoreArray(ctx->v,&pv);
90: return(0);
91: }
95: PetscErrorCode BVMultInPlaceTranspose_Svec(BV V,Mat Q,PetscInt s,PetscInt e)
96: {
98: BV_SVEC *ctx = (BV_SVEC*)V->data;
99: PetscScalar *pv,*q;
100: PetscInt ldq;
103: MatGetSize(Q,&ldq,NULL);
104: VecGetArray(ctx->v,&pv);
105: MatDenseGetArray(Q,&q);
106: 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);
107: MatDenseRestoreArray(Q,&q);
108: VecRestoreArray(ctx->v,&pv);
109: return(0);
110: }
114: PetscErrorCode BVDot_Svec(BV X,BV Y,Mat M)
115: {
116: PetscErrorCode ierr;
117: BV_SVEC *x = (BV_SVEC*)X->data,*y = (BV_SVEC*)Y->data;
118: const PetscScalar *px,*py;
119: PetscScalar *m;
120: PetscInt ldm;
123: MatGetSize(M,&ldm,NULL);
124: VecGetArrayRead(x->v,&px);
125: VecGetArrayRead(y->v,&py);
126: MatDenseGetArray(M,&m);
127: 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);
128: MatDenseRestoreArray(M,&m);
129: VecRestoreArrayRead(x->v,&px);
130: VecRestoreArrayRead(y->v,&py);
131: return(0);
132: }
136: PetscErrorCode BVDotVec_Svec(BV X,Vec y,PetscScalar *m)
137: {
138: PetscErrorCode ierr;
139: BV_SVEC *x = (BV_SVEC*)X->data;
140: const PetscScalar *px,*py;
141: Vec z = y;
144: if (X->matrix) {
145: BV_IPMatMult(X,y);
146: z = X->Bx;
147: }
148: VecGetArrayRead(x->v,&px);
149: VecGetArrayRead(z,&py);
150: BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,m,x->mpi);
151: VecRestoreArrayRead(z,&py);
152: VecRestoreArrayRead(x->v,&px);
153: return(0);
154: }
158: PetscErrorCode BVDotVec_Local_Svec(BV X,Vec y,PetscScalar *m)
159: {
161: BV_SVEC *x = (BV_SVEC*)X->data;
162: PetscScalar *px,*py;
163: Vec z = y;
166: if (X->matrix) {
167: BV_IPMatMult(X,y);
168: z = X->Bx;
169: }
170: VecGetArray(x->v,&px);
171: VecGetArray(z,&py);
172: BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,m,PETSC_FALSE);
173: VecRestoreArray(z,&py);
174: VecRestoreArray(x->v,&px);
175: return(0);
176: }
180: PetscErrorCode BVScale_Svec(BV bv,PetscInt j,PetscScalar alpha)
181: {
183: BV_SVEC *ctx = (BV_SVEC*)bv->data;
184: PetscScalar *array;
187: VecGetArray(ctx->v,&array);
188: if (j<0) {
189: BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->n,array+(bv->nc+bv->l)*bv->n,alpha);
190: } else {
191: BVScale_BLAS_Private(bv,bv->n,array+(bv->nc+j)*bv->n,alpha);
192: }
193: VecRestoreArray(ctx->v,&array);
194: return(0);
195: }
199: PetscErrorCode BVNorm_Svec(BV bv,PetscInt j,NormType type,PetscReal *val)
200: {
202: BV_SVEC *ctx = (BV_SVEC*)bv->data;
203: PetscScalar *array;
206: VecGetArray(ctx->v,&array);
207: if (j<0) {
208: BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,ctx->mpi);
209: } else {
210: BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,ctx->mpi);
211: }
212: VecRestoreArray(ctx->v,&array);
213: return(0);
214: }
218: PetscErrorCode BVNorm_Local_Svec(BV bv,PetscInt j,NormType type,PetscReal *val)
219: {
221: BV_SVEC *ctx = (BV_SVEC*)bv->data;
222: PetscScalar *array;
225: VecGetArray(ctx->v,&array);
226: if (j<0) {
227: BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,PETSC_FALSE);
228: } else {
229: BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,PETSC_FALSE);
230: }
231: VecRestoreArray(ctx->v,&array);
232: return(0);
233: }
237: PetscErrorCode BVOrthogonalize_Svec(BV V,Mat R)
238: {
240: BV_SVEC *ctx = (BV_SVEC*)V->data;
241: PetscScalar *pv,*r=NULL;
244: if (R) { MatDenseGetArray(R,&r); }
245: VecGetArray(ctx->v,&pv);
246: BVOrthogonalize_LAPACK_Private(V,V->n,V->k,pv+V->nc*V->n,r,ctx->mpi);
247: VecRestoreArray(ctx->v,&pv);
248: if (R) { MatDenseRestoreArray(R,&r); }
249: return(0);
250: }
254: PetscErrorCode BVMatMult_Svec(BV V,Mat A,BV W)
255: {
257: BV_SVEC *v = (BV_SVEC*)V->data,*w = (BV_SVEC*)W->data;
258: PetscScalar *pv,*pw,*pb,*pc;
259: PetscInt j,m;
260: PetscBool flg;
263: VecGetArray(v->v,&pv);
264: VecGetArray(w->v,&pw);
265: MatHasOperation(A,MATOP_MAT_MULT,&flg);
266: if (V->vmm && flg) {
267: m = V->k-V->l;
268: if (V->vmm==BV_MATMULT_MAT_SAVE) {
269: BV_AllocateMatMult(V,A,m);
270: MatDenseGetArray(V->B,&pb);
271: PetscMemcpy(pb,pv+(V->nc+V->l)*V->n,m*V->n*sizeof(PetscScalar));
272: MatDenseRestoreArray(V->B,&pb);
273: } else { /* BV_MATMULT_MAT */
274: MatCreateDense(PetscObjectComm((PetscObject)V),V->n,PETSC_DECIDE,V->N,m,pv+(V->nc+V->l)*V->n,&V->B);
275: }
276: if (!V->C) {
277: MatMatMultSymbolic(A,V->B,PETSC_DEFAULT,&V->C);
278: }
279: MatMatMultNumeric(A,V->B,V->C);
280: MatDenseGetArray(V->C,&pc);
281: PetscMemcpy(pw+(W->nc+W->l)*W->n,pc,m*V->n*sizeof(PetscScalar));
282: MatDenseRestoreArray(V->C,&pc);
283: if (V->vmm==BV_MATMULT_MAT) {
284: MatDestroy(&V->B);
285: MatDestroy(&V->C);
286: }
287: } else {
288: for (j=0;j<V->k-V->l;j++) {
289: VecPlaceArray(V->cv[1],pv+(V->nc+V->l+j)*V->n);
290: VecPlaceArray(W->cv[1],pw+(W->nc+W->l+j)*W->n);
291: MatMult(A,V->cv[1],W->cv[1]);
292: VecResetArray(V->cv[1]);
293: VecResetArray(W->cv[1]);
294: }
295: }
296: VecRestoreArray(v->v,&pv);
297: VecRestoreArray(w->v,&pw);
298: return(0);
299: }
303: PetscErrorCode BVCopy_Svec(BV V,BV W)
304: {
306: BV_SVEC *v = (BV_SVEC*)V->data,*w = (BV_SVEC*)W->data;
307: PetscScalar *pv,*pw,*pvc,*pwc;
310: VecGetArray(v->v,&pv);
311: VecGetArray(w->v,&pw);
312: pvc = pv+(V->nc+V->l)*V->n;
313: pwc = pw+(W->nc+W->l)*W->n;
314: PetscMemcpy(pwc,pvc,(V->k-V->l)*V->n*sizeof(PetscScalar));
315: VecRestoreArray(v->v,&pv);
316: VecRestoreArray(w->v,&pw);
317: return(0);
318: }
322: PetscErrorCode BVResize_Svec(BV bv,PetscInt m,PetscBool copy)
323: {
325: BV_SVEC *ctx = (BV_SVEC*)bv->data;
326: PetscScalar *pv,*pnew;
327: PetscInt bs;
328: Vec vnew;
329: char str[50];
332: VecGetBlockSize(bv->t,&bs);
333: VecCreate(PetscObjectComm((PetscObject)bv->t),&vnew);
334: VecSetType(vnew,((PetscObject)bv->t)->type_name);
335: VecSetSizes(vnew,m*bv->n,PETSC_DECIDE);
336: VecSetBlockSize(vnew,bs);
337: PetscLogObjectParent((PetscObject)bv,(PetscObject)vnew);
338: if (((PetscObject)bv)->name) {
339: PetscSNPrintf(str,50,"%s_0",((PetscObject)bv)->name);
340: PetscObjectSetName((PetscObject)vnew,str);
341: }
342: if (copy) {
343: VecGetArray(ctx->v,&pv);
344: VecGetArray(vnew,&pnew);
345: PetscMemcpy(pnew,pv,PetscMin(m,bv->m)*bv->n*sizeof(PetscScalar));
346: VecRestoreArray(ctx->v,&pv);
347: VecRestoreArray(vnew,&pnew);
348: }
349: VecDestroy(&ctx->v);
350: ctx->v = vnew;
351: return(0);
352: }
356: PetscErrorCode BVGetColumn_Svec(BV bv,PetscInt j,Vec *v)
357: {
359: BV_SVEC *ctx = (BV_SVEC*)bv->data;
360: PetscScalar *pv;
361: PetscInt l;
364: l = BVAvailableVec;
365: VecGetArray(ctx->v,&pv);
366: VecPlaceArray(bv->cv[l],pv+(bv->nc+j)*bv->n);
367: return(0);
368: }
372: PetscErrorCode BVRestoreColumn_Svec(BV bv,PetscInt j,Vec *v)
373: {
375: BV_SVEC *ctx = (BV_SVEC*)bv->data;
376: PetscInt l;
379: l = (j==bv->ci[0])? 0: 1;
380: VecResetArray(bv->cv[l]);
381: VecRestoreArray(ctx->v,NULL);
382: return(0);
383: }
387: PetscErrorCode BVGetArray_Svec(BV bv,PetscScalar **a)
388: {
390: BV_SVEC *ctx = (BV_SVEC*)bv->data;
393: VecGetArray(ctx->v,a);
394: return(0);
395: }
399: PetscErrorCode BVRestoreArray_Svec(BV bv,PetscScalar **a)
400: {
402: BV_SVEC *ctx = (BV_SVEC*)bv->data;
405: VecRestoreArray(ctx->v,a);
406: return(0);
407: }
411: PetscErrorCode BVGetArrayRead_Svec(BV bv,const PetscScalar **a)
412: {
414: BV_SVEC *ctx = (BV_SVEC*)bv->data;
417: VecGetArrayRead(ctx->v,a);
418: return(0);
419: }
423: PetscErrorCode BVRestoreArrayRead_Svec(BV bv,const PetscScalar **a)
424: {
426: BV_SVEC *ctx = (BV_SVEC*)bv->data;
429: VecRestoreArrayRead(ctx->v,a);
430: return(0);
431: }
435: PetscErrorCode BVView_Svec(BV bv,PetscViewer viewer)
436: {
437: PetscErrorCode ierr;
438: BV_SVEC *ctx = (BV_SVEC*)bv->data;
439: PetscViewerFormat format;
440: PetscBool isascii;
441: const char *bvname,*name;
444: VecView(ctx->v,viewer);
445: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
446: if (isascii) {
447: PetscViewerGetFormat(viewer,&format);
448: if (format == PETSC_VIEWER_ASCII_MATLAB) {
449: PetscObjectGetName((PetscObject)bv,&bvname);
450: PetscObjectGetName((PetscObject)ctx->v,&name);
451: PetscViewerASCIIPrintf(viewer,"%s=reshape(%s,%D,%D);clear %s\n",bvname,name,bv->N,bv->nc+bv->m,name);
452: if (bv->nc) {
453: PetscViewerASCIIPrintf(viewer,"%s=%s(:,%D:end);\n",bvname,bvname,bv->nc+1);
454: }
455: }
456: }
457: return(0);
458: }
462: PetscErrorCode BVDestroy_Svec(BV bv)
463: {
465: BV_SVEC *ctx = (BV_SVEC*)bv->data;
468: VecDestroy(&ctx->v);
469: VecDestroy(&bv->cv[0]);
470: VecDestroy(&bv->cv[1]);
471: PetscFree(bv->data);
472: return(0);
473: }
477: PETSC_EXTERN PetscErrorCode BVCreate_Svec(BV bv)
478: {
480: BV_SVEC *ctx;
481: PetscInt nloc,bs;
482: PetscBool seq;
483: char str[50];
486: PetscNewLog(bv,&ctx);
487: bv->data = (void*)ctx;
489: PetscObjectTypeCompare((PetscObject)bv->t,VECMPI,&ctx->mpi);
490: if (!ctx->mpi) {
491: PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
492: if (!seq) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot create a BVSVEC from a non-standard template vector");
493: }
495: VecGetLocalSize(bv->t,&nloc);
496: VecGetBlockSize(bv->t,&bs);
498: VecCreate(PetscObjectComm((PetscObject)bv->t),&ctx->v);
499: VecSetType(ctx->v,((PetscObject)bv->t)->type_name);
500: VecSetSizes(ctx->v,bv->m*nloc,PETSC_DECIDE);
501: VecSetBlockSize(ctx->v,bs);
502: PetscLogObjectParent((PetscObject)bv,(PetscObject)ctx->v);
503: if (((PetscObject)bv)->name) {
504: PetscSNPrintf(str,50,"%s_0",((PetscObject)bv)->name);
505: PetscObjectSetName((PetscObject)ctx->v,str);
506: }
508: if (ctx->mpi) {
509: VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,NULL,&bv->cv[0]);
510: VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,NULL,&bv->cv[1]);
511: } else {
512: VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,NULL,&bv->cv[0]);
513: VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,NULL,&bv->cv[1]);
514: }
516: bv->ops->mult = BVMult_Svec;
517: bv->ops->multvec = BVMultVec_Svec;
518: bv->ops->multinplace = BVMultInPlace_Svec;
519: bv->ops->multinplacetrans = BVMultInPlaceTranspose_Svec;
520: bv->ops->dot = BVDot_Svec;
521: bv->ops->dotvec = BVDotVec_Svec;
522: bv->ops->dotvec_local = BVDotVec_Local_Svec;
523: bv->ops->scale = BVScale_Svec;
524: bv->ops->norm = BVNorm_Svec;
525: bv->ops->norm_local = BVNorm_Local_Svec;
526: /*bv->ops->orthogonalize = BVOrthogonalize_Svec;*/
527: bv->ops->matmult = BVMatMult_Svec;
528: bv->ops->copy = BVCopy_Svec;
529: bv->ops->resize = BVResize_Svec;
530: bv->ops->getcolumn = BVGetColumn_Svec;
531: bv->ops->restorecolumn = BVRestoreColumn_Svec;
532: bv->ops->getarray = BVGetArray_Svec;
533: bv->ops->restorearray = BVRestoreArray_Svec;
534: bv->ops->getarrayread = BVGetArrayRead_Svec;
535: bv->ops->restorearrayread = BVRestoreArrayRead_Svec;
536: bv->ops->destroy = BVDestroy_Svec;
537: if (!ctx->mpi) bv->ops->view = BVView_Svec;
538: return(0);
539: }