Actual source code: bvfunc.c

slepc-3.15.0 2021-03-31
Report Typos and Errors
  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 (basis vectors) interface routines, callable by users
 12: */

 14: #include <slepc/private/bvimpl.h>

 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_Normalize = 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,MPIU_LAPY2;

 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);CHKERRMPI(ierr);
 41:   MPI_Op_free(&MPIU_LAPY2);CHKERRMPI(ierr);
 42:   BVPackageInitialized = PETSC_FALSE;
 43:   BVRegisterAllCalled  = PETSC_FALSE;
 44:   return(0);
 45: }

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

 52:    Level: developer

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

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

105: /*@C
106:    BVDestroy - Destroys BV context that was created with BVCreate().

108:    Collective on bv

110:    Input Parameter:
111: .  bv - the basis vectors context

113:    Level: beginner

115: .seealso: BVCreate()
116: @*/
117: PetscErrorCode BVDestroy(BV *bv)
118: {

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

145: /*@
146:    BVCreate - Creates a basis vectors context.

148:    Collective

150:    Input Parameter:
151: .  comm - MPI communicator

153:    Output Parameter:
154: .  bv - location to put the basis vectors context

156:    Level: beginner

158: .seealso: BVSetUp(), BVDestroy(), BV
159: @*/
160: PetscErrorCode BVCreate(MPI_Comm comm,BV *newbv)
161: {
163:   BV             bv;

167:   *newbv = 0;
168:   BVInitializePackage();
169:   SlepcHeaderCreate(bv,BV_CLASSID,"BV","Basis Vectors","BV",comm,BVDestroy,BVView);

171:   bv->t            = NULL;
172:   bv->n            = -1;
173:   bv->N            = -1;
174:   bv->m            = 0;
175:   bv->l            = 0;
176:   bv->k            = 0;
177:   bv->nc           = 0;
178:   bv->orthog_type  = BV_ORTHOG_CGS;
179:   bv->orthog_ref   = BV_ORTHOG_REFINE_IFNEEDED;
180:   bv->orthog_eta   = 0.7071;
181:   bv->orthog_block = BV_ORTHOG_BLOCK_GS;
182:   bv->matrix       = NULL;
183:   bv->indef        = PETSC_FALSE;
184:   bv->vmm          = BV_MATMULT_MAT;
185:   bv->rrandom      = PETSC_FALSE;
186:   bv->deftol       = 10*PETSC_MACHINE_EPSILON;

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

224:   *newbv = bv;
225:   return(0);
226: }

228: /*@
229:    BVCreateFromMat - Creates a basis vectors object from a dense Mat object.

231:    Collective on A

233:    Input Parameter:
234: .  A - a dense tall-skinny matrix

236:    Output Parameter:
237: .  bv - the new basis vectors context

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

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

245:    Level: intermediate

247: .seealso: BVCreate(), BVDestroy(), BVCreateMat()
248: @*/
249: PetscErrorCode BVCreateFromMat(Mat A,BV *bv)
250: {
252:   PetscInt       n,N,k;


258:   MatGetSize(A,&N,&k);
259:   MatGetLocalSize(A,&n,NULL);
260:   BVCreate(PetscObjectComm((PetscObject)A),bv);
261:   BVSetSizes(*bv,n,N,k);

263:   (*bv)->Acreate = A;
264:   PetscObjectReference((PetscObject)A);
265:   return(0);
266: }

268: /*@
269:    BVInsertVec - Insert a vector into the specified column.

271:    Collective on V

273:    Input Parameters:
274: +  V - basis vectors
275: .  j - the column of V to be overwritten
276: -  w - the vector to be copied

278:    Level: intermediate

280: .seealso: BVInsertVecs()
281: @*/
282: PetscErrorCode BVInsertVec(BV V,PetscInt j,Vec w)
283: {
285:   PetscInt       n,N;
286:   Vec            v;

293:   BVCheckSizes(V,1);

296:   VecGetSize(w,&N);
297:   VecGetLocalSize(w,&n);
298:   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);
299:   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);

301:   BVGetColumn(V,j,&v);
302:   VecCopy(w,v);
303:   BVRestoreColumn(V,j,&v);
304:   PetscObjectStateIncrease((PetscObject)V);
305:   return(0);
306: }

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

311:    Collective on V

313:    Input Parameters:
314: +  V - basis vectors
315: .  s - first column of V to be overwritten
316: .  W - set of vectors to be copied
317: -  orth - flag indicating if the vectors must be orthogonalized

319:    Input/Output Parameter:
320: .  m - number of input vectors, on output the number of linearly independent
321:        vectors

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

329:    Level: intermediate

331: .seealso: BVInsertVec(), BVOrthogonalizeColumn()
332: @*/
333: PetscErrorCode BVInsertVecs(BV V,PetscInt s,PetscInt *m,Vec *W,PetscBool orth)
334: {
336:   PetscInt       n,N,i,ndep;
337:   PetscBool      lindep;
338:   PetscReal      norm;
339:   Vec            v;

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

355:   VecGetSize(*W,&N);
356:   VecGetLocalSize(*W,&n);
357:   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);
358:   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);
359:   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);

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

381: /*@
382:    BVInsertConstraints - Insert a set of vectors as constraints.

384:    Collective on V

386:    Input Parameters:
387: +  V - basis vectors
388: -  C - set of vectors to be inserted as constraints

390:    Input/Output Parameter:
391: .  nc - number of input vectors, on output the number of linearly independent
392:        vectors

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

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

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

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

412:    Level: advanced

414: .seealso: BVInsertVecs(), BVOrthogonalizeColumn(), BVGetColumn(), BVGetNumConstraints()
415: @*/
416: PetscErrorCode BVInsertConstraints(BV V,PetscInt *nc,Vec *C)
417: {
419:   PetscInt       msave;

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

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

447: /*@C
448:    BVSetOptionsPrefix - Sets the prefix used for searching for all
449:    BV options in the database.

451:    Logically Collective on bv

453:    Input Parameters:
454: +  bv     - the basis vectors context
455: -  prefix - the prefix string to prepend to all BV option requests

457:    Notes:
458:    A hyphen (-) must NOT be given at the beginning of the prefix name.
459:    The first character of all runtime options is AUTOMATICALLY the
460:    hyphen.

462:    Level: advanced

464: .seealso: BVAppendOptionsPrefix(), BVGetOptionsPrefix()
465: @*/
466: PetscErrorCode BVSetOptionsPrefix(BV bv,const char *prefix)
467: {

472:   PetscObjectSetOptionsPrefix((PetscObject)bv,prefix);
473:   return(0);
474: }

476: /*@C
477:    BVAppendOptionsPrefix - Appends to the prefix used for searching for all
478:    BV options in the database.

480:    Logically Collective on bv

482:    Input Parameters:
483: +  bv     - the basis vectors context
484: -  prefix - the prefix string to prepend to all BV option requests

486:    Notes:
487:    A hyphen (-) must NOT be given at the beginning of the prefix name.
488:    The first character of all runtime options is AUTOMATICALLY the
489:    hyphen.

491:    Level: advanced

493: .seealso: BVSetOptionsPrefix(), BVGetOptionsPrefix()
494: @*/
495: PetscErrorCode BVAppendOptionsPrefix(BV bv,const char *prefix)
496: {

501:   PetscObjectAppendOptionsPrefix((PetscObject)bv,prefix);
502:   return(0);
503: }

505: /*@C
506:    BVGetOptionsPrefix - Gets the prefix used for searching for all
507:    BV options in the database.

509:    Not Collective

511:    Input Parameters:
512: .  bv - the basis vectors context

514:    Output Parameters:
515: .  prefix - pointer to the prefix string used, is returned

517:    Note:
518:    On the Fortran side, the user should pass in a string 'prefix' of
519:    sufficient length to hold the prefix.

521:    Level: advanced

523: .seealso: BVSetOptionsPrefix(), BVAppendOptionsPrefix()
524: @*/
525: PetscErrorCode BVGetOptionsPrefix(BV bv,const char *prefix[])
526: {

532:   PetscObjectGetOptionsPrefix((PetscObject)bv,prefix);
533:   return(0);
534: }

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

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

567: /*@C
568:    BVView - Prints the BV data structure.

570:    Collective on bv

572:    Input Parameters:
573: +  bv     - the BV context
574: -  viewer - optional visualization context

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

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

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

599:   if (!viewer) {
600:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)bv),&viewer);
601:   }

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

663: /*@C
664:    BVViewFromOptions - View from options

666:    Collective on BV

668:    Input Parameters:
669: +  bv   - the basis vectors context
670: .  obj  - optional object
671: -  name - command line option

673:    Level: intermediate

675: .seealso: BVView(), BVCreate()
676: @*/
677: PetscErrorCode BVViewFromOptions(BV bv,PetscObject obj,const char name[])
678: {

683:   PetscObjectViewFromOptions((PetscObject)bv,obj,name);
684:   return(0);
685: }

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

690:    Not collective

692:    Input Parameters:
693: +  name     - name of a new user-defined BV
694: -  function - routine to create context

696:    Notes:
697:    BVRegister() may be called multiple times to add several user-defined
698:    basis vectors.

700:    Level: advanced

702: .seealso: BVRegisterAll()
703: @*/
704: PetscErrorCode BVRegister(const char *name,PetscErrorCode (*function)(BV))
705: {

709:   BVInitializePackage();
710:   PetscFunctionListAdd(&BVList,name,function);
711:   return(0);
712: }

714: PetscErrorCode BVAllocateWork_Private(BV bv,PetscInt s)
715: {

719:   if (s>bv->lwork) {
720:     PetscFree(bv->work);
721:     PetscMalloc1(s,&bv->work);
722:     PetscLogObjectMemory((PetscObject)bv,(s-bv->lwork)*sizeof(PetscScalar));
723:     bv->lwork = s;
724:   }
725:   return(0);
726: }

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

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

743:   BVGetArray(bv,&array);
744:   BVGetSizes(bv,NULL,&N,&m);
745:   SlepcDebugViewMatrix(N,end-ini+1,array+ini*N,NULL,N,s,filename);
746:   BVRestoreArray(bv,&array);
747:   return(0);
748: }
749: #endif