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 as an array of Vecs sharing a contiguous array for elements
12 : */
13 :
14 : #include <slepc/private/bvimpl.h>
15 :
16 : typedef struct {
17 : Vec *V;
18 : PetscScalar *array;
19 : PetscBool mpi;
20 : } BV_CONTIGUOUS;
21 :
22 98 : static PetscErrorCode BVMult_Contiguous(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
23 : {
24 98 : BV_CONTIGUOUS *y = (BV_CONTIGUOUS*)Y->data,*x = (BV_CONTIGUOUS*)X->data;
25 98 : const PetscScalar *q;
26 98 : PetscInt ldq;
27 :
28 98 : PetscFunctionBegin;
29 98 : if (Q) {
30 78 : PetscCall(MatDenseGetLDA(Q,&ldq));
31 78 : PetscCall(MatDenseGetArrayRead(Q,&q));
32 78 : PetscCall(BVMult_BLAS_Private(Y,Y->n,Y->k-Y->l,X->k-X->l,alpha,x->array+(X->nc+X->l)*X->ld,X->ld,q+Y->l*ldq+X->l,ldq,beta,y->array+(Y->nc+Y->l)*Y->ld,Y->ld));
33 78 : PetscCall(MatDenseRestoreArrayRead(Q,&q));
34 20 : } else PetscCall(BVAXPY_BLAS_Private(Y,Y->n,Y->k-Y->l,alpha,x->array+(X->nc+X->l)*X->ld,X->ld,beta,y->array+(Y->nc+Y->l)*Y->ld,Y->ld));
35 98 : PetscFunctionReturn(PETSC_SUCCESS);
36 : }
37 :
38 343 : static PetscErrorCode BVMultVec_Contiguous(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
39 : {
40 343 : BV_CONTIGUOUS *x = (BV_CONTIGUOUS*)X->data;
41 343 : PetscScalar *py,*qq=q;
42 :
43 343 : PetscFunctionBegin;
44 343 : PetscCall(VecGetArray(y,&py));
45 343 : if (!q) PetscCall(VecGetArray(X->buffer,&qq));
46 343 : PetscCall(BVMultVec_BLAS_Private(X,X->n,X->k-X->l,alpha,x->array+(X->nc+X->l)*X->ld,X->ld,qq,beta,py));
47 343 : if (!q) PetscCall(VecRestoreArray(X->buffer,&qq));
48 343 : PetscCall(VecRestoreArray(y,&py));
49 343 : PetscFunctionReturn(PETSC_SUCCESS);
50 : }
51 :
52 46 : static PetscErrorCode BVMultInPlace_Contiguous(BV V,Mat Q,PetscInt s,PetscInt e)
53 : {
54 46 : BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)V->data;
55 46 : const PetscScalar *q;
56 46 : PetscInt ldq;
57 :
58 46 : PetscFunctionBegin;
59 46 : if (s>=e || !V->n) PetscFunctionReturn(PETSC_SUCCESS);
60 46 : PetscCall(MatDenseGetLDA(Q,&ldq));
61 46 : PetscCall(MatDenseGetArrayRead(Q,&q));
62 46 : PetscCall(BVMultInPlace_BLAS_Private(V,V->n,V->k-V->l,s-V->l,e-V->l,ctx->array+(V->nc+V->l)*V->ld,V->ld,q+V->l*ldq+V->l,ldq,PETSC_FALSE));
63 46 : PetscCall(MatDenseRestoreArrayRead(Q,&q));
64 46 : PetscFunctionReturn(PETSC_SUCCESS);
65 : }
66 :
67 1 : static PetscErrorCode BVMultInPlaceHermitianTranspose_Contiguous(BV V,Mat Q,PetscInt s,PetscInt e)
68 : {
69 1 : BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)V->data;
70 1 : const PetscScalar *q;
71 1 : PetscInt ldq;
72 :
73 1 : PetscFunctionBegin;
74 1 : if (s>=e || !V->n) PetscFunctionReturn(PETSC_SUCCESS);
75 1 : PetscCall(MatDenseGetLDA(Q,&ldq));
76 1 : PetscCall(MatDenseGetArrayRead(Q,&q));
77 1 : PetscCall(BVMultInPlace_BLAS_Private(V,V->n,V->k-V->l,s-V->l,e-V->l,ctx->array+(V->nc+V->l)*V->ld,V->ld,q+V->l*ldq+V->l,ldq,PETSC_TRUE));
78 1 : PetscCall(MatDenseRestoreArrayRead(Q,&q));
79 1 : PetscFunctionReturn(PETSC_SUCCESS);
80 : }
81 :
82 235 : static PetscErrorCode BVDot_Contiguous(BV X,BV Y,Mat M)
83 : {
84 235 : BV_CONTIGUOUS *x = (BV_CONTIGUOUS*)X->data,*y = (BV_CONTIGUOUS*)Y->data;
85 235 : PetscScalar *m;
86 235 : PetscInt ldm;
87 :
88 235 : PetscFunctionBegin;
89 235 : PetscCall(MatDenseGetLDA(M,&ldm));
90 235 : PetscCall(MatDenseGetArray(M,&m));
91 235 : PetscCall(BVDot_BLAS_Private(X,Y->k-Y->l,X->k-X->l,X->n,y->array+(Y->nc+Y->l)*Y->ld,Y->ld,x->array+(X->nc+X->l)*X->ld,X->ld,m+X->l*ldm+Y->l,ldm,x->mpi));
92 235 : PetscCall(MatDenseRestoreArray(M,&m));
93 235 : PetscFunctionReturn(PETSC_SUCCESS);
94 : }
95 :
96 351 : static PetscErrorCode BVDotVec_Contiguous(BV X,Vec y,PetscScalar *q)
97 : {
98 351 : BV_CONTIGUOUS *x = (BV_CONTIGUOUS*)X->data;
99 351 : const PetscScalar *py;
100 351 : PetscScalar *qq=q;
101 351 : Vec z = y;
102 :
103 351 : PetscFunctionBegin;
104 351 : if (PetscUnlikely(X->matrix)) {
105 57 : PetscCall(BV_IPMatMult(X,y));
106 57 : z = X->Bx;
107 : }
108 351 : PetscCall(VecGetArrayRead(z,&py));
109 351 : if (!q) PetscCall(VecGetArray(X->buffer,&qq));
110 351 : PetscCall(BVDotVec_BLAS_Private(X,X->n,X->k-X->l,x->array+(X->nc+X->l)*X->ld,X->ld,py,qq,x->mpi));
111 351 : if (!q) PetscCall(VecRestoreArray(X->buffer,&qq));
112 351 : PetscCall(VecRestoreArrayRead(z,&py));
113 351 : PetscFunctionReturn(PETSC_SUCCESS);
114 : }
115 :
116 4 : static PetscErrorCode BVDotVec_Local_Contiguous(BV X,Vec y,PetscScalar *m)
117 : {
118 4 : BV_CONTIGUOUS *x = (BV_CONTIGUOUS*)X->data;
119 4 : PetscScalar *py;
120 4 : Vec z = y;
121 :
122 4 : PetscFunctionBegin;
123 4 : if (PetscUnlikely(X->matrix)) {
124 0 : PetscCall(BV_IPMatMult(X,y));
125 0 : z = X->Bx;
126 : }
127 4 : PetscCall(VecGetArray(z,&py));
128 4 : PetscCall(BVDotVec_BLAS_Private(X,X->n,X->k-X->l,x->array+(X->nc+X->l)*X->ld,X->ld,py,m,PETSC_FALSE));
129 4 : PetscCall(VecRestoreArray(z,&py));
130 4 : PetscFunctionReturn(PETSC_SUCCESS);
131 : }
132 :
133 324 : static PetscErrorCode BVScale_Contiguous(BV bv,PetscInt j,PetscScalar alpha)
134 : {
135 324 : BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
136 :
137 324 : PetscFunctionBegin;
138 324 : if (!bv->n) PetscFunctionReturn(PETSC_SUCCESS);
139 324 : if (PetscUnlikely(j<0)) PetscCall(BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->ld,ctx->array+(bv->nc+bv->l)*bv->ld,alpha));
140 320 : else PetscCall(BVScale_BLAS_Private(bv,bv->n,ctx->array+(bv->nc+j)*bv->ld,alpha));
141 324 : PetscFunctionReturn(PETSC_SUCCESS);
142 : }
143 :
144 110 : static PetscErrorCode BVNorm_Contiguous(BV bv,PetscInt j,NormType type,PetscReal *val)
145 : {
146 110 : BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
147 :
148 110 : PetscFunctionBegin;
149 110 : if (PetscUnlikely(j<0)) PetscCall(BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,ctx->array+(bv->nc+bv->l)*bv->ld,bv->ld,type,val,ctx->mpi));
150 42 : else PetscCall(BVNorm_LAPACK_Private(bv,bv->n,1,ctx->array+(bv->nc+j)*bv->ld,bv->ld,type,val,ctx->mpi));
151 110 : PetscFunctionReturn(PETSC_SUCCESS);
152 : }
153 :
154 2 : static PetscErrorCode BVNorm_Local_Contiguous(BV bv,PetscInt j,NormType type,PetscReal *val)
155 : {
156 2 : BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
157 :
158 2 : PetscFunctionBegin;
159 2 : if (PetscUnlikely(j<0)) PetscCall(BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,ctx->array+(bv->nc+bv->l)*bv->ld,bv->ld,type,val,PETSC_FALSE));
160 2 : else PetscCall(BVNorm_LAPACK_Private(bv,bv->n,1,ctx->array+(bv->nc+j)*bv->ld,bv->ld,type,val,PETSC_FALSE));
161 2 : PetscFunctionReturn(PETSC_SUCCESS);
162 : }
163 :
164 6 : static PetscErrorCode BVNormalize_Contiguous(BV bv,PetscScalar *eigi)
165 : {
166 6 : BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
167 6 : PetscScalar *wi=NULL;
168 :
169 6 : PetscFunctionBegin;
170 6 : if (eigi) wi = eigi+bv->l;
171 6 : PetscCall(BVNormalize_LAPACK_Private(bv,bv->n,bv->k-bv->l,ctx->array+(bv->nc+bv->l)*bv->ld,bv->ld,wi,ctx->mpi));
172 6 : PetscFunctionReturn(PETSC_SUCCESS);
173 : }
174 :
175 107 : static PetscErrorCode BVMatMult_Contiguous(BV V,Mat A,BV W)
176 : {
177 107 : BV_CONTIGUOUS *v = (BV_CONTIGUOUS*)V->data,*w = (BV_CONTIGUOUS*)W->data;
178 107 : PetscInt j;
179 107 : Mat Vmat,Wmat;
180 :
181 107 : PetscFunctionBegin;
182 107 : if (V->vmm) {
183 100 : PetscCall(BVGetMat(V,&Vmat));
184 100 : PetscCall(BVGetMat(W,&Wmat));
185 100 : PetscCall(MatProductCreateWithMat(A,Vmat,NULL,Wmat));
186 100 : PetscCall(MatProductSetType(Wmat,MATPRODUCT_AB));
187 100 : PetscCall(MatProductSetFromOptions(Wmat));
188 100 : PetscCall(MatProductSymbolic(Wmat));
189 100 : PetscCall(MatProductNumeric(Wmat));
190 100 : PetscCall(MatProductClear(Wmat));
191 100 : PetscCall(BVRestoreMat(V,&Vmat));
192 100 : PetscCall(BVRestoreMat(W,&Wmat));
193 : } else {
194 50 : for (j=0;j<V->k-V->l;j++) PetscCall(MatMult(A,v->V[V->nc+V->l+j],w->V[W->nc+W->l+j]));
195 : }
196 107 : PetscFunctionReturn(PETSC_SUCCESS);
197 : }
198 :
199 97 : static PetscErrorCode BVCopy_Contiguous(BV V,BV W)
200 : {
201 97 : BV_CONTIGUOUS *v = (BV_CONTIGUOUS*)V->data,*w = (BV_CONTIGUOUS*)W->data;
202 97 : PetscInt j;
203 :
204 97 : PetscFunctionBegin;
205 726 : for (j=0;j<V->k-V->l;j++) PetscCall(PetscArraycpy(w->array+(W->nc+W->l+j)*W->ld,v->array+(V->nc+V->l+j)*V->ld,V->n));
206 97 : PetscFunctionReturn(PETSC_SUCCESS);
207 : }
208 :
209 2 : static PetscErrorCode BVCopyColumn_Contiguous(BV V,PetscInt j,PetscInt i)
210 : {
211 2 : BV_CONTIGUOUS *v = (BV_CONTIGUOUS*)V->data;
212 :
213 2 : PetscFunctionBegin;
214 2 : PetscCall(PetscArraycpy(v->array+(V->nc+i)*V->ld,v->array+(V->nc+j)*V->ld,V->n));
215 2 : PetscFunctionReturn(PETSC_SUCCESS);
216 : }
217 :
218 8 : static PetscErrorCode BVResize_Contiguous(BV bv,PetscInt m,PetscBool copy)
219 : {
220 8 : BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
221 8 : PetscInt j,bs;
222 8 : PetscScalar *newarray;
223 8 : Vec *newV;
224 8 : char str[50];
225 :
226 8 : PetscFunctionBegin;
227 8 : PetscCall(PetscLayoutGetBlockSize(bv->map,&bs));
228 8 : PetscCall(PetscCalloc1(m*bv->ld,&newarray));
229 8 : PetscCall(PetscMalloc1(m,&newV));
230 94 : for (j=0;j<m;j++) {
231 86 : if (ctx->mpi) PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv),bs,bv->n,PETSC_DECIDE,newarray+j*bv->ld,newV+j));
232 86 : else PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv),bs,bv->n,newarray+j*bv->ld,newV+j));
233 : }
234 8 : if (((PetscObject)bv)->name) {
235 94 : for (j=0;j<m;j++) {
236 86 : PetscCall(PetscSNPrintf(str,sizeof(str),"%s_%" PetscInt_FMT,((PetscObject)bv)->name,j));
237 86 : PetscCall(PetscObjectSetName((PetscObject)newV[j],str));
238 : }
239 : }
240 8 : if (copy) PetscCall(PetscArraycpy(newarray,ctx->array,PetscMin(m,bv->m)*bv->ld));
241 8 : PetscCall(VecDestroyVecs(bv->m,&ctx->V));
242 8 : ctx->V = newV;
243 8 : PetscCall(PetscFree(ctx->array));
244 8 : ctx->array = newarray;
245 8 : PetscFunctionReturn(PETSC_SUCCESS);
246 : }
247 :
248 2631 : static PetscErrorCode BVGetColumn_Contiguous(BV bv,PetscInt j,Vec *v)
249 : {
250 2631 : BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
251 2631 : PetscInt l;
252 :
253 2631 : PetscFunctionBegin;
254 2631 : l = BVAvailableVec;
255 2631 : bv->cv[l] = ctx->V[bv->nc+j];
256 2631 : PetscFunctionReturn(PETSC_SUCCESS);
257 : }
258 :
259 2631 : static PetscErrorCode BVRestoreColumn_Contiguous(BV bv,PetscInt j,Vec *v)
260 : {
261 2631 : PetscInt l;
262 :
263 2631 : PetscFunctionBegin;
264 2631 : l = (j==bv->ci[0])? 0: 1;
265 2631 : bv->cv[l] = NULL;
266 2631 : PetscFunctionReturn(PETSC_SUCCESS);
267 : }
268 :
269 251 : static PetscErrorCode BVGetArray_Contiguous(BV bv,PetscScalar **a)
270 : {
271 251 : BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
272 :
273 251 : PetscFunctionBegin;
274 251 : *a = ctx->array;
275 251 : PetscFunctionReturn(PETSC_SUCCESS);
276 : }
277 :
278 8 : static PetscErrorCode BVGetArrayRead_Contiguous(BV bv,const PetscScalar **a)
279 : {
280 8 : BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
281 :
282 8 : PetscFunctionBegin;
283 8 : *a = ctx->array;
284 8 : PetscFunctionReturn(PETSC_SUCCESS);
285 : }
286 :
287 339 : static PetscErrorCode BVDestroy_Contiguous(BV bv)
288 : {
289 339 : BV_CONTIGUOUS *ctx = (BV_CONTIGUOUS*)bv->data;
290 :
291 339 : PetscFunctionBegin;
292 339 : if (!bv->issplit) {
293 259 : PetscCall(VecDestroyVecs(bv->nc+bv->m,&ctx->V));
294 259 : PetscCall(PetscFree(ctx->array));
295 : }
296 339 : PetscCall(PetscFree(bv->data));
297 339 : PetscFunctionReturn(PETSC_SUCCESS);
298 : }
299 :
300 12 : static PetscErrorCode BVView_Contiguous(BV bv,PetscViewer viewer)
301 : {
302 12 : PetscInt j;
303 12 : Vec v;
304 12 : PetscViewerFormat format;
305 12 : PetscBool isascii,ismatlab=PETSC_FALSE;
306 12 : const char *bvname,*name;
307 :
308 12 : PetscFunctionBegin;
309 12 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
310 12 : if (isascii) {
311 12 : PetscCall(PetscViewerGetFormat(viewer,&format));
312 12 : if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(PETSC_SUCCESS);
313 10 : if (format == PETSC_VIEWER_ASCII_MATLAB) ismatlab = PETSC_TRUE;
314 : }
315 0 : if (ismatlab) {
316 0 : PetscCall(PetscObjectGetName((PetscObject)bv,&bvname));
317 0 : PetscCall(PetscViewerASCIIPrintf(viewer,"%s=[];\n",bvname));
318 : }
319 48 : for (j=0;j<bv->m;j++) {
320 38 : PetscCall(BVGetColumn(bv,j,&v));
321 38 : PetscCall(VecView(v,viewer));
322 38 : if (ismatlab) {
323 0 : PetscCall(PetscObjectGetName((PetscObject)v,&name));
324 0 : PetscCall(PetscViewerASCIIPrintf(viewer,"%s=[%s,%s];clear %s\n",bvname,bvname,name,name));
325 : }
326 38 : PetscCall(BVRestoreColumn(bv,j,&v));
327 : }
328 10 : PetscFunctionReturn(PETSC_SUCCESS);
329 : }
330 :
331 339 : SLEPC_EXTERN PetscErrorCode BVCreate_Contiguous(BV bv)
332 : {
333 339 : BV_CONTIGUOUS *ctx;
334 339 : PetscInt j,nloc,bs,lsplit,lda;
335 339 : PetscBool seq,isdense;
336 339 : PetscScalar *aa;
337 339 : char str[50];
338 339 : PetscScalar *array;
339 339 : BV parent;
340 339 : Vec *Vpar;
341 339 : MatType mtype;
342 :
343 339 : PetscFunctionBegin;
344 339 : PetscCall(PetscNew(&ctx));
345 339 : bv->data = (void*)ctx;
346 :
347 339 : PetscCall(PetscStrcmp(bv->vtype,VECMPI,&ctx->mpi));
348 339 : if (!ctx->mpi) {
349 61 : PetscCall(PetscStrcmp(bv->vtype,VECSEQ,&seq));
350 61 : PetscCheck(seq,PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Cannot create a contiguous BV from a non-standard vector type: %s",bv->vtype);
351 : }
352 :
353 339 : PetscCall(PetscLayoutGetLocalSize(bv->map,&nloc));
354 339 : PetscCall(PetscLayoutGetBlockSize(bv->map,&bs));
355 339 : PetscCall(BV_SetDefaultLD(bv,nloc));
356 :
357 339 : if (PetscUnlikely(bv->issplit)) {
358 80 : PetscCheck(bv->issplit>0,PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"BVCONTIGUOUS does not support BVGetSplitRows()");
359 : /* split BV: share memory and Vecs of the parent BV */
360 80 : parent = bv->splitparent;
361 80 : lsplit = parent->lsplit;
362 80 : Vpar = ((BV_CONTIGUOUS*)parent->data)->V;
363 80 : ctx->V = (bv->issplit==1)? Vpar: Vpar+lsplit;
364 80 : array = ((BV_CONTIGUOUS*)parent->data)->array;
365 80 : ctx->array = (bv->issplit==1)? array: array+lsplit*bv->ld;
366 : } else {
367 : /* regular BV: allocate memory and Vecs for the BV entries */
368 259 : PetscCall(PetscCalloc1(bv->m*bv->ld,&ctx->array));
369 259 : PetscCall(PetscMalloc1(bv->m,&ctx->V));
370 2233 : for (j=0;j<bv->m;j++) {
371 1974 : if (ctx->mpi) PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv),bs,nloc,PETSC_DECIDE,ctx->array+j*bv->ld,ctx->V+j));
372 1974 : else PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv),bs,nloc,ctx->array+j*bv->ld,ctx->V+j));
373 : }
374 : }
375 339 : if (((PetscObject)bv)->name) {
376 991 : for (j=0;j<bv->m;j++) {
377 878 : PetscCall(PetscSNPrintf(str,sizeof(str),"%s_%" PetscInt_FMT,((PetscObject)bv)->name,j));
378 878 : PetscCall(PetscObjectSetName((PetscObject)ctx->V[j],str));
379 : }
380 : }
381 :
382 339 : if (PetscUnlikely(bv->Acreate)) {
383 2 : PetscCall(MatGetType(bv->Acreate,&mtype));
384 2 : PetscCall(PetscStrcmpAny(mtype,&isdense,MATSEQDENSE,MATMPIDENSE,""));
385 2 : PetscCheck(isdense,PetscObjectComm((PetscObject)bv->Acreate),PETSC_ERR_SUP,"BVCONTIGUOUS requires a dense matrix in BVCreateFromMat()");
386 2 : PetscCall(MatDenseGetArray(bv->Acreate,&aa));
387 2 : PetscCall(MatDenseGetLDA(bv->Acreate,&lda));
388 18 : for (j=0;j<bv->m;j++) PetscCall(PetscArraycpy(ctx->array+j*bv->ld,aa+j*lda,bv->n));
389 2 : PetscCall(MatDenseRestoreArray(bv->Acreate,&aa));
390 2 : PetscCall(MatDestroy(&bv->Acreate));
391 : }
392 :
393 339 : bv->ops->mult = BVMult_Contiguous;
394 339 : bv->ops->multvec = BVMultVec_Contiguous;
395 339 : bv->ops->multinplace = BVMultInPlace_Contiguous;
396 339 : bv->ops->multinplacetrans = BVMultInPlaceHermitianTranspose_Contiguous;
397 339 : bv->ops->dot = BVDot_Contiguous;
398 339 : bv->ops->dotvec = BVDotVec_Contiguous;
399 339 : bv->ops->dotvec_local = BVDotVec_Local_Contiguous;
400 339 : bv->ops->scale = BVScale_Contiguous;
401 339 : bv->ops->norm = BVNorm_Contiguous;
402 339 : bv->ops->norm_local = BVNorm_Local_Contiguous;
403 339 : bv->ops->normalize = BVNormalize_Contiguous;
404 339 : bv->ops->matmult = BVMatMult_Contiguous;
405 339 : bv->ops->copy = BVCopy_Contiguous;
406 339 : bv->ops->copycolumn = BVCopyColumn_Contiguous;
407 339 : bv->ops->resize = BVResize_Contiguous;
408 339 : bv->ops->getcolumn = BVGetColumn_Contiguous;
409 339 : bv->ops->restorecolumn = BVRestoreColumn_Contiguous;
410 339 : bv->ops->getarray = BVGetArray_Contiguous;
411 339 : bv->ops->getarrayread = BVGetArrayRead_Contiguous;
412 339 : bv->ops->getmat = BVGetMat_Default;
413 339 : bv->ops->restoremat = BVRestoreMat_Default;
414 339 : bv->ops->destroy = BVDestroy_Contiguous;
415 339 : bv->ops->view = BVView_Contiguous;
416 339 : PetscFunctionReturn(PETSC_SUCCESS);
417 : }
|