Actual source code: vecs.c
slepc-3.5.2 2014-10-10
1: /*
2: BV implemented as an array of independent Vecs
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: } BV_VECS;
32: PetscErrorCode BVMult_Vecs(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
33: {
35: BV_VECS *y = (BV_VECS*)Y->data,*x = (BV_VECS*)X->data;
36: PetscScalar *q,*s=NULL;
37: PetscInt i,j,ldq;
40: MatGetSize(Q,&ldq,NULL);
41: if (alpha!=1.0) {
42: BVAllocateWork_Private(Y,X->k-X->l);
43: s = Y->work;
44: }
45: MatDenseGetArray(Q,&q);
46: for (j=Y->l;j<Y->k;j++) {
47: VecScale(y->V[Y->nc+j],beta);
48: if (alpha!=1.0) {
49: for (i=X->l;i<X->k;i++) s[i-X->l] = alpha*q[i+j*ldq];
50: } else s = q+j*ldq+X->l;
51: VecMAXPY(y->V[Y->nc+j],X->k-X->l,s,x->V+X->nc+X->l);
52: }
53: MatDenseRestoreArray(Q,&q);
54: return(0);
55: }
59: PetscErrorCode BVMultVec_Vecs(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
60: {
62: BV_VECS *x = (BV_VECS*)X->data;
63: PetscScalar *s=NULL;
64: PetscInt i;
67: if (alpha!=1.0) {
68: BVAllocateWork_Private(X,X->k-X->l);
69: s = X->work;
70: }
71: VecScale(y,beta);
72: if (alpha!=1.0) {
73: for (i=0;i<X->k-X->l;i++) s[i] = alpha*q[i];
74: } else s = q;
75: VecMAXPY(y,X->k-X->l,s,x->V+X->nc+X->l);
76: return(0);
77: }
81: /*
82: BVMultInPlace_Vecs - V(:,s:e-1) = V*Q(:,s:e-1) for regular vectors.
84: Writing V = [ V1 V2 V3 ] and Q(:,s:e-1) = [ Q1 Q2 Q3 ]', where V2
85: corresponds to the columns s:e-1, the computation is done as
86: V2 := V2*Q2 + V1*Q1 + V3*Q3
87: */
88: PetscErrorCode BVMultInPlace_Vecs(BV V,Mat Q,PetscInt s,PetscInt e)
89: {
91: BV_VECS *ctx = (BV_VECS*)V->data;
92: PetscScalar *q;
93: PetscInt i,ldq;
96: MatGetSize(Q,&ldq,NULL);
97: MatDenseGetArray(Q,&q);
98: /* V2 := V2*Q2 */
99: BVMultInPlace_Vecs_Private(V,V->n,e-s,V->k,ctx->V+V->nc+s,q+s*ldq+s,PETSC_FALSE);
100: /* V2 += V1*Q1 + V3*Q3 */
101: for (i=s;i<e;i++) {
102: if (s>V->l) {
103: VecMAXPY(ctx->V[V->nc+i],s-V->l,q+i*ldq+V->l,ctx->V+V->nc+V->l);
104: }
105: if (V->k>e) {
106: VecMAXPY(ctx->V[V->nc+i],V->k-e,q+i*ldq+e,ctx->V+V->nc+e);
107: }
108: }
109: MatDenseRestoreArray(Q,&q);
110: return(0);
111: }
115: /*
116: BVMultInPlaceTranspose_Vecs - V(:,s:e-1) = V*Q'(:,s:e-1) for regular vectors.
117: */
118: PetscErrorCode BVMultInPlaceTranspose_Vecs(BV V,Mat Q,PetscInt s,PetscInt e)
119: {
121: BV_VECS *ctx = (BV_VECS*)V->data;
122: PetscScalar *q;
123: PetscInt i,j,ldq,n;
126: MatGetSize(Q,&ldq,&n);
127: MatDenseGetArray(Q,&q);
128: /* V2 := V2*Q2' */
129: BVMultInPlace_Vecs_Private(V,V->n,e-s,ldq,ctx->V+V->nc+s,q+s*ldq+s,PETSC_TRUE);
130: /* V2 += V1*Q1' + V3*Q3' */
131: for (i=s;i<e;i++) {
132: for (j=V->l;j<s;j++) {
133: VecAXPY(ctx->V[V->nc+i],q[i+j*ldq],ctx->V[V->nc+j]);
134: }
135: for (j=e;j<n;j++) {
136: VecAXPY(ctx->V[V->nc+i],q[i+j*ldq],ctx->V[V->nc+j]);
137: }
138: }
139: MatDenseRestoreArray(Q,&q);
140: return(0);
141: }
145: PetscErrorCode BVAXPY_Vecs(BV Y,PetscScalar alpha,BV X)
146: {
148: BV_VECS *y = (BV_VECS*)Y->data,*x = (BV_VECS*)X->data;
149: PetscInt j;
152: for (j=0;j<Y->k-Y->l;j++) {
153: VecAXPY(y->V[Y->nc+Y->l+j],alpha,x->V[X->nc+X->l+j]);
154: }
155: return(0);
156: }
160: PetscErrorCode BVDot_Vecs(BV X,BV Y,Mat M)
161: {
163: BV_VECS *x = (BV_VECS*)X->data,*y = (BV_VECS*)Y->data;
164: PetscScalar *m;
165: PetscInt j,ldm;
168: MatGetSize(M,&ldm,NULL);
169: MatDenseGetArray(M,&m);
170: for (j=X->l;j<X->k;j++) {
171: VecMDot(x->V[X->nc+j],Y->k-Y->l,y->V+Y->nc+Y->l,m+j*ldm+Y->l);
172: }
173: MatDenseRestoreArray(M,&m);
174: return(0);
175: }
179: PetscErrorCode BVDotVec_Vecs(BV X,Vec y,PetscScalar *m)
180: {
182: BV_VECS *x = (BV_VECS*)X->data;
183: Vec z = y;
186: if (X->matrix) {
187: BV_IPMatMult(X,y);
188: z = X->Bx;
189: }
190: VecMDot(z,X->k-X->l,x->V+X->nc+X->l,m);
191: return(0);
192: }
196: PetscErrorCode BVScale_Vecs(BV bv,PetscInt j,PetscScalar alpha)
197: {
199: PetscInt i;
200: BV_VECS *ctx = (BV_VECS*)bv->data;
203: if (j<0) {
204: for (i=bv->l;i<bv->k;i++) {
205: VecScale(ctx->V[bv->nc+i],alpha);
206: }
207: } else {
208: VecScale(ctx->V[bv->nc+j],alpha);
209: }
210: return(0);
211: }
215: PetscErrorCode BVNorm_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val)
216: {
218: PetscInt i;
219: PetscReal nrm;
220: BV_VECS *ctx = (BV_VECS*)bv->data;
223: if (j<0) {
224: switch (type) {
225: case NORM_FROBENIUS:
226: *val = 0.0;
227: for (i=bv->l;i<bv->k;i++) {
228: VecNorm(ctx->V[bv->nc+i],NORM_2,&nrm);
229: *val += nrm*nrm;
230: }
231: *val = PetscSqrtReal(*val);
232: break;
233: default:
234: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Requested norm not implemented in BVVECS");
235: }
236: } else {
237: VecNorm(ctx->V[bv->nc+j],type,val);
238: }
239: return(0);
240: }
244: PetscErrorCode BVMatMult_Vecs(BV V,Mat A,BV W)
245: {
247: BV_VECS *v = (BV_VECS*)V->data,*w = (BV_VECS*)W->data;
248: PetscInt j;
251: for (j=0;j<V->k-V->l;j++) {
252: MatMult(A,v->V[V->nc+V->l+j],w->V[W->nc+W->l+j]);
253: }
254: return(0);
255: }
259: PetscErrorCode BVCopy_Vecs(BV V,BV W)
260: {
262: BV_VECS *v = (BV_VECS*)V->data,*w = (BV_VECS*)W->data;
263: PetscInt j;
266: for (j=0;j<V->k-V->l;j++) {
267: VecCopy(v->V[V->nc+V->l+j],w->V[W->nc+W->l+j]);
268: }
269: return(0);
270: }
274: PetscErrorCode BVResize_Vecs(BV bv,PetscInt m,PetscBool copy)
275: {
277: BV_VECS *ctx = (BV_VECS*)bv->data;
278: Vec *newV;
279: PetscInt j;
280: char str[50];
283: VecDuplicateVecs(bv->t,m,&newV);
284: PetscLogObjectParents(bv,m,newV);
285: if (((PetscObject)bv)->name) {
286: for (j=0;j<m;j++) {
287: PetscSNPrintf(str,50,"%s_%D",((PetscObject)bv)->name,j);
288: PetscObjectSetName((PetscObject)newV[j],str);
289: }
290: }
291: if (copy) {
292: for (j=0;j<PetscMin(m,bv->m);j++) {
293: VecCopy(ctx->V[j],newV[j]);
294: }
295: }
296: VecDestroyVecs(bv->m,&ctx->V);
297: ctx->V = newV;
298: return(0);
299: }
303: PetscErrorCode BVGetColumn_Vecs(BV bv,PetscInt j,Vec *v)
304: {
305: BV_VECS *ctx = (BV_VECS*)bv->data;
306: PetscInt l;
309: l = BVAvailableVec;
310: bv->cv[l] = ctx->V[bv->nc+j];
311: return(0);
312: }
316: PetscErrorCode BVGetArray_Vecs(BV bv,PetscScalar **a)
317: {
319: BV_VECS *ctx = (BV_VECS*)bv->data;
320: PetscInt j;
321: PetscScalar *p;
324: PetscMalloc((bv->nc+bv->m)*bv->n*sizeof(PetscScalar),a);
325: for (j=0;j<bv->nc+bv->m;j++) {
326: VecGetArray(ctx->V[j],&p);
327: PetscMemcpy(*a+j*bv->n,p,bv->n*sizeof(PetscScalar));
328: VecRestoreArray(ctx->V[j],&p);
329: }
330: return(0);
331: }
335: PetscErrorCode BVRestoreArray_Vecs(BV bv,PetscScalar **a)
336: {
338: BV_VECS *ctx = (BV_VECS*)bv->data;
339: PetscInt j;
340: PetscScalar *p;
343: for (j=0;j<bv->nc+bv->m;j++) {
344: VecGetArray(ctx->V[j],&p);
345: PetscMemcpy(p,*a+j*bv->n,bv->n*sizeof(PetscScalar));
346: VecRestoreArray(ctx->V[j],&p);
347: }
348: PetscFree(*a);
349: return(0);
350: }
354: PetscErrorCode BVView_Vecs(BV bv,PetscViewer viewer)
355: {
356: PetscErrorCode ierr;
357: BV_VECS *ctx = (BV_VECS*)bv->data;
358: PetscInt j;
359: PetscViewerFormat format;
360: PetscBool isascii,ismatlab=PETSC_FALSE;
363: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
364: if (isascii) {
365: PetscViewerGetFormat(viewer,&format);
366: if (format == PETSC_VIEWER_ASCII_MATLAB) ismatlab = PETSC_TRUE;
367: }
368: if (ismatlab) {
369: PetscViewerASCIIPrintf(viewer,"%s=[];\n",((PetscObject)bv)->name);
370: }
371: for (j=bv->nc;j<bv->nc+bv->m;j++) {
372: VecView(ctx->V[j],viewer);
373: if (ismatlab) {
374: PetscViewerASCIIPrintf(viewer,"%s=[%s,%s];clear %s\n",((PetscObject)bv)->name,((PetscObject)bv)->name,((PetscObject)ctx->V[j])->name,((PetscObject)ctx->V[j])->name);
375: }
376: }
377: return(0);
378: }
382: PetscErrorCode BVDestroy_Vecs(BV bv)
383: {
385: BV_VECS *ctx = (BV_VECS*)bv->data;
388: VecDestroyVecs(bv->nc+bv->m,&ctx->V);
389: PetscFree(bv->data);
390: return(0);
391: }
395: PETSC_EXTERN PetscErrorCode BVCreate_Vecs(BV bv)
396: {
398: BV_VECS *ctx;
399: PetscInt j;
400: char str[50];
403: PetscNewLog(bv,&ctx);
404: bv->data = (void*)ctx;
406: VecDuplicateVecs(bv->t,bv->m,&ctx->V);
407: PetscLogObjectParents(bv,bv->m,ctx->V);
408: if (((PetscObject)bv)->name) {
409: for (j=0;j<bv->m;j++) {
410: PetscSNPrintf(str,50,"%s_%D",((PetscObject)bv)->name,j);
411: PetscObjectSetName((PetscObject)ctx->V[j],str);
412: }
413: }
415: bv->ops->mult = BVMult_Vecs;
416: bv->ops->multvec = BVMultVec_Vecs;
417: bv->ops->multinplace = BVMultInPlace_Vecs;
418: bv->ops->multinplacetrans = BVMultInPlaceTranspose_Vecs;
419: bv->ops->axpy = BVAXPY_Vecs;
420: bv->ops->dot = BVDot_Vecs;
421: bv->ops->dotvec = BVDotVec_Vecs;
422: bv->ops->scale = BVScale_Vecs;
423: bv->ops->norm = BVNorm_Vecs;
424: bv->ops->matmult = BVMatMult_Vecs;
425: bv->ops->copy = BVCopy_Vecs;
426: bv->ops->resize = BVResize_Vecs;
427: bv->ops->getcolumn = BVGetColumn_Vecs;
428: bv->ops->getarray = BVGetArray_Vecs;
429: bv->ops->restorearray = BVRestoreArray_Vecs;
430: bv->ops->destroy = BVDestroy_Vecs;
431: bv->ops->view = BVView_Vecs;
432: return(0);
433: }