Actual source code: svec.c
slepc-3.16.0 2021-09-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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: PetscErrorCode ierr;
20: BV_SVEC *y = (BV_SVEC*)Y->data,*x = (BV_SVEC*)X->data;
21: const PetscScalar *px,*q;
22: PetscScalar *py;
23: PetscInt ldq;
26: VecGetArrayRead(x->v,&px);
27: VecGetArray(y->v,&py);
28: if (Q) {
29: MatGetSize(Q,&ldq,NULL);
30: MatDenseGetArrayRead(Q,&q);
31: 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);
32: MatDenseRestoreArrayRead(Q,&q);
33: } else {
34: 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);
35: }
36: VecRestoreArrayRead(x->v,&px);
37: VecRestoreArray(y->v,&py);
38: return(0);
39: }
41: PetscErrorCode BVMultVec_Svec(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
42: {
44: BV_SVEC *x = (BV_SVEC*)X->data;
45: PetscScalar *px,*py,*qq=q;
48: VecGetArray(x->v,&px);
49: VecGetArray(y,&py);
50: if (!q) { VecGetArray(X->buffer,&qq); }
51: BVMultVec_BLAS_Private(X,X->n,X->k-X->l,alpha,px+(X->nc+X->l)*X->n,qq,beta,py);
52: if (!q) { VecRestoreArray(X->buffer,&qq); }
53: VecRestoreArray(x->v,&px);
54: VecRestoreArray(y,&py);
55: return(0);
56: }
58: PetscErrorCode BVMultInPlace_Svec(BV V,Mat Q,PetscInt s,PetscInt e)
59: {
60: PetscErrorCode ierr;
61: BV_SVEC *ctx = (BV_SVEC*)V->data;
62: PetscScalar *pv;
63: const PetscScalar *q;
64: PetscInt ldq;
67: MatGetSize(Q,&ldq,NULL);
68: VecGetArray(ctx->v,&pv);
69: MatDenseGetArrayRead(Q,&q);
70: 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);
71: MatDenseRestoreArrayRead(Q,&q);
72: VecRestoreArray(ctx->v,&pv);
73: return(0);
74: }
76: PetscErrorCode BVMultInPlaceHermitianTranspose_Svec(BV V,Mat Q,PetscInt s,PetscInt e)
77: {
78: PetscErrorCode ierr;
79: BV_SVEC *ctx = (BV_SVEC*)V->data;
80: PetscScalar *pv;
81: const PetscScalar *q;
82: PetscInt ldq;
85: MatGetSize(Q,&ldq,NULL);
86: VecGetArray(ctx->v,&pv);
87: MatDenseGetArrayRead(Q,&q);
88: 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);
89: MatDenseRestoreArrayRead(Q,&q);
90: VecRestoreArray(ctx->v,&pv);
91: return(0);
92: }
94: PetscErrorCode BVDot_Svec(BV X,BV Y,Mat M)
95: {
96: PetscErrorCode ierr;
97: BV_SVEC *x = (BV_SVEC*)X->data,*y = (BV_SVEC*)Y->data;
98: const PetscScalar *px,*py;
99: PetscScalar *m;
100: PetscInt ldm;
103: MatGetSize(M,&ldm,NULL);
104: VecGetArrayRead(x->v,&px);
105: VecGetArrayRead(y->v,&py);
106: MatDenseGetArray(M,&m);
107: 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);
108: MatDenseRestoreArray(M,&m);
109: VecRestoreArrayRead(x->v,&px);
110: VecRestoreArrayRead(y->v,&py);
111: return(0);
112: }
114: PetscErrorCode BVDotVec_Svec(BV X,Vec y,PetscScalar *q)
115: {
116: PetscErrorCode ierr;
117: BV_SVEC *x = (BV_SVEC*)X->data;
118: const PetscScalar *px,*py;
119: PetscScalar *qq=q;
120: Vec z = y;
123: if (X->matrix) {
124: BV_IPMatMult(X,y);
125: z = X->Bx;
126: }
127: VecGetArrayRead(x->v,&px);
128: VecGetArrayRead(z,&py);
129: if (!q) { VecGetArray(X->buffer,&qq); }
130: BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,qq,x->mpi);
131: if (!q) { VecRestoreArray(X->buffer,&qq); }
132: VecRestoreArrayRead(z,&py);
133: VecRestoreArrayRead(x->v,&px);
134: return(0);
135: }
137: PetscErrorCode BVDotVec_Local_Svec(BV X,Vec y,PetscScalar *m)
138: {
140: BV_SVEC *x = (BV_SVEC*)X->data;
141: PetscScalar *px,*py;
142: Vec z = y;
145: if (X->matrix) {
146: BV_IPMatMult(X,y);
147: z = X->Bx;
148: }
149: VecGetArray(x->v,&px);
150: VecGetArray(z,&py);
151: BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,m,PETSC_FALSE);
152: VecRestoreArray(z,&py);
153: VecRestoreArray(x->v,&px);
154: return(0);
155: }
157: PetscErrorCode BVScale_Svec(BV bv,PetscInt j,PetscScalar alpha)
158: {
160: BV_SVEC *ctx = (BV_SVEC*)bv->data;
161: PetscScalar *array;
164: VecGetArray(ctx->v,&array);
165: if (j<0) {
166: BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->n,array+(bv->nc+bv->l)*bv->n,alpha);
167: } else {
168: BVScale_BLAS_Private(bv,bv->n,array+(bv->nc+j)*bv->n,alpha);
169: }
170: VecRestoreArray(ctx->v,&array);
171: return(0);
172: }
174: PetscErrorCode BVNorm_Svec(BV bv,PetscInt j,NormType type,PetscReal *val)
175: {
177: BV_SVEC *ctx = (BV_SVEC*)bv->data;
178: PetscScalar *array;
181: VecGetArray(ctx->v,&array);
182: if (j<0) {
183: BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,ctx->mpi);
184: } else {
185: BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,ctx->mpi);
186: }
187: VecRestoreArray(ctx->v,&array);
188: return(0);
189: }
191: PetscErrorCode BVNorm_Local_Svec(BV bv,PetscInt j,NormType type,PetscReal *val)
192: {
194: BV_SVEC *ctx = (BV_SVEC*)bv->data;
195: PetscScalar *array;
198: VecGetArray(ctx->v,&array);
199: if (j<0) {
200: BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,PETSC_FALSE);
201: } else {
202: BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,PETSC_FALSE);
203: }
204: VecRestoreArray(ctx->v,&array);
205: return(0);
206: }
208: PetscErrorCode BVNormalize_Svec(BV bv,PetscScalar *eigi)
209: {
211: BV_SVEC *ctx = (BV_SVEC*)bv->data;
212: PetscScalar *array,*wi=NULL;
215: VecGetArray(ctx->v,&array);
216: if (eigi) wi = eigi+bv->l;
217: BVNormalize_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,wi,ctx->mpi);
218: VecRestoreArray(ctx->v,&array);
219: return(0);
220: }
222: PetscErrorCode BVMatMult_Svec(BV V,Mat A,BV W)
223: {
225: PetscInt j;
226: Mat Vmat,Wmat;
227: Vec vv,ww;
230: if (V->vmm) {
231: BVGetMat(V,&Vmat);
232: BVGetMat(W,&Wmat);
233: MatProductCreateWithMat(A,Vmat,NULL,Wmat);
234: MatProductSetType(Wmat,MATPRODUCT_AB);
235: MatProductSetFromOptions(Wmat);
236: MatProductSymbolic(Wmat);
237: MatProductNumeric(Wmat);
238: MatProductClear(Wmat);
239: BVRestoreMat(V,&Vmat);
240: BVRestoreMat(W,&Wmat);
241: } else {
242: for (j=0;j<V->k-V->l;j++) {
243: BVGetColumn(V,V->l+j,&vv);
244: BVGetColumn(W,W->l+j,&ww);
245: MatMult(A,vv,ww);
246: BVRestoreColumn(V,V->l+j,&vv);
247: BVRestoreColumn(W,W->l+j,&ww);
248: }
249: }
250: return(0);
251: }
253: PetscErrorCode BVCopy_Svec(BV V,BV W)
254: {
256: BV_SVEC *v = (BV_SVEC*)V->data,*w = (BV_SVEC*)W->data;
257: PetscScalar *pv,*pw,*pvc,*pwc;
260: VecGetArray(v->v,&pv);
261: VecGetArray(w->v,&pw);
262: pvc = pv+(V->nc+V->l)*V->n;
263: pwc = pw+(W->nc+W->l)*W->n;
264: PetscArraycpy(pwc,pvc,(V->k-V->l)*V->n);
265: VecRestoreArray(v->v,&pv);
266: VecRestoreArray(w->v,&pw);
267: return(0);
268: }
270: PetscErrorCode BVCopyColumn_Svec(BV V,PetscInt j,PetscInt i)
271: {
273: BV_SVEC *v = (BV_SVEC*)V->data;
274: PetscScalar *pv;
277: VecGetArray(v->v,&pv);
278: PetscArraycpy(pv+(V->nc+i)*V->n,pv+(V->nc+j)*V->n,V->n);
279: VecRestoreArray(v->v,&pv);
280: return(0);
281: }
283: PetscErrorCode BVResize_Svec(BV bv,PetscInt m,PetscBool copy)
284: {
285: PetscErrorCode ierr;
286: BV_SVEC *ctx = (BV_SVEC*)bv->data;
287: PetscScalar *pnew;
288: const PetscScalar *pv;
289: PetscInt bs;
290: Vec vnew;
291: char str[50];
294: VecGetBlockSize(bv->t,&bs);
295: VecCreate(PetscObjectComm((PetscObject)bv->t),&vnew);
296: VecSetType(vnew,((PetscObject)bv->t)->type_name);
297: VecSetSizes(vnew,m*bv->n,PETSC_DECIDE);
298: VecSetBlockSize(vnew,bs);
299: PetscLogObjectParent((PetscObject)bv,(PetscObject)vnew);
300: if (((PetscObject)bv)->name) {
301: PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name);
302: PetscObjectSetName((PetscObject)vnew,str);
303: }
304: if (copy) {
305: VecGetArrayRead(ctx->v,&pv);
306: VecGetArray(vnew,&pnew);
307: PetscArraycpy(pnew,pv,PetscMin(m,bv->m)*bv->n);
308: VecRestoreArrayRead(ctx->v,&pv);
309: VecRestoreArray(vnew,&pnew);
310: }
311: VecDestroy(&ctx->v);
312: ctx->v = vnew;
313: return(0);
314: }
316: PetscErrorCode BVGetColumn_Svec(BV bv,PetscInt j,Vec *v)
317: {
319: BV_SVEC *ctx = (BV_SVEC*)bv->data;
320: PetscScalar *pv;
321: PetscInt l;
324: l = BVAvailableVec;
325: VecGetArray(ctx->v,&pv);
326: VecPlaceArray(bv->cv[l],pv+(bv->nc+j)*bv->n);
327: return(0);
328: }
330: PetscErrorCode BVRestoreColumn_Svec(BV bv,PetscInt j,Vec *v)
331: {
333: BV_SVEC *ctx = (BV_SVEC*)bv->data;
334: PetscInt l;
337: l = (j==bv->ci[0])? 0: 1;
338: VecResetArray(bv->cv[l]);
339: VecRestoreArray(ctx->v,NULL);
340: return(0);
341: }
343: PetscErrorCode BVGetArray_Svec(BV bv,PetscScalar **a)
344: {
346: BV_SVEC *ctx = (BV_SVEC*)bv->data;
349: VecGetArray(ctx->v,a);
350: return(0);
351: }
353: PetscErrorCode BVRestoreArray_Svec(BV bv,PetscScalar **a)
354: {
356: BV_SVEC *ctx = (BV_SVEC*)bv->data;
359: VecRestoreArray(ctx->v,a);
360: return(0);
361: }
363: PetscErrorCode BVGetArrayRead_Svec(BV bv,const PetscScalar **a)
364: {
366: BV_SVEC *ctx = (BV_SVEC*)bv->data;
369: VecGetArrayRead(ctx->v,a);
370: return(0);
371: }
373: PetscErrorCode BVRestoreArrayRead_Svec(BV bv,const PetscScalar **a)
374: {
376: BV_SVEC *ctx = (BV_SVEC*)bv->data;
379: VecRestoreArrayRead(ctx->v,a);
380: return(0);
381: }
383: PetscErrorCode BVView_Svec(BV bv,PetscViewer viewer)
384: {
385: PetscErrorCode ierr;
386: BV_SVEC *ctx = (BV_SVEC*)bv->data;
387: PetscViewerFormat format;
388: PetscBool isascii;
389: const char *bvname,*name;
392: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
393: if (isascii) {
394: PetscViewerGetFormat(viewer,&format);
395: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
396: VecView(ctx->v,viewer);
397: if (format == PETSC_VIEWER_ASCII_MATLAB) {
398: PetscObjectGetName((PetscObject)bv,&bvname);
399: PetscObjectGetName((PetscObject)ctx->v,&name);
400: PetscViewerASCIIPrintf(viewer,"%s=reshape(%s,%D,%D);clear %s\n",bvname,name,bv->N,bv->nc+bv->m,name);
401: if (bv->nc) {
402: PetscViewerASCIIPrintf(viewer,"%s=%s(:,%D:end);\n",bvname,bvname,bv->nc+1);
403: }
404: }
405: } else {
406: VecView(ctx->v,viewer);
407: }
408: return(0);
409: }
411: PetscErrorCode BVDestroy_Svec(BV bv)
412: {
414: BV_SVEC *ctx = (BV_SVEC*)bv->data;
417: VecDestroy(&ctx->v);
418: VecDestroy(&bv->cv[0]);
419: VecDestroy(&bv->cv[1]);
420: PetscFree(bv->data);
421: bv->cuda = PETSC_FALSE;
422: return(0);
423: }
425: SLEPC_EXTERN PetscErrorCode BVCreate_Svec(BV bv)
426: {
427: PetscErrorCode ierr;
428: BV_SVEC *ctx;
429: PetscInt nloc,N,bs,tglobal=0,tlocal,lsplit;
430: PetscBool seq;
431: PetscScalar *vv;
432: const PetscScalar *aa,*array,*ptr;
433: char str[50];
434: BV parent;
435: Vec vpar;
436: #if defined(PETSC_HAVE_CUDA)
437: PetscScalar *gpuarray,*gptr;
438: #endif
441: PetscNewLog(bv,&ctx);
442: bv->data = (void*)ctx;
444: PetscObjectTypeCompareAny((PetscObject)bv->t,&bv->cuda,VECSEQCUDA,VECMPICUDA,"");
445: PetscObjectTypeCompareAny((PetscObject)bv->t,&ctx->mpi,VECMPI,VECMPICUDA,"");
447: PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
448: if (!seq && !ctx->mpi && !bv->cuda) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"BVSVEC does not support the type of the provided template vector");
450: VecGetLocalSize(bv->t,&nloc);
451: VecGetSize(bv->t,&N);
452: VecGetBlockSize(bv->t,&bs);
453: tlocal = bv->m*nloc;
454: PetscIntMultError(bv->m,N,&tglobal);
456: if (bv->issplit) {
457: /* split BV: create Vec sharing the memory of the parent BV */
458: parent = bv->splitparent;
459: lsplit = parent->lsplit;
460: vpar = ((BV_SVEC*)parent->data)->v;
461: if (bv->cuda) {
462: #if defined(PETSC_HAVE_CUDA)
463: VecCUDAGetArray(vpar,&gpuarray);
464: gptr = (bv->issplit==1)? gpuarray: gpuarray+lsplit*nloc;
465: VecCUDARestoreArray(vpar,&gpuarray);
466: if (ctx->mpi) {
467: VecCreateMPICUDAWithArray(PetscObjectComm((PetscObject)bv->t),bs,tlocal,bv->m*N,NULL,&ctx->v);
468: } else {
469: VecCreateSeqCUDAWithArray(PetscObjectComm((PetscObject)bv->t),bs,tlocal,NULL,&ctx->v);
470: }
471: VecCUDAPlaceArray(ctx->v,gptr);
472: #endif
473: } else {
474: VecGetArrayRead(vpar,&array);
475: ptr = (bv->issplit==1)? array: array+lsplit*nloc;
476: VecRestoreArrayRead(vpar,&array);
477: if (ctx->mpi) {
478: VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,tlocal,bv->m*N,NULL,&ctx->v);
479: } else {
480: VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,tlocal,NULL,&ctx->v);
481: }
482: VecPlaceArray(ctx->v,ptr);
483: }
484: } else {
485: /* regular BV: create Vec to store the BV entries */
486: VecCreate(PetscObjectComm((PetscObject)bv->t),&ctx->v);
487: VecSetType(ctx->v,((PetscObject)bv->t)->type_name);
488: VecSetSizes(ctx->v,tlocal,tglobal);
489: VecSetBlockSize(ctx->v,bs);
490: }
491: PetscLogObjectParent((PetscObject)bv,(PetscObject)ctx->v);
492: if (((PetscObject)bv)->name) {
493: PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name);
494: PetscObjectSetName((PetscObject)ctx->v,str);
495: }
497: if (bv->Acreate) {
498: MatDenseGetArrayRead(bv->Acreate,&aa);
499: VecGetArray(ctx->v,&vv);
500: PetscArraycpy(vv,aa,tlocal);
501: VecRestoreArray(ctx->v,&vv);
502: MatDenseRestoreArrayRead(bv->Acreate,&aa);
503: MatDestroy(&bv->Acreate);
504: }
506: VecDuplicateEmpty(bv->t,&bv->cv[0]);
507: VecDuplicateEmpty(bv->t,&bv->cv[1]);
509: if (bv->cuda) {
510: #if defined(PETSC_HAVE_CUDA)
511: bv->ops->mult = BVMult_Svec_CUDA;
512: bv->ops->multvec = BVMultVec_Svec_CUDA;
513: bv->ops->multinplace = BVMultInPlace_Svec_CUDA;
514: bv->ops->multinplacetrans = BVMultInPlaceHermitianTranspose_Svec_CUDA;
515: bv->ops->dot = BVDot_Svec_CUDA;
516: bv->ops->dotvec = BVDotVec_Svec_CUDA;
517: bv->ops->dotvec_local = BVDotVec_Local_Svec_CUDA;
518: bv->ops->scale = BVScale_Svec_CUDA;
519: bv->ops->matmult = BVMatMult_Svec_CUDA;
520: bv->ops->copy = BVCopy_Svec_CUDA;
521: bv->ops->copycolumn = BVCopyColumn_Svec_CUDA;
522: bv->ops->resize = BVResize_Svec_CUDA;
523: bv->ops->getcolumn = BVGetColumn_Svec_CUDA;
524: bv->ops->restorecolumn = BVRestoreColumn_Svec_CUDA;
525: bv->ops->restoresplit = BVRestoreSplit_Svec_CUDA;
526: bv->ops->getmat = BVGetMat_Svec_CUDA;
527: bv->ops->restoremat = BVRestoreMat_Svec_CUDA;
528: #endif
529: } else {
530: bv->ops->mult = BVMult_Svec;
531: bv->ops->multvec = BVMultVec_Svec;
532: bv->ops->multinplace = BVMultInPlace_Svec;
533: bv->ops->multinplacetrans = BVMultInPlaceHermitianTranspose_Svec;
534: bv->ops->dot = BVDot_Svec;
535: bv->ops->dotvec = BVDotVec_Svec;
536: bv->ops->dotvec_local = BVDotVec_Local_Svec;
537: bv->ops->scale = BVScale_Svec;
538: bv->ops->matmult = BVMatMult_Svec;
539: bv->ops->copy = BVCopy_Svec;
540: bv->ops->copycolumn = BVCopyColumn_Svec;
541: bv->ops->resize = BVResize_Svec;
542: bv->ops->getcolumn = BVGetColumn_Svec;
543: bv->ops->restorecolumn = BVRestoreColumn_Svec;
544: }
545: bv->ops->norm = BVNorm_Svec;
546: bv->ops->norm_local = BVNorm_Local_Svec;
547: bv->ops->normalize = BVNormalize_Svec;
548: bv->ops->getarray = BVGetArray_Svec;
549: bv->ops->restorearray = BVRestoreArray_Svec;
550: bv->ops->getarrayread = BVGetArrayRead_Svec;
551: bv->ops->restorearrayread = BVRestoreArrayRead_Svec;
552: bv->ops->destroy = BVDestroy_Svec;
553: if (!ctx->mpi) bv->ops->view = BVView_Svec;
554: return(0);
555: }