1: /*
2: Basic BV routines.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2014, 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: }
158: return(0);
159: }
163: /*@
164: BVSetSizesFromVec - Sets the local and global sizes, and the number of columns.
165: Local and global sizes are specified indirectly by passing a template vector.
167: Collective on BV169: Input Parameters:
170: + bv - the basis vectors
171: . t - the template vectors
172: - m - the number of columns
174: Level: beginner
176: .seealso: BVSetSizes(), BVGetSizes(), BVResize()
177: @*/
178: PetscErrorCode BVSetSizesFromVec(BV bv,Vec t,PetscInt m)179: {
181: PetscInt ma;
188: if (m <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of columns %D must be positive",m);
189: if (bv->t) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Template vector was already set by a previous call to BVSetSizes/FromVec");
190: VecGetSize(t,&bv->N);
191: VecGetLocalSize(t,&bv->n);
192: if (bv->matrix) { /* check compatible dimensions of user-provided matrix */
193: MatGetLocalSize(bv->matrix,&ma,NULL);
194: 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);
195: }
196: bv->m = m;
197: bv->k = m;
198: bv->t = t;
199: PetscObjectReference((PetscObject)t);
200: if (bv->ops->create) {
201: (*bv->ops->create)(bv);
202: bv->ops->create = 0;
203: }
204: return(0);
205: }
209: /*@
210: BVGetSizes - Returns the local and global sizes, and the number of columns.
212: Not Collective
214: Input Parameter:
215: . bv - the basis vectors
217: Output Parameters:
218: + n - the local size
219: . N - the global size
220: - m - the number of columns
222: Note:
223: Normal usage requires that bv has already been given its sizes, otherwise
224: the call fails. However, this function can also be used to determine if
225: a BV object has been initialized completely (sizes and type). For this,
226: call with n=NULL and N=NULL, then a return value of m=0 indicates that
227: the BV object is not ready for use yet.
229: Level: beginner
231: .seealso: BVSetSizes(), BVSetSizesFromVec()
232: @*/
233: PetscErrorCode BVGetSizes(BV bv,PetscInt *n,PetscInt *N,PetscInt *m)234: {
236: if (!bv) {
237: if (m && !n && !N) *m = 0;
238: return(0);
239: }
241: if (n || N) BVCheckSizes(bv,1);
242: if (n) *n = bv->n;
243: if (N) *N = bv->N;
244: if (m) *m = bv->m;
245: if (m && !n && !N && !((PetscObject)bv)->type_name) *m = 0;
246: return(0);
247: }
251: /*@
252: BVSetNumConstraints - Set the number of constraints.
254: Logically Collective on BV256: Input Parameters:
257: + V - basis vectors
258: - nc - number of constraints
260: Notes:
261: This function sets the number of constraints to nc and marks all remaining
262: columns as regular. Normal user would call BVInsertConstraints() instead.
264: Level: developer
266: .seealso: BVInsertConstraints()
267: @*/
268: PetscErrorCode BVSetNumConstraints(BV V,PetscInt nc)269: {
271: PetscInt total;
276: if (!nc) return(0);
277: if (nc<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of constraints (given %D) cannot be negative",nc);
279: BVCheckSizes(V,1);
280: 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");
282: total = V->nc+V->m;
283: V->nc = nc;
284: V->m = total-nc;
285: if (V->m<=0) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Not enough columns for the given nc value");
286: V->ci[0] = -V->nc-1;
287: V->ci[1] = -V->nc-1;
288: V->l = 0;
289: V->k = V->m;
290: PetscObjectStateIncrease((PetscObject)V);
291: return(0);
292: }
296: /*@
297: BVGetNumConstraints - Returns the number of constraints.
299: Not Collective
301: Input Parameter:
302: . bv - the basis vectors
304: Output Parameters:
305: . nc - the number of constraints
307: Level: advanced
309: .seealso: BVGetSizes(), BVInsertConstraints()
310: @*/
311: PetscErrorCode BVGetNumConstraints(BV bv,PetscInt *nc)312: {
316: *nc = bv->nc;
317: return(0);
318: }
322: /*@
323: BVResize - Change the number of columns.
325: Collective on BV327: Input Parameters:
328: + bv - the basis vectors
329: . m - the new number of columns
330: - copy - a flag indicating whether current values should be kept
332: Note:
333: Internal storage is reallocated. If the copy flag is set to true, then
334: the contents are copied to the leading part of the new space.
336: Level: advanced
338: .seealso: BVSetSizes(), BVSetSizesFromVec()
339: @*/
340: PetscErrorCode BVResize(BV bv,PetscInt m,PetscBool copy)341: {
343: PetscReal *omega;
344: PetscInt i;
351: if (m <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of columns %D must be positive",m);
352: if (bv->nc) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot resize a BV with constraints");
353: if (bv->m == m) return(0);
355: PetscLogEventBegin(BV_Create,bv,0,0,0);
356: (*bv->ops->resize)(bv,m,copy);
357: PetscFree2(bv->h,bv->c);
358: if (bv->omega) {
359: PetscMalloc1(m,&omega);
360: PetscLogObjectMemory((PetscObject)bv,m*sizeof(PetscReal));
361: for (i=0;i<m;i++) omega[i] = 1.0;
362: if (copy) {
363: PetscMemcpy(omega,bv->omega,PetscMin(m,bv->m)*sizeof(PetscReal));
364: }
365: PetscFree(bv->omega);
366: bv->omega = omega;
367: }
368: bv->m = m;
369: bv->k = PetscMin(bv->k,m);
370: bv->l = PetscMin(bv->l,m);
371: PetscLogEventEnd(BV_Create,bv,0,0,0);
372: return(0);
373: }
377: /*@
378: BVSetActiveColumns - Specify the columns that will be involved in operations.
380: Logically Collective on BV382: Input Parameters:
383: + bv - the basis vectors context
384: . l - number of leading columns
385: - k - number of active columns
387: Notes:
388: In operations such as BVMult() or BVDot(), only the first k columns are
389: considered. This is useful when the BV is filled from left to right, so
390: the last m-k columns do not have relevant information.
392: Also in operations such as BVMult() or BVDot(), the first l columns are
393: normally not included in the computation. See the manpage of each
394: operation.
396: In orthogonalization operations, the first l columns are treated
397: differently: they participate in the orthogonalization but the computed
398: coefficients are not stored.
400: Level: intermediate
402: .seealso: BVGetActiveColumns(), BVSetSizes()
403: @*/
404: PetscErrorCode BVSetActiveColumns(BV bv,PetscInt l,PetscInt k)405: {
410: BVCheckSizes(bv,1);
411: if (k==PETSC_DECIDE || k==PETSC_DEFAULT) {
412: bv->k = bv->m;
413: } else {
414: if (k<0 || k>bv->m) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of k. Must be between 0 and m");
415: bv->k = k;
416: }
417: if (l==PETSC_DECIDE || l==PETSC_DEFAULT) {
418: bv->l = 0;
419: } else {
420: if (l<0 || l>bv->k) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of l. Must be between 0 and k");
421: bv->l = l;
422: }
423: return(0);
424: }
428: /*@
429: BVGetActiveColumns - Returns the current active dimensions.
431: Not Collective
433: Input Parameter:
434: . bv - the basis vectors context
436: Output Parameter:
437: + l - number of leading columns
438: - k - number of active columns
440: Level: intermediate
442: .seealso: BVSetActiveColumns()
443: @*/
444: PetscErrorCode BVGetActiveColumns(BV bv,PetscInt *l,PetscInt *k)445: {
448: if (l) *l = bv->l;
449: if (k) *k = bv->k;
450: return(0);
451: }
455: /*@
456: BVSetMatrix - Specifies the inner product to be used in orthogonalization.
458: Collective on BV460: Input Parameters:
461: + bv - the basis vectors context
462: . B - a symmetric matrix (may be NULL)
463: - indef - a flag indicating if the matrix is indefinite
465: Notes:
466: This is used to specify a non-standard inner product, whose matrix
467: representation is given by B. Then, all inner products required during
468: orthogonalization are computed as (x,y)_B=y^H*B*x rather than the
469: standard form (x,y)=y^H*x.
471: Matrix B must be real symmetric (or complex Hermitian). A genuine inner
472: product requires that B is also positive (semi-)definite. However, we
473: also allow for an indefinite B (setting indef=PETSC_TRUE), in which
474: case the orthogonalization uses an indefinite inner product.
476: This affects operations BVDot(), BVNorm(), BVOrthogonalize(), and variants.
478: Setting B=NULL has the same effect as if the identity matrix was passed.
480: Level: advanced
482: .seealso: BVGetMatrix(), BVDot(), BVNorm(), BVOrthogonalize()
483: @*/
484: PetscErrorCode BVSetMatrix(BV bv,Mat B,PetscBool indef)485: {
487: PetscInt m,n;
492: if (B) {
494: MatGetLocalSize(B,&m,&n);
495: if (m!=n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix must be square");
496: if (bv->m && bv->n!=n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension BV %D, Mat %D",bv->n,n);
497: }
498: MatDestroy(&bv->matrix);
499: if (B) PetscObjectReference((PetscObject)B);
500: bv->matrix = B;
501: bv->indef = indef;
502: if (B && !bv->Bx) {
503: MatGetVecs(B,&bv->Bx,NULL);
504: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->Bx);
505: }
506: return(0);
507: }
511: /*@C
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: }
538: /*@C
539: BVApplyMatrix - Multiplies a vector by the matrix representation of the
540: inner product.
542: Neighbor-wise Collective on BV and Vec
544: Input Parameter:
545: + bv - the basis vectors context
546: - x - the vector
548: Output Parameter:
549: . y - the result
551: Note:
552: If no matrix was specified this function copies the vector.
554: Level: advanced
556: .seealso: BVSetMatrix()
557: @*/
558: PetscErrorCode BVApplyMatrix(BV bv,Vec x,Vec y)559: {
566: if (bv->matrix) {
567: BV_IPMatMult(bv,x);
568: VecCopy(bv->Bx,y);
569: } else {
570: VecCopy(x,y);
571: }
572: return(0);
573: }
577: /*@
578: BVSetSignature - Sets the signature matrix to be used in orthogonalization.
580: Logically Collective on BV582: Input Parameter:
583: + bv - the basis vectors context
584: - omega - a vector representing the diagonal of the signature matrix
586: Note:
587: The signature matrix Omega = V'*B*V is relevant only for an indefinite B.
589: Level: developer
591: .seealso: BVSetMatrix(), BVGetSignature()
592: @*/
593: PetscErrorCode BVSetSignature(BV bv,Vec omega)594: {
596: PetscInt i,n;
597: PetscScalar *pomega;
601: BVCheckSizes(bv,1);
605: VecGetSize(omega,&n);
606: if (n!=bv->k) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Vec argument has %D elements, should be %D",n,bv->k);
607: BV_AllocateSignature(bv);
608: if (bv->indef) {
609: VecGetArray(omega,&pomega);
610: for (i=0;i<n;i++) bv->omega[bv->nc+i] = PetscRealPart(pomega[i]);
611: VecRestoreArray(omega,&pomega);
612: } else {
613: PetscInfo(bv,"Ignoring signature because BV is not indefinite\n");
614: }
615: return(0);
616: }
620: /*@
621: BVGetSignature - Retrieves the signature matrix from last orthogonalization.
623: Not collective
625: Input Parameter:
626: . bv - the basis vectors context
628: Output Parameter:
629: . omega - a vector representing the diagonal of the signature matrix
631: Note:
632: The signature matrix Omega = V'*B*V is relevant only for an indefinite B.
634: Level: developer
636: .seealso: BVSetMatrix(), BVSetSignature()
637: @*/
638: PetscErrorCode BVGetSignature(BV bv,Vec omega)639: {
641: PetscInt i,n;
642: PetscScalar *pomega;
646: BVCheckSizes(bv,1);
650: VecGetSize(omega,&n);
651: if (n!=bv->k) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Vec argument has %D elements, should be %D",n,bv->k);
652: if (bv->indef && bv->omega) {
653: VecGetArray(omega,&pomega);
654: for (i=0;i<n;i++) pomega[i] = bv->omega[bv->nc+i];
655: VecRestoreArray(omega,&pomega);
656: } else {
657: VecSet(omega,1.0);
658: }
659: return(0);
660: }
664: /*@
665: BVSetFromOptions - Sets BV options from the options database.
667: Collective on BV669: Input Parameter:
670: . bv - the basis vectors context
672: Level: beginner
673: @*/
674: PetscErrorCode BVSetFromOptions(BV bv)675: {
677: char type[256];
678: PetscBool flg;
679: PetscReal r;
680: PetscInt i,j;
681: const char *orth_list[2] = {"cgs","mgs"};
682: const char *ref_list[3] = {"ifneeded","never","always"};
686: if (!BVRegisterAllCalled) { BVRegisterAll(); }
687: PetscObjectOptionsBegin((PetscObject)bv);
688: PetscOptionsFList("-bv_type","Basis Vectors type","BVSetType",BVList,(char*)(((PetscObject)bv)->type_name?((PetscObject)bv)->type_name:BVSVEC),type,256,&flg);
689: if (flg) {
690: BVSetType(bv,type);
691: }
692: /*
693: Set the type if it was never set.
694: */
695: if (!((PetscObject)bv)->type_name) {
696: BVSetType(bv,BVSVEC);
697: }
699: i = bv->orthog_type;
700: PetscOptionsEList("-bv_orthog_type","Orthogonalization method","BVSetOrthogonalization",orth_list,2,orth_list[i],&i,NULL);
701: j = bv->orthog_ref;
702: PetscOptionsEList("-bv_orthog_refine","Iterative refinement mode during orthogonalization","BVSetOrthogonalization",ref_list,3,ref_list[j],&j,NULL);
703: r = bv->orthog_eta;
704: PetscOptionsReal("-bv_orthog_eta","Parameter of iterative refinement during orthogonalization","BVSetOrthogonalization",r,&r,NULL);
705: BVSetOrthogonalization(bv,(BVOrthogType)i,(BVOrthogRefineType)j,r);
707: if (bv->ops->setfromoptions) {
708: (*bv->ops->setfromoptions)(bv);
709: }
710: PetscObjectProcessOptionsHandlers((PetscObject)bv);
711: PetscOptionsEnd();
712: return(0);
713: }
717: /*@
718: BVSetOrthogonalization - Specifies the type of orthogonalization technique
719: to be used (classical or modified Gram-Schmidt with or without refinement).
721: Logically Collective on BV723: Input Parameters:
724: + bv - the basis vectors context
725: . type - the type of orthogonalization technique
726: . refine - type of refinement
727: - eta - parameter for selective refinement
729: Options Database Keys:
730: + -bv_orthog_type <type> - Where <type> is cgs for Classical Gram-Schmidt orthogonalization
731: (default) or mgs for Modified Gram-Schmidt orthogonalization
732: . -bv_orthog_refine <ref> - Where <ref> is one of never, ifneeded (default) or always
733: - -bv_orthog_eta <eta> - For setting the value of eta
735: Notes:
736: The default settings work well for most problems.
738: The parameter eta should be a real value between 0 and 1 (or PETSC_DEFAULT).
739: The value of eta is used only when the refinement type is "ifneeded".
741: When using several processors, MGS is likely to result in bad scalability.
743: Level: advanced
745: .seealso: BVOrthogonalizeColumn(), BVGetOrthogonalization(), BVOrthogType, BVOrthogRefineType746: @*/
747: PetscErrorCode BVSetOrthogonalization(BV bv,BVOrthogType type,BVOrthogRefineType refine,PetscReal eta)748: {
754: switch (type) {
755: case BV_ORTHOG_CGS:
756: case BV_ORTHOG_MGS:
757: bv->orthog_type = type;
758: break;
759: default:760: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
761: }
762: switch (refine) {
763: case BV_ORTHOG_REFINE_NEVER:
764: case BV_ORTHOG_REFINE_IFNEEDED:
765: case BV_ORTHOG_REFINE_ALWAYS:
766: bv->orthog_ref = refine;
767: break;
768: default:769: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown refinement type");
770: }
771: if (eta == PETSC_DEFAULT) {
772: bv->orthog_eta = 0.7071;
773: } else {
774: if (eta <= 0.0 || eta > 1.0) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Invalid eta value");
775: bv->orthog_eta = eta;
776: }
777: return(0);
778: }
782: /*@C
783: BVGetOrthogonalization - Gets the orthogonalization settings from the BV object.
785: Not Collective
787: Input Parameter:
788: . bv - basis vectors context
790: Output Parameter:
791: + type - type of orthogonalization technique
792: . refine - type of refinement
793: - eta - parameter for selective refinement
795: Level: advanced
797: .seealso: BVOrthogonalizeColumn(), BVSetOrthogonalization(), BVOrthogType, BVOrthogRefineType798: @*/
799: PetscErrorCode BVGetOrthogonalization(BV bv,BVOrthogType *type,BVOrthogRefineType *refine,PetscReal *eta)800: {
803: if (type) *type = bv->orthog_type;
804: if (refine) *refine = bv->orthog_ref;
805: if (eta) *eta = bv->orthog_eta;
806: return(0);
807: }
811: /*@
812: BVGetColumn - Returns a Vec object that contains the entries of the
813: requested column of the basis vectors object.
815: Logically Collective on BV817: Input Parameters:
818: + bv - the basis vectors context
819: - j - the index of the requested column
821: Output Parameter:
822: . v - vector containing the jth column
824: Notes:
825: The returned Vec must be seen as a reference (not a copy) of the BV826: column, that is, modifying the Vec will change the BV entries as well.
828: The returned Vec must not be destroyed. BVRestoreColumn() must be
829: called when it is no longer needed. At most, two columns can be fetched,
830: that is, this function can only be called twice before the corresponding
831: BVRestoreColumn() is invoked.
833: A negative index j selects the i-th constraint, where i=-j. Constraints
834: should not be modified.
836: Level: beginner
838: .seealso: BVRestoreColumn(), BVInsertConstraints()
839: @*/
840: PetscErrorCode BVGetColumn(BV bv,PetscInt j,Vec *v)841: {
843: PetscInt l;
848: BVCheckSizes(bv,1);
850: 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);
851: if (j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"You requested column %D but only %D are available",j,bv->m);
852: 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);
853: l = BVAvailableVec;
854: if (l==-1) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Too many requested columns; you must call BVReleaseColumn for one of the previously fetched columns");
855: (*bv->ops->getcolumn)(bv,j,v);
856: bv->ci[l] = j;
857: PetscObjectStateGet((PetscObject)bv->cv[l],&bv->st[l]);
858: PetscObjectGetId((PetscObject)bv->cv[l],&bv->id[l]);
859: *v = bv->cv[l];
860: return(0);
861: }
865: /*@
866: BVRestoreColumn - Restore a column obtained with BVGetColumn().
868: Logically Collective on BV870: Input Parameters:
871: + bv - the basis vectors context
872: . j - the index of the column
873: - v - vector obtained with BVGetColumn()
875: Note:
876: The arguments must match the corresponding call to BVGetColumn().
878: Level: beginner
880: .seealso: BVGetColumn()
881: @*/
882: PetscErrorCode BVRestoreColumn(BV bv,PetscInt j,Vec *v)883: {
884: PetscErrorCode ierr;
885: PetscObjectId id;
886: PetscObjectState st;
887: PetscInt l;
892: BVCheckSizes(bv,1);
896: 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);
897: if (j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"You requested column %D but only %D are available",j,bv->m);
898: 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);
899: l = (j==bv->ci[0])? 0: 1;
900: PetscObjectGetId((PetscObject)*v,&id);
901: 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");
902: PetscObjectStateGet((PetscObject)*v,&st);
903: if (st!=bv->st[l]) {
904: PetscObjectStateIncrease((PetscObject)bv);
905: }
906: if (bv->ops->restorecolumn) {
907: (*bv->ops->restorecolumn)(bv,j,v);
908: } else bv->cv[l] = NULL;
909: bv->ci[l] = -bv->nc-1;
910: bv->st[l] = -1;
911: bv->id[l] = 0;
912: *v = NULL;
913: return(0);
914: }
918: /*@C
919: BVGetArray - Returns a pointer to a contiguous array that contains this
920: processor's portion of the BV data.
922: Logically Collective on BV924: Input Parameters:
925: . bv - the basis vectors context
927: Output Parameter:
928: . a - location to put pointer to the array
930: Notes:
931: BVRestoreArray() must be called when access to the array is no longer needed.
932: This operation may imply a data copy, for BV types that do not store
933: data contiguously in memory.
935: The pointer will normally point to the first entry of the first column,
936: but if the BV has constraints then these go before the regular columns.
938: Level: advanced
940: .seealso: BVRestoreArray(), BVInsertConstraints()
941: @*/
942: PetscErrorCode BVGetArray(BV bv,PetscScalar **a)943: {
949: BVCheckSizes(bv,1);
950: (*bv->ops->getarray)(bv,a);
951: return(0);
952: }
956: /*@C
957: BVRestoreArray - Restore the BV object after BVGetArray() has been called.
959: Logically Collective on BV961: Input Parameters:
962: + bv - the basis vectors context
963: - a - location of pointer to array obtained from BVGetArray()
965: Note:
966: This operation may imply a data copy, for BV types that do not store
967: data contiguously in memory.
969: Level: advanced
971: .seealso: BVGetColumn()
972: @*/
973: PetscErrorCode BVRestoreArray(BV bv,PetscScalar **a)974: {
980: BVCheckSizes(bv,1);
981: if (bv->ops->restorearray) {
982: (*bv->ops->restorearray)(bv,a);
983: }
984: if (a) *a = NULL;
985: PetscObjectStateIncrease((PetscObject)bv);
986: return(0);
987: }
991: /*@
992: BVGetVec - Creates a new Vec object with the same type and dimensions
993: as the columns of the basis vectors object.
995: Collective on BV997: Input Parameter:
998: . bv - the basis vectors context
1000: Output Parameter:
1001: . v - the new vector
1003: Note:
1004: The user is responsible of destroying the returned vector.
1006: Level: beginner
1007: @*/
1008: PetscErrorCode BVGetVec(BV bv,Vec *v)1009: {
1014: BVCheckSizes(bv,1);
1016: VecDuplicate(bv->t,v);
1017: return(0);
1018: }
1022: /*@
1023: BVDuplicate - Creates a new basis vector object of the same type and
1024: dimensions as an existing one.
1026: Collective on BV1028: Input Parameter:
1029: . V - basis vectors context
1031: Output Parameter:
1032: . W - location to put the new BV1034: Notes:
1035: The new BV has the same type and dimensions as V, and it shares the same
1036: template vector. Also, the inner product matrix and orthogonalization
1037: options are copied.
1039: BVDuplicate() DOES NOT COPY the entries, but rather allocates storage
1040: for the new basis vectors. Use BVCopy() to copy the contents.
1042: Level: intermediate
1044: .seealso: BVCreate(), BVSetSizesFromVec(), BVCopy()
1045: @*/
1046: PetscErrorCode BVDuplicate(BV V,BV *W)1047: {
1053: BVCheckSizes(V,1);
1055: BVCreate(PetscObjectComm((PetscObject)V),W);
1056: BVSetSizesFromVec(*W,V->t,V->m);
1057: BVSetType(*W,((PetscObject)V)->type_name);
1058: BVSetMatrix(*W,V->matrix,V->indef);
1059: BVSetOrthogonalization(*W,V->orthog_type,V->orthog_ref,V->orthog_eta);
1060: PetscObjectStateIncrease((PetscObject)*W);
1061: return(0);
1062: }
1066: /*@
1067: BVDuplicateResize - Creates a new basis vector object of the same type and
1068: dimensions as an existing one, but with possibly different number of columns.
1070: Collective on BV1072: Input Parameter:
1073: + V - basis vectors context
1074: - m - the new number of columns
1076: Output Parameter:
1077: . W - location to put the new BV1079: Note:
1080: This is equivalent of a call to BVDuplicate() followed by BVResize(). The
1081: contents of V are not copied to W.
1083: Level: intermediate
1085: .seealso: BVDuplicate(), BVResize()
1086: @*/
1087: PetscErrorCode BVDuplicateResize(BV V,PetscInt m,BV *W)1088: {
1094: BVCheckSizes(V,1);
1097: BVCreate(PetscObjectComm((PetscObject)V),W);
1098: BVSetSizesFromVec(*W,V->t,m);
1099: BVSetType(*W,((PetscObject)V)->type_name);
1100: BVSetMatrix(*W,V->matrix,V->indef);
1101: BVSetOrthogonalization(*W,V->orthog_type,V->orthog_ref,V->orthog_eta);
1102: PetscObjectStateIncrease((PetscObject)*W);
1103: return(0);
1104: }
1108: /*@
1109: BVCopy - Copies a basis vector object into another one, W <- V.
1111: Logically Collective on BV1113: Input Parameter:
1114: . V - basis vectors context
1116: Output Parameter:
1117: . W - the copy
1119: Note:
1120: Both V and W must be distributed in the same manner; local copies are
1121: done. Only active columns (excluding the leading ones) are copied.
1122: In the destination W, columns are overwritten starting from the leading ones.
1123: Constraints are not copied.
1125: Level: beginner
1127: .seealso: BVCopyVec(), BVCopyColumn(), BVDuplicate(), BVSetActiveColumns()
1128: @*/
1129: PetscErrorCode BVCopy(BV V,BV W)1130: {
1136: BVCheckSizes(V,1);
1139: BVCheckSizes(W,2);
1141: if (V->n!=W->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension V %D, W %D",V->n,W->n);
1142: 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);
1143: if (!V->n) return(0);
1145: PetscLogEventBegin(BV_Copy,V,W,0,0);
1146: if (V->indef && V->matrix && V->indef==W->indef && V->matrix==W->matrix) {
1147: /* copy signature */
1148: BV_AllocateSignature(W);
1149: PetscMemcpy(W->omega+W->nc+W->l,V->omega+V->nc+V->l,(V->k-V->l)*sizeof(PetscReal));
1150: }
1151: (*V->ops->copy)(V,W);
1152: PetscLogEventEnd(BV_Copy,V,W,0,0);
1153: return(0);
1154: }
1158: /*@
1159: BVCopyVec - Copies one of the columns of a basis vectors object into a Vec.
1161: Logically Collective on BV1163: Input Parameter:
1164: + V - basis vectors context
1165: - j - the column number to be copied
1167: Output Parameter:
1168: . w - the copied column
1170: Note:
1171: Both V and w must be distributed in the same manner; local copies are done.
1173: Level: beginner
1175: .seealso: BVCopy(), BVCopyColumn(), BVInsertVec()
1176: @*/
1177: PetscErrorCode BVCopyVec(BV V,PetscInt j,Vec w)1178: {
1180: PetscInt n,N;
1181: Vec z;
1186: BVCheckSizes(V,1);
1191: VecGetSize(w,&N);
1192: VecGetLocalSize(w,&n);
1193: 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);
1195: PetscLogEventBegin(BV_Copy,V,w,0,0);
1196: BVGetColumn(V,j,&z);
1197: VecCopy(z,w);
1198: BVRestoreColumn(V,j,&z);
1199: PetscLogEventEnd(BV_Copy,V,w,0,0);
1200: return(0);
1201: }
1205: /*@
1206: BVCopyColumn - Copies the values from one of the columns to another one.
1208: Logically Collective on BV1210: Input Parameter:
1211: + V - basis vectors context
1212: . j - the number of the source column
1213: - i - the number of the destination column
1215: Level: beginner
1217: .seealso: BVCopy(), BVCopyVec()
1218: @*/
1219: PetscErrorCode BVCopyColumn(BV V,PetscInt j,PetscInt i)1220: {
1222: Vec z,w;
1227: BVCheckSizes(V,1);
1230: if (j==i) return(0);
1232: PetscLogEventBegin(BV_Copy,V,0,0,0);
1233: if (V->omega) V->omega[i] = V->omega[j];
1234: BVGetColumn(V,j,&z);
1235: BVGetColumn(V,i,&w);
1236: VecCopy(z,w);
1237: BVRestoreColumn(V,j,&z);
1238: BVRestoreColumn(V,i,&w);
1239: PetscLogEventEnd(BV_Copy,V,0,0,0);
1240: return(0);
1241: }