Actual source code: svec.c
slepc-3.6.1 2015-09-03
1: /*
2: BV implemented as a single Vec
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: 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: MatGetSize(Q,&ldq,NULL);
43: VecGetArrayRead(x->v,&px);
44: VecGetArray(y->v,&py);
45: MatDenseGetArray(Q,&q);
46: 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);
47: MatDenseRestoreArray(Q,&q);
48: VecRestoreArrayRead(x->v,&px);
49: VecRestoreArray(y->v,&py);
50: return(0);
51: }
55: PetscErrorCode BVMultVec_Svec(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
56: {
58: BV_SVEC *x = (BV_SVEC*)X->data;
59: PetscScalar *px,*py;
62: VecGetArray(x->v,&px);
63: VecGetArray(y,&py);
64: BVMultVec_BLAS_Private(X,X->n,X->k-X->l,alpha,px+(X->nc+X->l)*X->n,q,beta,py);
65: VecRestoreArray(x->v,&px);
66: VecRestoreArray(y,&py);
67: return(0);
68: }
72: PetscErrorCode BVMultInPlace_Svec(BV V,Mat Q,PetscInt s,PetscInt e)
73: {
75: BV_SVEC *ctx = (BV_SVEC*)V->data;
76: PetscScalar *pv,*q;
77: PetscInt ldq;
80: MatGetSize(Q,&ldq,NULL);
81: VecGetArray(ctx->v,&pv);
82: MatDenseGetArray(Q,&q);
83: 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);
84: MatDenseRestoreArray(Q,&q);
85: VecRestoreArray(ctx->v,&pv);
86: return(0);
87: }
91: PetscErrorCode BVMultInPlaceTranspose_Svec(BV V,Mat Q,PetscInt s,PetscInt e)
92: {
94: BV_SVEC *ctx = (BV_SVEC*)V->data;
95: PetscScalar *pv,*q;
96: PetscInt ldq;
99: MatGetSize(Q,&ldq,NULL);
100: VecGetArray(ctx->v,&pv);
101: MatDenseGetArray(Q,&q);
102: 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);
103: MatDenseRestoreArray(Q,&q);
104: VecRestoreArray(ctx->v,&pv);
105: return(0);
106: }
110: PetscErrorCode BVAXPY_Svec(BV Y,PetscScalar alpha,BV X)
111: {
113: BV_SVEC *x = (BV_SVEC*)X->data,*y = (BV_SVEC*)Y->data;
114: PetscScalar *px,*py;
117: VecGetArray(x->v,&px);
118: VecGetArray(y->v,&py);
119: 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);
120: VecRestoreArray(x->v,&px);
121: VecRestoreArray(y->v,&py);
122: return(0);
123: }
127: PetscErrorCode BVDot_Svec(BV X,BV Y,Mat M)
128: {
129: PetscErrorCode ierr;
130: BV_SVEC *x = (BV_SVEC*)X->data,*y = (BV_SVEC*)Y->data;
131: const PetscScalar *px,*py;
132: PetscScalar *m;
133: PetscInt ldm;
136: MatGetSize(M,&ldm,NULL);
137: VecGetArrayRead(x->v,&px);
138: VecGetArrayRead(y->v,&py);
139: MatDenseGetArray(M,&m);
140: 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);
141: MatDenseRestoreArray(M,&m);
142: VecRestoreArrayRead(x->v,&px);
143: VecRestoreArrayRead(y->v,&py);
144: return(0);
145: }
149: PetscErrorCode BVDotVec_Svec(BV X,Vec y,PetscScalar *m)
150: {
151: PetscErrorCode ierr;
152: BV_SVEC *x = (BV_SVEC*)X->data;
153: const PetscScalar *px,*py;
154: Vec z = y;
157: if (X->matrix) {
158: BV_IPMatMult(X,y);
159: z = X->Bx;
160: }
161: VecGetArrayRead(x->v,&px);
162: VecGetArrayRead(z,&py);
163: BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,m,x->mpi);
164: VecRestoreArrayRead(z,&py);
165: VecRestoreArrayRead(x->v,&px);
166: return(0);
167: }
171: PetscErrorCode BVDotVec_Local_Svec(BV X,Vec y,PetscScalar *m)
172: {
174: BV_SVEC *x = (BV_SVEC*)X->data;
175: PetscScalar *px,*py;
176: Vec z = y;
179: if (X->matrix) {
180: BV_IPMatMult(X,y);
181: z = X->Bx;
182: }
183: VecGetArray(x->v,&px);
184: VecGetArray(z,&py);
185: BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,m,PETSC_FALSE);
186: VecRestoreArray(z,&py);
187: VecRestoreArray(x->v,&px);
188: return(0);
189: }
193: PetscErrorCode BVScale_Svec(BV bv,PetscInt j,PetscScalar alpha)
194: {
196: BV_SVEC *ctx = (BV_SVEC*)bv->data;
197: PetscScalar *array;
200: VecGetArray(ctx->v,&array);
201: if (j<0) {
202: BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->n,array+(bv->nc+bv->l)*bv->n,alpha);
203: } else {
204: BVScale_BLAS_Private(bv,bv->n,array+(bv->nc+j)*bv->n,alpha);
205: }
206: VecRestoreArray(ctx->v,&array);
207: return(0);
208: }
212: PetscErrorCode BVNorm_Svec(BV bv,PetscInt j,NormType type,PetscReal *val)
213: {
215: BV_SVEC *ctx = (BV_SVEC*)bv->data;
216: PetscScalar *array;
219: VecGetArray(ctx->v,&array);
220: if (j<0) {
221: BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,ctx->mpi);
222: } else {
223: BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,ctx->mpi);
224: }
225: VecRestoreArray(ctx->v,&array);
226: return(0);
227: }
231: PetscErrorCode BVNorm_Local_Svec(BV bv,PetscInt j,NormType type,PetscReal *val)
232: {
234: BV_SVEC *ctx = (BV_SVEC*)bv->data;
235: PetscScalar *array;
238: VecGetArray(ctx->v,&array);
239: if (j<0) {
240: BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,PETSC_FALSE);
241: } else {
242: BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,PETSC_FALSE);
243: }
244: VecRestoreArray(ctx->v,&array);
245: return(0);
246: }
250: PetscErrorCode BVOrthogonalize_Svec(BV V,Mat R)
251: {
253: BV_SVEC *ctx = (BV_SVEC*)V->data;
254: PetscScalar *pv,*r=NULL;
257: if (R) { MatDenseGetArray(R,&r); }
258: VecGetArray(ctx->v,&pv);
259: BVOrthogonalize_LAPACK_Private(V,V->n,V->k,pv+V->nc*V->n,r,ctx->mpi);
260: VecRestoreArray(ctx->v,&pv);
261: if (R) { MatDenseRestoreArray(R,&r); }
262: return(0);
263: }
267: PetscErrorCode BVMatMult_Svec(BV V,Mat A,BV W)
268: {
270: BV_SVEC *v = (BV_SVEC*)V->data,*w = (BV_SVEC*)W->data;
271: PetscScalar *pv,*pw,*pb,*pc;
272: PetscInt j,m;
273: PetscBool flg;
276: VecGetArray(v->v,&pv);
277: VecGetArray(w->v,&pw);
278: MatHasOperation(A,MATOP_MAT_MULT,&flg);
279: if (V->vmm && flg) {
280: m = V->k-V->l;
281: if (V->vmm==BV_MATMULT_MAT_SAVE) {
282: BV_AllocateMatMult(V,A,m);
283: MatDenseGetArray(V->B,&pb);
284: PetscMemcpy(pb,pv+(V->nc+V->l)*V->n,m*V->n*sizeof(PetscScalar));
285: MatDenseRestoreArray(V->B,&pb);
286: } else { /* BV_MATMULT_MAT */
287: MatCreateDense(PetscObjectComm((PetscObject)V),V->n,PETSC_DECIDE,V->N,m,pv+(V->nc+V->l)*V->n,&V->B);
288: }
289: if (!V->C) {
290: MatMatMultSymbolic(A,V->B,PETSC_DEFAULT,&V->C);
291: }
292: MatMatMultNumeric(A,V->B,V->C);
293: MatDenseGetArray(V->C,&pc);
294: PetscMemcpy(pw+(W->nc+W->l)*W->n,pc,m*V->n*sizeof(PetscScalar));
295: MatDenseRestoreArray(V->C,&pc);
296: if (V->vmm==BV_MATMULT_MAT) {
297: MatDestroy(&V->B);
298: MatDestroy(&V->C);
299: }
300: } else {
301: for (j=0;j<V->k-V->l;j++) {
302: VecPlaceArray(V->cv[1],pv+(V->nc+V->l+j)*V->n);
303: VecPlaceArray(W->cv[1],pw+(W->nc+W->l+j)*W->n);
304: MatMult(A,V->cv[1],W->cv[1]);
305: VecResetArray(V->cv[1]);
306: VecResetArray(W->cv[1]);
307: }
308: }
309: VecRestoreArray(v->v,&pv);
310: VecRestoreArray(w->v,&pw);
311: return(0);
312: }
316: PetscErrorCode BVCopy_Svec(BV V,BV W)
317: {
319: BV_SVEC *v = (BV_SVEC*)V->data,*w = (BV_SVEC*)W->data;
320: PetscScalar *pv,*pw,*pvc,*pwc;
323: VecGetArray(v->v,&pv);
324: VecGetArray(w->v,&pw);
325: pvc = pv+(V->nc+V->l)*V->n;
326: pwc = pw+(W->nc+W->l)*W->n;
327: PetscMemcpy(pwc,pvc,(V->k-V->l)*V->n*sizeof(PetscScalar));
328: VecRestoreArray(v->v,&pv);
329: VecRestoreArray(w->v,&pw);
330: return(0);
331: }
335: PetscErrorCode BVResize_Svec(BV bv,PetscInt m,PetscBool copy)
336: {
338: BV_SVEC *ctx = (BV_SVEC*)bv->data;
339: PetscScalar *pv,*pnew;
340: PetscInt bs;
341: Vec vnew;
342: char str[50];
345: VecGetBlockSize(bv->t,&bs);
346: VecCreate(PetscObjectComm((PetscObject)bv->t),&vnew);
347: VecSetType(vnew,((PetscObject)bv->t)->type_name);
348: VecSetSizes(vnew,m*bv->n,PETSC_DECIDE);
349: VecSetBlockSize(vnew,bs);
350: PetscLogObjectParent((PetscObject)bv,(PetscObject)vnew);
351: if (((PetscObject)bv)->name) {
352: PetscSNPrintf(str,50,"%s_0",((PetscObject)bv)->name);
353: PetscObjectSetName((PetscObject)vnew,str);
354: }
355: if (copy) {
356: VecGetArray(ctx->v,&pv);
357: VecGetArray(vnew,&pnew);
358: PetscMemcpy(pnew,pv,PetscMin(m,bv->m)*bv->n*sizeof(PetscScalar));
359: VecRestoreArray(ctx->v,&pv);
360: VecRestoreArray(vnew,&pnew);
361: }
362: VecDestroy(&ctx->v);
363: ctx->v = vnew;
364: return(0);
365: }
369: PetscErrorCode BVGetColumn_Svec(BV bv,PetscInt j,Vec *v)
370: {
372: BV_SVEC *ctx = (BV_SVEC*)bv->data;
373: PetscScalar *pv;
374: PetscInt l;
377: l = BVAvailableVec;
378: VecGetArray(ctx->v,&pv);
379: VecPlaceArray(bv->cv[l],pv+(bv->nc+j)*bv->n);
380: return(0);
381: }
385: PetscErrorCode BVRestoreColumn_Svec(BV bv,PetscInt j,Vec *v)
386: {
388: BV_SVEC *ctx = (BV_SVEC*)bv->data;
389: PetscInt l;
392: l = (j==bv->ci[0])? 0: 1;
393: VecResetArray(bv->cv[l]);
394: VecRestoreArray(ctx->v,NULL);
395: return(0);
396: }
400: PetscErrorCode BVGetArray_Svec(BV bv,PetscScalar **a)
401: {
403: BV_SVEC *ctx = (BV_SVEC*)bv->data;
406: VecGetArray(ctx->v,a);
407: return(0);
408: }
412: PetscErrorCode BVRestoreArray_Svec(BV bv,PetscScalar **a)
413: {
415: BV_SVEC *ctx = (BV_SVEC*)bv->data;
418: VecRestoreArray(ctx->v,a);
419: return(0);
420: }
424: PetscErrorCode BVView_Svec(BV bv,PetscViewer viewer)
425: {
426: PetscErrorCode ierr;
427: BV_SVEC *ctx = (BV_SVEC*)bv->data;
428: PetscViewerFormat format;
429: PetscBool isascii;
430: const char *bvname,*name;
433: VecView(ctx->v,viewer);
434: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
435: if (isascii) {
436: PetscViewerGetFormat(viewer,&format);
437: if (format == PETSC_VIEWER_ASCII_MATLAB) {
438: PetscObjectGetName((PetscObject)bv,&bvname);
439: PetscObjectGetName((PetscObject)ctx->v,&name);
440: PetscViewerASCIIPrintf(viewer,"%s=reshape(%s,%D,%D);clear %s\n",bvname,name,bv->N,bv->nc+bv->m,name);
441: if (bv->nc) {
442: PetscViewerASCIIPrintf(viewer,"%s=%s(:,%D:end);\n",bvname,bvname,bv->nc+1);
443: }
444: }
445: }
446: return(0);
447: }
451: PetscErrorCode BVDestroy_Svec(BV bv)
452: {
454: BV_SVEC *ctx = (BV_SVEC*)bv->data;
457: VecDestroy(&ctx->v);
458: VecDestroy(&bv->cv[0]);
459: VecDestroy(&bv->cv[1]);
460: PetscFree(bv->data);
461: return(0);
462: }
466: PETSC_EXTERN PetscErrorCode BVCreate_Svec(BV bv)
467: {
469: BV_SVEC *ctx;
470: PetscInt nloc,bs;
471: PetscBool seq;
472: char str[50];
475: PetscNewLog(bv,&ctx);
476: bv->data = (void*)ctx;
478: PetscObjectTypeCompare((PetscObject)bv->t,VECMPI,&ctx->mpi);
479: if (!ctx->mpi) {
480: PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
481: if (!seq) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot create a BVSVEC from a non-standard template vector");
482: }
484: VecGetLocalSize(bv->t,&nloc);
485: VecGetBlockSize(bv->t,&bs);
487: VecCreate(PetscObjectComm((PetscObject)bv->t),&ctx->v);
488: VecSetType(ctx->v,((PetscObject)bv->t)->type_name);
489: VecSetSizes(ctx->v,bv->m*nloc,PETSC_DECIDE);
490: VecSetBlockSize(ctx->v,bs);
491: PetscLogObjectParent((PetscObject)bv,(PetscObject)ctx->v);
492: if (((PetscObject)bv)->name) {
493: PetscSNPrintf(str,50,"%s_0",((PetscObject)bv)->name);
494: PetscObjectSetName((PetscObject)ctx->v,str);
495: }
497: if (ctx->mpi) {
498: VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,NULL,&bv->cv[0]);
499: VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,PETSC_DECIDE,NULL,&bv->cv[1]);
500: } else {
501: VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,NULL,&bv->cv[0]);
502: VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,nloc,NULL,&bv->cv[1]);
503: }
505: bv->ops->mult = BVMult_Svec;
506: bv->ops->multvec = BVMultVec_Svec;
507: bv->ops->multinplace = BVMultInPlace_Svec;
508: bv->ops->multinplacetrans = BVMultInPlaceTranspose_Svec;
509: bv->ops->axpy = BVAXPY_Svec;
510: bv->ops->dot = BVDot_Svec;
511: bv->ops->dotvec = BVDotVec_Svec;
512: bv->ops->dotvec_local = BVDotVec_Local_Svec;
513: bv->ops->scale = BVScale_Svec;
514: bv->ops->norm = BVNorm_Svec;
515: bv->ops->norm_local = BVNorm_Local_Svec;
516: /*bv->ops->orthogonalize = BVOrthogonalize_Svec;*/
517: bv->ops->matmult = BVMatMult_Svec;
518: bv->ops->copy = BVCopy_Svec;
519: bv->ops->resize = BVResize_Svec;
520: bv->ops->getcolumn = BVGetColumn_Svec;
521: bv->ops->restorecolumn = BVRestoreColumn_Svec;
522: bv->ops->getarray = BVGetArray_Svec;
523: bv->ops->restorearray = BVRestoreArray_Svec;
524: bv->ops->destroy = BVDestroy_Svec;
525: if (!ctx->mpi) bv->ops->view = BVView_Svec;
526: return(0);
527: }