Actual source code: svec.c

slepc-3.14.0 2020-09-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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;
 22:   PetscScalar       *py,*q;
 23:   PetscInt          ldq;

 26:   VecGetArrayRead(x->v,&px);
 27:   VecGetArray(y->v,&py);
 28:   if (Q) {
 29:     MatGetSize(Q,&ldq,NULL);
 30:     MatDenseGetArray(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:     MatDenseRestoreArray(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: {
 61:   BV_SVEC        *ctx = (BV_SVEC*)V->data;
 62:   PetscScalar    *pv,*q;
 63:   PetscInt       ldq;

 66:   MatGetSize(Q,&ldq,NULL);
 67:   VecGetArray(ctx->v,&pv);
 68:   MatDenseGetArray(Q,&q);
 69:   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);
 70:   MatDenseRestoreArray(Q,&q);
 71:   VecRestoreArray(ctx->v,&pv);
 72:   return(0);
 73: }

 75: PetscErrorCode BVMultInPlaceTranspose_Svec(BV V,Mat Q,PetscInt s,PetscInt e)
 76: {
 78:   BV_SVEC        *ctx = (BV_SVEC*)V->data;
 79:   PetscScalar    *pv,*q;
 80:   PetscInt       ldq;

 83:   MatGetSize(Q,&ldq,NULL);
 84:   VecGetArray(ctx->v,&pv);
 85:   MatDenseGetArray(Q,&q);
 86:   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);
 87:   MatDenseRestoreArray(Q,&q);
 88:   VecRestoreArray(ctx->v,&pv);
 89:   return(0);
 90: }

 92: PetscErrorCode BVDot_Svec(BV X,BV Y,Mat M)
 93: {
 94:   PetscErrorCode    ierr;
 95:   BV_SVEC           *x = (BV_SVEC*)X->data,*y = (BV_SVEC*)Y->data;
 96:   const PetscScalar *px,*py;
 97:   PetscScalar       *m;
 98:   PetscInt          ldm;

101:   MatGetSize(M,&ldm,NULL);
102:   VecGetArrayRead(x->v,&px);
103:   VecGetArrayRead(y->v,&py);
104:   MatDenseGetArray(M,&m);
105:   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);
106:   MatDenseRestoreArray(M,&m);
107:   VecRestoreArrayRead(x->v,&px);
108:   VecRestoreArrayRead(y->v,&py);
109:   return(0);
110: }

112: PetscErrorCode BVDotVec_Svec(BV X,Vec y,PetscScalar *q)
113: {
114:   PetscErrorCode    ierr;
115:   BV_SVEC           *x = (BV_SVEC*)X->data;
116:   const PetscScalar *px,*py;
117:   PetscScalar       *qq=q;
118:   Vec               z = y;

121:   if (X->matrix) {
122:     BV_IPMatMult(X,y);
123:     z = X->Bx;
124:   }
125:   VecGetArrayRead(x->v,&px);
126:   VecGetArrayRead(z,&py);
127:   if (!q) { VecGetArray(X->buffer,&qq); }
128:   BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,qq,x->mpi);
129:   if (!q) { VecRestoreArray(X->buffer,&qq); }
130:   VecRestoreArrayRead(z,&py);
131:   VecRestoreArrayRead(x->v,&px);
132:   return(0);
133: }

135: PetscErrorCode BVDotVec_Local_Svec(BV X,Vec y,PetscScalar *m)
136: {
138:   BV_SVEC        *x = (BV_SVEC*)X->data;
139:   PetscScalar    *px,*py;
140:   Vec            z = y;

143:   if (X->matrix) {
144:     BV_IPMatMult(X,y);
145:     z = X->Bx;
146:   }
147:   VecGetArray(x->v,&px);
148:   VecGetArray(z,&py);
149:   BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->n,py,m,PETSC_FALSE);
150:   VecRestoreArray(z,&py);
151:   VecRestoreArray(x->v,&px);
152:   return(0);
153: }

155: PetscErrorCode BVScale_Svec(BV bv,PetscInt j,PetscScalar alpha)
156: {
158:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;
159:   PetscScalar    *array;

162:   VecGetArray(ctx->v,&array);
163:   if (j<0) {
164:     BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->n,array+(bv->nc+bv->l)*bv->n,alpha);
165:   } else {
166:     BVScale_BLAS_Private(bv,bv->n,array+(bv->nc+j)*bv->n,alpha);
167:   }
168:   VecRestoreArray(ctx->v,&array);
169:   return(0);
170: }

172: PetscErrorCode BVNorm_Svec(BV bv,PetscInt j,NormType type,PetscReal *val)
173: {
175:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;
176:   PetscScalar    *array;

179:   VecGetArray(ctx->v,&array);
180:   if (j<0) {
181:     BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,ctx->mpi);
182:   } else {
183:     BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,ctx->mpi);
184:   }
185:   VecRestoreArray(ctx->v,&array);
186:   return(0);
187: }

189: PetscErrorCode BVNorm_Local_Svec(BV bv,PetscInt j,NormType type,PetscReal *val)
190: {
192:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;
193:   PetscScalar    *array;

196:   VecGetArray(ctx->v,&array);
197:   if (j<0) {
198:     BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,type,val,PETSC_FALSE);
199:   } else {
200:     BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->n,type,val,PETSC_FALSE);
201:   }
202:   VecRestoreArray(ctx->v,&array);
203:   return(0);
204: }

206: PetscErrorCode BVNormalize_Svec(BV bv,PetscScalar *eigi)
207: {
209:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;
210:   PetscScalar    *array,*wi=NULL;

213:   VecGetArray(ctx->v,&array);
214:   if (eigi) wi = eigi+bv->l;
215:   BVNormalize_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->n,wi,ctx->mpi);
216:   VecRestoreArray(ctx->v,&array);
217:   return(0);
218: }

220: PetscErrorCode BVMatMult_Svec(BV V,Mat A,BV W)
221: {
223:   PetscInt       j;
224:   Mat            Vmat,Wmat;
225:   Vec            vv,ww;

228:   if (V->vmm) {
229:     BVGetMat(V,&Vmat);
230:     BVGetMat(W,&Wmat);
231:     MatProductCreateWithMat(A,Vmat,NULL,Wmat);
232:     MatProductSetType(Wmat,MATPRODUCT_AB);
233:     MatProductSetFromOptions(Wmat);
234:     MatProductSymbolic(Wmat);
235:     MatProductNumeric(Wmat);
236:     MatProductClear(Wmat);
237:     BVRestoreMat(V,&Vmat);
238:     BVRestoreMat(W,&Wmat);
239:   } else {
240:     for (j=0;j<V->k-V->l;j++) {
241:       BVGetColumn(V,V->l+j,&vv);
242:       BVGetColumn(W,W->l+j,&ww);
243:       MatMult(A,vv,ww);
244:       BVRestoreColumn(V,V->l+j,&vv);
245:       BVRestoreColumn(W,W->l+j,&ww);
246:     }
247:   }
248:   return(0);
249: }

251: PetscErrorCode BVCopy_Svec(BV V,BV W)
252: {
254:   BV_SVEC        *v = (BV_SVEC*)V->data,*w = (BV_SVEC*)W->data;
255:   PetscScalar    *pv,*pw,*pvc,*pwc;

258:   VecGetArray(v->v,&pv);
259:   VecGetArray(w->v,&pw);
260:   pvc = pv+(V->nc+V->l)*V->n;
261:   pwc = pw+(W->nc+W->l)*W->n;
262:   PetscArraycpy(pwc,pvc,(V->k-V->l)*V->n);
263:   VecRestoreArray(v->v,&pv);
264:   VecRestoreArray(w->v,&pw);
265:   return(0);
266: }

268: PetscErrorCode BVCopyColumn_Svec(BV V,PetscInt j,PetscInt i)
269: {
271:   BV_SVEC        *v = (BV_SVEC*)V->data;
272:   PetscScalar    *pv;

275:   VecGetArray(v->v,&pv);
276:   PetscArraycpy(pv+(V->nc+i)*V->n,pv+(V->nc+j)*V->n,V->n);
277:   VecRestoreArray(v->v,&pv);
278:   return(0);
279: }

281: PetscErrorCode BVResize_Svec(BV bv,PetscInt m,PetscBool copy)
282: {
283:   PetscErrorCode    ierr;
284:   BV_SVEC           *ctx = (BV_SVEC*)bv->data;
285:   PetscScalar       *pnew;
286:   const PetscScalar *pv;
287:   PetscInt          bs;
288:   Vec               vnew;
289:   char              str[50];

292:   VecGetBlockSize(bv->t,&bs);
293:   VecCreate(PetscObjectComm((PetscObject)bv->t),&vnew);
294:   VecSetType(vnew,((PetscObject)bv->t)->type_name);
295:   VecSetSizes(vnew,m*bv->n,PETSC_DECIDE);
296:   VecSetBlockSize(vnew,bs);
297:   PetscLogObjectParent((PetscObject)bv,(PetscObject)vnew);
298:   if (((PetscObject)bv)->name) {
299:     PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name);
300:     PetscObjectSetName((PetscObject)vnew,str);
301:   }
302:   if (copy) {
303:     VecGetArrayRead(ctx->v,&pv);
304:     VecGetArray(vnew,&pnew);
305:     PetscArraycpy(pnew,pv,PetscMin(m,bv->m)*bv->n);
306:     VecRestoreArrayRead(ctx->v,&pv);
307:     VecRestoreArray(vnew,&pnew);
308:   }
309:   VecDestroy(&ctx->v);
310:   ctx->v = vnew;
311:   return(0);
312: }

314: PetscErrorCode BVGetColumn_Svec(BV bv,PetscInt j,Vec *v)
315: {
317:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;
318:   PetscScalar    *pv;
319:   PetscInt       l;

322:   l = BVAvailableVec;
323:   VecGetArray(ctx->v,&pv);
324:   VecPlaceArray(bv->cv[l],pv+(bv->nc+j)*bv->n);
325:   return(0);
326: }

328: PetscErrorCode BVRestoreColumn_Svec(BV bv,PetscInt j,Vec *v)
329: {
331:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;
332:   PetscInt       l;

335:   l = (j==bv->ci[0])? 0: 1;
336:   VecResetArray(bv->cv[l]);
337:   VecRestoreArray(ctx->v,NULL);
338:   return(0);
339: }

341: PetscErrorCode BVGetArray_Svec(BV bv,PetscScalar **a)
342: {
344:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;

347:   VecGetArray(ctx->v,a);
348:   return(0);
349: }

351: PetscErrorCode BVRestoreArray_Svec(BV bv,PetscScalar **a)
352: {
354:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;

357:   VecRestoreArray(ctx->v,a);
358:   return(0);
359: }

361: PetscErrorCode BVGetArrayRead_Svec(BV bv,const PetscScalar **a)
362: {
364:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;

367:   VecGetArrayRead(ctx->v,a);
368:   return(0);
369: }

371: PetscErrorCode BVRestoreArrayRead_Svec(BV bv,const PetscScalar **a)
372: {
374:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;

377:   VecRestoreArrayRead(ctx->v,a);
378:   return(0);
379: }

381: PetscErrorCode BVView_Svec(BV bv,PetscViewer viewer)
382: {
383:   PetscErrorCode    ierr;
384:   BV_SVEC           *ctx = (BV_SVEC*)bv->data;
385:   PetscViewerFormat format;
386:   PetscBool         isascii;
387:   const char        *bvname,*name;

390:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
391:   if (isascii) {
392:     PetscViewerGetFormat(viewer,&format);
393:     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
394:     VecView(ctx->v,viewer);
395:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
396:       PetscObjectGetName((PetscObject)bv,&bvname);
397:       PetscObjectGetName((PetscObject)ctx->v,&name);
398:       PetscViewerASCIIPrintf(viewer,"%s=reshape(%s,%D,%D);clear %s\n",bvname,name,bv->N,bv->nc+bv->m,name);
399:       if (bv->nc) {
400:         PetscViewerASCIIPrintf(viewer,"%s=%s(:,%D:end);\n",bvname,bvname,bv->nc+1);
401:       }
402:     }
403:   } else {
404:     VecView(ctx->v,viewer);
405:   }
406:   return(0);
407: }

409: PetscErrorCode BVDestroy_Svec(BV bv)
410: {
412:   BV_SVEC        *ctx = (BV_SVEC*)bv->data;

415:   VecDestroy(&ctx->v);
416:   VecDestroy(&bv->cv[0]);
417:   VecDestroy(&bv->cv[1]);
418:   PetscFree(bv->data);
419:   bv->cuda = PETSC_FALSE;
420:   return(0);
421: }

423: SLEPC_EXTERN PetscErrorCode BVCreate_Svec(BV bv)
424: {
425:   PetscErrorCode    ierr;
426:   BV_SVEC           *ctx;
427:   PetscInt          nloc,N,bs,tglobal=0,tlocal,lsplit;
428:   PetscBool         seq;
429:   PetscScalar       *aa,*vv;
430:   const PetscScalar *array,*ptr;
431:   char              str[50];
432:   BV                parent;
433:   Vec               vpar;
434: #if defined(PETSC_HAVE_CUDA)
435:   PetscScalar       *gpuarray,*gptr;
436: #endif

439:   PetscNewLog(bv,&ctx);
440:   bv->data = (void*)ctx;

442:   PetscObjectTypeCompareAny((PetscObject)bv->t,&bv->cuda,VECSEQCUDA,VECMPICUDA,"");
443:   PetscObjectTypeCompareAny((PetscObject)bv->t,&ctx->mpi,VECMPI,VECMPICUDA,"");

445:   PetscObjectTypeCompare((PetscObject)bv->t,VECSEQ,&seq);
446:   if (!seq && !ctx->mpi && !bv->cuda) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"BVSVEC does not support the type of the provided template vector");

448:   VecGetLocalSize(bv->t,&nloc);
449:   VecGetSize(bv->t,&N);
450:   VecGetBlockSize(bv->t,&bs);
451:   tlocal  = bv->m*nloc;
452:   PetscIntMultError(bv->m,N,&tglobal);

454:   if (bv->issplit) {
455:     /* split BV: create Vec sharing the memory of the parent BV */
456:     parent = bv->splitparent;
457:     lsplit = parent->lsplit;
458:     vpar = ((BV_SVEC*)parent->data)->v;
459:     if (bv->cuda) {
460: #if defined(PETSC_HAVE_CUDA)
461:       VecCUDAGetArray(vpar,&gpuarray);
462:       gptr = (bv->issplit==1)? gpuarray: gpuarray+lsplit*nloc;
463:       VecCUDARestoreArray(vpar,&gpuarray);
464:       if (ctx->mpi) {
465:         VecCreateMPICUDAWithArray(PetscObjectComm((PetscObject)bv->t),bs,tlocal,bv->m*N,NULL,&ctx->v);
466:       } else {
467:         VecCreateSeqCUDAWithArray(PetscObjectComm((PetscObject)bv->t),bs,tlocal,NULL,&ctx->v);
468:       }
469:       VecCUDAPlaceArray(ctx->v,gptr);
470: #endif
471:     } else {
472:       VecGetArrayRead(vpar,&array);
473:       ptr = (bv->issplit==1)? array: array+lsplit*nloc;
474:       VecRestoreArrayRead(vpar,&array);
475:       if (ctx->mpi) {
476:         VecCreateMPIWithArray(PetscObjectComm((PetscObject)bv->t),bs,tlocal,bv->m*N,NULL,&ctx->v);
477:       } else {
478:         VecCreateSeqWithArray(PetscObjectComm((PetscObject)bv->t),bs,tlocal,NULL,&ctx->v);
479:       }
480:       VecPlaceArray(ctx->v,ptr);
481:     }
482:   } else {
483:     /* regular BV: create Vec to store the BV entries */
484:     VecCreate(PetscObjectComm((PetscObject)bv->t),&ctx->v);
485:     VecSetType(ctx->v,((PetscObject)bv->t)->type_name);
486:     VecSetSizes(ctx->v,tlocal,tglobal);
487:     VecSetBlockSize(ctx->v,bs);
488:   }
489:   PetscLogObjectParent((PetscObject)bv,(PetscObject)ctx->v);
490:   if (((PetscObject)bv)->name) {
491:     PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name);
492:     PetscObjectSetName((PetscObject)ctx->v,str);
493:   }

495:   if (bv->Acreate) {
496:     MatDenseGetArray(bv->Acreate,&aa);
497:     VecGetArray(ctx->v,&vv);
498:     PetscArraycpy(vv,aa,tlocal);
499:     VecRestoreArray(ctx->v,&vv);
500:     MatDenseRestoreArray(bv->Acreate,&aa);
501:     MatDestroy(&bv->Acreate);
502:   }

504:   VecDuplicateEmpty(bv->t,&bv->cv[0]);
505:   VecDuplicateEmpty(bv->t,&bv->cv[1]);

507:   if (bv->cuda) {
508: #if defined(PETSC_HAVE_CUDA)
509:     bv->ops->mult             = BVMult_Svec_CUDA;
510:     bv->ops->multvec          = BVMultVec_Svec_CUDA;
511:     bv->ops->multinplace      = BVMultInPlace_Svec_CUDA;
512:     bv->ops->multinplacetrans = BVMultInPlaceTranspose_Svec_CUDA;
513:     bv->ops->dot              = BVDot_Svec_CUDA;
514:     bv->ops->dotvec           = BVDotVec_Svec_CUDA;
515:     bv->ops->dotvec_local     = BVDotVec_Local_Svec_CUDA;
516:     bv->ops->scale            = BVScale_Svec_CUDA;
517:     bv->ops->matmult          = BVMatMult_Svec_CUDA;
518:     bv->ops->copy             = BVCopy_Svec_CUDA;
519:     bv->ops->copycolumn       = BVCopyColumn_Svec_CUDA;
520:     bv->ops->resize           = BVResize_Svec_CUDA;
521:     bv->ops->getcolumn        = BVGetColumn_Svec_CUDA;
522:     bv->ops->restorecolumn    = BVRestoreColumn_Svec_CUDA;
523:     bv->ops->restoresplit     = BVRestoreSplit_Svec_CUDA;
524:     bv->ops->getmat           = BVGetMat_Svec_CUDA;
525:     bv->ops->restoremat       = BVRestoreMat_Svec_CUDA;
526: #endif
527:   } else {
528:     bv->ops->mult             = BVMult_Svec;
529:     bv->ops->multvec          = BVMultVec_Svec;
530:     bv->ops->multinplace      = BVMultInPlace_Svec;
531:     bv->ops->multinplacetrans = BVMultInPlaceTranspose_Svec;
532:     bv->ops->dot              = BVDot_Svec;
533:     bv->ops->dotvec           = BVDotVec_Svec;
534:     bv->ops->dotvec_local     = BVDotVec_Local_Svec;
535:     bv->ops->scale            = BVScale_Svec;
536:     bv->ops->matmult          = BVMatMult_Svec;
537:     bv->ops->copy             = BVCopy_Svec;
538:     bv->ops->copycolumn       = BVCopyColumn_Svec;
539:     bv->ops->resize           = BVResize_Svec;
540:     bv->ops->getcolumn        = BVGetColumn_Svec;
541:     bv->ops->restorecolumn    = BVRestoreColumn_Svec;
542:   }
543:   bv->ops->norm             = BVNorm_Svec;
544:   bv->ops->norm_local       = BVNorm_Local_Svec;
545:   bv->ops->normalize        = BVNormalize_Svec;
546:   bv->ops->getarray         = BVGetArray_Svec;
547:   bv->ops->restorearray     = BVRestoreArray_Svec;
548:   bv->ops->getarrayread     = BVGetArrayRead_Svec;
549:   bv->ops->restorearrayread = BVRestoreArrayRead_Svec;
550:   bv->ops->destroy          = BVDestroy_Svec;
551:   if (!ctx->mpi) bv->ops->view = BVView_Svec;
552:   return(0);
553: }