Actual source code: bvmat.c
slepc-3.17.2 2022-08-09
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 with a dense Mat
12: */
14: #include <slepc/private/bvimpl.h>
16: typedef struct {
17: Mat A;
18: PetscBool mpi;
19: } BV_MAT;
21: PetscErrorCode BVMult_Mat(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
22: {
23: BV_MAT *y = (BV_MAT*)Y->data,*x = (BV_MAT*)X->data;
24: PetscScalar *py;
25: const PetscScalar *px,*q;
26: PetscInt ldq;
28: MatDenseGetArrayRead(x->A,&px);
29: MatDenseGetArray(y->A,&py);
30: if (Q) {
31: MatDenseGetLDA(Q,&ldq);
32: MatDenseGetArrayRead(Q,&q);
33: 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);
34: MatDenseRestoreArrayRead(Q,&q);
35: } else BVAXPY_BLAS_Private(Y,Y->n,Y->k-Y->l,alpha,px+(X->nc+X->l)*X->n,beta,py+(Y->nc+Y->l)*Y->n);
36: MatDenseRestoreArrayRead(x->A,&px);
37: MatDenseRestoreArray(y->A,&py);
38: PetscFunctionReturn(0);
39: }
41: PetscErrorCode BVMultVec_Mat(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
42: {
43: BV_MAT *x = (BV_MAT*)X->data;
44: PetscScalar *py,*qq=q;
45: const PetscScalar *px;
47: MatDenseGetArrayRead(x->A,&px);
48: VecGetArray(y,&py);
49: if (!q) VecGetArray(X->buffer,&qq);
50: BVMultVec_BLAS_Private(X,X->n,X->k-X->l,alpha,px+(X->nc+X->l)*X->n,qq,beta,py);
51: if (!q) VecRestoreArray(X->buffer,&qq);
52: MatDenseRestoreArrayRead(x->A,&px);
53: VecRestoreArray(y,&py);
54: PetscFunctionReturn(0);
55: }
57: PetscErrorCode BVMultInPlace_Mat(BV V,Mat Q,PetscInt s,PetscInt e)
58: {
59: BV_MAT *ctx = (BV_MAT*)V->data;
60: PetscScalar *pv;
61: const PetscScalar *q;
62: PetscInt ldq;
64: MatDenseGetLDA(Q,&ldq);
65: MatDenseGetArray(ctx->A,&pv);
66: MatDenseGetArrayRead(Q,&q);
67: 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);
68: MatDenseRestoreArrayRead(Q,&q);
69: MatDenseRestoreArray(ctx->A,&pv);
70: PetscFunctionReturn(0);
71: }
73: PetscErrorCode BVMultInPlaceHermitianTranspose_Mat(BV V,Mat Q,PetscInt s,PetscInt e)
74: {
75: BV_MAT *ctx = (BV_MAT*)V->data;
76: PetscScalar *pv;
77: const PetscScalar *q;
78: PetscInt ldq;
80: MatDenseGetLDA(Q,&ldq);
81: MatDenseGetArray(ctx->A,&pv);
82: MatDenseGetArrayRead(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_TRUE);
84: MatDenseRestoreArrayRead(Q,&q);
85: MatDenseRestoreArray(ctx->A,&pv);
86: PetscFunctionReturn(0);
87: }
89: PetscErrorCode BVDot_Mat(BV X,BV Y,Mat M)
90: {
91: BV_MAT *x = (BV_MAT*)X->data,*y = (BV_MAT*)Y->data;
92: PetscScalar *m;
93: const PetscScalar *px,*py;
94: PetscInt ldm;
96: MatDenseGetLDA(M,&ldm);
97: MatDenseGetArrayRead(x->A,&px);
98: MatDenseGetArrayRead(y->A,&py);
99: MatDenseGetArray(M,&m);
100: 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);
101: MatDenseRestoreArray(M,&m);
102: MatDenseRestoreArrayRead(x->A,&px);
103: MatDenseRestoreArrayRead(y->A,&py);
104: PetscFunctionReturn(0);
105: }
107: PetscErrorCode BVDotVec_Mat(BV X,Vec y,PetscScalar *q)
108: {
109: BV_MAT *x = (BV_MAT*)X->data;
110: PetscScalar *qq=q;
111: const PetscScalar *px,*py;
112: Vec z = y;
114: if (PetscUnlikely(X->matrix)) {
115: BV_IPMatMult(X,y);
116: z = X->Bx;
117: }
118: MatDenseGetArrayRead(x->A,&px);
119: VecGetArrayRead(z,&py);
120: if (!q) VecGetArray(X->buffer,&qq);
121: BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,qq,x->mpi);
122: if (!q) VecRestoreArray(X->buffer,&qq);
123: VecRestoreArrayRead(z,&py);
124: MatDenseRestoreArrayRead(x->A,&px);
125: PetscFunctionReturn(0);
126: }
128: PetscErrorCode BVDotVec_Local_Mat(BV X,Vec y,PetscScalar *m)
129: {
130: BV_MAT *x = (BV_MAT*)X->data;
131: const PetscScalar *px,*py;
132: Vec z = y;
134: if (PetscUnlikely(X->matrix)) {
135: BV_IPMatMult(X,y);
136: z = X->Bx;
137: }
138: MatDenseGetArrayRead(x->A,&px);
139: VecGetArrayRead(z,&py);
140: BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,m,PETSC_FALSE);
141: VecRestoreArrayRead(z,&py);
142: MatDenseRestoreArrayRead(x->A,&px);
143: PetscFunctionReturn(0);
144: }
146: PetscErrorCode BVScale_Mat(BV bv,PetscInt j,PetscScalar alpha)
147: {
148: BV_MAT *ctx = (BV_MAT*)bv->data;
149: PetscScalar *array;
151: MatDenseGetArray(ctx->A,&array);
152: if (PetscUnlikely(j<0)) BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->n,array+(bv->nc+bv->l)*bv->n,alpha);
153: else BVScale_BLAS_Private(bv,bv->n,array+(bv->nc+j)*bv->n,alpha);
154: MatDenseRestoreArray(ctx->A,&array);
155: PetscFunctionReturn(0);
156: }
158: PetscErrorCode BVNorm_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
159: {
160: BV_MAT *ctx = (BV_MAT*)bv->data;
161: const PetscScalar *array;
163: MatDenseGetArrayRead(ctx->A,&array);
164: if (PetscUnlikely(j<0)) BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,ctx->mpi);
165: else BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,ctx->mpi);
166: MatDenseRestoreArrayRead(ctx->A,&array);
167: PetscFunctionReturn(0);
168: }
170: PetscErrorCode BVNorm_Local_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
171: {
172: BV_MAT *ctx = (BV_MAT*)bv->data;
173: const PetscScalar *array;
175: MatDenseGetArrayRead(ctx->A,&array);
176: if (PetscUnlikely(j<0)) BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,PETSC_FALSE);
177: else BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,PETSC_FALSE);
178: MatDenseRestoreArrayRead(ctx->A,&array);
179: PetscFunctionReturn(0);
180: }
182: PetscErrorCode BVNormalize_Mat(BV bv,PetscScalar *eigi)
183: {
184: BV_MAT *ctx = (BV_MAT*)bv->data;
185: PetscScalar *array,*wi=NULL;
187: MatDenseGetArray(ctx->A,&array);
188: if (eigi) wi = eigi+bv->l;
189: BVNormalize_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,wi,ctx->mpi);
190: MatDenseRestoreArray(ctx->A,&array);
191: PetscFunctionReturn(0);
192: }
194: PetscErrorCode BVMatMult_Mat(BV V,Mat A,BV W)
195: {
196: PetscInt j;
197: Mat Vmat,Wmat;
198: Vec vv,ww;
200: if (V->vmm) {
201: BVGetMat(V,&Vmat);
202: BVGetMat(W,&Wmat);
203: MatProductCreateWithMat(A,Vmat,NULL,Wmat);
204: MatProductSetType(Wmat,MATPRODUCT_AB);
205: MatProductSetFromOptions(Wmat);
206: MatProductSymbolic(Wmat);
207: MatProductNumeric(Wmat);
208: MatProductClear(Wmat);
209: BVRestoreMat(V,&Vmat);
210: BVRestoreMat(W,&Wmat);
211: } else {
212: for (j=0;j<V->k-V->l;j++) {
213: BVGetColumn(V,V->l+j,&vv);
214: BVGetColumn(W,W->l+j,&ww);
215: MatMult(A,vv,ww);
216: BVRestoreColumn(V,V->l+j,&vv);
217: BVRestoreColumn(W,W->l+j,&ww);
218: }
219: }
220: PetscFunctionReturn(0);
221: }
223: PetscErrorCode BVCopy_Mat(BV V,BV W)
224: {
225: BV_MAT *v = (BV_MAT*)V->data,*w = (BV_MAT*)W->data;
226: PetscScalar *pw,*pwc;
227: const PetscScalar *pv,*pvc;
229: MatDenseGetArrayRead(v->A,&pv);
230: MatDenseGetArray(w->A,&pw);
231: pvc = pv+(V->nc+V->l)*V->n;
232: pwc = pw+(W->nc+W->l)*W->n;
233: PetscArraycpy(pwc,pvc,(V->k-V->l)*V->n);
234: MatDenseRestoreArrayRead(v->A,&pv);
235: MatDenseRestoreArray(w->A,&pw);
236: PetscFunctionReturn(0);
237: }
239: PetscErrorCode BVCopyColumn_Mat(BV V,PetscInt j,PetscInt i)
240: {
241: BV_MAT *v = (BV_MAT*)V->data;
242: PetscScalar *pv;
244: MatDenseGetArray(v->A,&pv);
245: PetscArraycpy(pv+(V->nc+i)*V->n,pv+(V->nc+j)*V->n,V->n);
246: MatDenseRestoreArray(v->A,&pv);
247: PetscFunctionReturn(0);
248: }
250: PetscErrorCode BVResize_Mat(BV bv,PetscInt m,PetscBool copy)
251: {
252: BV_MAT *ctx = (BV_MAT*)bv->data;
253: PetscScalar *pnew;
254: const PetscScalar *pA;
255: Mat A;
256: char str[50];
258: MatCreateDense(PetscObjectComm((PetscObject)bv->t),bv->n,PETSC_DECIDE,PETSC_DECIDE,m,NULL,&A);
259: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
260: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
261: PetscLogObjectParent((PetscObject)bv,(PetscObject)A);
262: if (((PetscObject)bv)->name) {
263: PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name);
264: PetscObjectSetName((PetscObject)A,str);
265: }
266: if (copy) {
267: MatDenseGetArrayRead(ctx->A,&pA);
268: MatDenseGetArrayWrite(A,&pnew);
269: PetscArraycpy(pnew,pA,PetscMin(m,bv->m)*bv->n);
270: MatDenseRestoreArrayRead(ctx->A,&pA);
271: MatDenseRestoreArrayWrite(A,&pnew);
272: }
273: MatDestroy(&ctx->A);
274: ctx->A = A;
275: PetscFunctionReturn(0);
276: }
278: PetscErrorCode BVGetColumn_Mat(BV bv,PetscInt j,Vec *v)
279: {
280: BV_MAT *ctx = (BV_MAT*)bv->data;
281: PetscScalar *pA;
282: PetscInt l;
284: l = BVAvailableVec;
285: MatDenseGetArray(ctx->A,&pA);
286: VecPlaceArray(bv->cv[l],pA+(bv->nc+j)*bv->n);
287: PetscFunctionReturn(0);
288: }
290: PetscErrorCode BVRestoreColumn_Mat(BV bv,PetscInt j,Vec *v)
291: {
292: BV_MAT *ctx = (BV_MAT*)bv->data;
293: PetscScalar *pA;
294: PetscInt l;
296: l = (j==bv->ci[0])? 0: 1;
297: VecResetArray(bv->cv[l]);
298: MatDenseRestoreArray(ctx->A,&pA);
299: PetscFunctionReturn(0);
300: }
302: PetscErrorCode BVGetArray_Mat(BV bv,PetscScalar **a)
303: {
304: BV_MAT *ctx = (BV_MAT*)bv->data;
306: MatDenseGetArray(ctx->A,a);
307: PetscFunctionReturn(0);
308: }
310: PetscErrorCode BVRestoreArray_Mat(BV bv,PetscScalar **a)
311: {
312: BV_MAT *ctx = (BV_MAT*)bv->data;
314: if (a) MatDenseRestoreArray(ctx->A,a);
315: PetscFunctionReturn(0);
316: }
318: PetscErrorCode BVGetArrayRead_Mat(BV bv,const PetscScalar **a)
319: {
320: BV_MAT *ctx = (BV_MAT*)bv->data;
322: MatDenseGetArrayRead(ctx->A,a);
323: PetscFunctionReturn(0);
324: }
326: PetscErrorCode BVRestoreArrayRead_Mat(BV bv,const PetscScalar **a)
327: {
328: BV_MAT *ctx = (BV_MAT*)bv->data;
330: if (a) MatDenseRestoreArrayRead(ctx->A,a);
331: PetscFunctionReturn(0);
332: }
334: PetscErrorCode BVView_Mat(BV bv,PetscViewer viewer)
335: {
336: BV_MAT *ctx = (BV_MAT*)bv->data;
337: PetscViewerFormat format;
338: PetscBool isascii;
339: const char *bvname,*name;
341: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
342: if (isascii) {
343: PetscViewerGetFormat(viewer,&format);
344: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0);
345: MatView(ctx->A,viewer);
346: if (format == PETSC_VIEWER_ASCII_MATLAB) {
347: PetscObjectGetName((PetscObject)bv,&bvname);
348: PetscObjectGetName((PetscObject)ctx->A,&name);
349: PetscViewerASCIIPrintf(viewer,"%s=%s;clear %s\n",bvname,name,name);
350: if (bv->nc) PetscViewerASCIIPrintf(viewer,"%s=%s(:,%" PetscInt_FMT ":end);\n",bvname,bvname,bv->nc+1);
351: }
352: } else MatView(ctx->A,viewer);
353: PetscFunctionReturn(0);
354: }
356: PetscErrorCode BVDestroy_Mat(BV bv)
357: {
358: BV_MAT *ctx = (BV_MAT*)bv->data;
360: MatDestroy(&ctx->A);
361: VecDestroy(&bv->cv[0]);
362: VecDestroy(&bv->cv[1]);
363: PetscFree(bv->data);
364: PetscFunctionReturn(0);
365: }
367: SLEPC_EXTERN PetscErrorCode BVCreate_Mat(BV bv)
368: {
369: BV_MAT *ctx;
370: PetscInt nloc,bs,lsplit;
371: PetscBool seq;
372: char str[50];
373: PetscScalar *array,*ptr;
374: BV parent;
376: PetscNewLog(bv,&ctx);
377: bv->data = (void*)ctx;
379: PetscObjectTypeCompare((PetscObject)bv->t,VECMPI,&ctx->mpi);
380: if (!ctx->mpi) {
381: PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
383: }
385: VecGetLocalSize(bv->t,&nloc);
386: VecGetBlockSize(bv->t,&bs);
388: if (PetscUnlikely(bv->issplit)) {
389: /* split BV: share the memory of the parent BV */
390: parent = bv->splitparent;
391: lsplit = parent->lsplit;
392: MatDenseGetArray(((BV_MAT*)parent->data)->A,&array);
393: ptr = (bv->issplit==1)? array: array+lsplit*nloc;
394: MatDenseRestoreArray(((BV_MAT*)parent->data)->A,&array);
395: } else {
396: /* regular BV: allocate memory for the BV entries */
397: ptr = NULL;
398: }
399: MatCreateDense(PetscObjectComm((PetscObject)bv->t),nloc,PETSC_DECIDE,PETSC_DECIDE,bv->m,ptr,&ctx->A);
400: MatAssemblyBegin(ctx->A,MAT_FINAL_ASSEMBLY);
401: MatAssemblyEnd(ctx->A,MAT_FINAL_ASSEMBLY);
402: PetscLogObjectParent((PetscObject)bv,(PetscObject)ctx->A);
403: if (((PetscObject)bv)->name) {
404: PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name);
405: PetscObjectSetName((PetscObject)ctx->A,str);
406: }
408: if (PetscUnlikely(bv->Acreate)) {
409: MatCopy(bv->Acreate,ctx->A,SAME_NONZERO_PATTERN);
410: MatDestroy(&bv->Acreate);
411: }
413: VecDuplicateEmpty(bv->t,&bv->cv[0]);
414: VecDuplicateEmpty(bv->t,&bv->cv[1]);
416: bv->ops->mult = BVMult_Mat;
417: bv->ops->multvec = BVMultVec_Mat;
418: bv->ops->multinplace = BVMultInPlace_Mat;
419: bv->ops->multinplacetrans = BVMultInPlaceHermitianTranspose_Mat;
420: bv->ops->dot = BVDot_Mat;
421: bv->ops->dotvec = BVDotVec_Mat;
422: bv->ops->dotvec_local = BVDotVec_Local_Mat;
423: bv->ops->scale = BVScale_Mat;
424: bv->ops->norm = BVNorm_Mat;
425: bv->ops->norm_local = BVNorm_Local_Mat;
426: bv->ops->normalize = BVNormalize_Mat;
427: bv->ops->matmult = BVMatMult_Mat;
428: bv->ops->copy = BVCopy_Mat;
429: bv->ops->copycolumn = BVCopyColumn_Mat;
430: bv->ops->resize = BVResize_Mat;
431: bv->ops->getcolumn = BVGetColumn_Mat;
432: bv->ops->restorecolumn = BVRestoreColumn_Mat;
433: bv->ops->getarray = BVGetArray_Mat;
434: bv->ops->restorearray = BVRestoreArray_Mat;
435: bv->ops->getarrayread = BVGetArrayRead_Mat;
436: bv->ops->restorearrayread = BVRestoreArrayRead_Mat;
437: bv->ops->destroy = BVDestroy_Mat;
438: if (!ctx->mpi) bv->ops->view = BVView_Mat;
439: PetscFunctionReturn(0);
440: }