Actual source code: vecs.c
slepc-3.17.0 2022-03-31
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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: {
25: BV_VECS *y = (BV_VECS*)Y->data,*x = (BV_VECS*)X->data;
26: PetscScalar *s=NULL;
27: const PetscScalar *q;
28: PetscInt i,j,ldq;
29: PetscBool trivial=(alpha==1.0)?PETSC_TRUE:PETSC_FALSE;
31: if (Q) {
32: MatDenseGetLDA(Q,&ldq);
33: if (!trivial) {
34: BVAllocateWork_Private(Y,X->k-X->l);
35: s = Y->work;
36: }
37: MatDenseGetArrayRead(Q,&q);
38: for (j=Y->l;j<Y->k;j++) {
39: VecScale(y->V[Y->nc+j],beta);
40: if (!trivial) {
41: for (i=X->l;i<X->k;i++) s[i-X->l] = alpha*q[i+j*ldq];
42: } else s = (PetscScalar*)(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: MatDenseRestoreArrayRead(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: PetscFunctionReturn(0);
53: }
55: PetscErrorCode BVMultVec_Vecs(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
56: {
57: BV_VECS *x = (BV_VECS*)X->data;
58: PetscScalar *s=NULL,*qq=q;
59: PetscInt i;
60: PetscBool trivial=(alpha==1.0)?PETSC_TRUE:PETSC_FALSE;
62: if (!trivial) {
63: BVAllocateWork_Private(X,X->k-X->l);
64: s = X->work;
65: }
66: if (!q) VecGetArray(X->buffer,&qq);
67: VecScale(y,beta);
68: if (!trivial) {
69: for (i=0;i<X->k-X->l;i++) s[i] = alpha*qq[i];
70: } else s = qq;
71: VecMAXPY(y,X->k-X->l,s,x->V+X->nc+X->l);
72: if (!q) VecRestoreArray(X->buffer,&qq);
73: PetscFunctionReturn(0);
74: }
76: /*
77: BVMultInPlace_Vecs_ME - V(:,s:e-1) = V*Q(:,s:e-1) for regular vectors.
79: Memory-efficient version, uses VecGetArray (default in CPU)
81: Writing V = [ V1 V2 V3 ] and Q(:,s:e-1) = [ Q1 Q2 Q3 ]', where V2
82: corresponds to the columns s:e-1, the computation is done as
83: V2 := V2*Q2 + V1*Q1 + V3*Q3
84: */
85: PetscErrorCode BVMultInPlace_Vecs_ME(BV V,Mat Q,PetscInt s,PetscInt e)
86: {
87: BV_VECS *ctx = (BV_VECS*)V->data;
88: const PetscScalar *q;
89: PetscInt i,ldq;
91: MatDenseGetLDA(Q,&ldq);
92: MatDenseGetArrayRead(Q,&q);
93: /* V2 := V2*Q2 */
94: BVMultInPlace_Vecs_Private(V,V->n,e-s,ldq,ctx->V+V->nc+s,q+s*ldq+s,PETSC_FALSE);
95: /* V2 += V1*Q1 + V3*Q3 */
96: for (i=s;i<e;i++) {
97: if (PetscUnlikely(s>V->l)) VecMAXPY(ctx->V[V->nc+i],s-V->l,q+i*ldq+V->l,ctx->V+V->nc+V->l);
98: if (V->k>e) VecMAXPY(ctx->V[V->nc+i],V->k-e,q+i*ldq+e,ctx->V+V->nc+e);
99: }
100: MatDenseRestoreArrayRead(Q,&q);
101: PetscFunctionReturn(0);
102: }
104: /*
105: BVMultInPlace_Vecs_Alloc - V(:,s:e-1) = V*Q(:,s:e-1) for regular vectors.
107: Version that allocates (e-s) work vectors in every call (default in GPU)
108: */
109: PetscErrorCode BVMultInPlace_Vecs_Alloc(BV V,Mat Q,PetscInt s,PetscInt e)
110: {
111: BV_VECS *ctx = (BV_VECS*)V->data;
112: const PetscScalar *q;
113: PetscInt i,ldq;
114: Vec *W;
116: MatDenseGetLDA(Q,&ldq);
117: MatDenseGetArrayRead(Q,&q);
118: VecDuplicateVecs(V->t,e-s,&W);
119: for (i=s;i<e;i++) VecMAXPY(W[i-s],V->k-V->l,q+i*ldq+V->l,ctx->V+V->nc+V->l);
120: for (i=s;i<e;i++) VecCopy(W[i-s],ctx->V[V->nc+i]);
121: VecDestroyVecs(e-s,&W);
122: MatDenseRestoreArrayRead(Q,&q);
123: PetscFunctionReturn(0);
124: }
126: /*
127: BVMultInPlaceHermitianTranspose_Vecs - V(:,s:e-1) = V*Q'(:,s:e-1) for regular vectors.
128: */
129: PetscErrorCode BVMultInPlaceHermitianTranspose_Vecs(BV V,Mat Q,PetscInt s,PetscInt e)
130: {
131: BV_VECS *ctx = (BV_VECS*)V->data;
132: const PetscScalar *q;
133: PetscInt i,j,ldq,n;
135: MatGetSize(Q,NULL,&n);
136: MatDenseGetLDA(Q,&ldq);
137: MatDenseGetArrayRead(Q,&q);
138: /* V2 := V2*Q2' */
139: BVMultInPlace_Vecs_Private(V,V->n,e-s,ldq,ctx->V+V->nc+s,q+s*ldq+s,PETSC_TRUE);
140: /* V2 += V1*Q1' + V3*Q3' */
141: for (i=s;i<e;i++) {
142: for (j=V->l;j<s;j++) VecAXPY(ctx->V[V->nc+i],q[i+j*ldq],ctx->V[V->nc+j]);
143: for (j=e;j<n;j++) VecAXPY(ctx->V[V->nc+i],q[i+j*ldq],ctx->V[V->nc+j]);
144: }
145: MatDenseRestoreArrayRead(Q,&q);
146: PetscFunctionReturn(0);
147: }
149: PetscErrorCode BVDot_Vecs(BV X,BV Y,Mat M)
150: {
151: BV_VECS *x = (BV_VECS*)X->data,*y = (BV_VECS*)Y->data;
152: PetscScalar *m;
153: PetscInt j,ldm;
155: MatDenseGetLDA(M,&ldm);
156: MatDenseGetArray(M,&m);
157: for (j=X->l;j<X->k;j++) VecMDot(x->V[X->nc+j],Y->k-Y->l,y->V+Y->nc+Y->l,m+j*ldm+Y->l);
158: MatDenseRestoreArray(M,&m);
159: PetscFunctionReturn(0);
160: }
162: PetscErrorCode BVDotVec_Vecs(BV X,Vec y,PetscScalar *q)
163: {
164: BV_VECS *x = (BV_VECS*)X->data;
165: Vec z = y;
166: PetscScalar *qq=q;
168: if (PetscUnlikely(X->matrix)) {
169: BV_IPMatMult(X,y);
170: z = X->Bx;
171: }
172: if (!q) VecGetArray(X->buffer,&qq);
173: VecMDot(z,X->k-X->l,x->V+X->nc+X->l,qq);
174: if (!q) VecRestoreArray(X->buffer,&qq);
175: PetscFunctionReturn(0);
176: }
178: PetscErrorCode BVDotVec_Begin_Vecs(BV X,Vec y,PetscScalar *m)
179: {
180: BV_VECS *x = (BV_VECS*)X->data;
181: Vec z = y;
183: if (PetscUnlikely(X->matrix)) {
184: BV_IPMatMult(X,y);
185: z = X->Bx;
186: }
187: VecMDotBegin(z,X->k-X->l,x->V+X->nc+X->l,m);
188: PetscFunctionReturn(0);
189: }
191: PetscErrorCode BVDotVec_End_Vecs(BV X,Vec y,PetscScalar *m)
192: {
193: BV_VECS *x = (BV_VECS*)X->data;
195: VecMDotEnd(y,X->k-X->l,x->V+X->nc+X->l,m);
196: PetscFunctionReturn(0);
197: }
199: PetscErrorCode BVScale_Vecs(BV bv,PetscInt j,PetscScalar alpha)
200: {
201: PetscInt i;
202: BV_VECS *ctx = (BV_VECS*)bv->data;
204: if (PetscUnlikely(j<0)) {
205: for (i=bv->l;i<bv->k;i++) VecScale(ctx->V[bv->nc+i],alpha);
206: } else VecScale(ctx->V[bv->nc+j],alpha);
207: PetscFunctionReturn(0);
208: }
210: PetscErrorCode BVNorm_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val)
211: {
212: PetscInt i;
213: PetscReal nrm;
214: BV_VECS *ctx = (BV_VECS*)bv->data;
216: if (PetscUnlikely(j<0)) {
218: *val = 0.0;
219: for (i=bv->l;i<bv->k;i++) {
220: VecNorm(ctx->V[bv->nc+i],NORM_2,&nrm);
221: *val += nrm*nrm;
222: }
223: *val = PetscSqrtReal(*val);
224: } else VecNorm(ctx->V[bv->nc+j],type,val);
225: PetscFunctionReturn(0);
226: }
228: PetscErrorCode BVNorm_Begin_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val)
229: {
230: BV_VECS *ctx = (BV_VECS*)bv->data;
233: VecNormBegin(ctx->V[bv->nc+j],type,val);
234: PetscFunctionReturn(0);
235: }
237: PetscErrorCode BVNorm_End_Vecs(BV bv,PetscInt j,NormType type,PetscReal *val)
238: {
239: BV_VECS *ctx = (BV_VECS*)bv->data;
242: VecNormEnd(ctx->V[bv->nc+j],type,val);
243: PetscFunctionReturn(0);
244: }
246: PetscErrorCode BVMatMult_Vecs(BV V,Mat A,BV W)
247: {
248: BV_VECS *v = (BV_VECS*)V->data,*w = (BV_VECS*)W->data;
249: PetscInt j;
250: Mat Vmat,Wmat;
252: if (V->vmm) {
253: BVGetMat(V,&Vmat);
254: BVGetMat(W,&Wmat);
255: MatProductCreateWithMat(A,Vmat,NULL,Wmat);
256: MatProductSetType(Wmat,MATPRODUCT_AB);
257: MatProductSetFromOptions(Wmat);
258: MatProductSymbolic(Wmat);
259: MatProductNumeric(Wmat);
260: MatProductClear(Wmat);
261: BVRestoreMat(V,&Vmat);
262: BVRestoreMat(W,&Wmat);
263: } else {
264: for (j=0;j<V->k-V->l;j++) MatMult(A,v->V[V->nc+V->l+j],w->V[W->nc+W->l+j]);
265: }
266: PetscFunctionReturn(0);
267: }
269: PetscErrorCode BVCopy_Vecs(BV V,BV W)
270: {
271: BV_VECS *v = (BV_VECS*)V->data,*w = (BV_VECS*)W->data;
272: PetscInt j;
274: for (j=0;j<V->k-V->l;j++) VecCopy(v->V[V->nc+V->l+j],w->V[W->nc+W->l+j]);
275: PetscFunctionReturn(0);
276: }
278: PetscErrorCode BVCopyColumn_Vecs(BV V,PetscInt j,PetscInt i)
279: {
280: BV_VECS *v = (BV_VECS*)V->data;
282: VecCopy(v->V[V->nc+j],v->V[V->nc+i]);
283: PetscFunctionReturn(0);
284: }
286: PetscErrorCode BVResize_Vecs(BV bv,PetscInt m,PetscBool copy)
287: {
288: BV_VECS *ctx = (BV_VECS*)bv->data;
289: Vec *newV;
290: PetscInt j;
291: char str[50];
293: VecDuplicateVecs(bv->t,m,&newV);
294: PetscLogObjectParents(bv,m,newV);
295: if (((PetscObject)bv)->name) {
296: for (j=0;j<m;j++) {
297: PetscSNPrintf(str,sizeof(str),"%s_%" PetscInt_FMT,((PetscObject)bv)->name,j);
298: PetscObjectSetName((PetscObject)newV[j],str);
299: }
300: }
301: if (copy) {
302: for (j=0;j<PetscMin(m,bv->m);j++) VecCopy(ctx->V[j],newV[j]);
303: }
304: VecDestroyVecs(bv->m,&ctx->V);
305: ctx->V = newV;
306: PetscFunctionReturn(0);
307: }
309: PetscErrorCode BVGetColumn_Vecs(BV bv,PetscInt j,Vec *v)
310: {
311: BV_VECS *ctx = (BV_VECS*)bv->data;
312: PetscInt l;
314: l = BVAvailableVec;
315: bv->cv[l] = ctx->V[bv->nc+j];
316: PetscFunctionReturn(0);
317: }
319: PetscErrorCode BVGetArray_Vecs(BV bv,PetscScalar **a)
320: {
321: BV_VECS *ctx = (BV_VECS*)bv->data;
322: PetscInt j;
323: const PetscScalar *p;
325: PetscMalloc1((bv->nc+bv->m)*bv->n,a);
326: for (j=0;j<bv->nc+bv->m;j++) {
327: VecGetArrayRead(ctx->V[j],&p);
328: PetscArraycpy(*a+j*bv->n,p,bv->n);
329: VecRestoreArrayRead(ctx->V[j],&p);
330: }
331: PetscFunctionReturn(0);
332: }
334: PetscErrorCode BVRestoreArray_Vecs(BV bv,PetscScalar **a)
335: {
336: BV_VECS *ctx = (BV_VECS*)bv->data;
337: PetscInt j;
338: PetscScalar *p;
340: for (j=0;j<bv->nc+bv->m;j++) {
341: VecGetArray(ctx->V[j],&p);
342: PetscArraycpy(p,*a+j*bv->n,bv->n);
343: VecRestoreArray(ctx->V[j],&p);
344: }
345: PetscFree(*a);
346: PetscFunctionReturn(0);
347: }
349: PetscErrorCode BVGetArrayRead_Vecs(BV bv,const PetscScalar **a)
350: {
351: BV_VECS *ctx = (BV_VECS*)bv->data;
352: PetscInt j;
353: const PetscScalar *p;
355: PetscMalloc1((bv->nc+bv->m)*bv->n,(PetscScalar**)a);
356: for (j=0;j<bv->nc+bv->m;j++) {
357: VecGetArrayRead(ctx->V[j],&p);
358: PetscArraycpy((PetscScalar*)*a+j*bv->n,p,bv->n);
359: VecRestoreArrayRead(ctx->V[j],&p);
360: }
361: PetscFunctionReturn(0);
362: }
364: PetscErrorCode BVRestoreArrayRead_Vecs(BV bv,const PetscScalar **a)
365: {
366: PetscFree(*a);
367: PetscFunctionReturn(0);
368: }
370: /*
371: Sets the value of vmip flag and resets ops->multinplace accordingly
372: */
373: static inline PetscErrorCode BVVecsSetVmip(BV bv,PetscInt vmip)
374: {
375: typedef PetscErrorCode (*fmultinplace)(BV,Mat,PetscInt,PetscInt);
376: fmultinplace multinplace[2] = {BVMultInPlace_Vecs_ME, BVMultInPlace_Vecs_Alloc};
377: BV_VECS *ctx = (BV_VECS*)bv->data;
379: ctx->vmip = vmip;
380: bv->ops->multinplace = multinplace[vmip];
381: PetscFunctionReturn(0);
382: }
384: PetscErrorCode BVSetFromOptions_Vecs(PetscOptionItems *PetscOptionsObject,BV bv)
385: {
386: BV_VECS *ctx = (BV_VECS*)bv->data;
388: PetscOptionsHead(PetscOptionsObject,"BV Vecs Options");
390: PetscOptionsRangeInt("-bv_vecs_vmip","Version of BVMultInPlace operation","",ctx->vmip,&ctx->vmip,NULL,0,1);
391: BVVecsSetVmip(bv,ctx->vmip);
393: PetscOptionsTail();
394: PetscFunctionReturn(0);
395: }
397: PetscErrorCode BVView_Vecs(BV bv,PetscViewer viewer)
398: {
399: BV_VECS *ctx = (BV_VECS*)bv->data;
400: PetscInt j;
401: PetscViewerFormat format;
402: PetscBool isascii,ismatlab=PETSC_FALSE;
403: const char *bvname,*name;
405: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
406: if (isascii) {
407: PetscViewerGetFormat(viewer,&format);
408: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0);
409: if (format == PETSC_VIEWER_ASCII_MATLAB) ismatlab = PETSC_TRUE;
410: }
411: if (ismatlab) {
412: PetscObjectGetName((PetscObject)bv,&bvname);
413: PetscViewerASCIIPrintf(viewer,"%s=[];\n",bvname);
414: }
415: for (j=bv->nc;j<bv->nc+bv->m;j++) {
416: VecView(ctx->V[j],viewer);
417: if (ismatlab) {
418: PetscObjectGetName((PetscObject)ctx->V[j],&name);
419: PetscViewerASCIIPrintf(viewer,"%s=[%s,%s];clear %s\n",bvname,bvname,name,name);
420: }
421: }
422: PetscFunctionReturn(0);
423: }
425: PetscErrorCode BVDestroy_Vecs(BV bv)
426: {
427: BV_VECS *ctx = (BV_VECS*)bv->data;
429: if (!bv->issplit) VecDestroyVecs(bv->nc+bv->m,&ctx->V);
430: PetscFree(bv->data);
431: PetscFunctionReturn(0);
432: }
434: PetscErrorCode BVDuplicate_Vecs(BV V,BV W)
435: {
436: BV_VECS *ctx = (BV_VECS*)V->data;
438: BVVecsSetVmip(W,ctx->vmip);
439: PetscFunctionReturn(0);
440: }
442: SLEPC_EXTERN PetscErrorCode BVCreate_Vecs(BV bv)
443: {
445: BV_VECS *ctx;
446: PetscInt j,lsplit;
447: PetscBool isgpu;
448: char str[50];
449: BV parent;
450: Vec *Vpar;
452: PetscNewLog(bv,&ctx);
453: bv->data = (void*)ctx;
455: if (PetscUnlikely(bv->issplit)) {
456: /* split BV: share the Vecs of the parent BV */
457: parent = bv->splitparent;
458: lsplit = parent->lsplit;
459: Vpar = ((BV_VECS*)parent->data)->V;
460: ctx->V = (bv->issplit==1)? Vpar: Vpar+lsplit;
461: } else {
462: /* regular BV: create array of Vecs to store the BV columns */
463: VecDuplicateVecs(bv->t,bv->m,&ctx->V);
464: PetscLogObjectParents(bv,bv->m,ctx->V);
465: if (((PetscObject)bv)->name) {
466: for (j=0;j<bv->m;j++) {
467: PetscSNPrintf(str,sizeof(str),"%s_%" PetscInt_FMT,((PetscObject)bv)->name,j);
468: PetscObjectSetName((PetscObject)ctx->V[j],str);
469: }
470: }
471: }
473: if (PetscUnlikely(bv->Acreate)) {
474: for (j=0;j<bv->m;j++) MatGetColumnVector(bv->Acreate,ctx->V[j],j);
475: MatDestroy(&bv->Acreate);
476: }
478: /* Default version of BVMultInPlace */
479: PetscObjectTypeCompareAny((PetscObject)bv->t,&isgpu,VECSEQCUDA,VECMPICUDA,"");
480: ctx->vmip = isgpu? 1: 0;
482: /* Default BVMatMult method */
483: bv->vmm = BV_MATMULT_VECS;
485: /* Deferred call to setfromoptions */
486: if (bv->defersfo) {
487: ierr = PetscObjectOptionsBegin((PetscObject)bv);
488: BVSetFromOptions_Vecs(PetscOptionsObject,bv);
489: ierr = PetscOptionsEnd();
490: }
491: BVVecsSetVmip(bv,ctx->vmip);
493: bv->ops->mult = BVMult_Vecs;
494: bv->ops->multvec = BVMultVec_Vecs;
495: bv->ops->multinplacetrans = BVMultInPlaceHermitianTranspose_Vecs;
496: bv->ops->dot = BVDot_Vecs;
497: bv->ops->dotvec = BVDotVec_Vecs;
498: bv->ops->dotvec_begin = BVDotVec_Begin_Vecs;
499: bv->ops->dotvec_end = BVDotVec_End_Vecs;
500: bv->ops->scale = BVScale_Vecs;
501: bv->ops->norm = BVNorm_Vecs;
502: bv->ops->norm_begin = BVNorm_Begin_Vecs;
503: bv->ops->norm_end = BVNorm_End_Vecs;
504: bv->ops->matmult = BVMatMult_Vecs;
505: bv->ops->copy = BVCopy_Vecs;
506: bv->ops->copycolumn = BVCopyColumn_Vecs;
507: bv->ops->resize = BVResize_Vecs;
508: bv->ops->getcolumn = BVGetColumn_Vecs;
509: bv->ops->getarray = BVGetArray_Vecs;
510: bv->ops->restorearray = BVRestoreArray_Vecs;
511: bv->ops->getarrayread = BVGetArrayRead_Vecs;
512: bv->ops->restorearrayread = BVRestoreArrayRead_Vecs;
513: bv->ops->destroy = BVDestroy_Vecs;
514: bv->ops->duplicate = BVDuplicate_Vecs;
515: bv->ops->setfromoptions = BVSetFromOptions_Vecs;
516: bv->ops->view = BVView_Vecs;
517: PetscFunctionReturn(0);
518: }