Actual source code: bvfunc.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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;

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

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

107:    Collective on BV

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

112:    Level: beginner

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

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

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

147:    Collective on MPI_Comm

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

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

155:    Level: beginner

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

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

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

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

220:   *newbv = bv;
221:   return(0);
222: }

224: /*@
225:    BVCreateFromMat - Creates a basis vectors object from a dense Mat object.

227:    Collective on Mat

229:    Input Parameter:
230: .  A - a dense tall-skinny matrix

232:    Output Parameter:
233: .  bv - the new basis vectors context

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

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

241:    Level: intermediate

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

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

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

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

266: /*@
267:    BVInsertVec - Insert a vector into the specified column.

269:    Collective on BV

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

276:    Level: intermediate

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

291:   BVCheckSizes(V,1);

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

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

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

309:    Collective on BV

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

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

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

327:    Level: intermediate

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

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

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

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

379: /*@
380:    BVInsertConstraints - Insert a set of vectors as constraints.

382:    Collective on BV

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

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

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

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

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

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

410:    Level: advanced

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

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

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

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

449:    Logically Collective on BV

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

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

460:    Level: advanced

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

470:   PetscObjectSetOptionsPrefix((PetscObject)bv,prefix);
471:   return(0);
472: }

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

478:    Logically Collective on BV

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

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

489:    Level: advanced

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

499:   PetscObjectAppendOptionsPrefix((PetscObject)bv,prefix);
500:   return(0);
501: }

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

507:    Not Collective

509:    Input Parameters:
510: .  bv - the basis vectors context

512:    Output Parameters:
513: .  prefix - pointer to the prefix string used, is returned

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

519:    Level: advanced

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

530:   PetscObjectGetOptionsPrefix((PetscObject)bv,prefix);
531:   return(0);
532: }

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

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

565: /*@C
566:    BVView - Prints the BV data structure.

568:    Collective on BV

570:    Input Parameters:
571: +  bv     - the BV context
572: -  viewer - optional visualization context

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

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

585:    Level: beginner

587: .seealso: PetscViewerASCIIOpen()
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,"  inner product matrix:\n");
631:         PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
632:         PetscViewerASCIIPushTab(viewer);
633:         MatView(bv->matrix,viewer);
634:         PetscViewerASCIIPopTab(viewer);
635:         PetscViewerPopFormat(viewer);
636:       }
637:       switch (bv->vmm) {
638:         case BV_MATMULT_VECS:
639:           PetscViewerASCIIPrintf(viewer,"  doing matmult as matrix-vector products\n");
640:           break;
641:         case BV_MATMULT_MAT:
642:           PetscViewerASCIIPrintf(viewer,"  doing matmult as a single matrix-matrix product\n");
643:           break;
644:         case BV_MATMULT_MAT_SAVE:
645:           PetscViewerASCIIPrintf(viewer,"  mat_save is deprecated, use mat\n");
646:           break;
647:       }
648:       if (bv->rrandom) {
649:         PetscViewerASCIIPrintf(viewer,"  generating random vectors independent of the number of processes\n");
650:       }
651:       if (bv->ops->view) { (*bv->ops->view)(bv,viewer); }
652:     } else {
653:       if (bv->ops->view) { (*bv->ops->view)(bv,viewer); }
654:       else { BVView_Default(bv,viewer); }
655:     }
656:   } else {
657:     (*bv->ops->view)(bv,viewer);
658:   }
659:   return(0);
660: }

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

665:    Not collective

667:    Input Parameters:
668: +  name     - name of a new user-defined BV
669: -  function - routine to create context

671:    Notes:
672:    BVRegister() may be called multiple times to add several user-defined
673:    basis vectors.

675:    Level: advanced

677: .seealso: BVRegisterAll()
678: @*/
679: PetscErrorCode BVRegister(const char *name,PetscErrorCode (*function)(BV))
680: {

684:   BVInitializePackage();
685:   PetscFunctionListAdd(&BVList,name,function);
686:   return(0);
687: }

689: PetscErrorCode BVAllocateWork_Private(BV bv,PetscInt s)
690: {

694:   if (s>bv->lwork) {
695:     PetscFree(bv->work);
696:     PetscMalloc1(s,&bv->work);
697:     PetscLogObjectMemory((PetscObject)bv,(s-bv->lwork)*sizeof(PetscScalar));
698:     bv->lwork = s;
699:   }
700:   return(0);
701: }

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

707:      ini, end: columns to be viewed
708:      s: name of Matlab variable
709:      filename: optionally write output to a file
710:  */
711: PETSC_UNUSED PetscErrorCode SlepcDebugBVView(BV bv,PetscInt ini,PetscInt end,const char *s,const char *filename)
712: {
714:   PetscInt       N,m;
715:   PetscScalar    *array;

718:   BVGetArray(bv,&array);
719:   BVGetSizes(bv,NULL,&N,&m);
720:   SlepcDebugViewMatrix(N,end-ini+1,array+ini*N,NULL,N,s,filename);
721:   BVRestoreArray(bv,&array);
722:   return(0);
723: }
724: #endif