Actual source code: bvbasic.c

slepc-3.14.2 2021-02-01
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:    Basic BV routines
 12: */

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

 16: PetscBool         BVRegisterAllCalled = PETSC_FALSE;
 17: PetscFunctionList BVList = 0;

 19: /*@C
 20:    BVSetType - Selects the type for the BV object.

 22:    Logically Collective on bv

 24:    Input Parameter:
 25: +  bv   - the basis vectors context
 26: -  type - a known type

 28:    Options Database Key:
 29: .  -bv_type <type> - Sets BV type

 31:    Level: intermediate

 33: .seealso: BVGetType()
 34: @*/
 35: PetscErrorCode BVSetType(BV bv,BVType type)
 36: {
 37:   PetscErrorCode ierr,(*r)(BV);
 38:   PetscBool      match;


 44:   PetscObjectTypeCompare((PetscObject)bv,type,&match);
 45:   if (match) return(0);
 46:   PetscStrcmp(type,BVTENSOR,&match);
 47:   if (match) SETERRQ(PetscObjectComm((PetscObject)bv),1,"Use BVCreateTensor() to create a BV of type tensor");

 49:    PetscFunctionListFind(BVList,type,&r);
 50:   if (!r) SETERRQ1(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested BV type %s",type);

 52:   if (bv->ops->destroy) { (*bv->ops->destroy)(bv); }
 53:   PetscMemzero(bv->ops,sizeof(struct _BVOps));

 55:   PetscObjectChangeTypeName((PetscObject)bv,type);
 56:   if (bv->n < 0 && bv->N < 0) {
 57:     bv->ops->create = r;
 58:   } else {
 59:     PetscLogEventBegin(BV_Create,bv,0,0,0);
 60:     (*r)(bv);
 61:     PetscLogEventEnd(BV_Create,bv,0,0,0);
 62:   }
 63:   return(0);
 64: }

 66: /*@C
 67:    BVGetType - Gets the BV type name (as a string) from the BV context.

 69:    Not Collective

 71:    Input Parameter:
 72: .  bv - the basis vectors context

 74:    Output Parameter:
 75: .  name - name of the type of basis vectors

 77:    Level: intermediate

 79: .seealso: BVSetType()
 80: @*/
 81: PetscErrorCode BVGetType(BV bv,BVType *type)
 82: {
 86:   *type = ((PetscObject)bv)->type_name;
 87:   return(0);
 88: }

 90: /*@
 91:    BVSetSizes - Sets the local and global sizes, and the number of columns.

 93:    Collective on bv

 95:    Input Parameters:
 96: +  bv - the basis vectors
 97: .  n  - the local size (or PETSC_DECIDE to have it set)
 98: .  N  - the global size (or PETSC_DECIDE)
 99: -  m  - the number of columns

101:    Notes:
102:    n and N cannot be both PETSC_DECIDE.
103:    If one processor calls this with N of PETSC_DECIDE then all processors must,
104:    otherwise the program will hang.

106:    Level: beginner

108: .seealso: BVSetSizesFromVec(), BVGetSizes(), BVResize()
109: @*/
110: PetscErrorCode BVSetSizes(BV bv,PetscInt n,PetscInt N,PetscInt m)
111: {
113:   PetscInt       ma;

119:   if (N >= 0 && n > N) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Local size %D cannot be larger than global size %D",n,N);
120:   if (m <= 0) SETERRQ1(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Number of columns %D must be positive",m);
121:   if ((bv->n >= 0 || bv->N >= 0) && (bv->n != n || bv->N != N)) SETERRQ4(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Cannot change/reset vector sizes to %D local %D global after previously setting them to %D local %D global",n,N,bv->n,bv->N);
122:   if (bv->m > 0 && bv->m != m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Cannot change the number of columns to %D after previously setting it to %D; use BVResize()",m,bv->m);
123:   bv->n = n;
124:   bv->N = N;
125:   bv->m = m;
126:   bv->k = m;
127:   if (!bv->t) {  /* create template vector and get actual dimensions */
128:     VecCreate(PetscObjectComm((PetscObject)bv),&bv->t);
129:     VecSetSizes(bv->t,bv->n,bv->N);
130:     VecSetFromOptions(bv->t);
131:     VecGetSize(bv->t,&bv->N);
132:     VecGetLocalSize(bv->t,&bv->n);
133:     if (bv->matrix) {  /* check compatible dimensions of user-provided matrix */
134:       MatGetLocalSize(bv->matrix,&ma,NULL);
135:       if (bv->n!=ma) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Local dimension %D does not match that of matrix given at BVSetMatrix %D",bv->n,ma);
136:     }
137:   }
138:   if (bv->ops->create) {
139:     PetscLogEventBegin(BV_Create,bv,0,0,0);
140:     (*bv->ops->create)(bv);
141:     PetscLogEventEnd(BV_Create,bv,0,0,0);
142:     bv->ops->create = 0;
143:     bv->defersfo = PETSC_FALSE;
144:   }
145:   return(0);
146: }

148: /*@
149:    BVSetSizesFromVec - Sets the local and global sizes, and the number of columns.
150:    Local and global sizes are specified indirectly by passing a template vector.

152:    Collective on bv

154:    Input Parameters:
155: +  bv - the basis vectors
156: .  t  - the template vector
157: -  m  - the number of columns

159:    Level: beginner

161: .seealso: BVSetSizes(), BVGetSizes(), BVResize()
162: @*/
163: PetscErrorCode BVSetSizesFromVec(BV bv,Vec t,PetscInt m)
164: {
166:   PetscInt       ma;

173:   if (m <= 0) SETERRQ1(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Number of columns %D must be positive",m);
174:   if (bv->t) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Template vector was already set by a previous call to BVSetSizes/FromVec");
175:   VecGetSize(t,&bv->N);
176:   VecGetLocalSize(t,&bv->n);
177:   if (bv->matrix) {  /* check compatible dimensions of user-provided matrix */
178:     MatGetLocalSize(bv->matrix,&ma,NULL);
179:     if (bv->n!=ma) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Local dimension %D does not match that of matrix given at BVSetMatrix %D",bv->n,ma);
180:   }
181:   bv->m = m;
182:   bv->k = m;
183:   bv->t = t;
184:   PetscObjectReference((PetscObject)t);
185:   if (bv->ops->create) {
186:     (*bv->ops->create)(bv);
187:     bv->ops->create = 0;
188:     bv->defersfo = PETSC_FALSE;
189:   }
190:   return(0);
191: }

193: /*@
194:    BVGetSizes - Returns the local and global sizes, and the number of columns.

196:    Not Collective

198:    Input Parameter:
199: .  bv - the basis vectors

201:    Output Parameters:
202: +  n  - the local size
203: .  N  - the global size
204: -  m  - the number of columns

206:    Note:
207:    Normal usage requires that bv has already been given its sizes, otherwise
208:    the call fails. However, this function can also be used to determine if
209:    a BV object has been initialized completely (sizes and type). For this,
210:    call with n=NULL and N=NULL, then a return value of m=0 indicates that
211:    the BV object is not ready for use yet.

213:    Level: beginner

215: .seealso: BVSetSizes(), BVSetSizesFromVec()
216: @*/
217: PetscErrorCode BVGetSizes(BV bv,PetscInt *n,PetscInt *N,PetscInt *m)
218: {
220:   if (!bv) {
221:     if (m && !n && !N) *m = 0;
222:     return(0);
223:   }
225:   if (n || N) BVCheckSizes(bv,1);
226:   if (n) *n = bv->n;
227:   if (N) *N = bv->N;
228:   if (m) *m = bv->m;
229:   if (m && !n && !N && !((PetscObject)bv)->type_name) *m = 0;
230:   return(0);
231: }

233: /*@
234:    BVSetNumConstraints - Set the number of constraints.

236:    Logically Collective on V

238:    Input Parameters:
239: +  V  - basis vectors
240: -  nc - number of constraints

242:    Notes:
243:    This function sets the number of constraints to nc and marks all remaining
244:    columns as regular. Normal user would call BVInsertConstraints() instead.

246:    If nc is smaller than the previously set value, then some of the constraints
247:    are discarded. In particular, using nc=0 removes all constraints preserving
248:    the content of regular columns.

250:    Level: developer

252: .seealso: BVInsertConstraints()
253: @*/
254: PetscErrorCode BVSetNumConstraints(BV V,PetscInt nc)
255: {
257:   PetscInt       total,diff,i;
258:   Vec            x,y;

263:   if (nc<0) SETERRQ1(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Number of constraints (given %D) cannot be negative",nc);
265:   BVCheckSizes(V,1);
266:   if (V->ci[0]!=-V->nc-1 || V->ci[1]!=-V->nc-1) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Cannot call BVSetNumConstraints after BVGetColumn");

268:   diff = nc-V->nc;
269:   if (!diff) return(0);
270:   total = V->nc+V->m;
271:   if (total-nc<=0) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Not enough columns for the given nc value");
272:   if (diff<0) {  /* lessen constraints, shift contents of BV */
273:     for (i=0;i<V->m;i++) {
274:       BVGetColumn(V,i,&x);
275:       BVGetColumn(V,i+diff,&y);
276:       VecCopy(x,y);
277:       BVRestoreColumn(V,i,&x);
278:       BVRestoreColumn(V,i+diff,&y);
279:     }
280:   }
281:   V->nc = nc;
282:   V->ci[0] = -V->nc-1;
283:   V->ci[1] = -V->nc-1;
284:   V->m = total-nc;
285:   V->l = PetscMin(V->l,V->m);
286:   V->k = PetscMin(V->k,V->m);
287:   PetscObjectStateIncrease((PetscObject)V);
288:   return(0);
289: }

291: /*@
292:    BVGetNumConstraints - Returns the number of constraints.

294:    Not Collective

296:    Input Parameter:
297: .  bv - the basis vectors

299:    Output Parameters:
300: .  nc - the number of constraints

302:    Level: advanced

304: .seealso: BVGetSizes(), BVInsertConstraints()
305: @*/
306: PetscErrorCode BVGetNumConstraints(BV bv,PetscInt *nc)
307: {
311:   *nc = bv->nc;
312:   return(0);
313: }

315: /*@
316:    BVResize - Change the number of columns.

318:    Collective on bv

320:    Input Parameters:
321: +  bv   - the basis vectors
322: .  m    - the new number of columns
323: -  copy - a flag indicating whether current values should be kept

325:    Note:
326:    Internal storage is reallocated. If the copy flag is set to true, then
327:    the contents are copied to the leading part of the new space.

329:    Level: advanced

331: .seealso: BVSetSizes(), BVSetSizesFromVec()
332: @*/
333: PetscErrorCode BVResize(BV bv,PetscInt m,PetscBool copy)
334: {
335:   PetscErrorCode    ierr;
336:   PetscScalar       *array;
337:   const PetscScalar *omega;
338:   Vec               v;

345:   if (m <= 0) SETERRQ1(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Number of columns %D must be positive",m);
346:   if (bv->nc && !bv->issplit) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"Cannot resize a BV with constraints");
347:   if (bv->m == m) return(0);
348:   BVCheckOp(bv,1,resize);

350:   PetscLogEventBegin(BV_Create,bv,0,0,0);
351:   (*bv->ops->resize)(bv,m,copy);
352:   VecDestroy(&bv->buffer);
353:   BVDestroy(&bv->cached);
354:   PetscFree2(bv->h,bv->c);
355:   if (bv->omega) {
356:     if (bv->cuda) {
357: #if defined(PETSC_HAVE_CUDA)
358:       VecCreateSeqCUDA(PETSC_COMM_SELF,m,&v);
359: #else
360:       SETERRQ(PetscObjectComm((PetscObject)bv),1,"Something wrong happened");
361: #endif
362:     } else {
363:       VecCreateSeq(PETSC_COMM_SELF,m,&v);
364:     }
365:     PetscLogObjectParent((PetscObject)bv,(PetscObject)v);
366:     if (copy) {
367:       VecGetArray(v,&array);
368:       VecGetArrayRead(bv->omega,&omega);
369:       PetscArraycpy(array,omega,PetscMin(m,bv->m));
370:       VecRestoreArrayRead(bv->omega,&omega);
371:       VecRestoreArray(v,&array);
372:     } else {
373:       VecSet(v,1.0);
374:     }
375:     VecDestroy(&bv->omega);
376:     bv->omega = v;
377:   }
378:   bv->m = m;
379:   bv->k = PetscMin(bv->k,m);
380:   bv->l = PetscMin(bv->l,m);
381:   PetscLogEventEnd(BV_Create,bv,0,0,0);
382:   PetscObjectStateIncrease((PetscObject)bv);
383:   return(0);
384: }

386: /*@
387:    BVSetActiveColumns - Specify the columns that will be involved in operations.

389:    Logically Collective on bv

391:    Input Parameters:
392: +  bv - the basis vectors context
393: .  l  - number of leading columns
394: -  k  - number of active columns

396:    Notes:
397:    In operations such as BVMult() or BVDot(), only the first k columns are
398:    considered. This is useful when the BV is filled from left to right, so
399:    the last m-k columns do not have relevant information.

401:    Also in operations such as BVMult() or BVDot(), the first l columns are
402:    normally not included in the computation. See the manpage of each
403:    operation.

405:    In orthogonalization operations, the first l columns are treated
406:    differently: they participate in the orthogonalization but the computed
407:    coefficients are not stored.

409:    Level: intermediate

411: .seealso: BVGetActiveColumns(), BVSetSizes()
412: @*/
413: PetscErrorCode BVSetActiveColumns(BV bv,PetscInt l,PetscInt k)
414: {
419:   BVCheckSizes(bv,1);
420:   if (k==PETSC_DECIDE || k==PETSC_DEFAULT) {
421:     bv->k = bv->m;
422:   } else {
423:     if (k<0 || k>bv->m) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of k. Must be between 0 and m");
424:     bv->k = k;
425:   }
426:   if (l==PETSC_DECIDE || l==PETSC_DEFAULT) {
427:     bv->l = 0;
428:   } else {
429:     if (l<0 || l>bv->k) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of l. Must be between 0 and k");
430:     bv->l = l;
431:   }
432:   return(0);
433: }

435: /*@
436:    BVGetActiveColumns - Returns the current active dimensions.

438:    Not Collective

440:    Input Parameter:
441: .  bv - the basis vectors context

443:    Output Parameter:
444: +  l  - number of leading columns
445: -  k  - number of active columns

447:    Level: intermediate

449: .seealso: BVSetActiveColumns()
450: @*/
451: PetscErrorCode BVGetActiveColumns(BV bv,PetscInt *l,PetscInt *k)
452: {
455:   if (l) *l = bv->l;
456:   if (k) *k = bv->k;
457:   return(0);
458: }

460: /*@
461:    BVSetMatrix - Specifies the inner product to be used in orthogonalization.

463:    Collective on bv

465:    Input Parameters:
466: +  bv    - the basis vectors context
467: .  B     - a symmetric matrix (may be NULL)
468: -  indef - a flag indicating if the matrix is indefinite

470:    Notes:
471:    This is used to specify a non-standard inner product, whose matrix
472:    representation is given by B. Then, all inner products required during
473:    orthogonalization are computed as (x,y)_B=y^H*B*x rather than the
474:    standard form (x,y)=y^H*x.

476:    Matrix B must be real symmetric (or complex Hermitian). A genuine inner
477:    product requires that B is also positive (semi-)definite. However, we
478:    also allow for an indefinite B (setting indef=PETSC_TRUE), in which
479:    case the orthogonalization uses an indefinite inner product.

481:    This affects operations BVDot(), BVNorm(), BVOrthogonalize(), and variants.

483:    Setting B=NULL has the same effect as if the identity matrix was passed.

485:    Level: advanced

487: .seealso: BVGetMatrix(), BVDot(), BVNorm(), BVOrthogonalize()
488: @*/
489: PetscErrorCode BVSetMatrix(BV bv,Mat B,PetscBool indef)
490: {
492:   PetscInt       m,n;

497:   if (B!=bv->matrix || (B && ((PetscObject)B)->id!=((PetscObject)bv->matrix)->id) || indef!=bv->indef) {
498:     if (B) {
500:       MatGetLocalSize(B,&m,&n);
501:       if (m!=n) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Matrix must be square");
502:       if (bv->m && bv->n!=n) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension BV %D, Mat %D",bv->n,n);
503:     }
504:     if (B) { PetscObjectReference((PetscObject)B); }
505:     MatDestroy(&bv->matrix);
506:     bv->matrix = B;
507:     bv->indef  = indef;
508:     PetscObjectStateIncrease((PetscObject)bv);
509:     if (bv->Bx) { PetscObjectStateIncrease((PetscObject)bv->Bx); }
510:     if (bv->cached) { PetscObjectStateIncrease((PetscObject)bv->cached); }
511:   }
512:   return(0);
513: }

515: /*@
516:    BVGetMatrix - Retrieves the matrix representation of the inner product.

518:    Not collective, though a parallel Mat may be returned

520:    Input Parameter:
521: .  bv    - the basis vectors context

523:    Output Parameter:
524: +  B     - the matrix of the inner product (may be NULL)
525: -  indef - the flag indicating if the matrix is indefinite

527:    Level: advanced

529: .seealso: BVSetMatrix()
530: @*/
531: PetscErrorCode BVGetMatrix(BV bv,Mat *B,PetscBool *indef)
532: {
535:   if (B)     *B     = bv->matrix;
536:   if (indef) *indef = bv->indef;
537:   return(0);
538: }

540: /*@
541:    BVApplyMatrix - Multiplies a vector by the matrix representation of the
542:    inner product.

544:    Neighbor-wise Collective on bv

546:    Input Parameter:
547: +  bv - the basis vectors context
548: -  x  - the vector

550:    Output Parameter:
551: .  y  - the result

553:    Note:
554:    If no matrix was specified this function copies the vector.

556:    Level: advanced

558: .seealso: BVSetMatrix(), BVApplyMatrixBV()
559: @*/
560: PetscErrorCode BVApplyMatrix(BV bv,Vec x,Vec y)
561: {

568:   if (bv->matrix) {
569:     BV_IPMatMult(bv,x);
570:     VecCopy(bv->Bx,y);
571:   } else {
572:     VecCopy(x,y);
573:   }
574:   return(0);
575: }

577: /*@
578:    BVApplyMatrixBV - Multiplies the BV vectors by the matrix representation
579:    of the inner product.

581:    Neighbor-wise Collective on X

583:    Input Parameter:
584: .  X - the basis vectors context

586:    Output Parameter:
587: .  Y - the basis vectors to store the result (optional)

589:    Note:
590:    This function computes Y = B*X, where B is the matrix given with
591:    BVSetMatrix(). This operation is computed as in BVMatMult().
592:    If no matrix was specified, then it just copies Y = X.

594:    If no Y is given, the result is stored internally in the cached BV.

596:    Level: developer

598: .seealso: BVSetMatrix(), BVApplyMatrix(), BVMatMult(), BVGetCachedBV()
599: @*/
600: PetscErrorCode BVApplyMatrixBV(BV X,BV Y)
601: {

606:   if (Y) {
608:     if (X->matrix) {
609:       BVMatMult(X,X->matrix,Y);
610:     } else {
611:       BVCopy(X,Y);
612:     }
613:   } else {
614:     BV_IPMatMultBV(X);
615:   }
616:   return(0);
617: }

619: /*@
620:    BVSetSignature - Sets the signature matrix to be used in orthogonalization.

622:    Logically Collective on bv

624:    Input Parameter:
625: +  bv    - the basis vectors context
626: -  omega - a vector representing the diagonal of the signature matrix

628:    Note:
629:    The signature matrix Omega = V'*B*V is relevant only for an indefinite B.

631:    Level: developer

633: .seealso: BVSetMatrix(), BVGetSignature()
634: @*/
635: PetscErrorCode BVSetSignature(BV bv,Vec omega)
636: {
637:   PetscErrorCode    ierr;
638:   PetscInt          i,n;
639:   const PetscScalar *pomega;
640:   PetscScalar       *intern;

644:   BVCheckSizes(bv,1);

648:   VecGetSize(omega,&n);
649:   if (n!=bv->k) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Vec argument has %D elements, should be %D",n,bv->k);
650:   BV_AllocateSignature(bv);
651:   if (bv->indef) {
652:     VecGetArrayRead(omega,&pomega);
653:     VecGetArray(bv->omega,&intern);
654:     for (i=0;i<n;i++) intern[bv->nc+i] = pomega[i];
655:     VecRestoreArray(bv->omega,&intern);
656:     VecRestoreArrayRead(omega,&pomega);
657:   } else {
658:     PetscInfo(bv,"Ignoring signature because BV is not indefinite\n");
659:   }
660:   PetscObjectStateIncrease((PetscObject)bv);
661:   return(0);
662: }

664: /*@
665:    BVGetSignature - Retrieves the signature matrix from last orthogonalization.

667:    Not collective

669:    Input Parameter:
670: .  bv    - the basis vectors context

672:    Output Parameter:
673: .  omega - a vector representing the diagonal of the signature matrix

675:    Note:
676:    The signature matrix Omega = V'*B*V is relevant only for an indefinite B.

678:    Level: developer

680: .seealso: BVSetMatrix(), BVSetSignature()
681: @*/
682: PetscErrorCode BVGetSignature(BV bv,Vec omega)
683: {
684:   PetscErrorCode    ierr;
685:   PetscInt          i,n;
686:   PetscScalar       *pomega;
687:   const PetscScalar *intern;

691:   BVCheckSizes(bv,1);

695:   VecGetSize(omega,&n);
696:   if (n!=bv->k) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Vec argument has %D elements, should be %D",n,bv->k);
697:   if (bv->indef && bv->omega) {
698:     VecGetArray(omega,&pomega);
699:     VecGetArrayRead(bv->omega,&intern);
700:     for (i=0;i<n;i++) pomega[i] = intern[bv->nc+i];
701:     VecRestoreArrayRead(bv->omega,&intern);
702:     VecRestoreArray(omega,&pomega);
703:   } else {
704:     VecSet(omega,1.0);
705:   }
706:   return(0);
707: }

709: /*@
710:    BVSetBufferVec - Attach a vector object to be used as buffer space for
711:    several operations.

713:    Collective on bv

715:    Input Parameters:
716: +  bv     - the basis vectors context)
717: -  buffer - the vector

719:    Notes:
720:    Use BVGetBufferVec() to retrieve the vector (for example, to free it
721:    at the end of the computations).

723:    The vector must be sequential of length (nc+m)*m, where m is the number
724:    of columns of bv and nc is the number of constraints.

726:    Level: developer

728: .seealso: BVGetBufferVec(), BVSetSizes(), BVGetNumConstraints()
729: @*/
730: PetscErrorCode BVSetBufferVec(BV bv,Vec buffer)
731: {
733:   PetscInt       ld,n;
734:   PetscMPIInt    size;

739:   BVCheckSizes(bv,1);
740:   VecGetSize(buffer,&n);
741:   ld = bv->m+bv->nc;
742:   if (n != ld*bv->m) SETERRQ1(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Buffer size must be %d",ld*bv->m);
743:   MPI_Comm_size(PetscObjectComm((PetscObject)buffer),&size);
744:   if (size>1) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Buffer must be a sequential vector");

746:   PetscObjectReference((PetscObject)buffer);
747:   VecDestroy(&bv->buffer);
748:   bv->buffer = buffer;
749:   PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->buffer);
750:   return(0);
751: }

753: /*@
754:    BVGetBufferVec - Obtain the buffer vector associated with the BV object.

756:    Not Collective, but Vec returned is parallel if BV is parallel

758:    Input Parameters:
759: .  bv - the basis vectors context

761:    Output Parameter:
762: .  buffer - vector

764:    Notes:
765:    The vector is created if not available previously. It is a sequential vector
766:    of length (nc+m)*m, where m is the number of columns of bv and nc is the number
767:    of constraints.

769:    Developer Notes:
770:    The buffer vector is viewed as a column-major matrix with leading dimension
771:    ld=nc+m, and m columns at most. In the most common usage, it has the structure
772: .vb
773:       | | C |
774:       |s|---|
775:       | | H |
776: .ve
777:    where H is an upper Hessenberg matrix of order m x (m-1), C contains coefficients
778:    related to orthogonalization against constraints (first nc rows), and s is the
779:    first column that contains scratch values computed during Gram-Schmidt
780:    orthogonalization. In particular, BVDotColumn() and BVMultColumn() use s to
781:    store the coefficients.

783:    Level: developer

785: .seealso: BVSetBufferVec(), BVSetSizes(), BVGetNumConstraints(), BVDotColumn(), BVMultColumn()
786: @*/
787: PetscErrorCode BVGetBufferVec(BV bv,Vec *buffer)
788: {
790:   PetscInt       ld;

795:   BVCheckSizes(bv,1);
796:   if (!bv->buffer) {
797:     ld = bv->m+bv->nc;
798:     VecCreate(PETSC_COMM_SELF,&bv->buffer);
799:     VecSetSizes(bv->buffer,PETSC_DECIDE,ld*bv->m);
800:     VecSetType(bv->buffer,((PetscObject)bv->t)->type_name);
801:     PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->buffer);
802:   }
803:   *buffer = bv->buffer;
804:   return(0);
805: }

807: /*@
808:    BVSetRandomContext - Sets the PetscRandom object associated with the BV,
809:    to be used in operations that need random numbers.

811:    Collective on bv

813:    Input Parameters:
814: +  bv   - the basis vectors context
815: -  rand - the random number generator context

817:    Level: advanced

819: .seealso: BVGetRandomContext(), BVSetRandom(), BVSetRandomColumn()
820: @*/
821: PetscErrorCode BVSetRandomContext(BV bv,PetscRandom rand)
822: {

829:   PetscObjectReference((PetscObject)rand);
830:   PetscRandomDestroy(&bv->rand);
831:   bv->rand = rand;
832:   PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->rand);
833:   return(0);
834: }

836: /*@
837:    BVGetRandomContext - Gets the PetscRandom object associated with the BV.

839:    Not Collective

841:    Input Parameter:
842: .  bv - the basis vectors context

844:    Output Parameter:
845: .  rand - the random number generator context

847:    Level: advanced

849: .seealso: BVSetRandomContext(), BVSetRandom(), BVSetRandomColumn()
850: @*/
851: PetscErrorCode BVGetRandomContext(BV bv,PetscRandom* rand)
852: {

858:   if (!bv->rand) {
859:     PetscRandomCreate(PetscObjectComm((PetscObject)bv),&bv->rand);
860:     PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->rand);
861:     if (bv->rrandom) {
862:       PetscRandomSetSeed(bv->rand,0x12345678);
863:       PetscRandomSeed(bv->rand);
864:     }
865:   }
866:   *rand = bv->rand;
867:   return(0);
868: }

870: /*@
871:    BVSetFromOptions - Sets BV options from the options database.

873:    Collective on bv

875:    Input Parameter:
876: .  bv - the basis vectors context

878:    Level: beginner
879: @*/
880: PetscErrorCode BVSetFromOptions(BV bv)
881: {
882:   PetscErrorCode     ierr;
883:   char               type[256];
884:   PetscBool          flg1,flg2,flg3,flg4;
885:   PetscReal          r;
886:   BVOrthogType       otype;
887:   BVOrthogRefineType orefine;
888:   BVOrthogBlockType  oblock;

892:   BVRegisterAll();
893:   PetscObjectOptionsBegin((PetscObject)bv);
894:     PetscOptionsFList("-bv_type","Basis Vectors type","BVSetType",BVList,(char*)(((PetscObject)bv)->type_name?((PetscObject)bv)->type_name:BVSVEC),type,sizeof(type),&flg1);
895:     if (flg1) {
896:       BVSetType(bv,type);
897:     } else if (!((PetscObject)bv)->type_name) {
898:       BVSetType(bv,BVSVEC);
899:     }

901:     otype = bv->orthog_type;
902:     PetscOptionsEnum("-bv_orthog_type","Orthogonalization method","BVSetOrthogonalization",BVOrthogTypes,(PetscEnum)otype,(PetscEnum*)&otype,&flg1);
903:     orefine = bv->orthog_ref;
904:     PetscOptionsEnum("-bv_orthog_refine","Iterative refinement mode during orthogonalization","BVSetOrthogonalization",BVOrthogRefineTypes,(PetscEnum)orefine,(PetscEnum*)&orefine,&flg2);
905:     r = bv->orthog_eta;
906:     PetscOptionsReal("-bv_orthog_eta","Parameter of iterative refinement during orthogonalization","BVSetOrthogonalization",r,&r,&flg3);
907:     oblock = bv->orthog_block;
908:     PetscOptionsEnum("-bv_orthog_block","Block orthogonalization method","BVSetOrthogonalization",BVOrthogBlockTypes,(PetscEnum)oblock,(PetscEnum*)&oblock,&flg4);
909:     if (flg1 || flg2 || flg3 || flg4) { BVSetOrthogonalization(bv,otype,orefine,r,oblock); }

911:     PetscOptionsEnum("-bv_matmult","Method for BVMatMult","BVSetMatMultMethod",BVMatMultTypes,(PetscEnum)bv->vmm,(PetscEnum*)&bv->vmm,NULL);

913:     /* undocumented option to generate random vectors that are independent of the number of processes */
914:     PetscOptionsGetBool(NULL,NULL,"-bv_reproducible_random",&bv->rrandom,NULL);

916:     if (bv->ops->create) bv->defersfo = PETSC_TRUE;   /* defer call to setfromoptions */
917:     else if (bv->ops->setfromoptions) {
918:       (*bv->ops->setfromoptions)(PetscOptionsObject,bv);
919:     }
920:     PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)bv);
921:   PetscOptionsEnd();

923:   if (!bv->rand) { BVGetRandomContext(bv,&bv->rand); }
924:   PetscRandomSetFromOptions(bv->rand);
925:   return(0);
926: }

928: /*@
929:    BVSetOrthogonalization - Specifies the method used for the orthogonalization of
930:    vectors (classical or modified Gram-Schmidt with or without refinement), and
931:    for the block-orthogonalization (simultaneous orthogonalization of a set of
932:    vectors).

934:    Logically Collective on bv

936:    Input Parameters:
937: +  bv     - the basis vectors context
938: .  type   - the method of vector orthogonalization
939: .  refine - type of refinement
940: .  eta    - parameter for selective refinement
941: -  block  - the method of block orthogonalization

943:    Options Database Keys:
944: +  -bv_orthog_type <type> - Where <type> is cgs for Classical Gram-Schmidt orthogonalization
945:                          (default) or mgs for Modified Gram-Schmidt orthogonalization
946: .  -bv_orthog_refine <ref> - Where <ref> is one of never, ifneeded (default) or always
947: .  -bv_orthog_eta <eta> -  For setting the value of eta
948: -  -bv_orthog_block <block> - Where <block> is the block-orthogonalization method

950:    Notes:
951:    The default settings work well for most problems.

953:    The parameter eta should be a real value between 0 and 1 (or PETSC_DEFAULT).
954:    The value of eta is used only when the refinement type is "ifneeded".

956:    When using several processors, MGS is likely to result in bad scalability.

958:    If the method set for block orthogonalization is GS, then the computation
959:    is done column by column with the vector orthogonalization.

961:    Level: advanced

963: .seealso: BVOrthogonalizeColumn(), BVGetOrthogonalization(), BVOrthogType, BVOrthogRefineType, BVOrthogBlockType
964: @*/
965: PetscErrorCode BVSetOrthogonalization(BV bv,BVOrthogType type,BVOrthogRefineType refine,PetscReal eta,BVOrthogBlockType block)
966: {
973:   switch (type) {
974:     case BV_ORTHOG_CGS:
975:     case BV_ORTHOG_MGS:
976:       bv->orthog_type = type;
977:       break;
978:     default:
979:       SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
980:   }
981:   switch (refine) {
982:     case BV_ORTHOG_REFINE_NEVER:
983:     case BV_ORTHOG_REFINE_IFNEEDED:
984:     case BV_ORTHOG_REFINE_ALWAYS:
985:       bv->orthog_ref = refine;
986:       break;
987:     default:
988:       SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown refinement type");
989:   }
990:   if (eta == PETSC_DEFAULT) {
991:     bv->orthog_eta = 0.7071;
992:   } else {
993:     if (eta <= 0.0 || eta > 1.0) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Invalid eta value");
994:     bv->orthog_eta = eta;
995:   }
996:   switch (block) {
997:     case BV_ORTHOG_BLOCK_GS:
998:     case BV_ORTHOG_BLOCK_CHOL:
999:     case BV_ORTHOG_BLOCK_TSQR:
1000:     case BV_ORTHOG_BLOCK_TSQRCHOL:
1001:     case BV_ORTHOG_BLOCK_SVQB:
1002:       bv->orthog_block = block;
1003:       break;
1004:     default:
1005:       SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown block orthogonalization type");
1006:   }
1007:   return(0);
1008: }

1010: /*@
1011:    BVGetOrthogonalization - Gets the orthogonalization settings from the BV object.

1013:    Not Collective

1015:    Input Parameter:
1016: .  bv - basis vectors context

1018:    Output Parameter:
1019: +  type   - the method of vector orthogonalization
1020: .  refine - type of refinement
1021: .  eta    - parameter for selective refinement
1022: -  block  - the method of block orthogonalization

1024:    Level: advanced

1026: .seealso: BVOrthogonalizeColumn(), BVSetOrthogonalization(), BVOrthogType, BVOrthogRefineType, BVOrthogBlockType
1027: @*/
1028: PetscErrorCode BVGetOrthogonalization(BV bv,BVOrthogType *type,BVOrthogRefineType *refine,PetscReal *eta,BVOrthogBlockType *block)
1029: {
1032:   if (type)   *type   = bv->orthog_type;
1033:   if (refine) *refine = bv->orthog_ref;
1034:   if (eta)    *eta    = bv->orthog_eta;
1035:   if (block)  *block  = bv->orthog_block;
1036:   return(0);
1037: }

1039: /*@
1040:    BVSetMatMultMethod - Specifies the method used for the BVMatMult() operation.

1042:    Logically Collective on bv

1044:    Input Parameters:
1045: +  bv     - the basis vectors context
1046: -  method - the method for the BVMatMult() operation

1048:    Options Database Keys:
1049: .  -bv_matmult <meth> - choose one of the methods: vecs, mat

1051:    Notes:
1052:    Allowed values are
1053: +  BV_MATMULT_VECS - perform a matrix-vector multiply per each column
1054: .  BV_MATMULT_MAT - carry out a Mat-Mat product with a dense matrix
1055: -  BV_MATMULT_MAT_SAVE - this case is deprecated

1057:    The default is BV_MATMULT_MAT except in the case of BVVECS.

1059:    Level: advanced

1061: .seealso: BVMatMult(), BVGetMatMultMethod(), BVMatMultType
1062: @*/
1063: PetscErrorCode BVSetMatMultMethod(BV bv,BVMatMultType method)
1064: {

1070:   switch (method) {
1071:     case BV_MATMULT_VECS:
1072:     case BV_MATMULT_MAT:
1073:       bv->vmm = method;
1074:       break;
1075:     case BV_MATMULT_MAT_SAVE:
1076:       PetscInfo(bv,"BV_MATMULT_MAT_SAVE is deprecated, using BV_MATMULT_MAT\n");
1077:       bv->vmm = BV_MATMULT_MAT;
1078:       break;
1079:     default:
1080:       SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown matmult method");
1081:   }
1082:   return(0);
1083: }

1085: /*@
1086:    BVGetMatMultMethod - Gets the method used for the BVMatMult() operation.

1088:    Not Collective

1090:    Input Parameter:
1091: .  bv - basis vectors context

1093:    Output Parameter:
1094: .  method - the method for the BVMatMult() operation

1096:    Level: advanced

1098: .seealso: BVMatMult(), BVSetMatMultMethod(), BVMatMultType
1099: @*/
1100: PetscErrorCode BVGetMatMultMethod(BV bv,BVMatMultType *method)
1101: {
1105:   *method = bv->vmm;
1106:   return(0);
1107: }

1109: /*@
1110:    BVGetColumn - Returns a Vec object that contains the entries of the
1111:    requested column of the basis vectors object.

1113:    Logically Collective on bv

1115:    Input Parameters:
1116: +  bv - the basis vectors context
1117: -  j  - the index of the requested column

1119:    Output Parameter:
1120: .  v  - vector containing the jth column

1122:    Notes:
1123:    The returned Vec must be seen as a reference (not a copy) of the BV
1124:    column, that is, modifying the Vec will change the BV entries as well.

1126:    The returned Vec must not be destroyed. BVRestoreColumn() must be
1127:    called when it is no longer needed. At most, two columns can be fetched,
1128:    that is, this function can only be called twice before the corresponding
1129:    BVRestoreColumn() is invoked.

1131:    A negative index j selects the i-th constraint, where i=-j. Constraints
1132:    should not be modified.

1134:    Level: beginner

1136: .seealso: BVRestoreColumn(), BVInsertConstraints()
1137: @*/
1138: PetscErrorCode BVGetColumn(BV bv,PetscInt j,Vec *v)
1139: {
1141:   PetscInt       l;

1146:   BVCheckSizes(bv,1);
1147:   BVCheckOp(bv,1,getcolumn);
1149:   if (j<0 && -j>bv->nc) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"You requested constraint %D but only %D are available",-j,bv->nc);
1150:   if (j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"You requested column %D but only %D are available",j,bv->m);
1151:   if (j==bv->ci[0] || j==bv->ci[1]) SETERRQ1(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Column %D already fetched in a previous call to BVGetColumn",j);
1152:   l = BVAvailableVec;
1153:   if (l==-1) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Too many requested columns; you must call BVRestoreColumn for one of the previously fetched columns");
1154:   (*bv->ops->getcolumn)(bv,j,v);
1155:   bv->ci[l] = j;
1156:   PetscObjectStateGet((PetscObject)bv->cv[l],&bv->st[l]);
1157:   PetscObjectGetId((PetscObject)bv->cv[l],&bv->id[l]);
1158:   *v = bv->cv[l];
1159:   return(0);
1160: }

1162: /*@
1163:    BVRestoreColumn - Restore a column obtained with BVGetColumn().

1165:    Logically Collective on bv

1167:    Input Parameters:
1168: +  bv - the basis vectors context
1169: .  j  - the index of the column
1170: -  v  - vector obtained with BVGetColumn()

1172:    Note:
1173:    The arguments must match the corresponding call to BVGetColumn().

1175:    Level: beginner

1177: .seealso: BVGetColumn()
1178: @*/
1179: PetscErrorCode BVRestoreColumn(BV bv,PetscInt j,Vec *v)
1180: {
1181:   PetscErrorCode   ierr;
1182:   PetscObjectId    id;
1183:   PetscObjectState st;
1184:   PetscInt         l;

1189:   BVCheckSizes(bv,1);
1193:   if (j<0 && -j>bv->nc) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"You requested constraint %D but only %D are available",-j,bv->nc);
1194:   if (j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"You requested column %D but only %D are available",j,bv->m);
1195:   if (j!=bv->ci[0] && j!=bv->ci[1]) SETERRQ1(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Column %D has not been fetched with a call to BVGetColumn",j);
1196:   l = (j==bv->ci[0])? 0: 1;
1197:   PetscObjectGetId((PetscObject)*v,&id);
1198:   if (id!=bv->id[l]) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Argument 3 is not the same Vec that was obtained with BVGetColumn");
1199:   PetscObjectStateGet((PetscObject)*v,&st);
1200:   if (st!=bv->st[l]) {
1201:     PetscObjectStateIncrease((PetscObject)bv);
1202:   }
1203:   if (bv->ops->restorecolumn) {
1204:     (*bv->ops->restorecolumn)(bv,j,v);
1205:   } else bv->cv[l] = NULL;
1206:   bv->ci[l] = -bv->nc-1;
1207:   bv->st[l] = -1;
1208:   bv->id[l] = 0;
1209:   *v = NULL;
1210:   return(0);
1211: }

1213: /*@C
1214:    BVGetArray - Returns a pointer to a contiguous array that contains this
1215:    processor's portion of the BV data.

1217:    Logically Collective on bv

1219:    Input Parameters:
1220: .  bv - the basis vectors context

1222:    Output Parameter:
1223: .  a  - location to put pointer to the array

1225:    Notes:
1226:    BVRestoreArray() must be called when access to the array is no longer needed.
1227:    This operation may imply a data copy, for BV types that do not store
1228:    data contiguously in memory.

1230:    The pointer will normally point to the first entry of the first column,
1231:    but if the BV has constraints then these go before the regular columns.

1233:    Level: advanced

1235: .seealso: BVRestoreArray(), BVInsertConstraints()
1236: @*/
1237: PetscErrorCode BVGetArray(BV bv,PetscScalar **a)
1238: {

1244:   BVCheckSizes(bv,1);
1245:   BVCheckOp(bv,1,getarray);
1246:   (*bv->ops->getarray)(bv,a);
1247:   return(0);
1248: }

1250: /*@C
1251:    BVRestoreArray - Restore the BV object after BVGetArray() has been called.

1253:    Logically Collective on bv

1255:    Input Parameters:
1256: +  bv - the basis vectors context
1257: -  a  - location of pointer to array obtained from BVGetArray()

1259:    Note:
1260:    This operation may imply a data copy, for BV types that do not store
1261:    data contiguously in memory.

1263:    Level: advanced

1265: .seealso: BVGetColumn()
1266: @*/
1267: PetscErrorCode BVRestoreArray(BV bv,PetscScalar **a)
1268: {

1274:   BVCheckSizes(bv,1);
1275:   if (bv->ops->restorearray) {
1276:     (*bv->ops->restorearray)(bv,a);
1277:   }
1278:   if (a) *a = NULL;
1279:   PetscObjectStateIncrease((PetscObject)bv);
1280:   return(0);
1281: }

1283: /*@C
1284:    BVGetArrayRead - Returns a read-only pointer to a contiguous array that
1285:    contains this processor's portion of the BV data.

1287:    Not Collective

1289:    Input Parameters:
1290: .  bv - the basis vectors context

1292:    Output Parameter:
1293: .  a  - location to put pointer to the array

1295:    Notes:
1296:    BVRestoreArrayRead() must be called when access to the array is no
1297:    longer needed. This operation may imply a data copy, for BV types that
1298:    do not store data contiguously in memory.

1300:    The pointer will normally point to the first entry of the first column,
1301:    but if the BV has constraints then these go before the regular columns.

1303:    Level: advanced

1305: .seealso: BVRestoreArray(), BVInsertConstraints()
1306: @*/
1307: PetscErrorCode BVGetArrayRead(BV bv,const PetscScalar **a)
1308: {

1314:   BVCheckSizes(bv,1);
1315:   BVCheckOp(bv,1,getarrayread);
1316:   (*bv->ops->getarrayread)(bv,a);
1317:   return(0);
1318: }

1320: /*@C
1321:    BVRestoreArrayRead - Restore the BV object after BVGetArrayRead() has
1322:    been called.

1324:    Not Collective

1326:    Input Parameters:
1327: +  bv - the basis vectors context
1328: -  a  - location of pointer to array obtained from BVGetArrayRead()

1330:    Level: advanced

1332: .seealso: BVGetColumn()
1333: @*/
1334: PetscErrorCode BVRestoreArrayRead(BV bv,const PetscScalar **a)
1335: {

1341:   BVCheckSizes(bv,1);
1342:   if (bv->ops->restorearrayread) {
1343:     (*bv->ops->restorearrayread)(bv,a);
1344:   }
1345:   if (a) *a = NULL;
1346:   return(0);
1347: }

1349: /*@
1350:    BVCreateVec - Creates a new Vec object with the same type and dimensions
1351:    as the columns of the basis vectors object.

1353:    Collective on bv

1355:    Input Parameter:
1356: .  bv - the basis vectors context

1358:    Output Parameter:
1359: .  v  - the new vector

1361:    Note:
1362:    The user is responsible of destroying the returned vector.

1364:    Level: beginner

1366: .seealso: BVCreateMat()
1367: @*/
1368: PetscErrorCode BVCreateVec(BV bv,Vec *v)
1369: {

1374:   BVCheckSizes(bv,1);
1376:   VecDuplicate(bv->t,v);
1377:   return(0);
1378: }

1380: /*@
1381:    BVCreateMat - Creates a new Mat object of dense type and copies the contents
1382:    of the BV object.

1384:    Collective on bv

1386:    Input Parameter:
1387: .  bv - the basis vectors context

1389:    Output Parameter:
1390: .  A  - the new matrix

1392:    Notes:
1393:    The user is responsible of destroying the returned matrix.

1395:    The matrix contains all columns of the BV, not just the active columns.

1397:    Level: intermediate

1399: .seealso: BVCreateFromMat(), BVCreateVec(), BVGetMat()
1400: @*/
1401: PetscErrorCode BVCreateMat(BV bv,Mat *A)
1402: {
1403:   PetscErrorCode    ierr;
1404:   PetscScalar       *aa;
1405:   const PetscScalar *vv;

1409:   BVCheckSizes(bv,1);

1412:   MatCreateDense(PetscObjectComm((PetscObject)bv->t),bv->n,PETSC_DECIDE,bv->N,bv->m,NULL,A);
1413:   MatDenseGetArray(*A,&aa);
1414:   BVGetArrayRead(bv,&vv);
1415:   PetscArraycpy(aa,vv,bv->m*bv->n);
1416:   BVRestoreArrayRead(bv,&vv);
1417:   MatDenseRestoreArray(*A,&aa);
1418:   MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
1419:   MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
1420:   return(0);
1421: }

1423: /*@
1424:    BVGetMat - Returns a Mat object of dense type that shares the memory of
1425:    the BV object.

1427:    Collective on bv

1429:    Input Parameter:
1430: .  bv - the basis vectors context

1432:    Output Parameter:
1433: .  A  - the matrix

1435:    Notes:
1436:    The returned matrix contains only the active columns. If the content of
1437:    the Mat is modified, these changes are also done in the BV object. The
1438:    user must call BVRestoreMat() when no longer needed.

1440:    This operation implies a call to BVGetArray(), which may result in data
1441:    copies.

1443:    Level: advanced

1445: .seealso: BVRestoreMat(), BVCreateMat(), BVGetArray()
1446: @*/
1447: PetscErrorCode BVGetMat(BV bv,Mat *A)
1448: {
1450:   PetscScalar    *vv,*aa;
1451:   PetscBool      create=PETSC_FALSE;
1452:   PetscInt       m,cols;

1456:   BVCheckSizes(bv,1);
1458:   if (bv->ops->getmat) {
1459:     (*bv->ops->getmat)(bv,A);
1460:   } else {
1461:     m = bv->k-bv->l;
1462:     if (!bv->Aget) create=PETSC_TRUE;
1463:     else {
1464:       MatDenseGetArray(bv->Aget,&aa);
1465:       if (aa) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"BVGetMat already called on this BV");
1466:       MatGetSize(bv->Aget,NULL,&cols);
1467:       if (cols!=m) {
1468:         MatDestroy(&bv->Aget);
1469:         create=PETSC_TRUE;
1470:       }
1471:     }
1472:     BVGetArray(bv,&vv);
1473:     if (create) {
1474:       MatCreateDense(PetscObjectComm((PetscObject)bv),bv->n,PETSC_DECIDE,bv->N,m,vv,&bv->Aget); /* pass a pointer to avoid allocation of storage */
1475:       MatDenseReplaceArray(bv->Aget,NULL);  /* replace with a null pointer, the value after BVRestoreMat */
1476:       PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->Aget);
1477:     }
1478:     MatDensePlaceArray(bv->Aget,vv+(bv->nc+bv->l)*bv->n);  /* set the actual pointer */
1479:     *A = bv->Aget;
1480:   }
1481:   return(0);
1482: }

1484: /*@
1485:    BVRestoreMat - Restores the Mat obtained with BVGetMat().

1487:    Logically Collective on bv

1489:    Input Parameters:
1490: +  bv - the basis vectors context
1491: -  A  - the fetched matrix

1493:    Note:
1494:    A call to this function must match a previous call of BVGetMat().
1495:    The effect is that the contents of the Mat are copied back to the
1496:    BV internal data structures.

1498:    Level: advanced

1500: .seealso: BVGetMat(), BVRestoreArray()
1501: @*/
1502: PetscErrorCode BVRestoreMat(BV bv,Mat *A)
1503: {
1505:   PetscScalar    *vv,*aa;

1509:   BVCheckSizes(bv,1);
1511:   if (!bv->Aget) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"BVRestoreMat must match a previous call to BVGetMat");
1512:   if (bv->Aget!=*A) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"Mat argument is not the same as the one obtained with BVGetMat");
1513:   if (bv->ops->restoremat) {
1514:     (*bv->ops->restoremat)(bv,A);
1515:   } else {
1516:     MatDenseGetArray(bv->Aget,&aa);
1517:     vv = aa-(bv->nc+bv->l)*bv->n;
1518:     MatDenseResetArray(bv->Aget);
1519:     BVRestoreArray(bv,&vv);
1520:     *A = NULL;
1521:   }
1522:   return(0);
1523: }

1525: /*
1526:    Copy all user-provided attributes of V to another BV object W
1527:  */
1528: PETSC_STATIC_INLINE PetscErrorCode BVDuplicate_Private(BV V,BV W)
1529: {

1533:   BVSetType(W,((PetscObject)V)->type_name);
1534:   W->orthog_type  = V->orthog_type;
1535:   W->orthog_ref   = V->orthog_ref;
1536:   W->orthog_eta   = V->orthog_eta;
1537:   W->orthog_block = V->orthog_block;
1538:   if (V->matrix) { PetscObjectReference((PetscObject)V->matrix); }
1539:   W->matrix       = V->matrix;
1540:   W->indef        = V->indef;
1541:   W->vmm          = V->vmm;
1542:   W->rrandom      = V->rrandom;
1543:   if (V->rand) { PetscObjectReference((PetscObject)V->rand); }
1544:   W->rand         = V->rand;
1545:   if (V->ops->duplicate) { (*V->ops->duplicate)(V,W); }
1546:   PetscObjectStateIncrease((PetscObject)W);
1547:   return(0);
1548: }

1550: /*@
1551:    BVDuplicate - Creates a new basis vector object of the same type and
1552:    dimensions as an existing one.

1554:    Collective on V

1556:    Input Parameter:
1557: .  V - basis vectors context

1559:    Output Parameter:
1560: .  W - location to put the new BV

1562:    Notes:
1563:    The new BV has the same type and dimensions as V, and it shares the same
1564:    template vector. Also, the inner product matrix and orthogonalization
1565:    options are copied.

1567:    BVDuplicate() DOES NOT COPY the entries, but rather allocates storage
1568:    for the new basis vectors. Use BVCopy() to copy the contents.

1570:    Level: intermediate

1572: .seealso: BVDuplicateResize(), BVCreate(), BVSetSizesFromVec(), BVCopy()
1573: @*/
1574: PetscErrorCode BVDuplicate(BV V,BV *W)
1575: {

1581:   BVCheckSizes(V,1);
1583:   BVCreate(PetscObjectComm((PetscObject)V),W);
1584:   BVSetSizesFromVec(*W,V->t,V->m);
1585:   BVDuplicate_Private(V,*W);
1586:   return(0);
1587: }

1589: /*@
1590:    BVDuplicateResize - Creates a new basis vector object of the same type and
1591:    dimensions as an existing one, but with possibly different number of columns.

1593:    Collective on V

1595:    Input Parameter:
1596: +  V - basis vectors context
1597: -  m - the new number of columns

1599:    Output Parameter:
1600: .  W - location to put the new BV

1602:    Note:
1603:    This is equivalent of a call to BVDuplicate() followed by BVResize(). The
1604:    contents of V are not copied to W.

1606:    Level: intermediate

1608: .seealso: BVDuplicate(), BVResize()
1609: @*/
1610: PetscErrorCode BVDuplicateResize(BV V,PetscInt m,BV *W)
1611: {

1617:   BVCheckSizes(V,1);
1620:   BVCreate(PetscObjectComm((PetscObject)V),W);
1621:   BVSetSizesFromVec(*W,V->t,m);
1622:   BVDuplicate_Private(V,*W);
1623:   return(0);
1624: }

1626: /*@
1627:    BVGetCachedBV - Returns a BV object stored internally that holds the
1628:    result of B*X after a call to BVApplyMatrixBV().

1630:    Not collective

1632:    Input Parameter:
1633: .  bv    - the basis vectors context

1635:    Output Parameter:
1636: .  cached - the cached BV

1638:    Note:
1639:    The cached BV is created if not available previously.

1641:    Level: developer

1643: .seealso: BVApplyMatrixBV()
1644: @*/
1645: PetscErrorCode BVGetCachedBV(BV bv,BV *cached)
1646: {

1652:   BVCheckSizes(bv,1);
1653:   if (!bv->cached) {
1654:     BVCreate(PetscObjectComm((PetscObject)bv),&bv->cached);
1655:     BVSetSizesFromVec(bv->cached,bv->t,bv->m);
1656:     BVDuplicate_Private(bv,bv->cached);
1657:     PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->cached);
1658:   }
1659:   *cached = bv->cached;
1660:   return(0);
1661: }

1663: /*@
1664:    BVCopy - Copies a basis vector object into another one, W <- V.

1666:    Logically Collective on V

1668:    Input Parameter:
1669: .  V - basis vectors context

1671:    Output Parameter:
1672: .  W - the copy

1674:    Note:
1675:    Both V and W must be distributed in the same manner; local copies are
1676:    done. Only active columns (excluding the leading ones) are copied.
1677:    In the destination W, columns are overwritten starting from the leading ones.
1678:    Constraints are not copied.

1680:    Level: beginner

1682: .seealso: BVCopyVec(), BVCopyColumn(), BVDuplicate(), BVSetActiveColumns()
1683: @*/
1684: PetscErrorCode BVCopy(BV V,BV W)
1685: {
1686:   PetscErrorCode    ierr;
1687:   PetscScalar       *womega;
1688:   const PetscScalar *vomega;

1693:   BVCheckSizes(V,1);
1694:   BVCheckOp(V,1,copy);
1697:   BVCheckSizes(W,2);
1699:   if (V->n!=W->n) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension V %D, W %D",V->n,W->n);
1700:   if (V->k-V->l>W->m-W->l) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"W has %D non-leading columns, not enough to store %D columns",W->m-W->l,V->k-V->l);
1701:   if (V==W || !V->n) return(0);

1703:   PetscLogEventBegin(BV_Copy,V,W,0,0);
1704:   if (V->indef && V->matrix && V->indef==W->indef && V->matrix==W->matrix) {
1705:     /* copy signature */
1706:     BV_AllocateSignature(W);
1707:     VecGetArrayRead(V->omega,&vomega);
1708:     VecGetArray(W->omega,&womega);
1709:     PetscArraycpy(womega+W->nc+W->l,vomega+V->nc+V->l,V->k-V->l);
1710:     VecRestoreArray(W->omega,&womega);
1711:     VecRestoreArrayRead(V->omega,&vomega);
1712:   }
1713:   (*V->ops->copy)(V,W);
1714:   PetscLogEventEnd(BV_Copy,V,W,0,0);
1715:   PetscObjectStateIncrease((PetscObject)W);
1716:   return(0);
1717: }

1719: /*@
1720:    BVCopyVec - Copies one of the columns of a basis vectors object into a Vec.

1722:    Logically Collective on V

1724:    Input Parameter:
1725: +  V - basis vectors context
1726: -  j - the column number to be copied

1728:    Output Parameter:
1729: .  w - the copied column

1731:    Note:
1732:    Both V and w must be distributed in the same manner; local copies are done.

1734:    Level: beginner

1736: .seealso: BVCopy(), BVCopyColumn(), BVInsertVec()
1737: @*/
1738: PetscErrorCode BVCopyVec(BV V,PetscInt j,Vec w)
1739: {
1741:   PetscInt       n,N;
1742:   Vec            z;

1747:   BVCheckSizes(V,1);

1752:   VecGetSize(w,&N);
1753:   VecGetLocalSize(w,&n);
1754:   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);

1756:   PetscLogEventBegin(BV_Copy,V,w,0,0);
1757:   BVGetColumn(V,j,&z);
1758:   VecCopy(z,w);
1759:   BVRestoreColumn(V,j,&z);
1760:   PetscLogEventEnd(BV_Copy,V,w,0,0);
1761:   return(0);
1762: }

1764: /*@
1765:    BVCopyColumn - Copies the values from one of the columns to another one.

1767:    Logically Collective on V

1769:    Input Parameter:
1770: +  V - basis vectors context
1771: .  j - the number of the source column
1772: -  i - the number of the destination column

1774:    Level: beginner

1776: .seealso: BVCopy(), BVCopyVec()
1777: @*/
1778: PetscErrorCode BVCopyColumn(BV V,PetscInt j,PetscInt i)
1779: {
1781:   Vec            z,w;
1782:   PetscScalar    *omega;

1787:   BVCheckSizes(V,1);
1790:   if (j==i) return(0);

1792:   PetscLogEventBegin(BV_Copy,V,0,0,0);
1793:   if (V->omega) {
1794:     VecGetArray(V->omega,&omega);
1795:     omega[i] = omega[j];
1796:     VecRestoreArray(V->omega,&omega);
1797:   }
1798:   if (V->ops->copycolumn) {
1799:     (*V->ops->copycolumn)(V,j,i);
1800:   } else {
1801:     BVGetColumn(V,j,&z);
1802:     BVGetColumn(V,i,&w);
1803:     VecCopy(z,w);
1804:     BVRestoreColumn(V,j,&z);
1805:     BVRestoreColumn(V,i,&w);
1806:   }
1807:   PetscLogEventEnd(BV_Copy,V,0,0,0);
1808:   PetscObjectStateIncrease((PetscObject)V);
1809:   return(0);
1810: }

1812: static PetscErrorCode BVGetSplit_Private(BV bv,PetscBool left,BV *split)
1813: {
1815:   PetscInt       ncols;

1818:   ncols = left? bv->nc+bv->l: bv->m-bv->l;
1819:   if (*split && ncols!=(*split)->m) { BVDestroy(split); }
1820:   if (!*split) {
1821:     BVCreate(PetscObjectComm((PetscObject)bv),split);
1822:     PetscLogObjectParent((PetscObject)bv,(PetscObject)*split);
1823:     (*split)->issplit = left? 1: 2;
1824:     (*split)->splitparent = bv;
1825:     BVSetSizesFromVec(*split,bv->t,ncols);
1826:     BVDuplicate_Private(bv,*split);
1827:   }
1828:   (*split)->l  = 0;
1829:   (*split)->k  = left? bv->l: bv->k-bv->l;
1830:   (*split)->nc = left? bv->nc: 0;
1831:   (*split)->m  = ncols-(*split)->nc;
1832:   if ((*split)->nc) {
1833:     (*split)->ci[0] = -(*split)->nc-1;
1834:     (*split)->ci[1] = -(*split)->nc-1;
1835:   }
1836:   if (left) {
1837:     PetscObjectStateGet((PetscObject)*split,&bv->lstate);
1838:   } else {
1839:     PetscObjectStateGet((PetscObject)*split,&bv->rstate);
1840:   }
1841:   return(0);
1842: }

1844: /*@
1845:    BVGetSplit - Splits the BV object into two BV objects that share the
1846:    internal data, one of them containing the leading columns and the other
1847:    one containing the remaining columns.

1849:    Logically Collective on bv

1851:    Input Parameters:
1852: .  bv - the basis vectors context

1854:    Output Parameters:
1855: +  L - left BV containing leading columns (can be NULL)
1856: -  R - right BV containing remaining columns (can be NULL)

1858:    Notes:
1859:    The columns are split in two sets. The leading columns (including the
1860:    constraints) are assigned to the left BV and the remaining columns
1861:    are assigned to the right BV. The number of leading columns, as
1862:    specified with BVSetActiveColumns(), must be between 1 and m-1 (to
1863:    guarantee that both L and R have at least one column).

1865:    The returned BV's must be seen as references (not copies) of the input
1866:    BV, that is, modifying them will change the entries of bv as well.
1867:    The returned BV's must not be destroyed. BVRestoreSplit() must be called
1868:    when they are no longer needed.

1870:    Pass NULL for any of the output BV's that is not needed.

1872:    Level: advanced

1874: .seealso: BVRestoreSplit(), BVSetActiveColumns(), BVSetNumConstraints()
1875: @*/
1876: PetscErrorCode BVGetSplit(BV bv,BV *L,BV *R)
1877: {

1883:   BVCheckSizes(bv,1);
1884:   if (!bv->l) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"Must indicate the number of leading columns with BVSetActiveColumns()");
1885:   if (bv->lsplit) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"Cannot get the split BV's twice before restoring them with BVRestoreSplit()");
1886:   bv->lsplit = bv->nc+bv->l;
1887:   BVGetSplit_Private(bv,PETSC_TRUE,&bv->L);
1888:   BVGetSplit_Private(bv,PETSC_FALSE,&bv->R);
1889:   if (L) *L = bv->L;
1890:   if (R) *R = bv->R;
1891:   return(0);
1892: }

1894: /*@
1895:    BVRestoreSplit - Restore the BV objects obtained with BVGetSplit().

1897:    Logically Collective on bv

1899:    Input Parameters:
1900: +  bv - the basis vectors context
1901: .  L  - left BV obtained with BVGetSplit()
1902: -  R  - right BV obtained with BVGetSplit()

1904:    Note:
1905:    The arguments must match the corresponding call to BVGetSplit().

1907:    Level: advanced

1909: .seealso: BVGetSplit()
1910: @*/
1911: PetscErrorCode BVRestoreSplit(BV bv,BV *L,BV *R)
1912: {

1918:   BVCheckSizes(bv,1);
1919:   if (!bv->lsplit) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"Must call BVGetSplit first");
1920:   if (L && *L!=bv->L) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Argument 2 is not the same BV that was obtained with BVGetSplit");
1921:   if (R && *R!=bv->R) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Argument 3 is not the same BV that was obtained with BVGetSplit");
1922:   if (L && ((*L)->ci[0]>(*L)->nc-1 || (*L)->ci[1]>(*L)->nc-1)) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"Argument 2 has unrestored columns, use BVRestoreColumn()");
1923:   if (R && ((*R)->ci[0]>(*R)->nc-1 || (*R)->ci[1]>(*R)->nc-1)) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"Argument 3 has unrestored columns, use BVRestoreColumn()");

1925:   if (bv->ops->restoresplit) {
1926:     (*bv->ops->restoresplit)(bv,L,R);
1927:   }
1928:   bv->lsplit = 0;
1929:   if (L) *L = NULL;
1930:   if (R) *R = NULL;
1931:   return(0);
1932: }