Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
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 : */
13 :
14 : #include <slepc/private/bvimpl.h>
15 : #include "bvmat.h"
16 :
17 6522 : static PetscErrorCode BVMult_Mat(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
18 : {
19 6522 : BV_MAT *y = (BV_MAT*)Y->data,*x = (BV_MAT*)X->data;
20 6522 : PetscScalar *py;
21 6522 : const PetscScalar *px,*q;
22 6522 : PetscInt ldq;
23 :
24 6522 : PetscFunctionBegin;
25 6522 : PetscCall(MatDenseGetArrayRead(x->A,&px));
26 6522 : PetscCall(MatDenseGetArray(y->A,&py));
27 6522 : if (Q) {
28 5437 : PetscCall(MatDenseGetLDA(Q,&ldq));
29 5437 : PetscCall(MatDenseGetArrayRead(Q,&q));
30 5437 : PetscCall(BVMult_BLAS_Private(Y,Y->n,Y->k-Y->l,X->k-X->l,alpha,px+(X->nc+X->l)*X->ld,X->ld,q+Y->l*ldq+X->l,ldq,beta,py+(Y->nc+Y->l)*Y->ld,Y->ld));
31 5437 : PetscCall(MatDenseRestoreArrayRead(Q,&q));
32 1085 : } else PetscCall(BVAXPY_BLAS_Private(Y,Y->n,Y->k-Y->l,alpha,px+(X->nc+X->l)*X->ld,X->ld,beta,py+(Y->nc+Y->l)*Y->ld,Y->ld));
33 6522 : PetscCall(MatDenseRestoreArrayRead(x->A,&px));
34 6522 : PetscCall(MatDenseRestoreArray(y->A,&py));
35 6522 : PetscFunctionReturn(PETSC_SUCCESS);
36 : }
37 :
38 604662 : static PetscErrorCode BVMultVec_Mat(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
39 : {
40 604662 : BV_MAT *x = (BV_MAT*)X->data;
41 604662 : PetscScalar *py,*qq=q;
42 604662 : const PetscScalar *px;
43 :
44 604662 : PetscFunctionBegin;
45 604662 : PetscCall(MatDenseGetArrayRead(x->A,&px));
46 604662 : PetscCall(VecGetArray(y,&py));
47 604662 : if (!q) PetscCall(VecGetArray(X->buffer,&qq));
48 604662 : PetscCall(BVMultVec_BLAS_Private(X,X->n,X->k-X->l,alpha,px+(X->nc+X->l)*X->ld,X->ld,qq,beta,py));
49 604662 : if (!q) PetscCall(VecRestoreArray(X->buffer,&qq));
50 604662 : PetscCall(MatDenseRestoreArrayRead(x->A,&px));
51 604662 : PetscCall(VecRestoreArray(y,&py));
52 604662 : PetscFunctionReturn(PETSC_SUCCESS);
53 : }
54 :
55 26979 : static PetscErrorCode BVMultInPlace_Mat(BV V,Mat Q,PetscInt s,PetscInt e)
56 : {
57 26979 : BV_MAT *ctx = (BV_MAT*)V->data;
58 26979 : PetscScalar *pv;
59 26979 : const PetscScalar *q;
60 26979 : PetscInt ldq;
61 :
62 26979 : PetscFunctionBegin;
63 26979 : if (s>=e || !V->n) PetscFunctionReturn(PETSC_SUCCESS);
64 25166 : PetscCall(MatDenseGetLDA(Q,&ldq));
65 25166 : PetscCall(MatDenseGetArray(ctx->A,&pv));
66 25166 : PetscCall(MatDenseGetArrayRead(Q,&q));
67 25166 : PetscCall(BVMultInPlace_BLAS_Private(V,V->n,V->k-V->l,s-V->l,e-V->l,pv+(V->nc+V->l)*V->ld,V->ld,q+V->l*ldq+V->l,ldq,PETSC_FALSE));
68 25166 : PetscCall(MatDenseRestoreArrayRead(Q,&q));
69 25166 : PetscCall(MatDenseRestoreArray(ctx->A,&pv));
70 25166 : PetscFunctionReturn(PETSC_SUCCESS);
71 : }
72 :
73 1 : static PetscErrorCode BVMultInPlaceHermitianTranspose_Mat(BV V,Mat Q,PetscInt s,PetscInt e)
74 : {
75 1 : BV_MAT *ctx = (BV_MAT*)V->data;
76 1 : PetscScalar *pv;
77 1 : const PetscScalar *q;
78 1 : PetscInt ldq;
79 :
80 1 : PetscFunctionBegin;
81 1 : if (s>=e || !V->n) PetscFunctionReturn(PETSC_SUCCESS);
82 1 : PetscCall(MatDenseGetLDA(Q,&ldq));
83 1 : PetscCall(MatDenseGetArray(ctx->A,&pv));
84 1 : PetscCall(MatDenseGetArrayRead(Q,&q));
85 1 : PetscCall(BVMultInPlace_BLAS_Private(V,V->n,V->k-V->l,s-V->l,e-V->l,pv+(V->nc+V->l)*V->ld,V->ld,q+V->l*ldq+V->l,ldq,PETSC_TRUE));
86 1 : PetscCall(MatDenseRestoreArrayRead(Q,&q));
87 1 : PetscCall(MatDenseRestoreArray(ctx->A,&pv));
88 1 : PetscFunctionReturn(PETSC_SUCCESS);
89 : }
90 :
91 26141 : static PetscErrorCode BVDot_Mat(BV X,BV Y,Mat M)
92 : {
93 26141 : BV_MAT *x = (BV_MAT*)X->data,*y = (BV_MAT*)Y->data;
94 26141 : PetscScalar *m;
95 26141 : const PetscScalar *px,*py;
96 26141 : PetscInt ldm;
97 :
98 26141 : PetscFunctionBegin;
99 26141 : PetscCall(MatDenseGetLDA(M,&ldm));
100 26141 : PetscCall(MatDenseGetArrayRead(x->A,&px));
101 26141 : PetscCall(MatDenseGetArrayRead(y->A,&py));
102 26141 : PetscCall(MatDenseGetArray(M,&m));
103 26141 : PetscCall(BVDot_BLAS_Private(X,Y->k-Y->l,X->k-X->l,X->n,py+(Y->nc+Y->l)*Y->ld,Y->ld,px+(X->nc+X->l)*X->ld,X->ld,m+X->l*ldm+Y->l,ldm,x->mpi));
104 26141 : PetscCall(MatDenseRestoreArray(M,&m));
105 26141 : PetscCall(MatDenseRestoreArrayRead(x->A,&px));
106 26141 : PetscCall(MatDenseRestoreArrayRead(y->A,&py));
107 26141 : PetscFunctionReturn(PETSC_SUCCESS);
108 : }
109 :
110 460825 : static PetscErrorCode BVDotVec_Mat(BV X,Vec y,PetscScalar *q)
111 : {
112 460825 : BV_MAT *x = (BV_MAT*)X->data;
113 460825 : PetscScalar *qq=q;
114 460825 : const PetscScalar *px,*py;
115 460825 : Vec z = y;
116 :
117 460825 : PetscFunctionBegin;
118 460825 : if (PetscUnlikely(X->matrix)) {
119 131307 : PetscCall(BV_IPMatMult(X,y));
120 131307 : z = X->Bx;
121 : }
122 460825 : PetscCall(MatDenseGetArrayRead(x->A,&px));
123 460825 : PetscCall(VecGetArrayRead(z,&py));
124 460825 : if (!q) PetscCall(VecGetArray(X->buffer,&qq));
125 460825 : PetscCall(BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->ld,X->ld,py,qq,x->mpi));
126 460825 : if (!q) PetscCall(VecRestoreArray(X->buffer,&qq));
127 460825 : PetscCall(VecRestoreArrayRead(z,&py));
128 460825 : PetscCall(MatDenseRestoreArrayRead(x->A,&px));
129 460825 : PetscFunctionReturn(PETSC_SUCCESS);
130 : }
131 :
132 2485 : static PetscErrorCode BVDotVec_Local_Mat(BV X,Vec y,PetscScalar *m)
133 : {
134 2485 : BV_MAT *x = (BV_MAT*)X->data;
135 2485 : const PetscScalar *px,*py;
136 2485 : Vec z = y;
137 :
138 2485 : PetscFunctionBegin;
139 2485 : if (PetscUnlikely(X->matrix)) {
140 0 : PetscCall(BV_IPMatMult(X,y));
141 0 : z = X->Bx;
142 : }
143 2485 : PetscCall(MatDenseGetArrayRead(x->A,&px));
144 2485 : PetscCall(VecGetArrayRead(z,&py));
145 2485 : PetscCall(BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->ld,X->ld,py,m,PETSC_FALSE));
146 2485 : PetscCall(VecRestoreArrayRead(z,&py));
147 2485 : PetscCall(MatDenseRestoreArrayRead(x->A,&px));
148 2485 : PetscFunctionReturn(PETSC_SUCCESS);
149 : }
150 :
151 229342 : static PetscErrorCode BVScale_Mat(BV bv,PetscInt j,PetscScalar alpha)
152 : {
153 229342 : BV_MAT *ctx = (BV_MAT*)bv->data;
154 229342 : PetscScalar *array;
155 :
156 229342 : PetscFunctionBegin;
157 229342 : if (!bv->n) PetscFunctionReturn(PETSC_SUCCESS);
158 229342 : PetscCall(MatDenseGetArray(ctx->A,&array));
159 229342 : if (PetscUnlikely(j<0)) PetscCall(BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->ld,array+(bv->nc+bv->l)*bv->ld,alpha));
160 229051 : else PetscCall(BVScale_BLAS_Private(bv,bv->n,array+(bv->nc+j)*bv->ld,alpha));
161 229342 : PetscCall(MatDenseRestoreArray(ctx->A,&array));
162 229342 : PetscFunctionReturn(PETSC_SUCCESS);
163 : }
164 :
165 5177 : static PetscErrorCode BVNorm_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
166 : {
167 5177 : BV_MAT *ctx = (BV_MAT*)bv->data;
168 5177 : const PetscScalar *array;
169 :
170 5177 : PetscFunctionBegin;
171 5177 : PetscCall(MatDenseGetArrayRead(ctx->A,&array));
172 5177 : if (PetscUnlikely(j<0)) PetscCall(BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->ld,bv->ld,type,val,ctx->mpi));
173 5107 : else PetscCall(BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->ld,bv->ld,type,val,ctx->mpi));
174 5177 : PetscCall(MatDenseRestoreArrayRead(ctx->A,&array));
175 5177 : PetscFunctionReturn(PETSC_SUCCESS);
176 : }
177 :
178 2308 : static PetscErrorCode BVNorm_Local_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
179 : {
180 2308 : BV_MAT *ctx = (BV_MAT*)bv->data;
181 2308 : const PetscScalar *array;
182 :
183 2308 : PetscFunctionBegin;
184 2308 : PetscCall(MatDenseGetArrayRead(ctx->A,&array));
185 2308 : if (PetscUnlikely(j<0)) PetscCall(BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->ld,bv->ld,type,val,PETSC_FALSE));
186 2308 : else PetscCall(BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->ld,bv->ld,type,val,PETSC_FALSE));
187 2308 : PetscCall(MatDenseRestoreArrayRead(ctx->A,&array));
188 2308 : PetscFunctionReturn(PETSC_SUCCESS);
189 : }
190 :
191 3498 : static PetscErrorCode BVNormalize_Mat(BV bv,PetscScalar *eigi)
192 : {
193 3498 : BV_MAT *ctx = (BV_MAT*)bv->data;
194 3498 : PetscScalar *array,*wi=NULL;
195 :
196 3498 : PetscFunctionBegin;
197 3498 : PetscCall(MatDenseGetArray(ctx->A,&array));
198 3498 : if (eigi) wi = eigi+bv->l;
199 3498 : PetscCall(BVNormalize_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->ld,bv->ld,wi,ctx->mpi));
200 3498 : PetscCall(MatDenseRestoreArray(ctx->A,&array));
201 3498 : PetscFunctionReturn(PETSC_SUCCESS);
202 : }
203 :
204 11272 : static PetscErrorCode BVMatMult_Mat(BV V,Mat A,BV W)
205 : {
206 11272 : PetscInt j;
207 11272 : Mat Vmat,Wmat;
208 11272 : Vec vv,ww;
209 :
210 11272 : PetscFunctionBegin;
211 11272 : if (V->vmm) {
212 11265 : PetscCall(BVGetMat(V,&Vmat));
213 11265 : PetscCall(BVGetMat(W,&Wmat));
214 11265 : PetscCall(MatProductCreateWithMat(A,Vmat,NULL,Wmat));
215 11265 : PetscCall(MatProductSetType(Wmat,MATPRODUCT_AB));
216 11265 : PetscCall(MatProductSetFromOptions(Wmat));
217 11265 : PetscCall(MatProductSymbolic(Wmat));
218 11265 : PetscCall(MatProductNumeric(Wmat));
219 11265 : PetscCall(MatProductClear(Wmat));
220 11265 : PetscCall(BVRestoreMat(V,&Vmat));
221 11265 : PetscCall(BVRestoreMat(W,&Wmat));
222 : } else {
223 50 : for (j=0;j<V->k-V->l;j++) {
224 43 : PetscCall(BVGetColumn(V,V->l+j,&vv));
225 43 : PetscCall(BVGetColumn(W,W->l+j,&ww));
226 43 : PetscCall(MatMult(A,vv,ww));
227 43 : PetscCall(BVRestoreColumn(V,V->l+j,&vv));
228 43 : PetscCall(BVRestoreColumn(W,W->l+j,&ww));
229 : }
230 : }
231 11272 : PetscFunctionReturn(PETSC_SUCCESS);
232 : }
233 :
234 13482 : static PetscErrorCode BVCopy_Mat(BV V,BV W)
235 : {
236 13482 : BV_MAT *v = (BV_MAT*)V->data,*w = (BV_MAT*)W->data;
237 13482 : const PetscScalar *pv;
238 13482 : PetscScalar *pw;
239 13482 : PetscInt j;
240 :
241 13482 : PetscFunctionBegin;
242 13482 : PetscCall(MatDenseGetArrayRead(v->A,&pv));
243 13482 : PetscCall(MatDenseGetArray(w->A,&pw));
244 80452 : for (j=0;j<V->k-V->l;j++) PetscCall(PetscArraycpy(pw+(W->nc+W->l+j)*W->ld,pv+(V->nc+V->l+j)*V->ld,V->n));
245 13482 : PetscCall(MatDenseRestoreArrayRead(v->A,&pv));
246 13482 : PetscCall(MatDenseRestoreArray(w->A,&pw));
247 13482 : PetscFunctionReturn(PETSC_SUCCESS);
248 : }
249 :
250 8379 : static PetscErrorCode BVCopyColumn_Mat(BV V,PetscInt j,PetscInt i)
251 : {
252 8379 : BV_MAT *v = (BV_MAT*)V->data;
253 8379 : PetscScalar *pv;
254 :
255 8379 : PetscFunctionBegin;
256 8379 : PetscCall(MatDenseGetArray(v->A,&pv));
257 8379 : PetscCall(PetscArraycpy(pv+(V->nc+i)*V->ld,pv+(V->nc+j)*V->ld,V->n));
258 8379 : PetscCall(MatDenseRestoreArray(v->A,&pv));
259 8379 : PetscFunctionReturn(PETSC_SUCCESS);
260 : }
261 :
262 98 : static PetscErrorCode BVResize_Mat(BV bv,PetscInt m,PetscBool copy)
263 : {
264 98 : BV_MAT *ctx = (BV_MAT*)bv->data;
265 98 : Mat A,Msrc,Mdst;
266 98 : char str[50];
267 :
268 98 : PetscFunctionBegin;
269 98 : PetscCall(MatCreateDenseFromVecType(PetscObjectComm((PetscObject)bv),bv->vtype,bv->n,PETSC_DECIDE,bv->N,m,bv->ld,NULL,&A));
270 98 : if (((PetscObject)bv)->name) {
271 18 : PetscCall(PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name));
272 18 : PetscCall(PetscObjectSetName((PetscObject)A,str));
273 : }
274 98 : if (copy) {
275 5 : PetscCall(MatDenseGetSubMatrix(ctx->A,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PetscMin(m,bv->m),&Msrc));
276 5 : PetscCall(MatDenseGetSubMatrix(A,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PetscMin(m,bv->m),&Mdst));
277 5 : PetscCall(MatCopy(Msrc,Mdst,SAME_NONZERO_PATTERN));
278 5 : PetscCall(MatDenseRestoreSubMatrix(ctx->A,&Msrc));
279 5 : PetscCall(MatDenseRestoreSubMatrix(A,&Mdst));
280 : }
281 98 : PetscCall(MatDestroy(&ctx->A));
282 98 : ctx->A = A;
283 98 : PetscFunctionReturn(PETSC_SUCCESS);
284 : }
285 :
286 1528826 : static PetscErrorCode BVGetColumn_Mat(BV bv,PetscInt j,Vec *v)
287 : {
288 1528826 : BV_MAT *ctx = (BV_MAT*)bv->data;
289 1528826 : PetscScalar *pA;
290 1528826 : PetscInt l;
291 :
292 1528826 : PetscFunctionBegin;
293 1528826 : l = BVAvailableVec;
294 1528826 : PetscCall(MatDenseGetArray(ctx->A,&pA));
295 1528826 : PetscCall(VecPlaceArray(bv->cv[l],pA+(bv->nc+j)*bv->ld));
296 1528826 : PetscFunctionReturn(PETSC_SUCCESS);
297 : }
298 :
299 1528826 : static PetscErrorCode BVRestoreColumn_Mat(BV bv,PetscInt j,Vec *v)
300 : {
301 1528826 : BV_MAT *ctx = (BV_MAT*)bv->data;
302 1528826 : PetscScalar *pA;
303 1528826 : PetscInt l;
304 :
305 1528826 : PetscFunctionBegin;
306 1528826 : l = (j==bv->ci[0])? 0: 1;
307 1528826 : PetscCall(VecResetArray(bv->cv[l]));
308 1528826 : PetscCall(MatDenseRestoreArray(ctx->A,&pA));
309 1528826 : PetscFunctionReturn(PETSC_SUCCESS);
310 : }
311 :
312 25541 : static PetscErrorCode BVGetArray_Mat(BV bv,PetscScalar **a)
313 : {
314 25541 : BV_MAT *ctx = (BV_MAT*)bv->data;
315 :
316 25541 : PetscFunctionBegin;
317 25541 : PetscCall(MatDenseGetArray(ctx->A,a));
318 25541 : PetscFunctionReturn(PETSC_SUCCESS);
319 : }
320 :
321 25541 : static PetscErrorCode BVRestoreArray_Mat(BV bv,PetscScalar **a)
322 : {
323 25541 : BV_MAT *ctx = (BV_MAT*)bv->data;
324 :
325 25541 : PetscFunctionBegin;
326 25541 : if (a) PetscCall(MatDenseRestoreArray(ctx->A,a));
327 25541 : PetscFunctionReturn(PETSC_SUCCESS);
328 : }
329 :
330 166 : static PetscErrorCode BVGetArrayRead_Mat(BV bv,const PetscScalar **a)
331 : {
332 166 : BV_MAT *ctx = (BV_MAT*)bv->data;
333 :
334 166 : PetscFunctionBegin;
335 166 : PetscCall(MatDenseGetArrayRead(ctx->A,a));
336 166 : PetscFunctionReturn(PETSC_SUCCESS);
337 : }
338 :
339 166 : static PetscErrorCode BVRestoreArrayRead_Mat(BV bv,const PetscScalar **a)
340 : {
341 166 : BV_MAT *ctx = (BV_MAT*)bv->data;
342 :
343 166 : PetscFunctionBegin;
344 166 : if (a) PetscCall(MatDenseRestoreArrayRead(ctx->A,a));
345 166 : PetscFunctionReturn(PETSC_SUCCESS);
346 : }
347 :
348 28 : static PetscErrorCode BVView_Mat(BV bv,PetscViewer viewer)
349 : {
350 28 : Mat A;
351 28 : PetscViewerFormat format;
352 28 : PetscBool isascii;
353 28 : const char *bvname,*name;
354 :
355 28 : PetscFunctionBegin;
356 28 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
357 28 : if (isascii) {
358 28 : PetscCall(PetscViewerGetFormat(viewer,&format));
359 28 : if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(PETSC_SUCCESS);
360 : }
361 10 : PetscCall(BVGetMat(bv,&A));
362 10 : PetscCall(MatView(A,viewer));
363 10 : if (format == PETSC_VIEWER_ASCII_MATLAB) {
364 0 : PetscCall(PetscObjectGetName((PetscObject)A,&name));
365 0 : PetscCall(PetscObjectGetName((PetscObject)bv,&bvname));
366 0 : PetscCall(PetscViewerASCIIPrintf(viewer,"%s=%s;clear %s\n",bvname,name,name));
367 0 : if (bv->nc) PetscCall(PetscViewerASCIIPrintf(viewer,"%s=%s(:,%" PetscInt_FMT ":end);\n",bvname,bvname,bv->nc+1));
368 : }
369 10 : PetscCall(BVRestoreMat(bv,&A));
370 10 : PetscFunctionReturn(PETSC_SUCCESS);
371 : }
372 :
373 13114 : static PetscErrorCode BVDestroy_Mat(BV bv)
374 : {
375 13114 : BV_MAT *ctx = (BV_MAT*)bv->data;
376 :
377 13114 : PetscFunctionBegin;
378 13114 : PetscCall(MatDestroy(&ctx->A));
379 13114 : PetscCall(VecDestroy(&bv->cv[0]));
380 13114 : PetscCall(VecDestroy(&bv->cv[1]));
381 13114 : PetscCall(PetscFree(bv->data));
382 13114 : PetscFunctionReturn(PETSC_SUCCESS);
383 : }
384 :
385 13114 : SLEPC_EXTERN PetscErrorCode BVCreate_Mat(BV bv)
386 : {
387 13114 : BV_MAT *ctx;
388 13114 : PetscInt nloc,lsplit;
389 13114 : PetscBool seq;
390 13114 : char str[50];
391 13114 : PetscScalar *array,*ptr=NULL;
392 13114 : BV parent;
393 13114 : Mat Apar;
394 :
395 13114 : PetscFunctionBegin;
396 13114 : PetscCall(PetscNew(&ctx));
397 13114 : bv->data = (void*)ctx;
398 :
399 13114 : PetscCall(PetscStrcmpAny(bv->vtype,&bv->cuda,VECSEQCUDA,VECMPICUDA,""));
400 13114 : PetscCall(PetscStrcmpAny(bv->vtype,&bv->hip,VECSEQHIP,VECMPIHIP,""));
401 13114 : PetscCall(PetscStrcmpAny(bv->vtype,&ctx->mpi,VECMPI,VECMPICUDA,VECMPIHIP,""));
402 :
403 13114 : PetscCall(PetscStrcmp(bv->vtype,VECSEQ,&seq));
404 13114 : PetscCheck(seq || ctx->mpi || bv->cuda || bv->hip,PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"BVMAT does not support the requested vector type: %s",bv->vtype);
405 :
406 13114 : PetscCall(PetscLayoutGetLocalSize(bv->map,&nloc));
407 13114 : PetscCall(BV_SetDefaultLD(bv,nloc));
408 :
409 13114 : if (PetscUnlikely(bv->issplit)) {
410 : /* split BV: share the memory of the parent BV */
411 120 : parent = bv->splitparent;
412 120 : lsplit = parent->lsplit;
413 120 : Apar = ((BV_MAT*)parent->data)->A;
414 120 : if (bv->cuda) {
415 : #if defined(PETSC_HAVE_CUDA)
416 : PetscCall(MatDenseCUDAGetArray(Apar,&array));
417 : if (bv->issplit>0) ptr = (bv->issplit==1)? array: array+lsplit*bv->ld;
418 : else ptr = (bv->issplit==1)? array: array-lsplit;
419 : PetscCall(MatDenseCUDARestoreArray(Apar,&array));
420 : #endif
421 120 : } else if (bv->hip) {
422 : #if defined(PETSC_HAVE_HIP)
423 : PetscCall(MatDenseHIPGetArray(Apar,&array));
424 : if (bv->issplit>0) ptr = (bv->issplit==1)? array: array+lsplit*bv->ld;
425 : else ptr = (bv->issplit==1)? array: array-lsplit;
426 : PetscCall(MatDenseHIPRestoreArray(Apar,&array));
427 : #endif
428 : } else {
429 120 : PetscCall(MatDenseGetArray(Apar,&array));
430 120 : if (bv->issplit>0) ptr = (bv->issplit==1)? array: array+lsplit*bv->ld;
431 40 : else ptr = (bv->issplit==-1)? array: array-lsplit;
432 120 : PetscCall(MatDenseRestoreArray(Apar,&array));
433 : }
434 : }
435 :
436 13114 : PetscCall(MatCreateDenseFromVecType(PetscObjectComm((PetscObject)bv),bv->vtype,nloc,PETSC_DECIDE,bv->N,bv->m,bv->ld,ptr,&ctx->A));
437 13114 : if (((PetscObject)bv)->name) {
438 126 : PetscCall(PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name));
439 126 : PetscCall(PetscObjectSetName((PetscObject)ctx->A,str));
440 : }
441 :
442 13114 : if (PetscUnlikely(bv->Acreate)) {
443 136 : PetscCall(MatConvert(bv->Acreate,bv->cuda?MATDENSECUDA:bv->hip?MATDENSEHIP:MATDENSE,MAT_REUSE_MATRIX,&ctx->A));
444 68 : PetscCall(MatDestroy(&bv->Acreate));
445 : }
446 :
447 13114 : PetscCall(BVCreateVecEmpty(bv,&bv->cv[0]));
448 13114 : PetscCall(BVCreateVecEmpty(bv,&bv->cv[1]));
449 :
450 13114 : if (bv->cuda) {
451 : #if defined(PETSC_HAVE_CUDA)
452 : bv->ops->mult = BVMult_Mat_CUDA;
453 : bv->ops->multvec = BVMultVec_Mat_CUDA;
454 : bv->ops->multinplace = BVMultInPlace_Mat_CUDA;
455 : bv->ops->multinplacetrans = BVMultInPlaceHermitianTranspose_Mat_CUDA;
456 : bv->ops->dot = BVDot_Mat_CUDA;
457 : bv->ops->dotvec = BVDotVec_Mat_CUDA;
458 : bv->ops->dotvec_local = BVDotVec_Local_Mat_CUDA;
459 : bv->ops->scale = BVScale_Mat_CUDA;
460 : bv->ops->norm = BVNorm_Mat_CUDA;
461 : bv->ops->norm_local = BVNorm_Local_Mat_CUDA;
462 : bv->ops->normalize = BVNormalize_Mat_CUDA;
463 : bv->ops->matmult = BVMatMult_Mat_CUDA;
464 : bv->ops->copy = BVCopy_Mat_CUDA;
465 : bv->ops->copycolumn = BVCopyColumn_Mat_CUDA;
466 : bv->ops->getcolumn = BVGetColumn_Mat_CUDA;
467 : bv->ops->restorecolumn = BVRestoreColumn_Mat_CUDA;
468 : bv->ops->restoresplit = BVRestoreSplit_Mat_CUDA;
469 : bv->ops->restoresplitrows = BVRestoreSplitRows_Mat_CUDA;
470 : bv->ops->getmat = BVGetMat_Mat_CUDA;
471 : bv->ops->restoremat = BVRestoreMat_Mat_CUDA;
472 : #endif
473 13114 : } else if (bv->hip) {
474 : #if defined(PETSC_HAVE_HIP)
475 : bv->ops->mult = BVMult_Mat_HIP;
476 : bv->ops->multvec = BVMultVec_Mat_HIP;
477 : bv->ops->multinplace = BVMultInPlace_Mat_HIP;
478 : bv->ops->multinplacetrans = BVMultInPlaceHermitianTranspose_Mat_HIP;
479 : bv->ops->dot = BVDot_Mat_HIP;
480 : bv->ops->dotvec = BVDotVec_Mat_HIP;
481 : bv->ops->dotvec_local = BVDotVec_Local_Mat_HIP;
482 : bv->ops->scale = BVScale_Mat_HIP;
483 : bv->ops->norm = BVNorm_Mat_HIP;
484 : bv->ops->norm_local = BVNorm_Local_Mat_HIP;
485 : bv->ops->normalize = BVNormalize_Mat_HIP;
486 : bv->ops->matmult = BVMatMult_Mat_HIP;
487 : bv->ops->copy = BVCopy_Mat_HIP;
488 : bv->ops->copycolumn = BVCopyColumn_Mat_HIP;
489 : bv->ops->getcolumn = BVGetColumn_Mat_HIP;
490 : bv->ops->restorecolumn = BVRestoreColumn_Mat_HIP;
491 : bv->ops->restoresplit = BVRestoreSplit_Mat_HIP;
492 : bv->ops->restoresplitrows = BVRestoreSplitRows_Mat_HIP;
493 : bv->ops->getmat = BVGetMat_Mat_HIP;
494 : bv->ops->restoremat = BVRestoreMat_Mat_HIP;
495 : #endif
496 : } else {
497 13114 : bv->ops->mult = BVMult_Mat;
498 13114 : bv->ops->multvec = BVMultVec_Mat;
499 13114 : bv->ops->multinplace = BVMultInPlace_Mat;
500 13114 : bv->ops->multinplacetrans = BVMultInPlaceHermitianTranspose_Mat;
501 13114 : bv->ops->dot = BVDot_Mat;
502 13114 : bv->ops->dotvec = BVDotVec_Mat;
503 13114 : bv->ops->dotvec_local = BVDotVec_Local_Mat;
504 13114 : bv->ops->scale = BVScale_Mat;
505 13114 : bv->ops->norm = BVNorm_Mat;
506 13114 : bv->ops->norm_local = BVNorm_Local_Mat;
507 13114 : bv->ops->normalize = BVNormalize_Mat;
508 13114 : bv->ops->matmult = BVMatMult_Mat;
509 13114 : bv->ops->copy = BVCopy_Mat;
510 13114 : bv->ops->copycolumn = BVCopyColumn_Mat;
511 13114 : bv->ops->getcolumn = BVGetColumn_Mat;
512 13114 : bv->ops->restorecolumn = BVRestoreColumn_Mat;
513 13114 : bv->ops->getmat = BVGetMat_Default;
514 13114 : bv->ops->restoremat = BVRestoreMat_Default;
515 : }
516 13114 : bv->ops->resize = BVResize_Mat;
517 13114 : bv->ops->getarray = BVGetArray_Mat;
518 13114 : bv->ops->restorearray = BVRestoreArray_Mat;
519 13114 : bv->ops->getarrayread = BVGetArrayRead_Mat;
520 13114 : bv->ops->restorearrayread = BVRestoreArrayRead_Mat;
521 13114 : bv->ops->destroy = BVDestroy_Mat;
522 13114 : bv->ops->view = BVView_Mat;
523 13114 : PetscFunctionReturn(PETSC_SUCCESS);
524 : }
|