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> /*I "slepcbv.h" I*/
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->l = 0;
285: V->m = total-nc;
286: 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) {
499: MatGetLocalSize(B,&m,&n);
500: if (m!=n) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Matrix must be square");
501: if (bv->m && bv->n!=n) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension BV %D, Mat %D",bv->n,n);
502: }
503: if (B) { PetscObjectReference((PetscObject)B); }
504: MatDestroy(&bv->matrix);
505: bv->matrix = B;
506: bv->indef = indef;
507: PetscObjectStateIncrease((PetscObject)bv);
508: return(0);
509: }
511: /*@
512: BVGetMatrix - Retrieves the matrix representation of the inner product.
514: Not collective, though a parallel Mat may be returned
516: Input Parameter:
517: . bv - the basis vectors context
519: Output Parameter:
520: + B - the matrix of the inner product (may be NULL)
521: - indef - the flag indicating if the matrix is indefinite
523: Level: advanced
525: .seealso: BVSetMatrix()
526: @*/
527: PetscErrorCode BVGetMatrix(BV bv,Mat *B,PetscBool *indef)528: {
531: if (B) *B = bv->matrix;
532: if (indef) *indef = bv->indef;
533: return(0);
534: }
536: /*@
537: BVApplyMatrix - Multiplies a vector by the matrix representation of the
538: inner product.
540: Neighbor-wise Collective on bv
542: Input Parameter:
543: + bv - the basis vectors context
544: - x - the vector
546: Output Parameter:
547: . y - the result
549: Note:
550: If no matrix was specified this function copies the vector.
552: Level: advanced
554: .seealso: BVSetMatrix(), BVApplyMatrixBV()
555: @*/
556: PetscErrorCode BVApplyMatrix(BV bv,Vec x,Vec y)557: {
564: if (bv->matrix) {
565: BV_IPMatMult(bv,x);
566: VecCopy(bv->Bx,y);
567: } else {
568: VecCopy(x,y);
569: }
570: return(0);
571: }
573: /*@
574: BVApplyMatrixBV - Multiplies the BV vectors by the matrix representation
575: of the inner product.
577: Neighbor-wise Collective on X
579: Input Parameter:
580: . X - the basis vectors context
582: Output Parameter:
583: . Y - the basis vectors to store the result (optional)
585: Note:
586: This function computes Y = B*X, where B is the matrix given with
587: BVSetMatrix(). This operation is computed as in BVMatMult().
588: If no matrix was specified, then it just copies Y = X.
590: If no Y is given, the result is stored internally in the cached BV.
592: Level: developer
594: .seealso: BVSetMatrix(), BVApplyMatrix(), BVMatMult(), BVGetCachedBV()
595: @*/
596: PetscErrorCode BVApplyMatrixBV(BV X,BV Y)597: {
602: if (Y) {
604: if (X->matrix) {
605: BVMatMult(X,X->matrix,Y);
606: } else {
607: BVCopy(X,Y);
608: }
609: } else {
610: BV_IPMatMultBV(X);
611: }
612: return(0);
613: }
615: /*@
616: BVSetSignature - Sets the signature matrix to be used in orthogonalization.
618: Logically Collective on bv
620: Input Parameter:
621: + bv - the basis vectors context
622: - omega - a vector representing the diagonal of the signature matrix
624: Note:
625: The signature matrix Omega = V'*B*V is relevant only for an indefinite B.
627: Level: developer
629: .seealso: BVSetMatrix(), BVGetSignature()
630: @*/
631: PetscErrorCode BVSetSignature(BV bv,Vec omega)632: {
633: PetscErrorCode ierr;
634: PetscInt i,n;
635: const PetscScalar *pomega;
636: PetscScalar *intern;
640: BVCheckSizes(bv,1);
644: VecGetSize(omega,&n);
645: if (n!=bv->k) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Vec argument has %D elements, should be %D",n,bv->k);
646: BV_AllocateSignature(bv);
647: if (bv->indef) {
648: VecGetArrayRead(omega,&pomega);
649: VecGetArray(bv->omega,&intern);
650: for (i=0;i<n;i++) intern[bv->nc+i] = pomega[i];
651: VecRestoreArray(bv->omega,&intern);
652: VecRestoreArrayRead(omega,&pomega);
653: } else {
654: PetscInfo(bv,"Ignoring signature because BV is not indefinite\n");
655: }
656: PetscObjectStateIncrease((PetscObject)bv);
657: return(0);
658: }
660: /*@
661: BVGetSignature - Retrieves the signature matrix from last orthogonalization.
663: Not collective
665: Input Parameter:
666: . bv - the basis vectors context
668: Output Parameter:
669: . omega - a vector representing the diagonal of the signature matrix
671: Note:
672: The signature matrix Omega = V'*B*V is relevant only for an indefinite B.
674: Level: developer
676: .seealso: BVSetMatrix(), BVSetSignature()
677: @*/
678: PetscErrorCode BVGetSignature(BV bv,Vec omega)679: {
680: PetscErrorCode ierr;
681: PetscInt i,n;
682: PetscScalar *pomega;
683: const PetscScalar *intern;
687: BVCheckSizes(bv,1);
691: VecGetSize(omega,&n);
692: if (n!=bv->k) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Vec argument has %D elements, should be %D",n,bv->k);
693: if (bv->indef && bv->omega) {
694: VecGetArray(omega,&pomega);
695: VecGetArrayRead(bv->omega,&intern);
696: for (i=0;i<n;i++) pomega[i] = intern[bv->nc+i];
697: VecRestoreArrayRead(bv->omega,&intern);
698: VecRestoreArray(omega,&pomega);
699: } else {
700: VecSet(omega,1.0);
701: }
702: return(0);
703: }
705: /*@
706: BVSetBufferVec - Attach a vector object to be used as buffer space for
707: several operations.
709: Collective on bv
711: Input Parameters:
712: + bv - the basis vectors context)
713: - buffer - the vector
715: Notes:
716: Use BVGetBufferVec() to retrieve the vector (for example, to free it
717: at the end of the computations).
719: The vector must be sequential of length (nc+m)*m, where m is the number
720: of columns of bv and nc is the number of constraints.
722: Level: developer
724: .seealso: BVGetBufferVec(), BVSetSizes(), BVGetNumConstraints()
725: @*/
726: PetscErrorCode BVSetBufferVec(BV bv,Vec buffer)727: {
729: PetscInt ld,n;
730: PetscMPIInt size;
735: BVCheckSizes(bv,1);
736: VecGetSize(buffer,&n);
737: ld = bv->m+bv->nc;
738: if (n != ld*bv->m) SETERRQ1(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Buffer size must be %d",ld*bv->m);
739: MPI_Comm_size(PetscObjectComm((PetscObject)buffer),&size);
740: if (size>1) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Buffer must be a sequential vector");
742: PetscObjectReference((PetscObject)buffer);
743: VecDestroy(&bv->buffer);
744: bv->buffer = buffer;
745: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->buffer);
746: return(0);
747: }
749: /*@
750: BVGetBufferVec - Obtain the buffer vector associated with the BV object.
752: Not Collective, but Vec returned is parallel if BV is parallel
754: Input Parameters:
755: . bv - the basis vectors context
757: Output Parameter:
758: . buffer - vector
760: Notes:
761: The vector is created if not available previously. It is a sequential vector
762: of length (nc+m)*m, where m is the number of columns of bv and nc is the number
763: of constraints.
765: Developer Notes:
766: The buffer vector is viewed as a column-major matrix with leading dimension
767: ld=nc+m, and m columns at most. In the most common usage, it has the structure
768: .vb
769: | | C |
770: |s|---|
771: | | H |
772: .ve
773: where H is an upper Hessenberg matrix of order m x (m-1), C contains coefficients
774: related to orthogonalization against constraints (first nc rows), and s is the
775: first column that contains scratch values computed during Gram-Schmidt
776: orthogonalization. In particular, BVDotColumn() and BVMultColumn() use s to
777: store the coefficients.
779: Level: developer
781: .seealso: BVSetBufferVec(), BVSetSizes(), BVGetNumConstraints(), BVDotColumn(), BVMultColumn()
782: @*/
783: PetscErrorCode BVGetBufferVec(BV bv,Vec *buffer)784: {
786: PetscInt ld;
791: BVCheckSizes(bv,1);
792: if (!bv->buffer) {
793: ld = bv->m+bv->nc;
794: VecCreate(PETSC_COMM_SELF,&bv->buffer);
795: VecSetSizes(bv->buffer,PETSC_DECIDE,ld*bv->m);
796: VecSetType(bv->buffer,((PetscObject)bv->t)->type_name);
797: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->buffer);
798: }
799: *buffer = bv->buffer;
800: return(0);
801: }
803: /*@
804: BVSetRandomContext - Sets the PetscRandom object associated with the BV,
805: to be used in operations that need random numbers.
807: Collective on bv
809: Input Parameters:
810: + bv - the basis vectors context
811: - rand - the random number generator context
813: Level: advanced
815: .seealso: BVGetRandomContext(), BVSetRandom(), BVSetRandomColumn()
816: @*/
817: PetscErrorCode BVSetRandomContext(BV bv,PetscRandom rand)818: {
825: PetscObjectReference((PetscObject)rand);
826: PetscRandomDestroy(&bv->rand);
827: bv->rand = rand;
828: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->rand);
829: return(0);
830: }
832: /*@
833: BVGetRandomContext - Gets the PetscRandom object associated with the BV.
835: Not Collective
837: Input Parameter:
838: . bv - the basis vectors context
840: Output Parameter:
841: . rand - the random number generator context
843: Level: advanced
845: .seealso: BVSetRandomContext(), BVSetRandom(), BVSetRandomColumn()
846: @*/
847: PetscErrorCode BVGetRandomContext(BV bv,PetscRandom* rand)848: {
854: if (!bv->rand) {
855: PetscRandomCreate(PetscObjectComm((PetscObject)bv),&bv->rand);
856: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->rand);
857: if (bv->rrandom) {
858: PetscRandomSetSeed(bv->rand,0x12345678);
859: PetscRandomSeed(bv->rand);
860: }
861: }
862: *rand = bv->rand;
863: return(0);
864: }
866: /*@
867: BVSetFromOptions - Sets BV options from the options database.
869: Collective on bv
871: Input Parameter:
872: . bv - the basis vectors context
874: Level: beginner
875: @*/
876: PetscErrorCode BVSetFromOptions(BV bv)877: {
878: PetscErrorCode ierr;
879: char type[256];
880: PetscBool flg1,flg2,flg3,flg4;
881: PetscReal r;
882: BVOrthogType otype;
883: BVOrthogRefineType orefine;
884: BVOrthogBlockType oblock;
888: BVRegisterAll();
889: PetscObjectOptionsBegin((PetscObject)bv);
890: PetscOptionsFList("-bv_type","Basis Vectors type","BVSetType",BVList,(char*)(((PetscObject)bv)->type_name?((PetscObject)bv)->type_name:BVSVEC),type,256,&flg1);
891: if (flg1) {
892: BVSetType(bv,type);
893: } else if (!((PetscObject)bv)->type_name) {
894: BVSetType(bv,BVSVEC);
895: }
897: otype = bv->orthog_type;
898: PetscOptionsEnum("-bv_orthog_type","Orthogonalization method","BVSetOrthogonalization",BVOrthogTypes,(PetscEnum)otype,(PetscEnum*)&otype,&flg1);
899: orefine = bv->orthog_ref;
900: PetscOptionsEnum("-bv_orthog_refine","Iterative refinement mode during orthogonalization","BVSetOrthogonalization",BVOrthogRefineTypes,(PetscEnum)orefine,(PetscEnum*)&orefine,&flg2);
901: r = bv->orthog_eta;
902: PetscOptionsReal("-bv_orthog_eta","Parameter of iterative refinement during orthogonalization","BVSetOrthogonalization",r,&r,&flg3);
903: oblock = bv->orthog_block;
904: PetscOptionsEnum("-bv_orthog_block","Block orthogonalization method","BVSetOrthogonalization",BVOrthogBlockTypes,(PetscEnum)oblock,(PetscEnum*)&oblock,&flg4);
905: if (flg1 || flg2 || flg3 || flg4) { BVSetOrthogonalization(bv,otype,orefine,r,oblock); }
907: PetscOptionsEnum("-bv_matmult","Method for BVMatMult","BVSetMatMultMethod",BVMatMultTypes,(PetscEnum)bv->vmm,(PetscEnum*)&bv->vmm,NULL);
909: /* undocumented option to generate random vectors that are independent of the number of processes */
910: PetscOptionsGetBool(NULL,NULL,"-bv_reproducible_random",&bv->rrandom,NULL);
912: if (bv->ops->create) bv->defersfo = PETSC_TRUE; /* defer call to setfromoptions */
913: else if (bv->ops->setfromoptions) {
914: (*bv->ops->setfromoptions)(PetscOptionsObject,bv);
915: }
916: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)bv);
917: PetscOptionsEnd();
919: if (!bv->rand) { BVGetRandomContext(bv,&bv->rand); }
920: PetscRandomSetFromOptions(bv->rand);
921: return(0);
922: }
924: /*@
925: BVSetOrthogonalization - Specifies the method used for the orthogonalization of
926: vectors (classical or modified Gram-Schmidt with or without refinement), and
927: for the block-orthogonalization (simultaneous orthogonalization of a set of
928: vectors).
930: Logically Collective on bv
932: Input Parameters:
933: + bv - the basis vectors context
934: . type - the method of vector orthogonalization
935: . refine - type of refinement
936: . eta - parameter for selective refinement
937: - block - the method of block orthogonalization
939: Options Database Keys:
940: + -bv_orthog_type <type> - Where <type> is cgs for Classical Gram-Schmidt orthogonalization
941: (default) or mgs for Modified Gram-Schmidt orthogonalization
942: . -bv_orthog_refine <ref> - Where <ref> is one of never, ifneeded (default) or always
943: . -bv_orthog_eta <eta> - For setting the value of eta
944: - -bv_orthog_block <block> - Where <block> is the block-orthogonalization method
946: Notes:
947: The default settings work well for most problems.
949: The parameter eta should be a real value between 0 and 1 (or PETSC_DEFAULT).
950: The value of eta is used only when the refinement type is "ifneeded".
952: When using several processors, MGS is likely to result in bad scalability.
954: If the method set for block orthogonalization is GS, then the computation
955: is done column by column with the vector orthogonalization.
957: Level: advanced
959: .seealso: BVOrthogonalizeColumn(), BVGetOrthogonalization(), BVOrthogType, BVOrthogRefineType, BVOrthogBlockType960: @*/
961: PetscErrorCode BVSetOrthogonalization(BV bv,BVOrthogType type,BVOrthogRefineType refine,PetscReal eta,BVOrthogBlockType block)962: {
969: switch (type) {
970: case BV_ORTHOG_CGS:
971: case BV_ORTHOG_MGS:
972: bv->orthog_type = type;
973: break;
974: default:975: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
976: }
977: switch (refine) {
978: case BV_ORTHOG_REFINE_NEVER:
979: case BV_ORTHOG_REFINE_IFNEEDED:
980: case BV_ORTHOG_REFINE_ALWAYS:
981: bv->orthog_ref = refine;
982: break;
983: default:984: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown refinement type");
985: }
986: if (eta == PETSC_DEFAULT) {
987: bv->orthog_eta = 0.7071;
988: } else {
989: if (eta <= 0.0 || eta > 1.0) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Invalid eta value");
990: bv->orthog_eta = eta;
991: }
992: switch (block) {
993: case BV_ORTHOG_BLOCK_GS:
994: case BV_ORTHOG_BLOCK_CHOL:
995: case BV_ORTHOG_BLOCK_TSQR:
996: case BV_ORTHOG_BLOCK_TSQRCHOL:
997: case BV_ORTHOG_BLOCK_SVQB:
998: bv->orthog_block = block;
999: break;
1000: default:1001: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown block orthogonalization type");
1002: }
1003: return(0);
1004: }
1006: /*@
1007: BVGetOrthogonalization - Gets the orthogonalization settings from the BV object.
1009: Not Collective
1011: Input Parameter:
1012: . bv - basis vectors context
1014: Output Parameter:
1015: + type - the method of vector orthogonalization
1016: . refine - type of refinement
1017: . eta - parameter for selective refinement
1018: - block - the method of block orthogonalization
1020: Level: advanced
1022: .seealso: BVOrthogonalizeColumn(), BVSetOrthogonalization(), BVOrthogType, BVOrthogRefineType, BVOrthogBlockType1023: @*/
1024: PetscErrorCode BVGetOrthogonalization(BV bv,BVOrthogType *type,BVOrthogRefineType *refine,PetscReal *eta,BVOrthogBlockType *block)1025: {
1028: if (type) *type = bv->orthog_type;
1029: if (refine) *refine = bv->orthog_ref;
1030: if (eta) *eta = bv->orthog_eta;
1031: if (block) *block = bv->orthog_block;
1032: return(0);
1033: }
1035: /*@
1036: BVSetMatMultMethod - Specifies the method used for the BVMatMult() operation.
1038: Logically Collective on bv
1040: Input Parameters:
1041: + bv - the basis vectors context
1042: - method - the method for the BVMatMult() operation
1044: Options Database Keys:
1045: . -bv_matmult <meth> - choose one of the methods: vecs, mat
1047: Notes:
1048: Allowed values are
1049: + BV_MATMULT_VECS - perform a matrix-vector multiply per each column
1050: . BV_MATMULT_MAT - carry out a Mat-Mat product with a dense matrix
1051: - BV_MATMULT_MAT_SAVE - this case is deprecated
1053: The default is BV_MATMULT_MAT except in the case of BVVECS.
1055: Level: advanced
1057: .seealso: BVMatMult(), BVGetMatMultMethod(), BVMatMultType1058: @*/
1059: PetscErrorCode BVSetMatMultMethod(BV bv,BVMatMultType method)1060: {
1066: switch (method) {
1067: case BV_MATMULT_VECS:
1068: case BV_MATMULT_MAT:
1069: bv->vmm = method;
1070: break;
1071: case BV_MATMULT_MAT_SAVE:
1072: PetscInfo(bv,"BV_MATMULT_MAT_SAVE is deprecated, using BV_MATMULT_MAT\n");
1073: bv->vmm = BV_MATMULT_MAT;
1074: break;
1075: default:1076: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown matmult method");
1077: }
1078: return(0);
1079: }
1081: /*@
1082: BVGetMatMultMethod - Gets the method used for the BVMatMult() operation.
1084: Not Collective
1086: Input Parameter:
1087: . bv - basis vectors context
1089: Output Parameter:
1090: . method - the method for the BVMatMult() operation
1092: Level: advanced
1094: .seealso: BVMatMult(), BVSetMatMultMethod(), BVMatMultType1095: @*/
1096: PetscErrorCode BVGetMatMultMethod(BV bv,BVMatMultType *method)1097: {
1101: *method = bv->vmm;
1102: return(0);
1103: }
1105: /*@
1106: BVGetColumn - Returns a Vec object that contains the entries of the
1107: requested column of the basis vectors object.
1109: Logically Collective on bv
1111: Input Parameters:
1112: + bv - the basis vectors context
1113: - j - the index of the requested column
1115: Output Parameter:
1116: . v - vector containing the jth column
1118: Notes:
1119: The returned Vec must be seen as a reference (not a copy) of the BV1120: column, that is, modifying the Vec will change the BV entries as well.
1122: The returned Vec must not be destroyed. BVRestoreColumn() must be
1123: called when it is no longer needed. At most, two columns can be fetched,
1124: that is, this function can only be called twice before the corresponding
1125: BVRestoreColumn() is invoked.
1127: A negative index j selects the i-th constraint, where i=-j. Constraints
1128: should not be modified.
1130: Level: beginner
1132: .seealso: BVRestoreColumn(), BVInsertConstraints()
1133: @*/
1134: PetscErrorCode BVGetColumn(BV bv,PetscInt j,Vec *v)1135: {
1137: PetscInt l;
1142: BVCheckSizes(bv,1);
1143: BVCheckOp(bv,1,getcolumn);
1145: 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);
1146: if (j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"You requested column %D but only %D are available",j,bv->m);
1147: 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);
1148: l = BVAvailableVec;
1149: 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");
1150: (*bv->ops->getcolumn)(bv,j,v);
1151: bv->ci[l] = j;
1152: PetscObjectStateGet((PetscObject)bv->cv[l],&bv->st[l]);
1153: PetscObjectGetId((PetscObject)bv->cv[l],&bv->id[l]);
1154: *v = bv->cv[l];
1155: return(0);
1156: }
1158: /*@
1159: BVRestoreColumn - Restore a column obtained with BVGetColumn().
1161: Logically Collective on bv
1163: Input Parameters:
1164: + bv - the basis vectors context
1165: . j - the index of the column
1166: - v - vector obtained with BVGetColumn()
1168: Note:
1169: The arguments must match the corresponding call to BVGetColumn().
1171: Level: beginner
1173: .seealso: BVGetColumn()
1174: @*/
1175: PetscErrorCode BVRestoreColumn(BV bv,PetscInt j,Vec *v)1176: {
1177: PetscErrorCode ierr;
1178: PetscObjectId id;
1179: PetscObjectState st;
1180: PetscInt l;
1185: BVCheckSizes(bv,1);
1189: 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);
1190: if (j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"You requested column %D but only %D are available",j,bv->m);
1191: 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);
1192: l = (j==bv->ci[0])? 0: 1;
1193: PetscObjectGetId((PetscObject)*v,&id);
1194: 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");
1195: PetscObjectStateGet((PetscObject)*v,&st);
1196: if (st!=bv->st[l]) {
1197: PetscObjectStateIncrease((PetscObject)bv);
1198: }
1199: if (bv->ops->restorecolumn) {
1200: (*bv->ops->restorecolumn)(bv,j,v);
1201: } else bv->cv[l] = NULL;
1202: bv->ci[l] = -bv->nc-1;
1203: bv->st[l] = -1;
1204: bv->id[l] = 0;
1205: *v = NULL;
1206: return(0);
1207: }
1209: /*@C
1210: BVGetArray - Returns a pointer to a contiguous array that contains this
1211: processor's portion of the BV data.
1213: Logically Collective on bv
1215: Input Parameters:
1216: . bv - the basis vectors context
1218: Output Parameter:
1219: . a - location to put pointer to the array
1221: Notes:
1222: BVRestoreArray() must be called when access to the array is no longer needed.
1223: This operation may imply a data copy, for BV types that do not store
1224: data contiguously in memory.
1226: The pointer will normally point to the first entry of the first column,
1227: but if the BV has constraints then these go before the regular columns.
1229: Level: advanced
1231: .seealso: BVRestoreArray(), BVInsertConstraints()
1232: @*/
1233: PetscErrorCode BVGetArray(BV bv,PetscScalar **a)1234: {
1240: BVCheckSizes(bv,1);
1241: BVCheckOp(bv,1,getarray);
1242: (*bv->ops->getarray)(bv,a);
1243: return(0);
1244: }
1246: /*@C
1247: BVRestoreArray - Restore the BV object after BVGetArray() has been called.
1249: Logically Collective on bv
1251: Input Parameters:
1252: + bv - the basis vectors context
1253: - a - location of pointer to array obtained from BVGetArray()
1255: Note:
1256: This operation may imply a data copy, for BV types that do not store
1257: data contiguously in memory.
1259: Level: advanced
1261: .seealso: BVGetColumn()
1262: @*/
1263: PetscErrorCode BVRestoreArray(BV bv,PetscScalar **a)1264: {
1270: BVCheckSizes(bv,1);
1271: if (bv->ops->restorearray) {
1272: (*bv->ops->restorearray)(bv,a);
1273: }
1274: if (a) *a = NULL;
1275: PetscObjectStateIncrease((PetscObject)bv);
1276: return(0);
1277: }
1279: /*@C
1280: BVGetArrayRead - Returns a read-only pointer to a contiguous array that
1281: contains this processor's portion of the BV data.
1283: Not Collective
1285: Input Parameters:
1286: . bv - the basis vectors context
1288: Output Parameter:
1289: . a - location to put pointer to the array
1291: Notes:
1292: BVRestoreArrayRead() must be called when access to the array is no
1293: longer needed. This operation may imply a data copy, for BV types that
1294: do not store data contiguously in memory.
1296: The pointer will normally point to the first entry of the first column,
1297: but if the BV has constraints then these go before the regular columns.
1299: Level: advanced
1301: .seealso: BVRestoreArray(), BVInsertConstraints()
1302: @*/
1303: PetscErrorCode BVGetArrayRead(BV bv,const PetscScalar **a)1304: {
1310: BVCheckSizes(bv,1);
1311: BVCheckOp(bv,1,getarrayread);
1312: (*bv->ops->getarrayread)(bv,a);
1313: return(0);
1314: }
1316: /*@C
1317: BVRestoreArrayRead - Restore the BV object after BVGetArrayRead() has
1318: been called.
1320: Not Collective
1322: Input Parameters:
1323: + bv - the basis vectors context
1324: - a - location of pointer to array obtained from BVGetArrayRead()
1326: Level: advanced
1328: .seealso: BVGetColumn()
1329: @*/
1330: PetscErrorCode BVRestoreArrayRead(BV bv,const PetscScalar **a)1331: {
1337: BVCheckSizes(bv,1);
1338: if (bv->ops->restorearrayread) {
1339: (*bv->ops->restorearrayread)(bv,a);
1340: }
1341: if (a) *a = NULL;
1342: return(0);
1343: }
1345: /*@
1346: BVCreateVec - Creates a new Vec object with the same type and dimensions
1347: as the columns of the basis vectors object.
1349: Collective on bv
1351: Input Parameter:
1352: . bv - the basis vectors context
1354: Output Parameter:
1355: . v - the new vector
1357: Note:
1358: The user is responsible of destroying the returned vector.
1360: Level: beginner
1362: .seealso: BVCreateMat()
1363: @*/
1364: PetscErrorCode BVCreateVec(BV bv,Vec *v)1365: {
1370: BVCheckSizes(bv,1);
1372: VecDuplicate(bv->t,v);
1373: return(0);
1374: }
1376: /*@
1377: BVCreateMat - Creates a new Mat object of dense type and copies the contents
1378: of the BV object.
1380: Collective on bv
1382: Input Parameter:
1383: . bv - the basis vectors context
1385: Output Parameter:
1386: . A - the new matrix
1388: Notes:
1389: The user is responsible of destroying the returned matrix.
1391: The matrix contains all columns of the BV, not just the active columns.
1393: Level: intermediate
1395: .seealso: BVCreateFromMat(), BVCreateVec(), BVGetMat()
1396: @*/
1397: PetscErrorCode BVCreateMat(BV bv,Mat *A)1398: {
1399: PetscErrorCode ierr;
1400: PetscScalar *aa;
1401: const PetscScalar *vv;
1405: BVCheckSizes(bv,1);
1408: MatCreateDense(PetscObjectComm((PetscObject)bv->t),bv->n,PETSC_DECIDE,bv->N,bv->m,NULL,A);
1409: MatDenseGetArray(*A,&aa);
1410: BVGetArrayRead(bv,&vv);
1411: PetscArraycpy(aa,vv,bv->m*bv->n);
1412: BVRestoreArrayRead(bv,&vv);
1413: MatDenseRestoreArray(*A,&aa);
1414: MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
1415: MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
1416: return(0);
1417: }
1419: /*@
1420: BVGetMat - Returns a Mat object of dense type that shares the memory of
1421: the BV object.
1423: Collective on bv
1425: Input Parameter:
1426: . bv - the basis vectors context
1428: Output Parameter:
1429: . A - the matrix
1431: Notes:
1432: The returned matrix contains only the active columns. If the content of
1433: the Mat is modified, these changes are also done in the BV object. The
1434: user must call BVRestoreMat() when no longer needed.
1436: This operation implies a call to BVGetArray(), which may result in data
1437: copies.
1439: Level: advanced
1441: .seealso: BVRestoreMat(), BVCreateMat(), BVGetArray()
1442: @*/
1443: PetscErrorCode BVGetMat(BV bv,Mat *A)1444: {
1446: PetscScalar *vv,*aa;
1447: PetscBool create=PETSC_FALSE;
1448: PetscInt m,cols;
1452: BVCheckSizes(bv,1);
1454: m = bv->k-bv->l;
1455: if (!bv->Aget) create=PETSC_TRUE;
1456: else {
1457: MatDenseGetArray(bv->Aget,&aa);
1458: if (aa) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"BVGetMat already called on this BV");
1459: MatGetSize(bv->Aget,NULL,&cols);
1460: if (cols!=m) {
1461: MatDestroy(&bv->Aget);
1462: create=PETSC_TRUE;
1463: }
1464: }
1465: BVGetArray(bv,&vv);
1466: if (create) {
1467: MatCreateDense(PetscObjectComm((PetscObject)bv),bv->n,PETSC_DECIDE,bv->N,m,vv,&bv->Aget); /* pass a pointer to avoid allocation of storage */
1468: MatDensePlaceArray(bv->Aget,NULL); /* replace with a null pointer, the value after BVRestoreMat */
1469: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->Aget);
1470: }
1471: MatDensePlaceArray(bv->Aget,vv+(bv->nc+bv->l)*bv->n); /* set the actual pointer */
1472: *A = bv->Aget;
1473: return(0);
1474: }
1476: /*@
1477: BVRestoreMat - Restores the Mat obtained with BVGetMat().
1479: Logically Collective on bv
1481: Input Parameters:
1482: + bv - the basis vectors context
1483: - A - the fetched matrix
1485: Note:
1486: A call to this function must match a previous call of BVGetMat().
1487: The effect is that the contents of the Mat are copied back to the
1488: BV internal data structures.
1490: Level: advanced
1492: .seealso: BVGetMat(), BVRestoreArray()
1493: @*/
1494: PetscErrorCode BVRestoreMat(BV bv,Mat *A)1495: {
1497: PetscScalar *vv,*aa;
1501: BVCheckSizes(bv,1);
1503: if (!bv->Aget) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"BVRestoreMat must match a previous call to BVGetMat");
1504: if (bv->Aget!=*A) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"Mat argument is not the same as the one obtained with BVGetMat");
1505: MatDenseGetArray(bv->Aget,&aa);
1506: vv = aa-(bv->nc+bv->l)*bv->n;
1507: MatDenseResetArray(bv->Aget);
1508: BVRestoreArray(bv,&vv);
1509: *A = NULL;
1510: return(0);
1511: }
1513: /*
1514: Copy all user-provided attributes of V to another BV object W
1515: */
1516: PETSC_STATIC_INLINE PetscErrorCode BVDuplicate_Private(BV V,BV W)1517: {
1521: BVSetType(W,((PetscObject)V)->type_name);
1522: W->orthog_type = V->orthog_type;
1523: W->orthog_ref = V->orthog_ref;
1524: W->orthog_eta = V->orthog_eta;
1525: W->orthog_block = V->orthog_block;
1526: if (V->matrix) { PetscObjectReference((PetscObject)V->matrix); }
1527: W->matrix = V->matrix;
1528: W->indef = V->indef;
1529: W->vmm = V->vmm;
1530: W->rrandom = V->rrandom;
1531: if (V->rand) { PetscObjectReference((PetscObject)V->rand); }
1532: W->rand = V->rand;
1533: if (V->ops->duplicate) { (*V->ops->duplicate)(V,W); }
1534: PetscObjectStateIncrease((PetscObject)W);
1535: return(0);
1536: }
1538: /*@
1539: BVDuplicate - Creates a new basis vector object of the same type and
1540: dimensions as an existing one.
1542: Collective on V
1544: Input Parameter:
1545: . V - basis vectors context
1547: Output Parameter:
1548: . W - location to put the new BV1550: Notes:
1551: The new BV has the same type and dimensions as V, and it shares the same
1552: template vector. Also, the inner product matrix and orthogonalization
1553: options are copied.
1555: BVDuplicate() DOES NOT COPY the entries, but rather allocates storage
1556: for the new basis vectors. Use BVCopy() to copy the contents.
1558: Level: intermediate
1560: .seealso: BVDuplicateResize(), BVCreate(), BVSetSizesFromVec(), BVCopy()
1561: @*/
1562: PetscErrorCode BVDuplicate(BV V,BV *W)1563: {
1569: BVCheckSizes(V,1);
1571: BVCreate(PetscObjectComm((PetscObject)V),W);
1572: BVSetSizesFromVec(*W,V->t,V->m);
1573: BVDuplicate_Private(V,*W);
1574: return(0);
1575: }
1577: /*@
1578: BVDuplicateResize - Creates a new basis vector object of the same type and
1579: dimensions as an existing one, but with possibly different number of columns.
1581: Collective on V
1583: Input Parameter:
1584: + V - basis vectors context
1585: - m - the new number of columns
1587: Output Parameter:
1588: . W - location to put the new BV1590: Note:
1591: This is equivalent of a call to BVDuplicate() followed by BVResize(). The
1592: contents of V are not copied to W.
1594: Level: intermediate
1596: .seealso: BVDuplicate(), BVResize()
1597: @*/
1598: PetscErrorCode BVDuplicateResize(BV V,PetscInt m,BV *W)1599: {
1605: BVCheckSizes(V,1);
1608: BVCreate(PetscObjectComm((PetscObject)V),W);
1609: BVSetSizesFromVec(*W,V->t,m);
1610: BVDuplicate_Private(V,*W);
1611: return(0);
1612: }
1614: /*@
1615: BVGetCachedBV - Returns a BV object stored internally that holds the
1616: result of B*X after a call to BVApplyMatrixBV().
1618: Not collective
1620: Input Parameter:
1621: . bv - the basis vectors context
1623: Output Parameter:
1624: . cached - the cached BV1626: Note:
1627: The cached BV is created if not available previously.
1629: Level: developer
1631: .seealso: BVApplyMatrixBV()
1632: @*/
1633: PetscErrorCode BVGetCachedBV(BV bv,BV *cached)1634: {
1640: BVCheckSizes(bv,1);
1641: if (!bv->cached) {
1642: BVCreate(PetscObjectComm((PetscObject)bv),&bv->cached);
1643: BVSetSizesFromVec(bv->cached,bv->t,bv->m);
1644: BVDuplicate_Private(bv,bv->cached);
1645: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->cached);
1646: }
1647: *cached = bv->cached;
1648: return(0);
1649: }
1651: /*@
1652: BVCopy - Copies a basis vector object into another one, W <- V.
1654: Logically Collective on V
1656: Input Parameter:
1657: . V - basis vectors context
1659: Output Parameter:
1660: . W - the copy
1662: Note:
1663: Both V and W must be distributed in the same manner; local copies are
1664: done. Only active columns (excluding the leading ones) are copied.
1665: In the destination W, columns are overwritten starting from the leading ones.
1666: Constraints are not copied.
1668: Level: beginner
1670: .seealso: BVCopyVec(), BVCopyColumn(), BVDuplicate(), BVSetActiveColumns()
1671: @*/
1672: PetscErrorCode BVCopy(BV V,BV W)1673: {
1674: PetscErrorCode ierr;
1675: PetscScalar *womega;
1676: const PetscScalar *vomega;
1681: BVCheckSizes(V,1);
1682: BVCheckOp(V,1,copy);
1685: BVCheckSizes(W,2);
1687: if (V->n!=W->n) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension V %D, W %D",V->n,W->n);
1688: 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);
1689: if (!V->n) return(0);
1691: PetscLogEventBegin(BV_Copy,V,W,0,0);
1692: if (V->indef && V->matrix && V->indef==W->indef && V->matrix==W->matrix) {
1693: /* copy signature */
1694: BV_AllocateSignature(W);
1695: VecGetArrayRead(V->omega,&vomega);
1696: VecGetArray(W->omega,&womega);
1697: PetscArraycpy(womega+W->nc+W->l,vomega+V->nc+V->l,V->k-V->l);
1698: VecRestoreArray(W->omega,&womega);
1699: VecRestoreArrayRead(V->omega,&vomega);
1700: }
1701: (*V->ops->copy)(V,W);
1702: PetscLogEventEnd(BV_Copy,V,W,0,0);
1703: PetscObjectStateIncrease((PetscObject)W);
1704: return(0);
1705: }
1707: /*@
1708: BVCopyVec - Copies one of the columns of a basis vectors object into a Vec.
1710: Logically Collective on V
1712: Input Parameter:
1713: + V - basis vectors context
1714: - j - the column number to be copied
1716: Output Parameter:
1717: . w - the copied column
1719: Note:
1720: Both V and w must be distributed in the same manner; local copies are done.
1722: Level: beginner
1724: .seealso: BVCopy(), BVCopyColumn(), BVInsertVec()
1725: @*/
1726: PetscErrorCode BVCopyVec(BV V,PetscInt j,Vec w)1727: {
1729: PetscInt n,N;
1730: Vec z;
1735: BVCheckSizes(V,1);
1740: VecGetSize(w,&N);
1741: VecGetLocalSize(w,&n);
1742: 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);
1744: PetscLogEventBegin(BV_Copy,V,w,0,0);
1745: BVGetColumn(V,j,&z);
1746: VecCopy(z,w);
1747: BVRestoreColumn(V,j,&z);
1748: PetscLogEventEnd(BV_Copy,V,w,0,0);
1749: return(0);
1750: }
1752: /*@
1753: BVCopyColumn - Copies the values from one of the columns to another one.
1755: Logically Collective on V
1757: Input Parameter:
1758: + V - basis vectors context
1759: . j - the number of the source column
1760: - i - the number of the destination column
1762: Level: beginner
1764: .seealso: BVCopy(), BVCopyVec()
1765: @*/
1766: PetscErrorCode BVCopyColumn(BV V,PetscInt j,PetscInt i)1767: {
1769: Vec z,w;
1770: PetscScalar *omega;
1775: BVCheckSizes(V,1);
1778: if (j==i) return(0);
1780: PetscLogEventBegin(BV_Copy,V,0,0,0);
1781: if (V->omega) {
1782: VecGetArray(V->omega,&omega);
1783: omega[i] = omega[j];
1784: VecRestoreArray(V->omega,&omega);
1785: }
1786: if (V->ops->copycolumn) {
1787: (*V->ops->copycolumn)(V,j,i);
1788: } else {
1789: BVGetColumn(V,j,&z);
1790: BVGetColumn(V,i,&w);
1791: VecCopy(z,w);
1792: BVRestoreColumn(V,j,&z);
1793: BVRestoreColumn(V,i,&w);
1794: }
1795: PetscLogEventEnd(BV_Copy,V,0,0,0);
1796: PetscObjectStateIncrease((PetscObject)V);
1797: return(0);
1798: }
1800: static PetscErrorCode BVGetSplit_Private(BV bv,PetscBool left,BV *split)1801: {
1803: PetscInt ncols;
1806: ncols = left? bv->nc+bv->l: bv->m-bv->l;
1807: if (!*split) {
1808: BVCreate(PetscObjectComm((PetscObject)bv),split);
1809: PetscLogObjectParent((PetscObject)bv,(PetscObject)*split);
1810: (*split)->issplit = left? 1: 2;
1811: (*split)->splitparent = bv;
1812: BVSetSizesFromVec(*split,bv->t,ncols);
1813: BVDuplicate_Private(bv,*split);
1814: } else {
1815: BVResize(*split,ncols,PETSC_FALSE);
1816: }
1817: (*split)->l = 0;
1818: (*split)->k = left? bv->l: bv->k-bv->l;
1819: (*split)->nc = left? bv->nc: 0;
1820: (*split)->m = ncols-(*split)->nc;
1821: if ((*split)->nc) {
1822: (*split)->ci[0] = -(*split)->nc-1;
1823: (*split)->ci[1] = -(*split)->nc-1;
1824: }
1825: if (left) {
1826: PetscObjectStateGet((PetscObject)*split,&bv->lstate);
1827: } else {
1828: PetscObjectStateGet((PetscObject)*split,&bv->rstate);
1829: }
1830: return(0);
1831: }
1833: /*@
1834: BVGetSplit - Splits the BV object into two BV objects that share the
1835: internal data, one of them containing the leading columns and the other
1836: one containing the remaining columns.
1838: Logically Collective on bv
1840: Input Parameters:
1841: . bv - the basis vectors context
1843: Output Parameters:
1844: + L - left BV containing leading columns (can be NULL)
1845: - R - right BV containing remaining columns (can be NULL)
1847: Notes:
1848: The columns are split in two sets. The leading columns (including the
1849: constraints) are assigned to the left BV and the remaining columns
1850: are assigned to the right BV. The number of leading columns, as
1851: specified with BVSetActiveColumns(), must be between 1 and m-1 (to
1852: guarantee that both L and R have at least one column).
1854: The returned BV's must be seen as references (not copies) of the input
1855: BV, that is, modifying them will change the entries of bv as well.
1856: The returned BV's must not be destroyed. BVRestoreSplit() must be called
1857: when they are no longer needed.
1859: Pass NULL for any of the output BV's that is not needed.
1861: Level: advanced
1863: .seealso: BVRestoreSplit(), BVSetActiveColumns(), BVSetNumConstraints()
1864: @*/
1865: PetscErrorCode BVGetSplit(BV bv,BV *L,BV *R)1866: {
1872: BVCheckSizes(bv,1);
1873: if (!bv->l) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"Must indicate the number of leading columns with BVSetActiveColumns()");
1874: if (bv->lsplit) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"Cannot get the split BV's twice before restoring them with BVRestoreSplit()");
1875: bv->lsplit = bv->nc+bv->l;
1876: BVGetSplit_Private(bv,PETSC_TRUE,&bv->L);
1877: BVGetSplit_Private(bv,PETSC_FALSE,&bv->R);
1878: if (L) *L = bv->L;
1879: if (R) *R = bv->R;
1880: return(0);
1881: }
1883: /*@
1884: BVRestoreSplit - Restore the BV objects obtained with BVGetSplit().
1886: Logically Collective on bv
1888: Input Parameters:
1889: + bv - the basis vectors context
1890: . L - left BV obtained with BVGetSplit()
1891: - R - right BV obtained with BVGetSplit()
1893: Note:
1894: The arguments must match the corresponding call to BVGetSplit().
1896: Level: advanced
1898: .seealso: BVGetSplit()
1899: @*/
1900: PetscErrorCode BVRestoreSplit(BV bv,BV *L,BV *R)1901: {
1907: BVCheckSizes(bv,1);
1908: if (!bv->lsplit) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"Must call BVGetSplit first");
1909: 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");
1910: 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");
1911: 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()");
1912: 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()");
1914: if (bv->ops->restoresplit) {
1915: (*bv->ops->restoresplit)(bv,L,R);
1916: }
1917: bv->lsplit = 0;
1918: if (L) *L = NULL;
1919: if (R) *R = NULL;
1920: return(0);
1921: }