Actual source code: bvfunc.c

slepc-3.13.0 2020-03-31
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 (basis vectors) interface routines, callable by users
 12: */

 14: #include <slepc/private/bvimpl.h>            /*I "slepcbv.h" I*/

 16: PetscClassId     BV_CLASSID = 0;
 17: PetscLogEvent    BV_Create = 0,BV_Copy = 0,BV_Mult = 0,BV_MultVec = 0,BV_MultInPlace = 0,BV_Dot = 0,BV_DotVec = 0,BV_Orthogonalize = 0,BV_OrthogonalizeVec = 0,BV_Scale = 0,BV_Norm = 0,BV_NormVec = 0,BV_SetRandom = 0,BV_MatMult = 0,BV_MatMultVec = 0,BV_MatProject = 0;
 18: static PetscBool BVPackageInitialized = PETSC_FALSE;
 19: MPI_Op MPIU_TSQR = 0;

 21: const char *BVOrthogTypes[] = {"CGS","MGS","BVOrthogType","BV_ORTHOG_",0};
 22: const char *BVOrthogRefineTypes[] = {"IFNEEDED","NEVER","ALWAYS","BVOrthogRefineType","BV_ORTHOG_REFINE_",0};
 23: const char *BVOrthogBlockTypes[] = {"GS","CHOL","TSQR","TSQRCHOL","SVQB","BVOrthogBlockType","BV_ORTHOG_BLOCK_",0};
 24: const char *BVMatMultTypes[] = {"VECS","MAT","MAT_SAVE","BVMatMultType","BV_MATMULT_",0};

 26: /*@C
 27:    BVFinalizePackage - This function destroys everything in the Slepc interface
 28:    to the BV package. It is called from SlepcFinalize().

 30:    Level: developer

 32: .seealso: SlepcFinalize()
 33: @*/
 34: PetscErrorCode BVFinalizePackage(void)
 35: {

 39:   PetscFunctionListDestroy(&BVList);
 40:   MPI_Op_free(&MPIU_TSQR);
 41:   BVPackageInitialized = PETSC_FALSE;
 42:   BVRegisterAllCalled  = PETSC_FALSE;
 43:   return(0);
 44: }

 46: /*@C
 47:    BVInitializePackage - This function initializes everything in the BV package.
 48:    It is called from PetscDLLibraryRegister() when using dynamic libraries, and
 49:    on the first call to BVCreate() when using static libraries.

 51:    Level: developer

 53: .seealso: SlepcInitialize()
 54: @*/
 55: PetscErrorCode BVInitializePackage(void)
 56: {
 57:   char           logList[256];
 58:   PetscBool      opt,pkg;
 59:   PetscClassId   classids[1];

 63:   if (BVPackageInitialized) return(0);
 64:   BVPackageInitialized = PETSC_TRUE;
 65:   /* Register Classes */
 66:   PetscClassIdRegister("Basis Vectors",&BV_CLASSID);
 67:   /* Register Constructors */
 68:   BVRegisterAll();
 69:   /* Register Events */
 70:   PetscLogEventRegister("BVCreate",BV_CLASSID,&BV_Create);
 71:   PetscLogEventRegister("BVCopy",BV_CLASSID,&BV_Copy);
 72:   PetscLogEventRegister("BVMult",BV_CLASSID,&BV_Mult);
 73:   PetscLogEventRegister("BVMultVec",BV_CLASSID,&BV_MultVec);
 74:   PetscLogEventRegister("BVMultInPlace",BV_CLASSID,&BV_MultInPlace);
 75:   PetscLogEventRegister("BVDot",BV_CLASSID,&BV_Dot);
 76:   PetscLogEventRegister("BVDotVec",BV_CLASSID,&BV_DotVec);
 77:   PetscLogEventRegister("BVOrthogonalize",BV_CLASSID,&BV_Orthogonalize);
 78:   PetscLogEventRegister("BVOrthogonalizeV",BV_CLASSID,&BV_OrthogonalizeVec);
 79:   PetscLogEventRegister("BVScale",BV_CLASSID,&BV_Scale);
 80:   PetscLogEventRegister("BVNorm",BV_CLASSID,&BV_Norm);
 81:   PetscLogEventRegister("BVNormVec",BV_CLASSID,&BV_NormVec);
 82:   PetscLogEventRegister("BVSetRandom",BV_CLASSID,&BV_SetRandom);
 83:   PetscLogEventRegister("BVMatMult",BV_CLASSID,&BV_MatMult);
 84:   PetscLogEventRegister("BVMatMultVec",BV_CLASSID,&BV_MatMultVec);
 85:   PetscLogEventRegister("BVMatProject",BV_CLASSID,&BV_MatProject);
 86:   /* MPI reduction operation used in BVOrthogonalize */
 87:   MPI_Op_create(SlepcGivensPacked,PETSC_TRUE,&MPIU_TSQR);
 88:   /* Process Info */
 89:   classids[0] = BV_CLASSID;
 90:   PetscInfoProcessClass("bv",1,&classids[0]);
 91:   /* Process summary exclusions */
 92:   PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt);
 93:   if (opt) {
 94:     PetscStrInList("bv",logList,',',&pkg);
 95:     if (pkg) { PetscLogEventDeactivateClass(BV_CLASSID); }
 96:   }
 97:   /* Register package finalizer */
 98:   PetscRegisterFinalize(BVFinalizePackage);
 99:   return(0);
100: }

102: /*@
103:    BVDestroy - Destroys BV context that was created with BVCreate().

105:    Collective on bv

107:    Input Parameter:
108: .  bv - the basis vectors context

110:    Level: beginner

112: .seealso: BVCreate()
113: @*/
114: PetscErrorCode BVDestroy(BV *bv)
115: {

119:   if (!*bv) return(0);
121:   if ((*bv)->lsplit) SETERRQ(PetscObjectComm((PetscObject)(*bv)),PETSC_ERR_ARG_WRONGSTATE,"Must call BVRestoreSplit before destroying the BV");
122:   if (--((PetscObject)(*bv))->refct > 0) { *bv = 0; return(0); }
123:   if ((*bv)->ops->destroy) { (*(*bv)->ops->destroy)(*bv); }
124:   VecDestroy(&(*bv)->t);
125:   MatDestroy(&(*bv)->matrix);
126:   VecDestroy(&(*bv)->Bx);
127:   VecDestroy(&(*bv)->buffer);
128:   BVDestroy(&(*bv)->cached);
129:   BVDestroy(&(*bv)->L);
130:   BVDestroy(&(*bv)->R);
131:   PetscFree((*bv)->work);
132:   PetscFree2((*bv)->h,(*bv)->c);
133:   VecDestroy(&(*bv)->omega);
134:   MatDestroy(&(*bv)->Acreate);
135:   MatDestroy(&(*bv)->Aget);
136:   MatDestroy(&(*bv)->Abuffer);
137:   PetscRandomDestroy(&(*bv)->rand);
138:   PetscHeaderDestroy(bv);
139:   return(0);
140: }

142: /*@
143:    BVCreate - Creates a basis vectors context.

145:    Collective

147:    Input Parameter:
148: .  comm - MPI communicator

150:    Output Parameter:
151: .  bv - location to put the basis vectors context

153:    Level: beginner

155: .seealso: BVSetUp(), BVDestroy(), BV
156: @*/
157: PetscErrorCode BVCreate(MPI_Comm comm,BV *newbv)
158: {
160:   BV             bv;

164:   *newbv = 0;
165:   BVInitializePackage();
166:   SlepcHeaderCreate(bv,BV_CLASSID,"BV","Basis Vectors","BV",comm,BVDestroy,BVView);

168:   bv->t            = NULL;
169:   bv->n            = -1;
170:   bv->N            = -1;
171:   bv->m            = 0;
172:   bv->l            = 0;
173:   bv->k            = 0;
174:   bv->nc           = 0;
175:   bv->orthog_type  = BV_ORTHOG_CGS;
176:   bv->orthog_ref   = BV_ORTHOG_REFINE_IFNEEDED;
177:   bv->orthog_eta   = 0.7071;
178:   bv->orthog_block = BV_ORTHOG_BLOCK_GS;
179:   bv->matrix       = NULL;
180:   bv->indef        = PETSC_FALSE;
181:   bv->vmm          = BV_MATMULT_MAT;

183:   bv->Bx           = NULL;
184:   bv->buffer       = NULL;
185:   bv->Abuffer      = NULL;
186:   bv->xid          = 0;
187:   bv->xstate       = 0;
188:   bv->cv[0]        = NULL;
189:   bv->cv[1]        = NULL;
190:   bv->ci[0]        = -1;
191:   bv->ci[1]        = -1;
192:   bv->st[0]        = -1;
193:   bv->st[1]        = -1;
194:   bv->id[0]        = 0;
195:   bv->id[1]        = 0;
196:   bv->h            = NULL;
197:   bv->c            = NULL;
198:   bv->omega        = NULL;
199:   bv->defersfo     = PETSC_FALSE;
200:   bv->cached       = NULL;
201:   bv->bvstate      = 0;
202:   bv->L            = NULL;
203:   bv->R            = NULL;
204:   bv->lstate       = 0;
205:   bv->rstate       = 0;
206:   bv->lsplit       = 0;
207:   bv->issplit      = 0;
208:   bv->splitparent  = NULL;
209:   bv->rand         = NULL;
210:   bv->rrandom      = PETSC_FALSE;
211:   bv->Acreate      = NULL;
212:   bv->Aget         = NULL;
213:   bv->cuda         = PETSC_FALSE;
214:   bv->work         = NULL;
215:   bv->lwork        = 0;
216:   bv->data         = NULL;

218:   *newbv = bv;
219:   return(0);
220: }

222: /*@
223:    BVCreateFromMat - Creates a basis vectors object from a dense Mat object.

225:    Collective on A

227:    Input Parameter:
228: .  A - a dense tall-skinny matrix

230:    Output Parameter:
231: .  bv - the new basis vectors context

233:    Notes:
234:    The matrix values are copied to the BV data storage, memory is not shared.

236:    The communicator of the BV object will be the same as A, and so will be
237:    the dimensions.

239:    Level: intermediate

241: .seealso: BVCreate(), BVDestroy(), BVCreateMat()
242: @*/
243: PetscErrorCode BVCreateFromMat(Mat A,BV *bv)
244: {
246:   PetscBool      match;
247:   PetscInt       n,N,k;

251:   PetscObjectTypeCompareAny((PetscObject)A,&match,MATSEQDENSE,MATMPIDENSE,"");
252:   if (!match) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Matrix must be of type dense");

254:   MatGetSize(A,&N,&k);
255:   MatGetLocalSize(A,&n,NULL);
256:   BVCreate(PetscObjectComm((PetscObject)A),bv);
257:   BVSetSizes(*bv,n,N,k);

259:   (*bv)->Acreate = A;
260:   PetscObjectReference((PetscObject)A);
261:   return(0);
262: }

264: /*@
265:    BVInsertVec - Insert a vector into the specified column.

267:    Collective on V

269:    Input Parameters:
270: +  V - basis vectors
271: .  j - the column of V to be overwritten
272: -  w - the vector to be copied

274:    Level: intermediate

276: .seealso: BVInsertVecs()
277: @*/
278: PetscErrorCode BVInsertVec(BV V,PetscInt j,Vec w)
279: {
281:   PetscInt       n,N;
282:   Vec            v;

289:   BVCheckSizes(V,1);

292:   VecGetSize(w,&N);
293:   VecGetLocalSize(w,&n);
294:   if (N!=V->N || n!=V->n) SETERRQ4(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Vec sizes (global %D, local %D) do not match BV sizes (global %D, local %D)",N,n,V->N,V->n);
295:   if (j<-V->nc || j>=V->m) SETERRQ3(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument j has wrong value %D, should be between %D and %D",j,-V->nc,V->m-1);

297:   BVGetColumn(V,j,&v);
298:   VecCopy(w,v);
299:   BVRestoreColumn(V,j,&v);
300:   PetscObjectStateIncrease((PetscObject)V);
301:   return(0);
302: }

304: /*@
305:    BVInsertVecs - Insert a set of vectors into the specified columns.

307:    Collective on V

309:    Input Parameters:
310: +  V - basis vectors
311: .  s - first column of V to be overwritten
312: .  W - set of vectors to be copied
313: -  orth - flag indicating if the vectors must be orthogonalized

315:    Input/Output Parameter:
316: .  m - number of input vectors, on output the number of linearly independent
317:        vectors

319:    Notes:
320:    Copies the contents of vectors W to V(:,s:s+n). If the orthogonalization
321:    flag is set, then the vectors are copied one by one and then orthogonalized
322:    against the previous ones. If any of them is linearly dependent then it
323:    is discarded and the value of m is decreased.

325:    Level: intermediate

327: .seealso: BVInsertVec(), BVOrthogonalizeColumn()
328: @*/
329: PetscErrorCode BVInsertVecs(BV V,PetscInt s,PetscInt *m,Vec *W,PetscBool orth)
330: {
332:   PetscInt       n,N,i,ndep;
333:   PetscBool      lindep;
334:   PetscReal      norm;
335:   Vec            v;

342:   if (!*m) return(0);
343:   if (*m<0) SETERRQ1(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Number of vectors (given %D) cannot be negative",*m);
348:   BVCheckSizes(V,1);

351:   VecGetSize(*W,&N);
352:   VecGetLocalSize(*W,&n);
353:   if (N!=V->N || n!=V->n) SETERRQ4(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Vec sizes (global %D, local %D) do not match BV sizes (global %D, local %D)",N,n,V->N,V->n);
354:   if (s<0 || s>=V->m) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument s has wrong value %D, should be between 0 and %D",s,V->m-1);
355:   if (s+(*m)>V->m) SETERRQ1(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Too many vectors provided, there is only room for %D",V->m);

357:   ndep = 0;
358:   for (i=0;i<*m;i++) {
359:     BVGetColumn(V,s+i-ndep,&v);
360:     VecCopy(W[i],v);
361:     BVRestoreColumn(V,s+i-ndep,&v);
362:     if (orth) {
363:       BVOrthogonalizeColumn(V,s+i-ndep,NULL,&norm,&lindep);
364:       if (norm==0.0 || lindep) {
365:         PetscInfo1(V,"Removing linearly dependent vector %D\n",i);
366:         ndep++;
367:       } else {
368:         BVScaleColumn(V,s+i-ndep,1.0/norm);
369:       }
370:     }
371:   }
372:   *m -= ndep;
373:   PetscObjectStateIncrease((PetscObject)V);
374:   return(0);
375: }

377: /*@
378:    BVInsertConstraints - Insert a set of vectors as constraints.

380:    Collective on V

382:    Input Parameters:
383: +  V - basis vectors
384: -  C - set of vectors to be inserted as constraints

386:    Input/Output Parameter:
387: .  nc - number of input vectors, on output the number of linearly independent
388:        vectors

390:    Notes:
391:    The constraints are relevant only during orthogonalization. Constraint
392:    vectors span a subspace that is deflated in every orthogonalization
393:    operation, so they are intended for removing those directions from the
394:    orthogonal basis computed in regular BV columns.

396:    Constraints are not stored in regular BV colums, but in a special part of
397:    the storage. They can be accessed with negative indices in BVGetColumn().

399:    This operation is DESTRUCTIVE, meaning that all data contained in the
400:    columns of V is lost. This is typically invoked just after creating the BV.
401:    Once a set of constraints has been set, it is not allowed to call this
402:    function again.

404:    The vectors are copied one by one and then orthogonalized against the
405:    previous ones. If any of them is linearly dependent then it is discarded
406:    and the value of nc is decreased. The behaviour is similar to BVInsertVecs().

408:    Level: advanced

410: .seealso: BVInsertVecs(), BVOrthogonalizeColumn(), BVGetColumn(), BVGetNumConstraints()
411: @*/
412: PetscErrorCode BVInsertConstraints(BV V,PetscInt *nc,Vec *C)
413: {
415:   PetscInt       msave;

421:   if (!*nc) return(0);
422:   if (*nc<0) SETERRQ1(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Number of constraints (given %D) cannot be negative",*nc);
426:   BVCheckSizes(V,1);
428:   if (V->issplit) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_WRONGSTATE,"Operation not permitted for a BV obtained from BVGetSplit");
429:   if (V->nc) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_WRONGSTATE,"Constraints already present in this BV object");
430:   if (V->ci[0]!=-1 || V->ci[1]!=-1) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Cannot call BVInsertConstraints after BVGetColumn");

432:   msave = V->m;
433:   BVResize(V,*nc+V->m,PETSC_FALSE);
434:   BVInsertVecs(V,0,nc,C,PETSC_TRUE);
435:   V->nc = *nc;
436:   V->m  = msave;
437:   V->ci[0] = -V->nc-1;
438:   V->ci[1] = -V->nc-1;
439:   PetscObjectStateIncrease((PetscObject)V);
440:   return(0);
441: }

443: /*@C
444:    BVSetOptionsPrefix - Sets the prefix used for searching for all
445:    BV options in the database.

447:    Logically Collective on bv

449:    Input Parameters:
450: +  bv     - the basis vectors context
451: -  prefix - the prefix string to prepend to all BV option requests

453:    Notes:
454:    A hyphen (-) must NOT be given at the beginning of the prefix name.
455:    The first character of all runtime options is AUTOMATICALLY the
456:    hyphen.

458:    Level: advanced

460: .seealso: BVAppendOptionsPrefix(), BVGetOptionsPrefix()
461: @*/
462: PetscErrorCode BVSetOptionsPrefix(BV bv,const char *prefix)
463: {

468:   PetscObjectSetOptionsPrefix((PetscObject)bv,prefix);
469:   return(0);
470: }

472: /*@C
473:    BVAppendOptionsPrefix - Appends to the prefix used for searching for all
474:    BV options in the database.

476:    Logically Collective on bv

478:    Input Parameters:
479: +  bv     - the basis vectors context
480: -  prefix - the prefix string to prepend to all BV option requests

482:    Notes:
483:    A hyphen (-) must NOT be given at the beginning of the prefix name.
484:    The first character of all runtime options is AUTOMATICALLY the
485:    hyphen.

487:    Level: advanced

489: .seealso: BVSetOptionsPrefix(), BVGetOptionsPrefix()
490: @*/
491: PetscErrorCode BVAppendOptionsPrefix(BV bv,const char *prefix)
492: {

497:   PetscObjectAppendOptionsPrefix((PetscObject)bv,prefix);
498:   return(0);
499: }

501: /*@C
502:    BVGetOptionsPrefix - Gets the prefix used for searching for all
503:    BV options in the database.

505:    Not Collective

507:    Input Parameters:
508: .  bv - the basis vectors context

510:    Output Parameters:
511: .  prefix - pointer to the prefix string used, is returned

513:    Note:
514:    On the Fortran side, the user should pass in a string 'prefix' of
515:    sufficient length to hold the prefix.

517:    Level: advanced

519: .seealso: BVSetOptionsPrefix(), BVAppendOptionsPrefix()
520: @*/
521: PetscErrorCode BVGetOptionsPrefix(BV bv,const char *prefix[])
522: {

528:   PetscObjectGetOptionsPrefix((PetscObject)bv,prefix);
529:   return(0);
530: }

532: static PetscErrorCode BVView_Default(BV bv,PetscViewer viewer)
533: {
534:   PetscErrorCode    ierr;
535:   PetscInt          j;
536:   Vec               v;
537:   PetscViewerFormat format;
538:   PetscBool         isascii,ismatlab=PETSC_FALSE;
539:   const char        *bvname,*name;

542:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
543:   if (isascii) {
544:     PetscViewerGetFormat(viewer,&format);
545:     if (format == PETSC_VIEWER_ASCII_MATLAB) ismatlab = PETSC_TRUE;
546:   }
547:   if (ismatlab) {
548:     PetscObjectGetName((PetscObject)bv,&bvname);
549:     PetscViewerASCIIPrintf(viewer,"%s=[];\n",bvname);
550:   }
551:   for (j=-bv->nc;j<bv->m;j++) {
552:     BVGetColumn(bv,j,&v);
553:     VecView(v,viewer);
554:     if (ismatlab) {
555:       PetscObjectGetName((PetscObject)v,&name);
556:       PetscViewerASCIIPrintf(viewer,"%s=[%s,%s];clear %s\n",bvname,bvname,name,name);
557:     }
558:     BVRestoreColumn(bv,j,&v);
559:   }
560:   return(0);
561: }

563: /*@C
564:    BVView - Prints the BV data structure.

566:    Collective on bv

568:    Input Parameters:
569: +  bv     - the BV context
570: -  viewer - optional visualization context

572:    Note:
573:    The available visualization contexts include
574: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
575: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
576:          output where only the first processor opens
577:          the file.  All other processors send their
578:          data to the first processor to print.

580:    The user can open an alternative visualization contexts with
581:    PetscViewerASCIIOpen() (output to a specified file).

583:    Level: beginner
584: @*/
585: PetscErrorCode BVView(BV bv,PetscViewer viewer)
586: {
587:   PetscErrorCode    ierr;
588:   PetscBool         isascii;
589:   PetscViewerFormat format;
590:   const char        *orthname[2] = {"classical","modified"};
591:   const char        *refname[3] = {"if needed","never","always"};

595:   if (!viewer) {
596:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)bv),&viewer);
597:   }

600:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
601:   if (isascii) {
602:     PetscObjectPrintClassNamePrefixType((PetscObject)bv,viewer);
603:     PetscViewerGetFormat(viewer,&format);
604:     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
605:       PetscViewerASCIIPrintf(viewer,"  %D columns of global length %D%s\n",bv->m,bv->N,bv->cuda?" (CUDA)":"");
606:       if (bv->nc>0) {
607:         PetscViewerASCIIPrintf(viewer,"  number of constraints: %D\n",bv->nc);
608:       }
609:       PetscViewerASCIIPrintf(viewer,"  vector orthogonalization method: %s Gram-Schmidt\n",orthname[bv->orthog_type]);
610:       switch (bv->orthog_ref) {
611:         case BV_ORTHOG_REFINE_IFNEEDED:
612:           PetscViewerASCIIPrintf(viewer,"  orthogonalization refinement: %s (eta: %g)\n",refname[bv->orthog_ref],(double)bv->orthog_eta);
613:           break;
614:         case BV_ORTHOG_REFINE_NEVER:
615:         case BV_ORTHOG_REFINE_ALWAYS:
616:           PetscViewerASCIIPrintf(viewer,"  orthogonalization refinement: %s\n",refname[bv->orthog_ref]);
617:           break;
618:       }
619:       PetscViewerASCIIPrintf(viewer,"  block orthogonalization method: %s\n",BVOrthogBlockTypes[bv->orthog_block]);
620:       if (bv->matrix) {
621:         if (bv->indef) {
622:           PetscViewerASCIIPrintf(viewer,"  indefinite inner product\n");
623:         } else {
624:           PetscViewerASCIIPrintf(viewer,"  non-standard inner product\n");
625:         }
626:         PetscViewerASCIIPrintf(viewer,"  inner product matrix:\n");
627:         PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
628:         PetscViewerASCIIPushTab(viewer);
629:         MatView(bv->matrix,viewer);
630:         PetscViewerASCIIPopTab(viewer);
631:         PetscViewerPopFormat(viewer);
632:       }
633:       switch (bv->vmm) {
634:         case BV_MATMULT_VECS:
635:           PetscViewerASCIIPrintf(viewer,"  doing matmult as matrix-vector products\n");
636:           break;
637:         case BV_MATMULT_MAT:
638:           PetscViewerASCIIPrintf(viewer,"  doing matmult as a single matrix-matrix product\n");
639:           break;
640:         case BV_MATMULT_MAT_SAVE:
641:           PetscViewerASCIIPrintf(viewer,"  mat_save is deprecated, use mat\n");
642:           break;
643:       }
644:       if (bv->rrandom) {
645:         PetscViewerASCIIPrintf(viewer,"  generating random vectors independent of the number of processes\n");
646:       }
647:       if (bv->ops->view) { (*bv->ops->view)(bv,viewer); }
648:     } else {
649:       if (bv->ops->view) { (*bv->ops->view)(bv,viewer); }
650:       else { BVView_Default(bv,viewer); }
651:     }
652:   } else {
653:     (*bv->ops->view)(bv,viewer);
654:   }
655:   return(0);
656: }

658: /*@C
659:    BVViewFromOptions - View from options

661:    Collective on BV

663:    Input Parameters:
664: +  bv   - the basis vectors context
665: .  obj  - optional object
666: -  name - command line option

668:    Level: intermediate

670: .seealso: BVView(), BVCreate()
671: @*/
672: PetscErrorCode BVViewFromOptions(BV bv,PetscObject obj,const char name[])
673: {

678:   PetscObjectViewFromOptions((PetscObject)bv,obj,name);
679:   return(0);
680: }

682: /*@C
683:    BVRegister - Adds a new storage format to the BV package.

685:    Not collective

687:    Input Parameters:
688: +  name     - name of a new user-defined BV
689: -  function - routine to create context

691:    Notes:
692:    BVRegister() may be called multiple times to add several user-defined
693:    basis vectors.

695:    Level: advanced

697: .seealso: BVRegisterAll()
698: @*/
699: PetscErrorCode BVRegister(const char *name,PetscErrorCode (*function)(BV))
700: {

704:   BVInitializePackage();
705:   PetscFunctionListAdd(&BVList,name,function);
706:   return(0);
707: }

709: PetscErrorCode BVAllocateWork_Private(BV bv,PetscInt s)
710: {

714:   if (s>bv->lwork) {
715:     PetscFree(bv->work);
716:     PetscMalloc1(s,&bv->work);
717:     PetscLogObjectMemory((PetscObject)bv,(s-bv->lwork)*sizeof(PetscScalar));
718:     bv->lwork = s;
719:   }
720:   return(0);
721: }

723: #if defined(PETSC_USE_DEBUG)
724: /*
725:    SlepcDebugBVView - partially view a BV object, to be used from within a debugger.

727:      ini, end: columns to be viewed
728:      s: name of Matlab variable
729:      filename: optionally write output to a file
730:  */
731: PETSC_UNUSED PetscErrorCode SlepcDebugBVView(BV bv,PetscInt ini,PetscInt end,const char *s,const char *filename)
732: {
734:   PetscInt       N,m;
735:   PetscScalar    *array;

738:   BVGetArray(bv,&array);
739:   BVGetSizes(bv,NULL,&N,&m);
740:   SlepcDebugViewMatrix(N,end-ini+1,array+ini*N,NULL,N,s,filename);
741:   BVRestoreArray(bv,&array);
742:   return(0);
743: }
744: #endif