1: /*
2: Basic BV routines.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc/private/bvimpl.h> /*I "slepcbv.h" I*/
26: PetscBool BVRegisterAllCalled = PETSC_FALSE;
27: PetscFunctionList BVList = 0;
31: /*@C
32: BVSetType - Selects the type for the BV object.
34: Logically Collective on BV 36: Input Parameter:
37: + bv - the basis vectors context
38: - type - a known type
40: Options Database Key:
41: . -bv_type <type> - Sets BV type
43: Level: intermediate
45: .seealso: BVGetType()
46: @*/
47: PetscErrorCode BVSetType(BV bv,BVType type) 48: {
49: PetscErrorCode ierr,(*r)(BV);
50: PetscBool match;
56: PetscObjectTypeCompare((PetscObject)bv,type,&match);
57: if (match) return(0);
59: PetscFunctionListFind(BVList,type,&r);
60: if (!r) SETERRQ1(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested BV type %s",type);
62: if (bv->ops->destroy) { (*bv->ops->destroy)(bv); }
63: PetscMemzero(bv->ops,sizeof(struct _BVOps));
65: PetscObjectChangeTypeName((PetscObject)bv,type);
66: if (bv->n < 0 && bv->N < 0) {
67: bv->ops->create = r;
68: } else {
69: PetscLogEventBegin(BV_Create,bv,0,0,0);
70: (*r)(bv);
71: PetscLogEventEnd(BV_Create,bv,0,0,0);
72: }
73: return(0);
74: }
78: /*@C
79: BVGetType - Gets the BV type name (as a string) from the BV context.
81: Not Collective
83: Input Parameter:
84: . bv - the basis vectors context
86: Output Parameter:
87: . name - name of the type of basis vectors
89: Level: intermediate
91: .seealso: BVSetType()
92: @*/
93: PetscErrorCode BVGetType(BV bv,BVType *type) 94: {
98: *type = ((PetscObject)bv)->type_name;
99: return(0);
100: }
104: /*@
105: BVSetSizes - Sets the local and global sizes, and the number of columns.
107: Collective on BV109: Input Parameters:
110: + bv - the basis vectors
111: . n - the local size (or PETSC_DECIDE to have it set)
112: . N - the global size (or PETSC_DECIDE)
113: - m - the number of columns
115: Notes:
116: n and N cannot be both PETSC_DECIDE.
117: If one processor calls this with N of PETSC_DECIDE then all processors must,
118: otherwise the program will hang.
120: Level: beginner
122: .seealso: BVSetSizesFromVec(), BVGetSizes(), BVResize()
123: @*/
124: PetscErrorCode BVSetSizes(BV bv,PetscInt n,PetscInt N,PetscInt m)125: {
127: PetscInt ma;
133: if (N >= 0 && n > N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local size %D cannot be larger than global size %D",n,N);
134: if (m <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of columns %D must be positive",m);
135: if ((bv->n >= 0 || bv->N >= 0) && (bv->n != n || bv->N != N)) SETERRQ4(PETSC_COMM_SELF,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);
136: if (bv->m > 0 && bv->m != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change the number of columns to %D after previously setting it to %D; use BVResize()",m,bv->m);
137: bv->n = n;
138: bv->N = N;
139: bv->m = m;
140: bv->k = m;
141: if (!bv->t) { /* create template vector and get actual dimensions */
142: VecCreate(PetscObjectComm((PetscObject)bv),&bv->t);
143: VecSetSizes(bv->t,bv->n,bv->N);
144: VecSetFromOptions(bv->t);
145: VecGetSize(bv->t,&bv->N);
146: VecGetLocalSize(bv->t,&bv->n);
147: if (bv->matrix) { /* check compatible dimensions of user-provided matrix */
148: MatGetLocalSize(bv->matrix,&ma,NULL);
149: if (bv->n!=ma) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local dimension %D does not match that of matrix given at BVSetMatrix %D",bv->n,ma);
150: }
151: }
152: if (bv->ops->create) {
153: PetscLogEventBegin(BV_Create,bv,0,0,0);
154: (*bv->ops->create)(bv);
155: PetscLogEventEnd(BV_Create,bv,0,0,0);
156: bv->ops->create = 0;
157: bv->defersfo = PETSC_FALSE;
158: }
159: return(0);
160: }
164: /*@
165: BVSetSizesFromVec - Sets the local and global sizes, and the number of columns.
166: Local and global sizes are specified indirectly by passing a template vector.
168: Collective on BV170: Input Parameters:
171: + bv - the basis vectors
172: . t - the template vectors
173: - m - the number of columns
175: Level: beginner
177: .seealso: BVSetSizes(), BVGetSizes(), BVResize()
178: @*/
179: PetscErrorCode BVSetSizesFromVec(BV bv,Vec t,PetscInt m)180: {
182: PetscInt ma;
189: if (m <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of columns %D must be positive",m);
190: if (bv->t) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Template vector was already set by a previous call to BVSetSizes/FromVec");
191: VecGetSize(t,&bv->N);
192: VecGetLocalSize(t,&bv->n);
193: if (bv->matrix) { /* check compatible dimensions of user-provided matrix */
194: MatGetLocalSize(bv->matrix,&ma,NULL);
195: if (bv->n!=ma) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local dimension %D does not match that of matrix given at BVSetMatrix %D",bv->n,ma);
196: }
197: bv->m = m;
198: bv->k = m;
199: bv->t = t;
200: PetscObjectReference((PetscObject)t);
201: if (bv->ops->create) {
202: (*bv->ops->create)(bv);
203: bv->ops->create = 0;
204: }
205: return(0);
206: }
210: /*@
211: BVGetSizes - Returns the local and global sizes, and the number of columns.
213: Not Collective
215: Input Parameter:
216: . bv - the basis vectors
218: Output Parameters:
219: + n - the local size
220: . N - the global size
221: - m - the number of columns
223: Note:
224: Normal usage requires that bv has already been given its sizes, otherwise
225: the call fails. However, this function can also be used to determine if
226: a BV object has been initialized completely (sizes and type). For this,
227: call with n=NULL and N=NULL, then a return value of m=0 indicates that
228: the BV object is not ready for use yet.
230: Level: beginner
232: .seealso: BVSetSizes(), BVSetSizesFromVec()
233: @*/
234: PetscErrorCode BVGetSizes(BV bv,PetscInt *n,PetscInt *N,PetscInt *m)235: {
237: if (!bv) {
238: if (m && !n && !N) *m = 0;
239: return(0);
240: }
242: if (n || N) BVCheckSizes(bv,1);
243: if (n) *n = bv->n;
244: if (N) *N = bv->N;
245: if (m) *m = bv->m;
246: if (m && !n && !N && !((PetscObject)bv)->type_name) *m = 0;
247: return(0);
248: }
252: /*@
253: BVSetNumConstraints - Set the number of constraints.
255: Logically Collective on BV257: Input Parameters:
258: + V - basis vectors
259: - nc - number of constraints
261: Notes:
262: This function sets the number of constraints to nc and marks all remaining
263: columns as regular. Normal user would call BVInsertConstraints() instead.
265: If nc is smaller than the previously set value, then some of the constraints
266: are discarded. In particular, using nc=0 removes all constraints preserving
267: the content of regular columns.
269: Level: developer
271: .seealso: BVInsertConstraints()
272: @*/
273: PetscErrorCode BVSetNumConstraints(BV V,PetscInt nc)274: {
276: PetscInt total,diff,i;
277: Vec x,y;
282: if (nc<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of constraints (given %D) cannot be negative",nc);
284: BVCheckSizes(V,1);
285: 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");
287: diff = nc-V->nc;
288: if (!diff) return(0);
289: total = V->nc+V->m;
290: if (total-nc<=0) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Not enough columns for the given nc value");
291: if (diff<0) { /* lessen constraints, shift contents of BV */
292: for (i=0;i<V->m;i++) {
293: BVGetColumn(V,i,&x);
294: BVGetColumn(V,i+diff,&y);
295: VecCopy(x,y);
296: BVRestoreColumn(V,i,&x);
297: BVRestoreColumn(V,i+diff,&y);
298: }
299: }
300: V->nc = nc;
301: V->ci[0] = -V->nc-1;
302: V->ci[1] = -V->nc-1;
303: V->l = 0;
304: V->m = total-nc;
305: V->k = V->m;
306: PetscObjectStateIncrease((PetscObject)V);
307: return(0);
308: }
312: /*@
313: BVGetNumConstraints - Returns the number of constraints.
315: Not Collective
317: Input Parameter:
318: . bv - the basis vectors
320: Output Parameters:
321: . nc - the number of constraints
323: Level: advanced
325: .seealso: BVGetSizes(), BVInsertConstraints()
326: @*/
327: PetscErrorCode BVGetNumConstraints(BV bv,PetscInt *nc)328: {
332: *nc = bv->nc;
333: return(0);
334: }
338: /*@
339: BVResize - Change the number of columns.
341: Collective on BV343: Input Parameters:
344: + bv - the basis vectors
345: . m - the new number of columns
346: - copy - a flag indicating whether current values should be kept
348: Note:
349: Internal storage is reallocated. If the copy flag is set to true, then
350: the contents are copied to the leading part of the new space.
352: Level: advanced
354: .seealso: BVSetSizes(), BVSetSizesFromVec()
355: @*/
356: PetscErrorCode BVResize(BV bv,PetscInt m,PetscBool copy)357: {
359: PetscReal *omega;
360: PetscInt i;
367: if (m <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of columns %D must be positive",m);
368: if (bv->nc) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot resize a BV with constraints");
369: if (bv->m == m) return(0);
371: PetscLogEventBegin(BV_Create,bv,0,0,0);
372: (*bv->ops->resize)(bv,m,copy);
373: PetscFree2(bv->h,bv->c);
374: if (bv->omega) {
375: PetscMalloc1(m,&omega);
376: PetscLogObjectMemory((PetscObject)bv,m*sizeof(PetscReal));
377: for (i=0;i<m;i++) omega[i] = 1.0;
378: if (copy) {
379: PetscMemcpy(omega,bv->omega,PetscMin(m,bv->m)*sizeof(PetscReal));
380: }
381: PetscFree(bv->omega);
382: bv->omega = omega;
383: }
384: bv->m = m;
385: bv->k = PetscMin(bv->k,m);
386: bv->l = PetscMin(bv->l,m);
387: PetscLogEventEnd(BV_Create,bv,0,0,0);
388: return(0);
389: }
393: /*@
394: BVSetActiveColumns - Specify the columns that will be involved in operations.
396: Logically Collective on BV398: Input Parameters:
399: + bv - the basis vectors context
400: . l - number of leading columns
401: - k - number of active columns
403: Notes:
404: In operations such as BVMult() or BVDot(), only the first k columns are
405: considered. This is useful when the BV is filled from left to right, so
406: the last m-k columns do not have relevant information.
408: Also in operations such as BVMult() or BVDot(), the first l columns are
409: normally not included in the computation. See the manpage of each
410: operation.
412: In orthogonalization operations, the first l columns are treated
413: differently: they participate in the orthogonalization but the computed
414: coefficients are not stored.
416: Level: intermediate
418: .seealso: BVGetActiveColumns(), BVSetSizes()
419: @*/
420: PetscErrorCode BVSetActiveColumns(BV bv,PetscInt l,PetscInt k)421: {
426: BVCheckSizes(bv,1);
427: if (k==PETSC_DECIDE || k==PETSC_DEFAULT) {
428: bv->k = bv->m;
429: } else {
430: if (k<0 || k>bv->m) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of k. Must be between 0 and m");
431: bv->k = k;
432: }
433: if (l==PETSC_DECIDE || l==PETSC_DEFAULT) {
434: bv->l = 0;
435: } else {
436: if (l<0 || l>bv->k) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of l. Must be between 0 and k");
437: bv->l = l;
438: }
439: return(0);
440: }
444: /*@
445: BVGetActiveColumns - Returns the current active dimensions.
447: Not Collective
449: Input Parameter:
450: . bv - the basis vectors context
452: Output Parameter:
453: + l - number of leading columns
454: - k - number of active columns
456: Level: intermediate
458: .seealso: BVSetActiveColumns()
459: @*/
460: PetscErrorCode BVGetActiveColumns(BV bv,PetscInt *l,PetscInt *k)461: {
464: if (l) *l = bv->l;
465: if (k) *k = bv->k;
466: return(0);
467: }
471: /*@
472: BVSetMatrix - Specifies the inner product to be used in orthogonalization.
474: Collective on BV476: Input Parameters:
477: + bv - the basis vectors context
478: . B - a symmetric matrix (may be NULL)
479: - indef - a flag indicating if the matrix is indefinite
481: Notes:
482: This is used to specify a non-standard inner product, whose matrix
483: representation is given by B. Then, all inner products required during
484: orthogonalization are computed as (x,y)_B=y^H*B*x rather than the
485: standard form (x,y)=y^H*x.
487: Matrix B must be real symmetric (or complex Hermitian). A genuine inner
488: product requires that B is also positive (semi-)definite. However, we
489: also allow for an indefinite B (setting indef=PETSC_TRUE), in which
490: case the orthogonalization uses an indefinite inner product.
492: This affects operations BVDot(), BVNorm(), BVOrthogonalize(), and variants.
494: Setting B=NULL has the same effect as if the identity matrix was passed.
496: Level: advanced
498: .seealso: BVGetMatrix(), BVDot(), BVNorm(), BVOrthogonalize()
499: @*/
500: PetscErrorCode BVSetMatrix(BV bv,Mat B,PetscBool indef)501: {
503: PetscInt m,n;
508: if (B) {
510: MatGetLocalSize(B,&m,&n);
511: if (m!=n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix must be square");
512: if (bv->m && bv->n!=n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension BV %D, Mat %D",bv->n,n);
513: }
514: MatDestroy(&bv->matrix);
515: if (B) PetscObjectReference((PetscObject)B);
516: bv->matrix = B;
517: bv->indef = indef;
518: if (B && !bv->Bx) {
519: MatCreateVecs(B,&bv->Bx,NULL);
520: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->Bx);
521: }
522: return(0);
523: }
527: /*@
528: BVGetMatrix - Retrieves the matrix representation of the inner product.
530: Not collective, though a parallel Mat may be returned
532: Input Parameter:
533: . bv - the basis vectors context
535: Output Parameter:
536: + B - the matrix of the inner product (may be NULL)
537: - indef - the flag indicating if the matrix is indefinite
539: Level: advanced
541: .seealso: BVSetMatrix()
542: @*/
543: PetscErrorCode BVGetMatrix(BV bv,Mat *B,PetscBool *indef)544: {
547: if (B) *B = bv->matrix;
548: if (indef) *indef = bv->indef;
549: return(0);
550: }
554: /*@
555: BVApplyMatrix - Multiplies a vector by the matrix representation of the
556: inner product.
558: Neighbor-wise Collective on BV and Vec
560: Input Parameter:
561: + bv - the basis vectors context
562: - x - the vector
564: Output Parameter:
565: . y - the result
567: Note:
568: If no matrix was specified this function copies the vector.
570: Level: advanced
572: .seealso: BVSetMatrix(), BVApplyMatrixBV()
573: @*/
574: PetscErrorCode BVApplyMatrix(BV bv,Vec x,Vec y)575: {
582: if (bv->matrix) {
583: BV_IPMatMult(bv,x);
584: VecCopy(bv->Bx,y);
585: } else {
586: VecCopy(x,y);
587: }
588: return(0);
589: }
593: /*@
594: BVApplyMatrixBV - Multiplies the BV vectors by the matrix representation
595: of the inner product.
597: Neighbor-wise Collective on BV599: Input Parameter:
600: + X - the basis vectors context
602: Output Parameter:
603: . Y - the basis vectors to store the result (optional)
605: Note:
606: This function computes Y = B*X, where B is the matrix given with
607: BVSetMatrix(). This operation is computed as in BVMatMult().
608: If no matrix was specified, then it just copies Y = X.
610: If no Y is given, the result is stored internally in the cached BV.
612: Level: developer
614: .seealso: BVSetMatrix(), BVApplyMatrix(), BVMatMult(), BVGetCachedBV()
615: @*/
616: PetscErrorCode BVApplyMatrixBV(BV X,BV Y)617: {
622: if (Y) {
624: if (X->matrix) {
625: BVMatMult(X,X->matrix,Y);
626: } else {
627: BVCopy(X,Y);
628: }
629: } else {
630: BV_IPMatMultBV(X);
631: }
632: return(0);
633: }
637: /*@
638: BVGetCachedBV - Returns a BV object stored internally that holds the
639: result of B*X after a call to BVApplyMatrixBV().
641: Not collective
643: Input Parameter:
644: . bv - the basis vectors context
646: Output Parameter:
647: . cached - the cached BV649: Note:
650: The function will return a NULL if BVApplyMatrixBV() was not called yet.
652: Level: developer
654: .seealso: BVApplyMatrixBV()
655: @*/
656: PetscErrorCode BVGetCachedBV(BV bv,BV *cached)657: {
661: *cached = bv->cached;
662: return(0);
663: }
667: /*@
668: BVSetSignature - Sets the signature matrix to be used in orthogonalization.
670: Logically Collective on BV672: Input Parameter:
673: + bv - the basis vectors context
674: - omega - a vector representing the diagonal of the signature matrix
676: Note:
677: The signature matrix Omega = V'*B*V is relevant only for an indefinite B.
679: Level: developer
681: .seealso: BVSetMatrix(), BVGetSignature()
682: @*/
683: PetscErrorCode BVSetSignature(BV bv,Vec omega)684: {
685: PetscErrorCode ierr;
686: PetscInt i,n;
687: const PetscScalar *pomega;
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: BV_AllocateSignature(bv);
698: if (bv->indef) {
699: VecGetArrayRead(omega,&pomega);
700: for (i=0;i<n;i++) bv->omega[bv->nc+i] = PetscRealPart(pomega[i]);
701: VecRestoreArrayRead(omega,&pomega);
702: } else {
703: PetscInfo(bv,"Ignoring signature because BV is not indefinite\n");
704: }
705: return(0);
706: }
710: /*@
711: BVGetSignature - Retrieves the signature matrix from last orthogonalization.
713: Not collective
715: Input Parameter:
716: . bv - the basis vectors context
718: Output Parameter:
719: . omega - a vector representing the diagonal of the signature matrix
721: Note:
722: The signature matrix Omega = V'*B*V is relevant only for an indefinite B.
724: Level: developer
726: .seealso: BVSetMatrix(), BVSetSignature()
727: @*/
728: PetscErrorCode BVGetSignature(BV bv,Vec omega)729: {
731: PetscInt i,n;
732: PetscScalar *pomega;
736: BVCheckSizes(bv,1);
740: VecGetSize(omega,&n);
741: if (n!=bv->k) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Vec argument has %D elements, should be %D",n,bv->k);
742: if (bv->indef && bv->omega) {
743: VecGetArray(omega,&pomega);
744: for (i=0;i<n;i++) pomega[i] = bv->omega[bv->nc+i];
745: VecRestoreArray(omega,&pomega);
746: } else {
747: VecSet(omega,1.0);
748: }
749: return(0);
750: }
754: /*@
755: BVSetFromOptions - Sets BV options from the options database.
757: Collective on BV759: Input Parameter:
760: . bv - the basis vectors context
762: Level: beginner
763: @*/
764: PetscErrorCode BVSetFromOptions(BV bv)765: {
767: char type[256];
768: PetscBool flg;
769: PetscReal r;
770: PetscInt i,j,k;
771: const char *orth_list[2] = {"cgs","mgs"};
772: const char *ref_list[3] = {"ifneeded","never","always"};
773: const char *borth_list[2] = {"gs","chol"};
777: BVRegisterAll();
778: PetscObjectOptionsBegin((PetscObject)bv);
779: PetscOptionsFList("-bv_type","Basis Vectors type","BVSetType",BVList,(char*)(((PetscObject)bv)->type_name?((PetscObject)bv)->type_name:BVSVEC),type,256,&flg);
780: if (flg) {
781: BVSetType(bv,type);
782: }
783: /*
784: Set the type if it was never set.
785: */
786: if (!((PetscObject)bv)->type_name) {
787: BVSetType(bv,BVSVEC);
788: }
790: i = bv->orthog_type;
791: PetscOptionsEList("-bv_orthog_type","Orthogonalization method","BVSetOrthogonalization",orth_list,2,orth_list[i],&i,NULL);
792: j = bv->orthog_ref;
793: PetscOptionsEList("-bv_orthog_refine","Iterative refinement mode during orthogonalization","BVSetOrthogonalization",ref_list,3,ref_list[j],&j,NULL);
794: r = bv->orthog_eta;
795: PetscOptionsReal("-bv_orthog_eta","Parameter of iterative refinement during orthogonalization","BVSetOrthogonalization",r,&r,NULL);
796: k = bv->orthog_block;
797: PetscOptionsEList("-bv_orthog_block","Block orthogonalization method","BVSetOrthogonalization",borth_list,2,borth_list[k],&k,NULL);
798: BVSetOrthogonalization(bv,(BVOrthogType)i,(BVOrthogRefineType)j,r,(BVOrthogBlockType)k);
800: PetscOptionsBoolGroupBegin("-bv_matmult_vecs","Do matmult as matrix-vector products","BVSetMatMultMethod",&flg);
801: if (flg) { BVSetMatMultMethod(bv,BV_MATMULT_VECS); }
802: PetscOptionsBoolGroup("-bv_matmult_mat","Do matmult as a single matrix-matrix product","BVSetMatMultMethod",&flg);
803: if (flg) { BVSetMatMultMethod(bv,BV_MATMULT_MAT); }
804: PetscOptionsBoolGroupEnd("-bv_matmult_mat_save","Do matmult as a single matrix-matrix product and save auxiliary matrices","BVSetMatMultMethod",&flg);
805: if (flg) { BVSetMatMultMethod(bv,BV_MATMULT_MAT_SAVE); }
807: if (bv->ops->create) bv->defersfo = PETSC_TRUE; /* defer call to setfromoptions */
808: else if (bv->ops->setfromoptions) {
809: (*bv->ops->setfromoptions)(PetscOptionsObject,bv);
810: }
811: PetscObjectProcessOptionsHandlers((PetscObject)bv);
812: PetscOptionsEnd();
813: return(0);
814: }
818: /*@
819: BVSetOrthogonalization - Specifies the method used for the orthogonalization of
820: vectors (classical or modified Gram-Schmidt with or without refinement), and
821: for the block-orthogonalization (simultaneous orthogonalization of a set of
822: vectors).
824: Logically Collective on BV826: Input Parameters:
827: + bv - the basis vectors context
828: . type - the method of vector orthogonalization
829: . refine - type of refinement
830: . eta - parameter for selective refinement
831: - block - the method of block orthogonalization
833: Options Database Keys:
834: + -bv_orthog_type <type> - Where <type> is cgs for Classical Gram-Schmidt orthogonalization
835: (default) or mgs for Modified Gram-Schmidt orthogonalization
836: . -bv_orthog_refine <ref> - Where <ref> is one of never, ifneeded (default) or always
837: . -bv_orthog_eta <eta> - For setting the value of eta
838: - -bv_orthog_block <block> - Where <block> is the block-orthogonalization method
840: Notes:
841: The default settings work well for most problems.
843: The parameter eta should be a real value between 0 and 1 (or PETSC_DEFAULT).
844: The value of eta is used only when the refinement type is "ifneeded".
846: When using several processors, MGS is likely to result in bad scalability.
848: If the method set for block orthogonalization is GS, then the computation
849: is done column by column with the vector orthogonalization.
851: Level: advanced
853: .seealso: BVOrthogonalizeColumn(), BVGetOrthogonalization(), BVOrthogType, BVOrthogRefineType, BVOrthogBlockType854: @*/
855: PetscErrorCode BVSetOrthogonalization(BV bv,BVOrthogType type,BVOrthogRefineType refine,PetscReal eta,BVOrthogBlockType block)856: {
863: switch (type) {
864: case BV_ORTHOG_CGS:
865: case BV_ORTHOG_MGS:
866: bv->orthog_type = type;
867: break;
868: default:869: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
870: }
871: switch (refine) {
872: case BV_ORTHOG_REFINE_NEVER:
873: case BV_ORTHOG_REFINE_IFNEEDED:
874: case BV_ORTHOG_REFINE_ALWAYS:
875: bv->orthog_ref = refine;
876: break;
877: default:878: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown refinement type");
879: }
880: if (eta == PETSC_DEFAULT) {
881: bv->orthog_eta = 0.7071;
882: } else {
883: if (eta <= 0.0 || eta > 1.0) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Invalid eta value");
884: bv->orthog_eta = eta;
885: }
886: switch (block) {
887: case BV_ORTHOG_BLOCK_GS:
888: case BV_ORTHOG_BLOCK_CHOL:
889: bv->orthog_block = block;
890: break;
891: default:892: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown block orthogonalization type");
893: }
894: return(0);
895: }
899: /*@
900: BVGetOrthogonalization - Gets the orthogonalization settings from the BV object.
902: Not Collective
904: Input Parameter:
905: . bv - basis vectors context
907: Output Parameter:
908: + type - the method of vector orthogonalization
909: . refine - type of refinement
910: . eta - parameter for selective refinement
911: - block - the method of block orthogonalization
913: Level: advanced
915: .seealso: BVOrthogonalizeColumn(), BVSetOrthogonalization(), BVOrthogType, BVOrthogRefineType, BVOrthogBlockType916: @*/
917: PetscErrorCode BVGetOrthogonalization(BV bv,BVOrthogType *type,BVOrthogRefineType *refine,PetscReal *eta,BVOrthogBlockType *block)918: {
921: if (type) *type = bv->orthog_type;
922: if (refine) *refine = bv->orthog_ref;
923: if (eta) *eta = bv->orthog_eta;
924: if (block) *block = bv->orthog_block;
925: return(0);
926: }
930: /*@
931: BVSetMatMultMethod - Specifies the method used for the BVMatMult() operation.
933: Logically Collective on BV935: Input Parameters:
936: + bv - the basis vectors context
937: - method - the method for the BVMatMult() operation
939: Options Database Keys:
940: + -bv_matmult_vecs - perform a matrix-vector multiply per each column
941: . -bv_matmult_mat - carry out a MatMatMult() product with a dense matrix
942: - -bv_matmult_mat_save - call MatMatMult() and keep auxiliary matrices
944: Note:
945: The default is BV_MATMULT_MAT.
947: Level: advanced
949: .seealso: BVGetMatMultMethod(), BVMatMultType950: @*/
951: PetscErrorCode BVSetMatMultMethod(BV bv,BVMatMultType method)952: {
956: switch (method) {
957: case BV_MATMULT_VECS:
958: case BV_MATMULT_MAT:
959: case BV_MATMULT_MAT_SAVE:
960: bv->vmm = method;
961: break;
962: default:963: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown matmult method");
964: }
965: return(0);
966: }
970: /*@
971: BVGetMatMultMethod - Gets the method used for the BVMatMult() operation.
973: Not Collective
975: Input Parameter:
976: . bv - basis vectors context
978: Output Parameter:
979: . method - the method for the BVMatMult() operation
981: Level: advanced
983: .seealso: BVSetMatMultMethod(), BVMatMultType984: @*/
985: PetscErrorCode BVGetMatMultMethod(BV bv,BVMatMultType *method)986: {
990: *method = bv->vmm;
991: return(0);
992: }
996: /*@
997: BVGetColumn - Returns a Vec object that contains the entries of the
998: requested column of the basis vectors object.
1000: Logically Collective on BV1002: Input Parameters:
1003: + bv - the basis vectors context
1004: - j - the index of the requested column
1006: Output Parameter:
1007: . v - vector containing the jth column
1009: Notes:
1010: The returned Vec must be seen as a reference (not a copy) of the BV1011: column, that is, modifying the Vec will change the BV entries as well.
1013: The returned Vec must not be destroyed. BVRestoreColumn() must be
1014: called when it is no longer needed. At most, two columns can be fetched,
1015: that is, this function can only be called twice before the corresponding
1016: BVRestoreColumn() is invoked.
1018: A negative index j selects the i-th constraint, where i=-j. Constraints
1019: should not be modified.
1021: Level: beginner
1023: .seealso: BVRestoreColumn(), BVInsertConstraints()
1024: @*/
1025: PetscErrorCode BVGetColumn(BV bv,PetscInt j,Vec *v)1026: {
1028: PetscInt l;
1033: BVCheckSizes(bv,1);
1035: 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);
1036: if (j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"You requested column %D but only %D are available",j,bv->m);
1037: 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);
1038: l = BVAvailableVec;
1039: 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");
1040: (*bv->ops->getcolumn)(bv,j,v);
1041: bv->ci[l] = j;
1042: PetscObjectStateGet((PetscObject)bv->cv[l],&bv->st[l]);
1043: PetscObjectGetId((PetscObject)bv->cv[l],&bv->id[l]);
1044: *v = bv->cv[l];
1045: return(0);
1046: }
1050: /*@
1051: BVRestoreColumn - Restore a column obtained with BVGetColumn().
1053: Logically Collective on BV1055: Input Parameters:
1056: + bv - the basis vectors context
1057: . j - the index of the column
1058: - v - vector obtained with BVGetColumn()
1060: Note:
1061: The arguments must match the corresponding call to BVGetColumn().
1063: Level: beginner
1065: .seealso: BVGetColumn()
1066: @*/
1067: PetscErrorCode BVRestoreColumn(BV bv,PetscInt j,Vec *v)1068: {
1069: PetscErrorCode ierr;
1070: PetscObjectId id;
1071: PetscObjectState st;
1072: PetscInt l;
1077: BVCheckSizes(bv,1);
1081: 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);
1082: if (j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"You requested column %D but only %D are available",j,bv->m);
1083: 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);
1084: l = (j==bv->ci[0])? 0: 1;
1085: PetscObjectGetId((PetscObject)*v,&id);
1086: 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");
1087: PetscObjectStateGet((PetscObject)*v,&st);
1088: if (st!=bv->st[l]) {
1089: PetscObjectStateIncrease((PetscObject)bv);
1090: }
1091: if (bv->ops->restorecolumn) {
1092: (*bv->ops->restorecolumn)(bv,j,v);
1093: } else bv->cv[l] = NULL;
1094: bv->ci[l] = -bv->nc-1;
1095: bv->st[l] = -1;
1096: bv->id[l] = 0;
1097: *v = NULL;
1098: return(0);
1099: }
1103: /*@C
1104: BVGetArray - Returns a pointer to a contiguous array that contains this
1105: processor's portion of the BV data.
1107: Logically Collective on BV1109: Input Parameters:
1110: . bv - the basis vectors context
1112: Output Parameter:
1113: . a - location to put pointer to the array
1115: Notes:
1116: BVRestoreArray() must be called when access to the array is no longer needed.
1117: This operation may imply a data copy, for BV types that do not store
1118: data contiguously in memory.
1120: The pointer will normally point to the first entry of the first column,
1121: but if the BV has constraints then these go before the regular columns.
1123: Level: advanced
1125: .seealso: BVRestoreArray(), BVInsertConstraints()
1126: @*/
1127: PetscErrorCode BVGetArray(BV bv,PetscScalar **a)1128: {
1134: BVCheckSizes(bv,1);
1135: (*bv->ops->getarray)(bv,a);
1136: return(0);
1137: }
1141: /*@C
1142: BVRestoreArray - Restore the BV object after BVGetArray() has been called.
1144: Logically Collective on BV1146: Input Parameters:
1147: + bv - the basis vectors context
1148: - a - location of pointer to array obtained from BVGetArray()
1150: Note:
1151: This operation may imply a data copy, for BV types that do not store
1152: data contiguously in memory.
1154: Level: advanced
1156: .seealso: BVGetColumn()
1157: @*/
1158: PetscErrorCode BVRestoreArray(BV bv,PetscScalar **a)1159: {
1165: BVCheckSizes(bv,1);
1166: if (bv->ops->restorearray) {
1167: (*bv->ops->restorearray)(bv,a);
1168: }
1169: if (a) *a = NULL;
1170: PetscObjectStateIncrease((PetscObject)bv);
1171: return(0);
1172: }
1176: /*@
1177: BVCreateVec - Creates a new Vec object with the same type and dimensions
1178: as the columns of the basis vectors object.
1180: Collective on BV1182: Input Parameter:
1183: . bv - the basis vectors context
1185: Output Parameter:
1186: . v - the new vector
1188: Note:
1189: The user is responsible of destroying the returned vector.
1191: Level: beginner
1192: @*/
1193: PetscErrorCode BVCreateVec(BV bv,Vec *v)1194: {
1199: BVCheckSizes(bv,1);
1201: VecDuplicate(bv->t,v);
1202: return(0);
1203: }
1207: /*@
1208: BVDuplicate - Creates a new basis vector object of the same type and
1209: dimensions as an existing one.
1211: Collective on BV1213: Input Parameter:
1214: . V - basis vectors context
1216: Output Parameter:
1217: . W - location to put the new BV1219: Notes:
1220: The new BV has the same type and dimensions as V, and it shares the same
1221: template vector. Also, the inner product matrix and orthogonalization
1222: options are copied.
1224: BVDuplicate() DOES NOT COPY the entries, but rather allocates storage
1225: for the new basis vectors. Use BVCopy() to copy the contents.
1227: Level: intermediate
1229: .seealso: BVDuplicateResize(), BVCreate(), BVSetSizesFromVec(), BVCopy()
1230: @*/
1231: PetscErrorCode BVDuplicate(BV V,BV *W)1232: {
1238: BVCheckSizes(V,1);
1240: BVCreate(PetscObjectComm((PetscObject)V),W);
1241: BVSetSizesFromVec(*W,V->t,V->m);
1242: BVSetType(*W,((PetscObject)V)->type_name);
1243: BVSetMatrix(*W,V->matrix,V->indef);
1244: BVSetOrthogonalization(*W,V->orthog_type,V->orthog_ref,V->orthog_eta,V->orthog_block);
1245: PetscObjectStateIncrease((PetscObject)*W);
1246: return(0);
1247: }
1251: /*@
1252: BVDuplicateResize - Creates a new basis vector object of the same type and
1253: dimensions as an existing one, but with possibly different number of columns.
1255: Collective on BV1257: Input Parameter:
1258: + V - basis vectors context
1259: - m - the new number of columns
1261: Output Parameter:
1262: . W - location to put the new BV1264: Note:
1265: This is equivalent of a call to BVDuplicate() followed by BVResize(). The
1266: contents of V are not copied to W.
1268: Level: intermediate
1270: .seealso: BVDuplicate(), BVResize()
1271: @*/
1272: PetscErrorCode BVDuplicateResize(BV V,PetscInt m,BV *W)1273: {
1279: BVCheckSizes(V,1);
1282: BVCreate(PetscObjectComm((PetscObject)V),W);
1283: BVSetSizesFromVec(*W,V->t,m);
1284: BVSetType(*W,((PetscObject)V)->type_name);
1285: BVSetMatrix(*W,V->matrix,V->indef);
1286: BVSetOrthogonalization(*W,V->orthog_type,V->orthog_ref,V->orthog_eta,V->orthog_block);
1287: PetscObjectStateIncrease((PetscObject)*W);
1288: return(0);
1289: }
1293: /*@
1294: BVCopy - Copies a basis vector object into another one, W <- V.
1296: Logically Collective on BV1298: Input Parameter:
1299: . V - basis vectors context
1301: Output Parameter:
1302: . W - the copy
1304: Note:
1305: Both V and W must be distributed in the same manner; local copies are
1306: done. Only active columns (excluding the leading ones) are copied.
1307: In the destination W, columns are overwritten starting from the leading ones.
1308: Constraints are not copied.
1310: Level: beginner
1312: .seealso: BVCopyVec(), BVCopyColumn(), BVDuplicate(), BVSetActiveColumns()
1313: @*/
1314: PetscErrorCode BVCopy(BV V,BV W)1315: {
1321: BVCheckSizes(V,1);
1324: BVCheckSizes(W,2);
1326: if (V->n!=W->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension V %D, W %D",V->n,W->n);
1327: if (V->k-V->l>W->m-W->l) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"W has %D non-leading columns, not enough to store %D columns",W->m-W->l,V->k-V->l);
1328: if (!V->n) return(0);
1330: PetscLogEventBegin(BV_Copy,V,W,0,0);
1331: if (V->indef && V->matrix && V->indef==W->indef && V->matrix==W->matrix) {
1332: /* copy signature */
1333: BV_AllocateSignature(W);
1334: PetscMemcpy(W->omega+W->nc+W->l,V->omega+V->nc+V->l,(V->k-V->l)*sizeof(PetscReal));
1335: }
1336: (*V->ops->copy)(V,W);
1337: PetscLogEventEnd(BV_Copy,V,W,0,0);
1338: return(0);
1339: }
1343: /*@
1344: BVCopyVec - Copies one of the columns of a basis vectors object into a Vec.
1346: Logically Collective on BV1348: Input Parameter:
1349: + V - basis vectors context
1350: - j - the column number to be copied
1352: Output Parameter:
1353: . w - the copied column
1355: Note:
1356: Both V and w must be distributed in the same manner; local copies are done.
1358: Level: beginner
1360: .seealso: BVCopy(), BVCopyColumn(), BVInsertVec()
1361: @*/
1362: PetscErrorCode BVCopyVec(BV V,PetscInt j,Vec w)1363: {
1365: PetscInt n,N;
1366: Vec z;
1371: BVCheckSizes(V,1);
1376: VecGetSize(w,&N);
1377: VecGetLocalSize(w,&n);
1378: 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);
1380: PetscLogEventBegin(BV_Copy,V,w,0,0);
1381: BVGetColumn(V,j,&z);
1382: VecCopy(z,w);
1383: BVRestoreColumn(V,j,&z);
1384: PetscLogEventEnd(BV_Copy,V,w,0,0);
1385: return(0);
1386: }
1390: /*@
1391: BVCopyColumn - Copies the values from one of the columns to another one.
1393: Logically Collective on BV1395: Input Parameter:
1396: + V - basis vectors context
1397: . j - the number of the source column
1398: - i - the number of the destination column
1400: Level: beginner
1402: .seealso: BVCopy(), BVCopyVec()
1403: @*/
1404: PetscErrorCode BVCopyColumn(BV V,PetscInt j,PetscInt i)1405: {
1407: Vec z,w;
1412: BVCheckSizes(V,1);
1415: if (j==i) return(0);
1417: PetscLogEventBegin(BV_Copy,V,0,0,0);
1418: if (V->omega) V->omega[i] = V->omega[j];
1419: BVGetColumn(V,j,&z);
1420: BVGetColumn(V,i,&w);
1421: VecCopy(z,w);
1422: BVRestoreColumn(V,j,&z);
1423: BVRestoreColumn(V,i,&w);
1424: PetscLogEventEnd(BV_Copy,V,0,0,0);
1425: return(0);
1426: }