1: /*
2: Basic BV routines.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2016, 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 vector
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: BVSetRandomContext - Sets the PetscRandom object associated with the BV,
756: to be used in operations that need random numbers.
758: Collective on BV760: Input Parameters:
761: + bv - the basis vectors context
762: - rand - the random number generator context
764: Level: advanced
766: .seealso: BVGetRandomContext(), BVSetRandom(), BVSetRandomColumn()
767: @*/
768: PetscErrorCode BVSetRandomContext(BV bv,PetscRandom rand)769: {
776: PetscObjectReference((PetscObject)rand);
777: PetscRandomDestroy(&bv->rand);
778: bv->rand = rand;
779: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->rand);
780: return(0);
781: }
785: /*@
786: BVGetRandomContext - Gets the PetscRandom object associated with the BV.
788: Not Collective
790: Input Parameter:
791: . bv - the basis vectors context
793: Output Parameter:
794: . rand - the random number generator context
796: Level: advanced
798: .seealso: BVSetRandomContext(), BVSetRandom(), BVSetRandomColumn()
799: @*/
800: PetscErrorCode BVGetRandomContext(BV bv,PetscRandom* rand)801: {
807: if (!bv->rand) {
808: PetscRandomCreate(PetscObjectComm((PetscObject)bv),&bv->rand);
809: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->rand);
810: if (bv->rrandom) {
811: PetscRandomSetSeed(bv->rand,0x12345678);
812: PetscRandomSeed(bv->rand);
813: }
814: }
815: *rand = bv->rand;
816: return(0);
817: }
821: /*@
822: BVSetFromOptions - Sets BV options from the options database.
824: Collective on BV826: Input Parameter:
827: . bv - the basis vectors context
829: Level: beginner
830: @*/
831: PetscErrorCode BVSetFromOptions(BV bv)832: {
834: char type[256];
835: PetscBool flg;
836: PetscReal r;
840: BVRegisterAll();
841: PetscObjectOptionsBegin((PetscObject)bv);
842: PetscOptionsFList("-bv_type","Basis Vectors type","BVSetType",BVList,(char*)(((PetscObject)bv)->type_name?((PetscObject)bv)->type_name:BVSVEC),type,256,&flg);
843: if (flg) {
844: BVSetType(bv,type);
845: }
846: /*
847: Set the type if it was never set.
848: */
849: if (!((PetscObject)bv)->type_name) {
850: BVSetType(bv,BVSVEC);
851: }
853: PetscOptionsEnum("-bv_orthog_type","Orthogonalization method","BVSetOrthogonalization",BVOrthogTypes,(PetscEnum)bv->orthog_type,(PetscEnum*)&bv->orthog_type,NULL);
854: PetscOptionsEnum("-bv_orthog_refine","Iterative refinement mode during orthogonalization","BVSetOrthogonalization",BVOrthogRefineTypes,(PetscEnum)bv->orthog_ref,(PetscEnum*)&bv->orthog_ref,NULL);
855: PetscOptionsEnum("-bv_orthog_block","Block orthogonalization method","BVSetOrthogonalization",BVOrthogBlockTypes,(PetscEnum)bv->orthog_block,(PetscEnum*)&bv->orthog_block,NULL);
856: r = bv->orthog_eta;
857: PetscOptionsReal("-bv_orthog_eta","Parameter of iterative refinement during orthogonalization","BVSetOrthogonalization",r,&r,NULL);
858: BVSetOrthogonalization(bv,bv->orthog_type,bv->orthog_ref,r,bv->orthog_block);
860: PetscOptionsEnum("-bv_matmult","Method for BVMatMult","BVSetMatMultMethod",BVMatMultTypes,(PetscEnum)bv->vmm,(PetscEnum*)&bv->vmm,NULL);
862: /* undocumented option to generate random vectors that are independent of the number of processes */
863: PetscOptionsGetBool(NULL,NULL,"-bv_reproducible_random",&bv->rrandom,NULL);
865: if (!bv->rand) { BVGetRandomContext(bv,&bv->rand); }
866: PetscRandomSetFromOptions(bv->rand);
868: if (bv->ops->create) bv->defersfo = PETSC_TRUE; /* defer call to setfromoptions */
869: else if (bv->ops->setfromoptions) {
870: (*bv->ops->setfromoptions)(PetscOptionsObject,bv);
871: }
872: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)bv);
873: PetscOptionsEnd();
874: return(0);
875: }
879: /*@
880: BVSetOrthogonalization - Specifies the method used for the orthogonalization of
881: vectors (classical or modified Gram-Schmidt with or without refinement), and
882: for the block-orthogonalization (simultaneous orthogonalization of a set of
883: vectors).
885: Logically Collective on BV887: Input Parameters:
888: + bv - the basis vectors context
889: . type - the method of vector orthogonalization
890: . refine - type of refinement
891: . eta - parameter for selective refinement
892: - block - the method of block orthogonalization
894: Options Database Keys:
895: + -bv_orthog_type <type> - Where <type> is cgs for Classical Gram-Schmidt orthogonalization
896: (default) or mgs for Modified Gram-Schmidt orthogonalization
897: . -bv_orthog_refine <ref> - Where <ref> is one of never, ifneeded (default) or always
898: . -bv_orthog_eta <eta> - For setting the value of eta
899: - -bv_orthog_block <block> - Where <block> is the block-orthogonalization method
901: Notes:
902: The default settings work well for most problems.
904: The parameter eta should be a real value between 0 and 1 (or PETSC_DEFAULT).
905: The value of eta is used only when the refinement type is "ifneeded".
907: When using several processors, MGS is likely to result in bad scalability.
909: If the method set for block orthogonalization is GS, then the computation
910: is done column by column with the vector orthogonalization.
912: Level: advanced
914: .seealso: BVOrthogonalizeColumn(), BVGetOrthogonalization(), BVOrthogType, BVOrthogRefineType, BVOrthogBlockType915: @*/
916: PetscErrorCode BVSetOrthogonalization(BV bv,BVOrthogType type,BVOrthogRefineType refine,PetscReal eta,BVOrthogBlockType block)917: {
924: switch (type) {
925: case BV_ORTHOG_CGS:
926: case BV_ORTHOG_MGS:
927: bv->orthog_type = type;
928: break;
929: default:930: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
931: }
932: switch (refine) {
933: case BV_ORTHOG_REFINE_NEVER:
934: case BV_ORTHOG_REFINE_IFNEEDED:
935: case BV_ORTHOG_REFINE_ALWAYS:
936: bv->orthog_ref = refine;
937: break;
938: default:939: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown refinement type");
940: }
941: if (eta == PETSC_DEFAULT) {
942: bv->orthog_eta = 0.7071;
943: } else {
944: if (eta <= 0.0 || eta > 1.0) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Invalid eta value");
945: bv->orthog_eta = eta;
946: }
947: switch (block) {
948: case BV_ORTHOG_BLOCK_GS:
949: case BV_ORTHOG_BLOCK_CHOL:
950: bv->orthog_block = block;
951: break;
952: default:953: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown block orthogonalization type");
954: }
955: return(0);
956: }
960: /*@
961: BVGetOrthogonalization - Gets the orthogonalization settings from the BV object.
963: Not Collective
965: Input Parameter:
966: . bv - basis vectors context
968: Output Parameter:
969: + type - the method of vector orthogonalization
970: . refine - type of refinement
971: . eta - parameter for selective refinement
972: - block - the method of block orthogonalization
974: Level: advanced
976: .seealso: BVOrthogonalizeColumn(), BVSetOrthogonalization(), BVOrthogType, BVOrthogRefineType, BVOrthogBlockType977: @*/
978: PetscErrorCode BVGetOrthogonalization(BV bv,BVOrthogType *type,BVOrthogRefineType *refine,PetscReal *eta,BVOrthogBlockType *block)979: {
982: if (type) *type = bv->orthog_type;
983: if (refine) *refine = bv->orthog_ref;
984: if (eta) *eta = bv->orthog_eta;
985: if (block) *block = bv->orthog_block;
986: return(0);
987: }
991: /*@
992: BVSetMatMultMethod - Specifies the method used for the BVMatMult() operation.
994: Logically Collective on BV996: Input Parameters:
997: + bv - the basis vectors context
998: - method - the method for the BVMatMult() operation
1000: Options Database Keys:
1001: . -bv_matmult <meth> - choose one of the methods: vecs, mat, mat_save
1003: Note:
1004: Allowed values are:
1005: + BV_MATMULT_VECS - perform a matrix-vector multiply per each column
1006: . BV_MATMULT_MAT - carry out a MatMatMult() product with a dense matrix
1007: - BV_MATMULT_MAT_SAVE - call MatMatMult() and keep auxiliary matrices
1008: The default is BV_MATMULT_MAT.
1010: Level: advanced
1012: .seealso: BVMatMult(), BVGetMatMultMethod(), BVMatMultType1013: @*/
1014: PetscErrorCode BVSetMatMultMethod(BV bv,BVMatMultType method)1015: {
1019: switch (method) {
1020: case BV_MATMULT_VECS:
1021: case BV_MATMULT_MAT:
1022: case BV_MATMULT_MAT_SAVE:
1023: bv->vmm = method;
1024: break;
1025: default:1026: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown matmult method");
1027: }
1028: return(0);
1029: }
1033: /*@
1034: BVGetMatMultMethod - Gets the method used for the BVMatMult() operation.
1036: Not Collective
1038: Input Parameter:
1039: . bv - basis vectors context
1041: Output Parameter:
1042: . method - the method for the BVMatMult() operation
1044: Level: advanced
1046: .seealso: BVMatMult(), BVSetMatMultMethod(), BVMatMultType1047: @*/
1048: PetscErrorCode BVGetMatMultMethod(BV bv,BVMatMultType *method)1049: {
1053: *method = bv->vmm;
1054: return(0);
1055: }
1059: /*@
1060: BVGetColumn - Returns a Vec object that contains the entries of the
1061: requested column of the basis vectors object.
1063: Logically Collective on BV1065: Input Parameters:
1066: + bv - the basis vectors context
1067: - j - the index of the requested column
1069: Output Parameter:
1070: . v - vector containing the jth column
1072: Notes:
1073: The returned Vec must be seen as a reference (not a copy) of the BV1074: column, that is, modifying the Vec will change the BV entries as well.
1076: The returned Vec must not be destroyed. BVRestoreColumn() must be
1077: called when it is no longer needed. At most, two columns can be fetched,
1078: that is, this function can only be called twice before the corresponding
1079: BVRestoreColumn() is invoked.
1081: A negative index j selects the i-th constraint, where i=-j. Constraints
1082: should not be modified.
1084: Level: beginner
1086: .seealso: BVRestoreColumn(), BVInsertConstraints()
1087: @*/
1088: PetscErrorCode BVGetColumn(BV bv,PetscInt j,Vec *v)1089: {
1091: PetscInt l;
1096: BVCheckSizes(bv,1);
1098: 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);
1099: if (j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"You requested column %D but only %D are available",j,bv->m);
1100: 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);
1101: l = BVAvailableVec;
1102: 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");
1103: (*bv->ops->getcolumn)(bv,j,v);
1104: bv->ci[l] = j;
1105: PetscObjectStateGet((PetscObject)bv->cv[l],&bv->st[l]);
1106: PetscObjectGetId((PetscObject)bv->cv[l],&bv->id[l]);
1107: *v = bv->cv[l];
1108: return(0);
1109: }
1113: /*@
1114: BVRestoreColumn - Restore a column obtained with BVGetColumn().
1116: Logically Collective on BV1118: Input Parameters:
1119: + bv - the basis vectors context
1120: . j - the index of the column
1121: - v - vector obtained with BVGetColumn()
1123: Note:
1124: The arguments must match the corresponding call to BVGetColumn().
1126: Level: beginner
1128: .seealso: BVGetColumn()
1129: @*/
1130: PetscErrorCode BVRestoreColumn(BV bv,PetscInt j,Vec *v)1131: {
1132: PetscErrorCode ierr;
1133: PetscObjectId id;
1134: PetscObjectState st;
1135: PetscInt l;
1140: BVCheckSizes(bv,1);
1144: 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);
1145: if (j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"You requested column %D but only %D are available",j,bv->m);
1146: 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);
1147: l = (j==bv->ci[0])? 0: 1;
1148: PetscObjectGetId((PetscObject)*v,&id);
1149: 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");
1150: PetscObjectStateGet((PetscObject)*v,&st);
1151: if (st!=bv->st[l]) {
1152: PetscObjectStateIncrease((PetscObject)bv);
1153: }
1154: if (bv->ops->restorecolumn) {
1155: (*bv->ops->restorecolumn)(bv,j,v);
1156: } else bv->cv[l] = NULL;
1157: bv->ci[l] = -bv->nc-1;
1158: bv->st[l] = -1;
1159: bv->id[l] = 0;
1160: *v = NULL;
1161: return(0);
1162: }
1166: /*@C
1167: BVGetArray - Returns a pointer to a contiguous array that contains this
1168: processor's portion of the BV data.
1170: Logically Collective on BV1172: Input Parameters:
1173: . bv - the basis vectors context
1175: Output Parameter:
1176: . a - location to put pointer to the array
1178: Notes:
1179: BVRestoreArray() must be called when access to the array is no longer needed.
1180: This operation may imply a data copy, for BV types that do not store
1181: data contiguously in memory.
1183: The pointer will normally point to the first entry of the first column,
1184: but if the BV has constraints then these go before the regular columns.
1186: Level: advanced
1188: .seealso: BVRestoreArray(), BVInsertConstraints()
1189: @*/
1190: PetscErrorCode BVGetArray(BV bv,PetscScalar **a)1191: {
1197: BVCheckSizes(bv,1);
1198: (*bv->ops->getarray)(bv,a);
1199: return(0);
1200: }
1204: /*@C
1205: BVRestoreArray - Restore the BV object after BVGetArray() has been called.
1207: Logically Collective on BV1209: Input Parameters:
1210: + bv - the basis vectors context
1211: - a - location of pointer to array obtained from BVGetArray()
1213: Note:
1214: This operation may imply a data copy, for BV types that do not store
1215: data contiguously in memory.
1217: Level: advanced
1219: .seealso: BVGetColumn()
1220: @*/
1221: PetscErrorCode BVRestoreArray(BV bv,PetscScalar **a)1222: {
1228: BVCheckSizes(bv,1);
1229: if (bv->ops->restorearray) {
1230: (*bv->ops->restorearray)(bv,a);
1231: }
1232: if (a) *a = NULL;
1233: PetscObjectStateIncrease((PetscObject)bv);
1234: return(0);
1235: }
1239: /*@C
1240: BVGetArrayRead - Returns a read-only pointer to a contiguous array that
1241: contains this processor's portion of the BV data.
1243: Not Collective
1245: Input Parameters:
1246: . bv - the basis vectors context
1248: Output Parameter:
1249: . a - location to put pointer to the array
1251: Notes:
1252: BVRestoreArrayRead() must be called when access to the array is no
1253: longer needed. This operation may imply a data copy, for BV types that
1254: do not store data contiguously in memory.
1256: The pointer will normally point to the first entry of the first column,
1257: but if the BV has constraints then these go before the regular columns.
1259: Level: advanced
1261: .seealso: BVRestoreArray(), BVInsertConstraints()
1262: @*/
1263: PetscErrorCode BVGetArrayRead(BV bv,const PetscScalar **a)1264: {
1270: BVCheckSizes(bv,1);
1271: (*bv->ops->getarrayread)(bv,a);
1272: return(0);
1273: }
1277: /*@C
1278: BVRestoreArrayRead - Restore the BV object after BVGetArrayRead() has
1279: been called.
1281: Logically Collective on BV1283: Input Parameters:
1284: + bv - the basis vectors context
1285: - a - location of pointer to array obtained from BVGetArrayRead()
1287: Level: advanced
1289: .seealso: BVGetColumn()
1290: @*/
1291: PetscErrorCode BVRestoreArrayRead(BV bv,const PetscScalar **a)1292: {
1298: BVCheckSizes(bv,1);
1299: if (bv->ops->restorearrayread) {
1300: (*bv->ops->restorearrayread)(bv,a);
1301: }
1302: if (a) *a = NULL;
1303: return(0);
1304: }
1308: /*@
1309: BVCreateVec - Creates a new Vec object with the same type and dimensions
1310: as the columns of the basis vectors object.
1312: Collective on BV1314: Input Parameter:
1315: . bv - the basis vectors context
1317: Output Parameter:
1318: . v - the new vector
1320: Note:
1321: The user is responsible of destroying the returned vector.
1323: Level: beginner
1324: @*/
1325: PetscErrorCode BVCreateVec(BV bv,Vec *v)1326: {
1331: BVCheckSizes(bv,1);
1333: VecDuplicate(bv->t,v);
1334: return(0);
1335: }
1339: PETSC_STATIC_INLINE PetscErrorCode BVDuplicate_Private(BV V,PetscInt m,BV *W)1340: {
1344: BVCreate(PetscObjectComm((PetscObject)V),W);
1345: BVSetSizesFromVec(*W,V->t,m);
1346: BVSetType(*W,((PetscObject)V)->type_name);
1347: BVSetMatrix(*W,V->matrix,V->indef);
1348: BVSetOrthogonalization(*W,V->orthog_type,V->orthog_ref,V->orthog_eta,V->orthog_block);
1349: (*W)->vmm = V->vmm;
1350: (*W)->rrandom = V->rrandom;
1351: if (V->ops->duplicate) { (*V->ops->duplicate)(V,W); }
1352: PetscObjectStateIncrease((PetscObject)*W);
1353: return(0);
1354: }
1358: /*@
1359: BVDuplicate - Creates a new basis vector object of the same type and
1360: dimensions as an existing one.
1362: Collective on BV1364: Input Parameter:
1365: . V - basis vectors context
1367: Output Parameter:
1368: . W - location to put the new BV1370: Notes:
1371: The new BV has the same type and dimensions as V, and it shares the same
1372: template vector. Also, the inner product matrix and orthogonalization
1373: options are copied.
1375: BVDuplicate() DOES NOT COPY the entries, but rather allocates storage
1376: for the new basis vectors. Use BVCopy() to copy the contents.
1378: Level: intermediate
1380: .seealso: BVDuplicateResize(), BVCreate(), BVSetSizesFromVec(), BVCopy()
1381: @*/
1382: PetscErrorCode BVDuplicate(BV V,BV *W)1383: {
1389: BVCheckSizes(V,1);
1391: BVDuplicate_Private(V,V->m,W);
1392: return(0);
1393: }
1397: /*@
1398: BVDuplicateResize - Creates a new basis vector object of the same type and
1399: dimensions as an existing one, but with possibly different number of columns.
1401: Collective on BV1403: Input Parameter:
1404: + V - basis vectors context
1405: - m - the new number of columns
1407: Output Parameter:
1408: . W - location to put the new BV1410: Note:
1411: This is equivalent of a call to BVDuplicate() followed by BVResize(). The
1412: contents of V are not copied to W.
1414: Level: intermediate
1416: .seealso: BVDuplicate(), BVResize()
1417: @*/
1418: PetscErrorCode BVDuplicateResize(BV V,PetscInt m,BV *W)1419: {
1425: BVCheckSizes(V,1);
1428: BVDuplicate_Private(V,m,W);
1429: return(0);
1430: }
1434: /*@
1435: BVCopy - Copies a basis vector object into another one, W <- V.
1437: Logically Collective on BV1439: Input Parameter:
1440: . V - basis vectors context
1442: Output Parameter:
1443: . W - the copy
1445: Note:
1446: Both V and W must be distributed in the same manner; local copies are
1447: done. Only active columns (excluding the leading ones) are copied.
1448: In the destination W, columns are overwritten starting from the leading ones.
1449: Constraints are not copied.
1451: Level: beginner
1453: .seealso: BVCopyVec(), BVCopyColumn(), BVDuplicate(), BVSetActiveColumns()
1454: @*/
1455: PetscErrorCode BVCopy(BV V,BV W)1456: {
1462: BVCheckSizes(V,1);
1465: BVCheckSizes(W,2);
1467: if (V->n!=W->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension V %D, W %D",V->n,W->n);
1468: 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);
1469: if (!V->n) return(0);
1471: PetscLogEventBegin(BV_Copy,V,W,0,0);
1472: if (V->indef && V->matrix && V->indef==W->indef && V->matrix==W->matrix) {
1473: /* copy signature */
1474: BV_AllocateSignature(W);
1475: PetscMemcpy(W->omega+W->nc+W->l,V->omega+V->nc+V->l,(V->k-V->l)*sizeof(PetscReal));
1476: }
1477: (*V->ops->copy)(V,W);
1478: PetscLogEventEnd(BV_Copy,V,W,0,0);
1479: return(0);
1480: }
1484: /*@
1485: BVCopyVec - Copies one of the columns of a basis vectors object into a Vec.
1487: Logically Collective on BV1489: Input Parameter:
1490: + V - basis vectors context
1491: - j - the column number to be copied
1493: Output Parameter:
1494: . w - the copied column
1496: Note:
1497: Both V and w must be distributed in the same manner; local copies are done.
1499: Level: beginner
1501: .seealso: BVCopy(), BVCopyColumn(), BVInsertVec()
1502: @*/
1503: PetscErrorCode BVCopyVec(BV V,PetscInt j,Vec w)1504: {
1506: PetscInt n,N;
1507: Vec z;
1512: BVCheckSizes(V,1);
1517: VecGetSize(w,&N);
1518: VecGetLocalSize(w,&n);
1519: 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);
1521: PetscLogEventBegin(BV_Copy,V,w,0,0);
1522: BVGetColumn(V,j,&z);
1523: VecCopy(z,w);
1524: BVRestoreColumn(V,j,&z);
1525: PetscLogEventEnd(BV_Copy,V,w,0,0);
1526: return(0);
1527: }
1531: /*@
1532: BVCopyColumn - Copies the values from one of the columns to another one.
1534: Logically Collective on BV1536: Input Parameter:
1537: + V - basis vectors context
1538: . j - the number of the source column
1539: - i - the number of the destination column
1541: Level: beginner
1543: .seealso: BVCopy(), BVCopyVec()
1544: @*/
1545: PetscErrorCode BVCopyColumn(BV V,PetscInt j,PetscInt i)1546: {
1548: Vec z,w;
1553: BVCheckSizes(V,1);
1556: if (j==i) return(0);
1558: PetscLogEventBegin(BV_Copy,V,0,0,0);
1559: if (V->omega) V->omega[i] = V->omega[j];
1560: BVGetColumn(V,j,&z);
1561: BVGetColumn(V,i,&w);
1562: VecCopy(z,w);
1563: BVRestoreColumn(V,j,&z);
1564: BVRestoreColumn(V,i,&w);
1565: PetscLogEventEnd(BV_Copy,V,0,0,0);
1566: return(0);
1567: }