Actual source code: svec.c
slepc-3.18.2 2023-01-26
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 a single Vec
12: */
14: #include <slepc/private/bvimpl.h>
15: #include "svec.h"
17: PetscErrorCode BVMult_Svec(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
18: {
19: BV_SVEC *y = (BV_SVEC*)Y->data,*x = (BV_SVEC*)X->data;
20: const PetscScalar *px,*q;
21: PetscScalar *py;
22: PetscInt ldq;
24: VecGetArrayRead(x->v,&px);
25: VecGetArray(y->v,&py);
26: if (Q) {
27: MatDenseGetLDA(Q,&ldq);
28: MatDenseGetArrayRead(Q,&q);
29: 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);
30: MatDenseRestoreArrayRead(Q,&q);
31: } 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);
32: VecRestoreArrayRead(x->v,&px);
33: VecRestoreArray(y->v,&py);
34: return 0;
35: }
37: PetscErrorCode BVMultVec_Svec(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
38: {
39: BV_SVEC *x = (BV_SVEC*)X->data;
40: PetscScalar *px,*py,*qq=q;
42: VecGetArray(x->v,&px);
43: VecGetArray(y,&py);
44: if (!q) VecGetArray(X->buffer,&qq);
45: BVMultVec_BLAS_Private(X,X->n,X->k-X->l,alpha,px+(X->nc+X->l)*X->n,qq,beta,py);
46: if (!q) VecRestoreArray(X->buffer,&qq);
47: VecRestoreArray(x->v,&px);
48: VecRestoreArray(y,&py);
49: return 0;
50: }
52: PetscErrorCode BVMultInPlace_Svec(BV V,Mat Q,PetscInt s,PetscInt e)
53: {
54: BV_SVEC *ctx = (BV_SVEC*)V->data;
55: PetscScalar *pv;
56: const PetscScalar *q;
57: PetscInt ldq;
59: MatDenseGetLDA(Q,&ldq);
60: VecGetArray(ctx->v,&pv);
61: MatDenseGetArrayRead(Q,&q);
62: 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);
63: MatDenseRestoreArrayRead(Q,&q);
64: VecRestoreArray(ctx->v,&pv);
65: return 0;
66: }
68: PetscErrorCode BVMultInPlaceHermitianTranspose_Svec(BV V,Mat Q,PetscInt s,PetscInt e)
69: {
70: BV_SVEC *ctx = (BV_SVEC*)V->data;
71: PetscScalar *pv;
72: const PetscScalar *q;
73: PetscInt ldq;
75: MatDenseGetLDA(Q,&ldq);
76: VecGetArray(ctx->v,&pv);
77: MatDenseGetArrayRead(Q,&q);
78: 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);
79: MatDenseRestoreArrayRead(Q,&q);
80: VecRestoreArray(ctx->v,&pv);
81: return 0;
82: }
84: PetscErrorCode BVDot_Svec(BV X,BV Y,Mat M)
85: {
86: BV_SVEC *x = (BV_SVEC*)X->data,*y = (BV_SVEC*)Y->data;
87: const PetscScalar *px,*py;
88: PetscScalar *m;
89: PetscInt ldm;
91: MatDenseGetLDA(M,&ldm);
92: VecGetArrayRead(x->v,&px);
93: VecGetArrayRead(y->v,&py);
94: MatDenseGetArray(M,&m);
95: 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);
96: MatDenseRestoreArray(M,&m);
97: VecRestoreArrayRead(x->v,&px);
98: VecRestoreArrayRead(y->v,&py);
99: return 0;
100: }
102: PetscErrorCode BVDotVec_Svec(BV X,Vec y,PetscScalar *q)
103: {
104: BV_SVEC *x = (BV_SVEC*)X->data;
105: const PetscScalar *px,*py;
106: PetscScalar *qq=q;
107: Vec z = y;
109: if (PetscUnlikely(X->matrix)) {
110: BV_IPMatMult(X,y);
111: z = X->Bx;
112: }
113: VecGetArrayRead(x->v,&px);
114: VecGetArrayRead(z,&py);
115: if (!q) VecGetArray(X->buffer,&qq);
116: BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,qq,x->mpi);
117: if (!q) VecRestoreArray(X->buffer,&qq);
118: VecRestoreArrayRead(z,&py);
119: VecRestoreArrayRead(x->v,&px);
120: return 0;
121: }
123: PetscErrorCode BVDotVec_Local_Svec(BV X,Vec y,PetscScalar *m)
124: {
125: BV_SVEC *x = (BV_SVEC*)X->data;
126: PetscScalar *px,*py;
127: Vec z = y;
129: if (PetscUnlikely(X->matrix)) {
130: BV_IPMatMult(X,y);
131: z = X->Bx;
132: }
133: VecGetArray(x->v,&px);
134: VecGetArray(z,&py);
135: BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,m,PETSC_FALSE);
136: VecRestoreArray(z,&py);
137: VecRestoreArray(x->v,&px);
138: return 0;
139: }
141: PetscErrorCode BVScale_Svec(BV bv,PetscInt j,PetscScalar alpha)
142: {
143: BV_SVEC *ctx = (BV_SVEC*)bv->data;
144: PetscScalar *array;
146: VecGetArray(ctx->v,&array);
147: if (PetscUnlikely(j<0)) BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->n,array+(bv->nc+bv->l)*bv->n,alpha);
148: else BVScale_BLAS_Private(bv,bv->n,array+(bv->nc+j)*bv->n,alpha);
149: VecRestoreArray(ctx->v,&array);
150: return 0;
151: }
153: PetscErrorCode BVNorm_Svec(BV bv,PetscInt j,NormType type,PetscReal *val)
154: {
155: BV_SVEC *ctx = (BV_SVEC*)bv->data;
156: PetscScalar *array;
158: VecGetArray(ctx->v,&array);
159: 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);
160: else BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,ctx->mpi);
161: VecRestoreArray(ctx->v,&array);
162: return 0;
163: }
165: PetscErrorCode BVNorm_Local_Svec(BV bv,PetscInt j,NormType type,PetscReal *val)
166: {
167: BV_SVEC *ctx = (BV_SVEC*)bv->data;
168: PetscScalar *array;
170: VecGetArray(ctx->v,&array);
171: 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);
172: else BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,PETSC_FALSE);
173: VecRestoreArray(ctx->v,&array);
174: return 0;
175: }
177: PetscErrorCode BVNormalize_Svec(BV bv,PetscScalar *eigi)
178: {
179: BV_SVEC *ctx = (BV_SVEC*)bv->data;
180: PetscScalar *array,*wi=NULL;
182: VecGetArray(ctx->v,&array);
183: if (eigi) wi = eigi+bv->l;
184: BVNormalize_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,wi,ctx->mpi);
185: VecRestoreArray(ctx->v,&array);
186: return 0;
187: }
189: PetscErrorCode BVMatMult_Svec(BV V,Mat A,BV W)
190: {
191: PetscInt j;
192: Mat Vmat,Wmat;
193: Vec vv,ww;
195: if (V->vmm) {
196: BVGetMat(V,&Vmat);
197: BVGetMat(W,&Wmat);
198: MatProductCreateWithMat(A,Vmat,NULL,Wmat);
199: MatProductSetType(Wmat,MATPRODUCT_AB);
200: MatProductSetFromOptions(Wmat);
201: MatProductSymbolic(Wmat);
202: MatProductNumeric(Wmat);
203: MatProductClear(Wmat);
204: BVRestoreMat(V,&Vmat);
205: BVRestoreMat(W,&Wmat);
206: } else {
207: for (j=0;j<V->k-V->l;j++) {
208: BVGetColumn(V,V->l+j,&vv);
209: BVGetColumn(W,W->l+j,&ww);
210: MatMult(A,vv,ww);
211: BVRestoreColumn(V,V->l+j,&vv);
212: BVRestoreColumn(W,W->l+j,&ww);
213: }
214: }
215: return 0;
216: }
218: PetscErrorCode BVCopy_Svec(BV V,BV W)
219: {
220: BV_SVEC *v = (BV_SVEC*)V->data,*w = (BV_SVEC*)W->data;
221: PetscScalar *pv,*pw,*pvc,*pwc;
223: VecGetArray(v->v,&pv);
224: VecGetArray(w->v,&pw);
225: pvc = pv+(V->nc+V->l)*V->n;
226: pwc = pw+(W->nc+W->l)*W->n;
227: PetscArraycpy(pwc,pvc,(V->k-V->l)*V->n);
228: VecRestoreArray(v->v,&pv);
229: VecRestoreArray(w->v,&pw);
230: return 0;
231: }
233: PetscErrorCode BVCopyColumn_Svec(BV V,PetscInt j,PetscInt i)
234: {
235: BV_SVEC *v = (BV_SVEC*)V->data;
236: PetscScalar *pv;
238: VecGetArray(v->v,&pv);
239: PetscArraycpy(pv+(V->nc+i)*V->n,pv+(V->nc+j)*V->n,V->n);
240: VecRestoreArray(v->v,&pv);
241: return 0;
242: }
244: PetscErrorCode BVResize_Svec(BV bv,PetscInt m,PetscBool copy)
245: {
246: BV_SVEC *ctx = (BV_SVEC*)bv->data;
247: PetscScalar *pnew;
248: const PetscScalar *pv;
249: PetscInt bs;
250: Vec vnew;
251: char str[50];
253: VecGetBlockSize(bv->t,&bs);
254: VecCreate(PetscObjectComm((PetscObject)bv->t),&vnew);
255: VecSetType(vnew,((PetscObject)bv->t)->type_name);
256: VecSetSizes(vnew,m*bv->n,PETSC_DECIDE);
257: VecSetBlockSize(vnew,bs);
258: if (((PetscObject)bv)->name) {
259: PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name);
260: PetscObjectSetName((PetscObject)vnew,str);
261: }
262: if (copy) {
263: VecGetArrayRead(ctx->v,&pv);
264: VecGetArray(vnew,&pnew);
265: PetscArraycpy(pnew,pv,PetscMin(m,bv->m)*bv->n);
266: VecRestoreArrayRead(ctx->v,&pv);
267: VecRestoreArray(vnew,&pnew);
268: }
269: VecDestroy(&ctx->v);
270: ctx->v = vnew;
271: return 0;
272: }
274: PetscErrorCode BVGetColumn_Svec(BV bv,PetscInt j,Vec *v)
275: {
276: BV_SVEC *ctx = (BV_SVEC*)bv->data;
277: PetscScalar *pv;
278: PetscInt l;
280: l = BVAvailableVec;
281: VecGetArray(ctx->v,&pv);
282: VecPlaceArray(bv->cv[l],pv+(bv->nc+j)*bv->n);
283: return 0;
284: }
286: PetscErrorCode BVRestoreColumn_Svec(BV bv,PetscInt j,Vec *v)
287: {
288: BV_SVEC *ctx = (BV_SVEC*)bv->data;
289: PetscInt l;
291: l = (j==bv->ci[0])? 0: 1;
292: VecResetArray(bv->cv[l]);
293: VecRestoreArray(ctx->v,NULL);
294: return 0;
295: }
297: PetscErrorCode BVGetArray_Svec(BV bv,PetscScalar **a)
298: {
299: BV_SVEC *ctx = (BV_SVEC*)bv->data;
301: VecGetArray(ctx->v,a);
302: return 0;
303: }
305: PetscErrorCode BVRestoreArray_Svec(BV bv,PetscScalar **a)
306: {
307: BV_SVEC *ctx = (BV_SVEC*)bv->data;
309: VecRestoreArray(ctx->v,a);
310: return 0;
311: }
313: PetscErrorCode BVGetArrayRead_Svec(BV bv,const PetscScalar **a)
314: {
315: BV_SVEC *ctx = (BV_SVEC*)bv->data;
317: VecGetArrayRead(ctx->v,a);
318: return 0;
319: }
321: PetscErrorCode BVRestoreArrayRead_Svec(BV bv,const PetscScalar **a)
322: {
323: BV_SVEC *ctx = (BV_SVEC*)bv->data;
325: VecRestoreArrayRead(ctx->v,a);
326: return 0;
327: }
329: PetscErrorCode BVView_Svec(BV bv,PetscViewer viewer)
330: {
331: BV_SVEC *ctx = (BV_SVEC*)bv->data;
332: PetscViewerFormat format;
333: PetscBool isascii;
334: const char *bvname,*name;
336: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
337: if (isascii) {
338: PetscViewerGetFormat(viewer,&format);
339: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return 0;
340: VecView(ctx->v,viewer);
341: if (format == PETSC_VIEWER_ASCII_MATLAB) {
342: PetscObjectGetName((PetscObject)bv,&bvname);
343: PetscObjectGetName((PetscObject)ctx->v,&name);
344: PetscViewerASCIIPrintf(viewer,"%s=reshape(%s,%" PetscInt_FMT ",%" PetscInt_FMT ");clear %s\n",bvname,name,bv->N,bv->nc+bv->m,name);
345: if (bv->nc) PetscViewerASCIIPrintf(viewer,"%s=%s(:,%" PetscInt_FMT ":end);\n",bvname,bvname,bv->nc+1);
346: }
347: } else VecView(ctx->v,viewer);
348: return 0;
349: }
351: PetscErrorCode BVDestroy_Svec(BV bv)
352: {
353: BV_SVEC *ctx = (BV_SVEC*)bv->data;
355: VecDestroy(&ctx->v);
356: VecDestroy(&bv->cv[0]);
357: VecDestroy(&bv->cv[1]);
358: PetscFree(bv->data);
359: bv->cuda = PETSC_FALSE;
360: return 0;
361: }
363: SLEPC_EXTERN PetscErrorCode BVCreate_Svec(BV bv)
364: {
365: BV_SVEC *ctx;
366: PetscInt nloc,N,bs,tglobal=0,tlocal,lsplit;
367: PetscBool seq;
368: PetscScalar *vv;
369: const PetscScalar *aa,*array,*ptr;
370: char str[50];
371: BV parent;
372: Vec vpar;
373: #if defined(PETSC_HAVE_CUDA)
374: PetscScalar *gpuarray,*gptr;
375: #endif
377: PetscNew(&ctx);
378: bv->data = (void*)ctx;
380: PetscObjectTypeCompareAny((PetscObject)bv->t,&bv->cuda,VECSEQCUDA,VECMPICUDA,"");
381: PetscObjectTypeCompareAny((PetscObject)bv->t,&ctx->mpi,VECMPI,VECMPICUDA,"");
383: PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
386: VecGetLocalSize(bv->t,&nloc);
387: VecGetSize(bv->t,&N);
388: VecGetBlockSize(bv->t,&bs);
389: tlocal = bv->m*nloc;
390: PetscIntMultError(bv->m,N,&tglobal);
392: if (PetscUnlikely(bv->issplit)) {
393: /* split BV: create Vec sharing the memory of the parent BV */
394: parent = bv->splitparent;
395: lsplit = parent->lsplit;
396: vpar = ((BV_SVEC*)parent->data)->v;
397: if (bv->cuda) {
398: #if defined(PETSC_HAVE_CUDA)
399: VecCUDAGetArray(vpar,&gpuarray);
400: gptr = (bv->issplit==1)? gpuarray: gpuarray+lsplit*nloc;
401: VecCUDARestoreArray(vpar,&gpuarray);
402: if (ctx->mpi) VecCreateMPICUDAWithArray(PetscObjectComm((PetscObject)bv->t),bs,tlocal,bv->m*N,NULL,&ctx->v);
403: else VecCreateSeqCUDAWithArray(PetscObjectComm((PetscObject)bv->t),bs,tlocal,NULL,&ctx->v);
404: VecCUDAPlaceArray(ctx->v,gptr);
405: #endif
406: } else {
407: VecGetArrayRead(vpar,&array);
408: ptr = (bv->issplit==1)? array: array+lsplit*nloc;
409: VecRestoreArrayRead(vpar,&array);
410: if (ctx->mpi) VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,tlocal,bv->m*N,NULL,&ctx->v);
411: else VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,tlocal,NULL,&ctx->v);
412: VecPlaceArray(ctx->v,ptr);
413: }
414: } else {
415: /* regular BV: create Vec to store the BV entries */
416: VecCreate(PetscObjectComm((PetscObject)bv->t),&ctx->v);
417: VecSetType(ctx->v,((PetscObject)bv->t)->type_name);
418: VecSetSizes(ctx->v,tlocal,tglobal);
419: VecSetBlockSize(ctx->v,bs);
420: }
421: if (((PetscObject)bv)->name) {
422: PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name);
423: PetscObjectSetName((PetscObject)ctx->v,str);
424: }
426: if (PetscUnlikely(bv->Acreate)) {
427: MatDenseGetArrayRead(bv->Acreate,&aa);
428: VecGetArray(ctx->v,&vv);
429: PetscArraycpy(vv,aa,tlocal);
430: VecRestoreArray(ctx->v,&vv);
431: MatDenseRestoreArrayRead(bv->Acreate,&aa);
432: MatDestroy(&bv->Acreate);
433: }
435: VecDuplicateEmpty(bv->t,&bv->cv[0]);
436: VecDuplicateEmpty(bv->t,&bv->cv[1]);
438: if (bv->cuda) {
439: #if defined(PETSC_HAVE_CUDA)
440: bv->ops->mult = BVMult_Svec_CUDA;
441: bv->ops->multvec = BVMultVec_Svec_CUDA;
442: bv->ops->multinplace = BVMultInPlace_Svec_CUDA;
443: bv->ops->multinplacetrans = BVMultInPlaceHermitianTranspose_Svec_CUDA;
444: bv->ops->dot = BVDot_Svec_CUDA;
445: bv->ops->dotvec = BVDotVec_Svec_CUDA;
446: bv->ops->dotvec_local = BVDotVec_Local_Svec_CUDA;
447: bv->ops->scale = BVScale_Svec_CUDA;
448: bv->ops->matmult = BVMatMult_Svec_CUDA;
449: bv->ops->copy = BVCopy_Svec_CUDA;
450: bv->ops->copycolumn = BVCopyColumn_Svec_CUDA;
451: bv->ops->resize = BVResize_Svec_CUDA;
452: bv->ops->getcolumn = BVGetColumn_Svec_CUDA;
453: bv->ops->restorecolumn = BVRestoreColumn_Svec_CUDA;
454: bv->ops->restoresplit = BVRestoreSplit_Svec_CUDA;
455: bv->ops->getmat = BVGetMat_Svec_CUDA;
456: bv->ops->restoremat = BVRestoreMat_Svec_CUDA;
457: #endif
458: } else {
459: bv->ops->mult = BVMult_Svec;
460: bv->ops->multvec = BVMultVec_Svec;
461: bv->ops->multinplace = BVMultInPlace_Svec;
462: bv->ops->multinplacetrans = BVMultInPlaceHermitianTranspose_Svec;
463: bv->ops->dot = BVDot_Svec;
464: bv->ops->dotvec = BVDotVec_Svec;
465: bv->ops->dotvec_local = BVDotVec_Local_Svec;
466: bv->ops->scale = BVScale_Svec;
467: bv->ops->matmult = BVMatMult_Svec;
468: bv->ops->copy = BVCopy_Svec;
469: bv->ops->copycolumn = BVCopyColumn_Svec;
470: bv->ops->resize = BVResize_Svec;
471: bv->ops->getcolumn = BVGetColumn_Svec;
472: bv->ops->restorecolumn = BVRestoreColumn_Svec;
473: bv->ops->getmat = BVGetMat_Default;
474: bv->ops->restoremat = BVRestoreMat_Default;
475: }
476: bv->ops->norm = BVNorm_Svec;
477: bv->ops->norm_local = BVNorm_Local_Svec;
478: bv->ops->normalize = BVNormalize_Svec;
479: bv->ops->getarray = BVGetArray_Svec;
480: bv->ops->restorearray = BVRestoreArray_Svec;
481: bv->ops->getarrayread = BVGetArrayRead_Svec;
482: bv->ops->restorearrayread = BVRestoreArrayRead_Svec;
483: bv->ops->destroy = BVDestroy_Svec;
484: if (!ctx->mpi) bv->ops->view = BVView_Svec;
485: return 0;
486: }