Actual source code: vecs.c
slepc-3.8.2 2017-12-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2017, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: BV implemented as an array of independent Vecs
12: */
14: #include <slepc/private/bvimpl.h>
16: typedef struct {
17: Vec *V;
18: PetscInt vmip; /* Version of BVMultInPlace:
19: 0: memory-efficient version, uses VecGetArray (default in CPU)
20: 1: version that allocates (e-s) work vectors in every call (default in GPU) */
21: } BV_VECS;
23: PetscErrorCode BVMult_Vecs(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
24: {
26: BV_VECS *y = (BV_VECS*)Y->data,*x = (BV_VECS*)X->data;
27: PetscScalar *q,*s=NULL;
28: PetscInt i,j,ldq;
31: if (Q) {
32: MatGetSize(Q,&ldq,NULL);
33: if (alpha!=1.0) {
34: BVAllocateWork_Private(Y,X->k-X->l);
35: s = Y->work;
36: }
37: MatDenseGetArray(Q,&q);
38: for (j=Y->l;j<Y->k;j++) {
39: VecScale(y->V[Y->nc+j],beta);
40: if (alpha!=1.0) {
41: for (i=X->l;i<X->k;i++) s[i-X->l] = alpha*q[i+j*ldq];
42: } else s = q+j*ldq+X->l;
43: VecMAXPY(y->V[Y->nc+j],X->k-X->l,s,x->V+X->nc+X->l);
44: }
45: MatDenseRestoreArray(Q,&q);
46: } else {
47: for (j=0;j<Y->k-Y->l;j++) {
48: VecScale(y->V[Y->nc+Y->l+j],beta);
49: VecAXPY(y->V[Y->nc+Y->l+j],alpha,x->V[X->nc+X->l+j]);
50: }
51: }
52: return(0);
53: }
55: PetscErrorCode BVMultVec_Vecs(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
56: {
58: BV_VECS *x = (BV_VECS*)X->data;
59: PetscScalar *s=NULL,*qq=q;
60: PetscInt i;
63: if (alpha!=1.0) {
64: BVAllocateWork_Private(X,X->k-X->l);
65: s = X->work;
66: }
67: if (!q) { VecGetArray(X->buffer,&qq); }
68: VecScale(y,beta);
69: if (alpha!=1.0) {
70: for (i=0;i<X->k-X->l;i++) s[i] = alpha*qq[i];
71: } else s = qq;
72: VecMAXPY(y,X->k-X->l,s,x->V+X->nc+X->l);
73: if (!q) { VecRestoreArray(X->buffer,&qq); }
74: return(0);
75: }
77: /*
78: BVMultInPlace_Vecs_ME - V(:,s:e-1) = V*Q(:,s:e-1) for regular vectors.
80: Memory-efficient version, uses VecGetArray (default in CPU)
82: Writing V = [ V1 V2 V3 ] and Q(:,s:e-1) = [ Q1 Q2 Q3 ]', where V2
83: corresponds to the columns s:e-1, the computation is done as
84: V2 := V2*Q2 + V1*Q1 + V3*Q3
85: */
86: PetscErrorCode BVMultInPlace_Vecs_ME(BV V,Mat Q,PetscInt s,PetscInt e)
87: {
89: BV_VECS *ctx = (BV_VECS*)V->data;
90: PetscScalar *q;
91: PetscInt i,ldq;
94: MatGetSize(Q,&ldq,NULL);
95: MatDenseGetArray(Q,&q);
96: /* V2 := V2*Q2 */
97: BVMultInPlace_Vecs_Private(V,V->n,e-s,ldq,ctx->V+V->nc+s,q+s*ldq+s,PETSC_FALSE);
98: /* V2 += V1*Q1 + V3*Q3 */
99: for (i=s;i<e;i++) {
100: if (s>V->l) {
101: VecMAXPY(ctx->V[V->nc+i],s-V->l,q+i*ldq+V->l,ctx->V+V->nc+V->l);
102: }
103: if (V->k>e) {
104: VecMAXPY(ctx->V[V->nc+i],V->k-e,q+i*ldq+e,ctx->V+V->nc+e);
105: }
106: }
107: MatDenseRestoreArray(Q,&q);
108: return(0);
109: }
111: /*
112: BVMultInPlace_Vecs_Alloc - V(:,s:e-1) = V*Q(:,s:e-1) for regular vectors.
114: Version that allocates (e-s) work vectors in every call (default in GPU)
115: */
116: PetscErrorCode BVMultInPlace_Vecs_Alloc(BV V,Mat Q,PetscInt s,PetscInt e)
117: {
119: BV_VECS *ctx = (BV_VECS*)V->data;
120: PetscScalar *q;
121: PetscInt i,ldq;
122: Vec *W;
125: MatGetSize(Q,&ldq,NULL);
126: MatDenseGetArray(Q,&q);
127: VecDuplicateVecs(V->t,e-s,&W);
128: for (i=s;i<e;i++) {
129: VecMAXPY(W[i-s],V->k-V->l,q+i*ldq+V->l,ctx->V+V->nc+V->l);
130: }
131: for (i=s;i<e;i++) {
132: VecCopy(W[i-s],ctx->V[V->nc+i]);
133: }
134: VecDestroyVecs(e-s,&W);
135: MatDenseRestoreArray(Q,&q);
136: return(0);
137: }
139: /*
140: BVMultInPlaceTranspose_Vecs - V(:,s:e-1) = V*Q'(:,s:e-1) for regular vectors.
141: */
142: PetscErrorCode BVMultInPlaceTranspose_Vecs(BV V,Mat Q,PetscInt s,PetscInt e)
143: {
145: BV_VECS *ctx = (BV_VECS*)V->data;
146: PetscScalar *q;
147: PetscInt i,j,ldq,n;
150: MatGetSize(Q,&ldq,&n);
151: MatDenseGetArray(Q,&q);
152: /* V2 := V2*Q2' */
153: BVMultInPlace_Vecs_Private(V,V->n,e-s,ldq,ctx->V+V->nc+s,q+s*ldq+s,PETSC_TRUE);
154: /* V2 += V1*Q1' + V3*Q3' */
155: for (i=s;i<e;i++) {
156: for (j=V->l;j<s;j++) {
157: VecAXPY(ctx->V[V->nc+i],q[i+j*ldq],ctx->V[V->nc+j]);
158: }
159: for (j=e;j<n;j++) {
160: VecAXPY(ctx->V[V->nc+i],q[i+j*ldq],ctx->V[V->nc+j]);
161: }
162: }
163: MatDenseRestoreArray(Q,&q);
164: return(0);
165: }
167: PetscErrorCode BVDot_Vecs(BV X,BV Y,Mat M)
168: {
170: BV_VECS *x = (BV_VECS*)X->data,*y = (BV_VECS*)Y->data;
171: PetscScalar *m;
172: PetscInt j,ldm;
175: MatGetSize(M,&ldm,NULL);
176: MatDenseGetArray(M,&m);
177: for (j=X->l;j<X->k;j++) {
178: VecMDot(x->V[X->nc+j],Y->k-Y->l,y->V+Y->nc+Y->l,m+j*ldm+Y->l);
179: }
180: MatDenseRestoreArray(M,&m);
181: return(0);
182: }
184: PetscErrorCode BVDotVec_Vecs(BV X,Vec y,PetscScalar *q)
185: {
187: BV_VECS *x = (BV_VECS*)X->data;
188: Vec z = y;
189: PetscScalar *qq=q;
192: if (X->matrix) {
193: BV_IPMatMult(X,y);
194: z = X->Bx;
195: }
196: if (!q) { VecGetArray(X->buffer,&qq); }
197: VecMDot(z,X->k-X->l,x->V+X->nc+X->l,qq);
198: if (!q) { VecRestoreArray(X->buffer,&qq); }
199: return(0);
200: }
202: PetscErrorCode BVDotVec_Begin_Vecs(BV X,Vec y,PetscScalar *m)
203: {
205: BV_VECS *x = (BV_VECS*)X->data;
206: Vec z = y;
209: if (X->matrix) {
210: BV_IPMatMult(X,y);
211: z = X->Bx;
212: }
213: VecMDotBegin(z,X->k-X->l,x->V+X->nc+X->l,m);
214: return(0);
215: }
217: PetscErrorCode BVDotVec_End_Vecs(BV X,Vec y,PetscScalar *m)
218: {
220: BV_VECS *x = (BV_VECS*)X->data;
223: VecMDotEnd(y,X->k-X->l,x->V+X->nc+X->l,m);
224: return(0);
225: }
227: PetscErrorCode BVScale_Vecs(BV bv,PetscInt j,PetscScalar alpha)
228: {
230: PetscInt i;
231: BV_VECS *ctx = (BV_VECS*)bv->data;
234: if (j<0) {
235: for (i=bv->l;i<bv->k;i++) {
236: VecScale(ctx->V[bv->nc+i],alpha);
237: }
238: } else {
239: VecScale(ctx->V[bv->nc+j],alpha);
240: }
241: return(0);
242: }
244: PetscErrorCode BVNorm_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val)
245: {
247: PetscInt i;
248: PetscReal nrm;
249: BV_VECS *ctx = (BV_VECS*)bv->data;
252: if (j<0) {
253: if (type==NORM_FROBENIUS) {
254: *val = 0.0;
255: for (i=bv->l;i<bv->k;i++) {
256: VecNorm(ctx->V[bv->nc+i],NORM_2,&nrm);
257: *val += nrm*nrm;
258: }
259: *val = PetscSqrtReal(*val);
260: } else SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not implemented in BVVECS");
261: } else {
262: VecNorm(ctx->V[bv->nc+j],type,val);
263: }
264: return(0);
265: }
267: PetscErrorCode BVNorm_Begin_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val)
268: {
270: BV_VECS *ctx = (BV_VECS*)bv->data;
273: if (j<0) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not implemented in BVVECS");
274: else {
275: VecNormBegin(ctx->V[bv->nc+j],type,val);
276: }
277: return(0);
278: }
280: PetscErrorCode BVNorm_End_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val)
281: {
283: BV_VECS *ctx = (BV_VECS*)bv->data;
286: if (j<0) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not implemented in BVVECS");
287: else {
288: VecNormEnd(ctx->V[bv->nc+j],type,val);
289: }
290: return(0);
291: }
293: PetscErrorCode BVMatMult_Vecs(BV V,Mat A,BV W)
294: {
296: BV_VECS *v = (BV_VECS*)V->data,*w = (BV_VECS*)W->data;
297: PetscInt j;
300: if (V->vmm) { PetscInfo(V,"BVMatMult_Vecs: ignoring method\n"); }
301: for (j=0;j<V->k-V->l;j++) {
302: MatMult(A,v->V[V->nc+V->l+j],w->V[W->nc+W->l+j]);
303: }
304: return(0);
305: }
307: PetscErrorCode BVCopy_Vecs(BV V,BV W)
308: {
310: BV_VECS *v = (BV_VECS*)V->data,*w = (BV_VECS*)W->data;
311: PetscInt j;
314: for (j=0;j<V->k-V->l;j++) {
315: VecCopy(v->V[V->nc+V->l+j],w->V[W->nc+W->l+j]);
316: }
317: return(0);
318: }
320: PetscErrorCode BVResize_Vecs(BV bv,PetscInt m,PetscBool copy)
321: {
323: BV_VECS *ctx = (BV_VECS*)bv->data;
324: Vec *newV;
325: PetscInt j;
326: char str[50];
329: VecDuplicateVecs(bv->t,m,&newV);
330: PetscLogObjectParents(bv,m,newV);
331: if (((PetscObject)bv)->name) {
332: for (j=0;j<m;j++) {
333: PetscSNPrintf(str,50,"%s_%D",((PetscObject)bv)->name,j);
334: PetscObjectSetName((PetscObject)newV[j],str);
335: }
336: }
337: if (copy) {
338: for (j=0;j<PetscMin(m,bv->m);j++) {
339: VecCopy(ctx->V[j],newV[j]);
340: }
341: }
342: VecDestroyVecs(bv->m,&ctx->V);
343: ctx->V = newV;
344: return(0);
345: }
347: PetscErrorCode BVGetColumn_Vecs(BV bv,PetscInt j,Vec *v)
348: {
349: BV_VECS *ctx = (BV_VECS*)bv->data;
350: PetscInt l;
353: l = BVAvailableVec;
354: bv->cv[l] = ctx->V[bv->nc+j];
355: return(0);
356: }
358: PetscErrorCode BVGetArray_Vecs(BV bv,PetscScalar **a)
359: {
360: PetscErrorCode ierr;
361: BV_VECS *ctx = (BV_VECS*)bv->data;
362: PetscInt j;
363: const PetscScalar *p;
366: PetscMalloc1((bv->nc+bv->m)*bv->n,a);
367: for (j=0;j<bv->nc+bv->m;j++) {
368: VecGetArrayRead(ctx->V[j],&p);
369: PetscMemcpy(*a+j*bv->n,p,bv->n*sizeof(PetscScalar));
370: VecRestoreArrayRead(ctx->V[j],&p);
371: }
372: return(0);
373: }
375: PetscErrorCode BVRestoreArray_Vecs(BV bv,PetscScalar **a)
376: {
378: BV_VECS *ctx = (BV_VECS*)bv->data;
379: PetscInt j;
380: PetscScalar *p;
383: for (j=0;j<bv->nc+bv->m;j++) {
384: VecGetArray(ctx->V[j],&p);
385: PetscMemcpy(p,*a+j*bv->n,bv->n*sizeof(PetscScalar));
386: VecRestoreArray(ctx->V[j],&p);
387: }
388: PetscFree(*a);
389: return(0);
390: }
392: PetscErrorCode BVGetArrayRead_Vecs(BV bv,const PetscScalar **a)
393: {
394: PetscErrorCode ierr;
395: BV_VECS *ctx = (BV_VECS*)bv->data;
396: PetscInt j;
397: const PetscScalar *p;
400: PetscMalloc1((bv->nc+bv->m)*bv->n,a);
401: for (j=0;j<bv->nc+bv->m;j++) {
402: VecGetArrayRead(ctx->V[j],&p);
403: PetscMemcpy((PetscScalar*)*a+j*bv->n,p,bv->n*sizeof(PetscScalar));
404: VecRestoreArrayRead(ctx->V[j],&p);
405: }
406: return(0);
407: }
409: PetscErrorCode BVRestoreArrayRead_Vecs(BV bv,const PetscScalar **a)
410: {
414: PetscFree(*a);
415: return(0);
416: }
418: /*
419: Sets the value of vmip flag and resets ops->multinplace accordingly
420: */
421: PETSC_STATIC_INLINE PetscErrorCode BVVecsSetVmip(BV bv,PetscInt vmip)
422: {
423: typedef PetscErrorCode (*fmultinplace)(BV,Mat,PetscInt,PetscInt);
424: fmultinplace multinplace[2] = {BVMultInPlace_Vecs_ME, BVMultInPlace_Vecs_Alloc};
425: BV_VECS *ctx = (BV_VECS*)bv->data;
428: ctx->vmip = vmip;
429: bv->ops->multinplace = multinplace[vmip];
430: return(0);
431: }
433: PetscErrorCode BVSetFromOptions_Vecs(PetscOptionItems *PetscOptionsObject,BV bv)
434: {
436: BV_VECS *ctx = (BV_VECS*)bv->data;
439: PetscOptionsHead(PetscOptionsObject,"BV Vecs Options");
441: PetscOptionsInt("-bv_vecs_vmip","Version of BVMultInPlace operation","",ctx->vmip,&ctx->vmip,NULL);
442: if (ctx->vmip<0 || ctx->vmip>1) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Wrong version of BVMultInPlace");
443: BVVecsSetVmip(bv,ctx->vmip);
445: PetscOptionsTail();
446: return(0);
447: }
449: PetscErrorCode BVView_Vecs(BV bv,PetscViewer viewer)
450: {
451: PetscErrorCode ierr;
452: BV_VECS *ctx = (BV_VECS*)bv->data;
453: PetscInt j;
454: PetscViewerFormat format;
455: PetscBool isascii,ismatlab=PETSC_FALSE;
456: const char *bvname,*name;
459: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
460: if (isascii) {
461: PetscViewerGetFormat(viewer,&format);
462: if (format == PETSC_VIEWER_ASCII_MATLAB) ismatlab = PETSC_TRUE;
463: }
464: if (ismatlab) {
465: PetscObjectGetName((PetscObject)bv,&bvname);
466: PetscViewerASCIIPrintf(viewer,"%s=[];\n",bvname);
467: }
468: for (j=bv->nc;j<bv->nc+bv->m;j++) {
469: VecView(ctx->V[j],viewer);
470: if (ismatlab) {
471: PetscObjectGetName((PetscObject)ctx->V[j],&name);
472: PetscViewerASCIIPrintf(viewer,"%s=[%s,%s];clear %s\n",bvname,bvname,name,name);
473: }
474: }
475: return(0);
476: }
478: PetscErrorCode BVDestroy_Vecs(BV bv)
479: {
481: BV_VECS *ctx = (BV_VECS*)bv->data;
484: VecDestroyVecs(bv->nc+bv->m,&ctx->V);
485: PetscFree(bv->data);
486: return(0);
487: }
489: PetscErrorCode BVDuplicate_Vecs(BV V,BV *W)
490: {
492: BV_VECS *ctx = (BV_VECS*)V->data;
495: BVVecsSetVmip(*W,ctx->vmip);
496: return(0);
497: }
499: PETSC_EXTERN PetscErrorCode BVCreate_Vecs(BV bv)
500: {
502: BV_VECS *ctx;
503: PetscInt j;
504: PetscBool isgpu;
505: char str[50];
508: PetscNewLog(bv,&ctx);
509: bv->data = (void*)ctx;
511: VecDuplicateVecs(bv->t,bv->m,&ctx->V);
512: PetscLogObjectParents(bv,bv->m,ctx->V);
513: if (((PetscObject)bv)->name) {
514: for (j=0;j<bv->m;j++) {
515: PetscSNPrintf(str,50,"%s_%D",((PetscObject)bv)->name,j);
516: PetscObjectSetName((PetscObject)ctx->V[j],str);
517: }
518: }
520: if (bv->Acreate) {
521: for (j=0;j<bv->m;j++) {
522: MatGetColumnVector(bv->Acreate,ctx->V[j],j);
523: }
524: MatDestroy(&bv->Acreate);
525: }
527: /* Default version of BVMultInPlace */
528: PetscObjectTypeCompareAny((PetscObject)bv->t,&isgpu,VECSEQCUDA,VECMPICUDA,VECSEQCUSP,VECMPICUSP,"");
529: ctx->vmip = isgpu? 1: 0;
531: /* Deferred call to setfromoptions */
532: if (bv->defersfo) {
533: PetscObjectOptionsBegin((PetscObject)bv);
534: BVSetFromOptions_Vecs(PetscOptionsObject,bv);
535: PetscOptionsEnd();
536: }
537: BVVecsSetVmip(bv,ctx->vmip);
539: bv->ops->mult = BVMult_Vecs;
540: bv->ops->multvec = BVMultVec_Vecs;
541: bv->ops->multinplacetrans = BVMultInPlaceTranspose_Vecs;
542: bv->ops->dot = BVDot_Vecs;
543: bv->ops->dotvec = BVDotVec_Vecs;
544: bv->ops->dotvec_begin = BVDotVec_Begin_Vecs;
545: bv->ops->dotvec_end = BVDotVec_End_Vecs;
546: bv->ops->scale = BVScale_Vecs;
547: bv->ops->norm = BVNorm_Vecs;
548: bv->ops->norm_begin = BVNorm_Begin_Vecs;
549: bv->ops->norm_end = BVNorm_End_Vecs;
550: bv->ops->matmult = BVMatMult_Vecs;
551: bv->ops->copy = BVCopy_Vecs;
552: bv->ops->resize = BVResize_Vecs;
553: bv->ops->getcolumn = BVGetColumn_Vecs;
554: bv->ops->getarray = BVGetArray_Vecs;
555: bv->ops->restorearray = BVRestoreArray_Vecs;
556: bv->ops->getarrayread = BVGetArrayRead_Vecs;
557: bv->ops->restorearrayread = BVRestoreArrayRead_Vecs;
558: bv->ops->destroy = BVDestroy_Vecs;
559: bv->ops->duplicate = BVDuplicate_Vecs;
560: bv->ops->setfromoptions = BVSetFromOptions_Vecs;
561: bv->ops->view = BVView_Vecs;
562: return(0);
563: }