Actual source code: vecs.c
slepc-3.7.4 2017-05-17
1: /*
2: BV implemented as an array of independent Vecs
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: PetscInt vmip; /* Version of BVMultInPlace:
29: 0: memory-efficient version, uses VecGetArray (default in CPU)
30: 1: version that allocates (e-s) work vectors in every call (default in GPU) */
31: } BV_VECS;
35: PetscErrorCode BVMult_Vecs(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
36: {
38: BV_VECS *y = (BV_VECS*)Y->data,*x = (BV_VECS*)X->data;
39: PetscScalar *q,*s=NULL;
40: PetscInt i,j,ldq;
43: if (Q) {
44: MatGetSize(Q,&ldq,NULL);
45: if (alpha!=1.0) {
46: BVAllocateWork_Private(Y,X->k-X->l);
47: s = Y->work;
48: }
49: MatDenseGetArray(Q,&q);
50: for (j=Y->l;j<Y->k;j++) {
51: VecScale(y->V[Y->nc+j],beta);
52: if (alpha!=1.0) {
53: for (i=X->l;i<X->k;i++) s[i-X->l] = alpha*q[i+j*ldq];
54: } else s = q+j*ldq+X->l;
55: VecMAXPY(y->V[Y->nc+j],X->k-X->l,s,x->V+X->nc+X->l);
56: }
57: MatDenseRestoreArray(Q,&q);
58: } else {
59: for (j=0;j<Y->k-Y->l;j++) {
60: VecScale(y->V[Y->nc+Y->l+j],beta);
61: VecAXPY(y->V[Y->nc+Y->l+j],alpha,x->V[X->nc+X->l+j]);
62: }
63: }
64: return(0);
65: }
69: PetscErrorCode BVMultVec_Vecs(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
70: {
72: BV_VECS *x = (BV_VECS*)X->data;
73: PetscScalar *s=NULL;
74: PetscInt i;
77: if (alpha!=1.0) {
78: BVAllocateWork_Private(X,X->k-X->l);
79: s = X->work;
80: }
81: VecScale(y,beta);
82: if (alpha!=1.0) {
83: for (i=0;i<X->k-X->l;i++) s[i] = alpha*q[i];
84: } else s = q;
85: VecMAXPY(y,X->k-X->l,s,x->V+X->nc+X->l);
86: return(0);
87: }
91: /*
92: BVMultInPlace_Vecs_ME - V(:,s:e-1) = V*Q(:,s:e-1) for regular vectors.
94: Memory-efficient version, uses VecGetArray (default in CPU)
96: Writing V = [ V1 V2 V3 ] and Q(:,s:e-1) = [ Q1 Q2 Q3 ]', where V2
97: corresponds to the columns s:e-1, the computation is done as
98: V2 := V2*Q2 + V1*Q1 + V3*Q3
99: */
100: PetscErrorCode BVMultInPlace_Vecs_ME(BV V,Mat Q,PetscInt s,PetscInt e)
101: {
103: BV_VECS *ctx = (BV_VECS*)V->data;
104: PetscScalar *q;
105: PetscInt i,ldq;
108: MatGetSize(Q,&ldq,NULL);
109: MatDenseGetArray(Q,&q);
110: /* V2 := V2*Q2 */
111: BVMultInPlace_Vecs_Private(V,V->n,e-s,ldq,ctx->V+V->nc+s,q+s*ldq+s,PETSC_FALSE);
112: /* V2 += V1*Q1 + V3*Q3 */
113: for (i=s;i<e;i++) {
114: if (s>V->l) {
115: VecMAXPY(ctx->V[V->nc+i],s-V->l,q+i*ldq+V->l,ctx->V+V->nc+V->l);
116: }
117: if (V->k>e) {
118: VecMAXPY(ctx->V[V->nc+i],V->k-e,q+i*ldq+e,ctx->V+V->nc+e);
119: }
120: }
121: MatDenseRestoreArray(Q,&q);
122: return(0);
123: }
127: /*
128: BVMultInPlace_Vecs_Alloc - V(:,s:e-1) = V*Q(:,s:e-1) for regular vectors.
130: Version that allocates (e-s) work vectors in every call (default in GPU)
131: */
132: PetscErrorCode BVMultInPlace_Vecs_Alloc(BV V,Mat Q,PetscInt s,PetscInt e)
133: {
135: BV_VECS *ctx = (BV_VECS*)V->data;
136: PetscScalar *q;
137: PetscInt i,ldq;
138: Vec *W;
141: MatGetSize(Q,&ldq,NULL);
142: MatDenseGetArray(Q,&q);
143: VecDuplicateVecs(V->t,e-s,&W);
144: for (i=s;i<e;i++) {
145: VecMAXPY(W[i-s],V->k-V->l,q+i*ldq+V->l,ctx->V+V->nc+V->l);
146: }
147: for (i=s;i<e;i++) {
148: VecCopy(W[i-s],ctx->V[V->nc+i]);
149: }
150: VecDestroyVecs(e-s,&W);
151: MatDenseRestoreArray(Q,&q);
152: return(0);
153: }
157: /*
158: BVMultInPlaceTranspose_Vecs - V(:,s:e-1) = V*Q'(:,s:e-1) for regular vectors.
159: */
160: PetscErrorCode BVMultInPlaceTranspose_Vecs(BV V,Mat Q,PetscInt s,PetscInt e)
161: {
163: BV_VECS *ctx = (BV_VECS*)V->data;
164: PetscScalar *q;
165: PetscInt i,j,ldq,n;
168: MatGetSize(Q,&ldq,&n);
169: MatDenseGetArray(Q,&q);
170: /* V2 := V2*Q2' */
171: BVMultInPlace_Vecs_Private(V,V->n,e-s,ldq,ctx->V+V->nc+s,q+s*ldq+s,PETSC_TRUE);
172: /* V2 += V1*Q1' + V3*Q3' */
173: for (i=s;i<e;i++) {
174: for (j=V->l;j<s;j++) {
175: VecAXPY(ctx->V[V->nc+i],q[i+j*ldq],ctx->V[V->nc+j]);
176: }
177: for (j=e;j<n;j++) {
178: VecAXPY(ctx->V[V->nc+i],q[i+j*ldq],ctx->V[V->nc+j]);
179: }
180: }
181: MatDenseRestoreArray(Q,&q);
182: return(0);
183: }
187: PetscErrorCode BVDot_Vecs(BV X,BV Y,Mat M)
188: {
190: BV_VECS *x = (BV_VECS*)X->data,*y = (BV_VECS*)Y->data;
191: PetscScalar *m;
192: PetscInt j,ldm;
195: MatGetSize(M,&ldm,NULL);
196: MatDenseGetArray(M,&m);
197: for (j=X->l;j<X->k;j++) {
198: VecMDot(x->V[X->nc+j],Y->k-Y->l,y->V+Y->nc+Y->l,m+j*ldm+Y->l);
199: }
200: MatDenseRestoreArray(M,&m);
201: return(0);
202: }
206: PetscErrorCode BVDotVec_Vecs(BV X,Vec y,PetscScalar *m)
207: {
209: BV_VECS *x = (BV_VECS*)X->data;
210: Vec z = y;
213: if (X->matrix) {
214: BV_IPMatMult(X,y);
215: z = X->Bx;
216: }
217: VecMDot(z,X->k-X->l,x->V+X->nc+X->l,m);
218: return(0);
219: }
223: PetscErrorCode BVDotVec_Begin_Vecs(BV X,Vec y,PetscScalar *m)
224: {
226: BV_VECS *x = (BV_VECS*)X->data;
227: Vec z = y;
230: if (X->matrix) {
231: BV_IPMatMult(X,y);
232: z = X->Bx;
233: }
234: VecMDotBegin(z,X->k-X->l,x->V+X->nc+X->l,m);
235: return(0);
236: }
240: PetscErrorCode BVDotVec_End_Vecs(BV X,Vec y,PetscScalar *m)
241: {
243: BV_VECS *x = (BV_VECS*)X->data;
246: VecMDotEnd(y,X->k-X->l,x->V+X->nc+X->l,m);
247: return(0);
248: }
252: PetscErrorCode BVScale_Vecs(BV bv,PetscInt j,PetscScalar alpha)
253: {
255: PetscInt i;
256: BV_VECS *ctx = (BV_VECS*)bv->data;
259: if (j<0) {
260: for (i=bv->l;i<bv->k;i++) {
261: VecScale(ctx->V[bv->nc+i],alpha);
262: }
263: } else {
264: VecScale(ctx->V[bv->nc+j],alpha);
265: }
266: return(0);
267: }
271: PetscErrorCode BVNorm_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val)
272: {
274: PetscInt i;
275: PetscReal nrm;
276: BV_VECS *ctx = (BV_VECS*)bv->data;
279: if (j<0) {
280: switch (type) {
281: case NORM_FROBENIUS:
282: *val = 0.0;
283: for (i=bv->l;i<bv->k;i++) {
284: VecNorm(ctx->V[bv->nc+i],NORM_2,&nrm);
285: *val += nrm*nrm;
286: }
287: *val = PetscSqrtReal(*val);
288: break;
289: default:
290: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Requested norm not implemented in BVVECS");
291: }
292: } else {
293: VecNorm(ctx->V[bv->nc+j],type,val);
294: }
295: return(0);
296: }
300: PetscErrorCode BVNorm_Begin_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val)
301: {
303: BV_VECS *ctx = (BV_VECS*)bv->data;
306: if (j<0) {
307: switch (type) {
308: default:
309: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Requested norm not implemented in BVVECS");
310: }
311: } else {
312: VecNormBegin(ctx->V[bv->nc+j],type,val);
313: }
314: return(0);
315: }
319: PetscErrorCode BVNorm_End_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val)
320: {
322: BV_VECS *ctx = (BV_VECS*)bv->data;
325: if (j<0) {
326: switch (type) {
327: default:
328: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Requested norm not implemented in BVVECS");
329: }
330: } else {
331: VecNormEnd(ctx->V[bv->nc+j],type,val);
332: }
333: return(0);
334: }
338: PetscErrorCode BVMatMult_Vecs(BV V,Mat A,BV W)
339: {
341: BV_VECS *v = (BV_VECS*)V->data,*w = (BV_VECS*)W->data;
342: PetscInt j;
345: if (V->vmm) { PetscInfo(V,"BVMatMult_Vecs: ignoring method\n"); }
346: for (j=0;j<V->k-V->l;j++) {
347: MatMult(A,v->V[V->nc+V->l+j],w->V[W->nc+W->l+j]);
348: }
349: return(0);
350: }
354: PetscErrorCode BVCopy_Vecs(BV V,BV W)
355: {
357: BV_VECS *v = (BV_VECS*)V->data,*w = (BV_VECS*)W->data;
358: PetscInt j;
361: for (j=0;j<V->k-V->l;j++) {
362: VecCopy(v->V[V->nc+V->l+j],w->V[W->nc+W->l+j]);
363: }
364: return(0);
365: }
369: PetscErrorCode BVResize_Vecs(BV bv,PetscInt m,PetscBool copy)
370: {
372: BV_VECS *ctx = (BV_VECS*)bv->data;
373: Vec *newV;
374: PetscInt j;
375: char str[50];
378: VecDuplicateVecs(bv->t,m,&newV);
379: PetscLogObjectParents(bv,m,newV);
380: if (((PetscObject)bv)->name) {
381: for (j=0;j<m;j++) {
382: PetscSNPrintf(str,50,"%s_%D",((PetscObject)bv)->name,j);
383: PetscObjectSetName((PetscObject)newV[j],str);
384: }
385: }
386: if (copy) {
387: for (j=0;j<PetscMin(m,bv->m);j++) {
388: VecCopy(ctx->V[j],newV[j]);
389: }
390: }
391: VecDestroyVecs(bv->m,&ctx->V);
392: ctx->V = newV;
393: return(0);
394: }
398: PetscErrorCode BVGetColumn_Vecs(BV bv,PetscInt j,Vec *v)
399: {
400: BV_VECS *ctx = (BV_VECS*)bv->data;
401: PetscInt l;
404: l = BVAvailableVec;
405: bv->cv[l] = ctx->V[bv->nc+j];
406: return(0);
407: }
411: PetscErrorCode BVGetArray_Vecs(BV bv,PetscScalar **a)
412: {
413: PetscErrorCode ierr;
414: BV_VECS *ctx = (BV_VECS*)bv->data;
415: PetscInt j;
416: const PetscScalar *p;
419: PetscMalloc1((bv->nc+bv->m)*bv->n,a);
420: for (j=0;j<bv->nc+bv->m;j++) {
421: VecGetArrayRead(ctx->V[j],&p);
422: PetscMemcpy(*a+j*bv->n,p,bv->n*sizeof(PetscScalar));
423: VecRestoreArrayRead(ctx->V[j],&p);
424: }
425: return(0);
426: }
430: PetscErrorCode BVRestoreArray_Vecs(BV bv,PetscScalar **a)
431: {
433: BV_VECS *ctx = (BV_VECS*)bv->data;
434: PetscInt j;
435: PetscScalar *p;
438: for (j=0;j<bv->nc+bv->m;j++) {
439: VecGetArray(ctx->V[j],&p);
440: PetscMemcpy(p,*a+j*bv->n,bv->n*sizeof(PetscScalar));
441: VecRestoreArray(ctx->V[j],&p);
442: }
443: PetscFree(*a);
444: return(0);
445: }
449: PetscErrorCode BVGetArrayRead_Vecs(BV bv,const PetscScalar **a)
450: {
451: PetscErrorCode ierr;
452: BV_VECS *ctx = (BV_VECS*)bv->data;
453: PetscInt j;
454: const PetscScalar *p;
457: PetscMalloc1((bv->nc+bv->m)*bv->n,a);
458: for (j=0;j<bv->nc+bv->m;j++) {
459: VecGetArrayRead(ctx->V[j],&p);
460: PetscMemcpy((PetscScalar*)*a+j*bv->n,p,bv->n*sizeof(PetscScalar));
461: VecRestoreArrayRead(ctx->V[j],&p);
462: }
463: return(0);
464: }
468: PetscErrorCode BVRestoreArrayRead_Vecs(BV bv,const PetscScalar **a)
469: {
473: PetscFree(*a);
474: return(0);
475: }
479: /*
480: Sets the value of vmip flag and resets ops->multinplace accordingly
481: */
482: PETSC_STATIC_INLINE PetscErrorCode BVVecsSetVmip(BV bv,PetscInt vmip)
483: {
484: typedef PetscErrorCode (*fmultinplace)(BV,Mat,PetscInt,PetscInt);
485: fmultinplace multinplace[2] = {BVMultInPlace_Vecs_ME, BVMultInPlace_Vecs_Alloc};
486: BV_VECS *ctx = (BV_VECS*)bv->data;
489: ctx->vmip = vmip;
490: bv->ops->multinplace = multinplace[vmip];
491: return(0);
492: }
496: PetscErrorCode BVSetFromOptions_Vecs(PetscOptionItems *PetscOptionsObject,BV bv)
497: {
499: BV_VECS *ctx = (BV_VECS*)bv->data;
502: PetscOptionsHead(PetscOptionsObject,"BV Vecs Options");
503: PetscOptionsInt("-bv_vecs_vmip","Version of BVMultInPlace operation","",ctx->vmip,&ctx->vmip,NULL);
504: if (ctx->vmip<0 || ctx->vmip>1) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Wrong version of BVMultInPlace");
505: BVVecsSetVmip(bv,ctx->vmip);
506: PetscOptionsTail();
507: return(0);
508: }
512: PetscErrorCode BVView_Vecs(BV bv,PetscViewer viewer)
513: {
514: PetscErrorCode ierr;
515: BV_VECS *ctx = (BV_VECS*)bv->data;
516: PetscInt j;
517: PetscViewerFormat format;
518: PetscBool isascii,ismatlab=PETSC_FALSE;
519: const char *bvname,*name;
522: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
523: if (isascii) {
524: PetscViewerGetFormat(viewer,&format);
525: if (format == PETSC_VIEWER_ASCII_MATLAB) ismatlab = PETSC_TRUE;
526: }
527: if (ismatlab) {
528: PetscObjectGetName((PetscObject)bv,&bvname);
529: PetscViewerASCIIPrintf(viewer,"%s=[];\n",bvname);
530: }
531: for (j=bv->nc;j<bv->nc+bv->m;j++) {
532: VecView(ctx->V[j],viewer);
533: if (ismatlab) {
534: PetscObjectGetName((PetscObject)ctx->V[j],&name);
535: PetscViewerASCIIPrintf(viewer,"%s=[%s,%s];clear %s\n",bvname,bvname,name,name);
536: }
537: }
538: return(0);
539: }
543: PetscErrorCode BVDestroy_Vecs(BV bv)
544: {
546: BV_VECS *ctx = (BV_VECS*)bv->data;
549: VecDestroyVecs(bv->nc+bv->m,&ctx->V);
550: PetscFree(bv->data);
551: return(0);
552: }
556: PetscErrorCode BVDuplicate_Vecs(BV V,BV *W)
557: {
559: BV_VECS *ctx = (BV_VECS*)V->data;
562: BVVecsSetVmip(*W,ctx->vmip);
563: return(0);
564: }
568: PETSC_EXTERN PetscErrorCode BVCreate_Vecs(BV bv)
569: {
571: BV_VECS *ctx;
572: PetscInt j;
573: PetscBool iscusp;
574: char str[50];
577: PetscNewLog(bv,&ctx);
578: bv->data = (void*)ctx;
580: VecDuplicateVecs(bv->t,bv->m,&ctx->V);
581: PetscLogObjectParents(bv,bv->m,ctx->V);
582: if (((PetscObject)bv)->name) {
583: for (j=0;j<bv->m;j++) {
584: PetscSNPrintf(str,50,"%s_%D",((PetscObject)bv)->name,j);
585: PetscObjectSetName((PetscObject)ctx->V[j],str);
586: }
587: }
589: /* Default version of BVMultInPlace */
590: PetscObjectTypeCompareAny((PetscObject)bv->t,&iscusp,VECSEQCUSP,VECMPICUSP,"");
591: ctx->vmip = iscusp? 1: 0;
593: /* Deferred call to setfromoptions */
594: if (bv->defersfo) {
595: PetscObjectOptionsBegin((PetscObject)bv);
596: BVSetFromOptions_Vecs(PetscOptionsObject,bv);
597: PetscOptionsEnd();
598: }
599: BVVecsSetVmip(bv,ctx->vmip);
601: bv->ops->mult = BVMult_Vecs;
602: bv->ops->multvec = BVMultVec_Vecs;
603: bv->ops->multinplacetrans = BVMultInPlaceTranspose_Vecs;
604: bv->ops->dot = BVDot_Vecs;
605: bv->ops->dotvec = BVDotVec_Vecs;
606: bv->ops->dotvec_begin = BVDotVec_Begin_Vecs;
607: bv->ops->dotvec_end = BVDotVec_End_Vecs;
608: bv->ops->scale = BVScale_Vecs;
609: bv->ops->norm = BVNorm_Vecs;
610: bv->ops->norm_begin = BVNorm_Begin_Vecs;
611: bv->ops->norm_end = BVNorm_End_Vecs;
612: bv->ops->matmult = BVMatMult_Vecs;
613: bv->ops->copy = BVCopy_Vecs;
614: bv->ops->resize = BVResize_Vecs;
615: bv->ops->getcolumn = BVGetColumn_Vecs;
616: bv->ops->getarray = BVGetArray_Vecs;
617: bv->ops->restorearray = BVRestoreArray_Vecs;
618: bv->ops->getarrayread = BVGetArrayRead_Vecs;
619: bv->ops->restorearrayread = BVRestoreArrayRead_Vecs;
620: bv->ops->destroy = BVDestroy_Vecs;
621: bv->ops->duplicate = BVDuplicate_Vecs;
622: bv->ops->setfromoptions = BVSetFromOptions_Vecs;
623: bv->ops->view = BVView_Vecs;
624: return(0);
625: }