1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2017, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: Basic BV routines
12: */
14: #include <slepc/private/bvimpl.h> /*I "slepcbv.h" I*/
16: PetscBool BVRegisterAllCalled = PETSC_FALSE;
17: PetscFunctionList BVList = 0;
19: /*@C
20: BVSetType - Selects the type for the BV object.
22: Logically Collective on BV 24: Input Parameter:
25: + bv - the basis vectors context
26: - type - a known type
28: Options Database Key:
29: . -bv_type <type> - Sets BV type
31: Level: intermediate
33: .seealso: BVGetType()
34: @*/
35: PetscErrorCode BVSetType(BV bv,BVType type) 36: {
37: PetscErrorCode ierr,(*r)(BV);
38: PetscBool match;
44: PetscObjectTypeCompare((PetscObject)bv,type,&match);
45: if (match) return(0);
47: PetscFunctionListFind(BVList,type,&r);
48: if (!r) SETERRQ1(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested BV type %s",type);
50: if (bv->ops->destroy) { (*bv->ops->destroy)(bv); }
51: PetscMemzero(bv->ops,sizeof(struct _BVOps));
53: PetscObjectChangeTypeName((PetscObject)bv,type);
54: if (bv->n < 0 && bv->N < 0) {
55: bv->ops->create = r;
56: } else {
57: PetscLogEventBegin(BV_Create,bv,0,0,0);
58: (*r)(bv);
59: PetscLogEventEnd(BV_Create,bv,0,0,0);
60: }
61: return(0);
62: }
64: /*@C
65: BVGetType - Gets the BV type name (as a string) from the BV context.
67: Not Collective
69: Input Parameter:
70: . bv - the basis vectors context
72: Output Parameter:
73: . name - name of the type of basis vectors
75: Level: intermediate
77: .seealso: BVSetType()
78: @*/
79: PetscErrorCode BVGetType(BV bv,BVType *type) 80: {
84: *type = ((PetscObject)bv)->type_name;
85: return(0);
86: }
88: /*@
89: BVSetSizes - Sets the local and global sizes, and the number of columns.
91: Collective on BV 93: Input Parameters:
94: + bv - the basis vectors
95: . n - the local size (or PETSC_DECIDE to have it set)
96: . N - the global size (or PETSC_DECIDE)
97: - m - the number of columns
99: Notes:
100: n and N cannot be both PETSC_DECIDE.
101: If one processor calls this with N of PETSC_DECIDE then all processors must,
102: otherwise the program will hang.
104: Level: beginner
106: .seealso: BVSetSizesFromVec(), BVGetSizes(), BVResize()
107: @*/
108: PetscErrorCode BVSetSizes(BV bv,PetscInt n,PetscInt N,PetscInt m)109: {
111: PetscInt ma;
117: if (N >= 0 && n > N) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Local size %D cannot be larger than global size %D",n,N);
118: if (m <= 0) SETERRQ1(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Number of columns %D must be positive",m);
119: if ((bv->n >= 0 || bv->N >= 0) && (bv->n != n || bv->N != N)) SETERRQ4(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Cannot change/reset vector sizes to %D local %D global after previously setting them to %D local %D global",n,N,bv->n,bv->N);
120: if (bv->m > 0 && bv->m != m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Cannot change the number of columns to %D after previously setting it to %D; use BVResize()",m,bv->m);
121: bv->n = n;
122: bv->N = N;
123: bv->m = m;
124: bv->k = m;
125: if (!bv->t) { /* create template vector and get actual dimensions */
126: VecCreate(PetscObjectComm((PetscObject)bv),&bv->t);
127: VecSetSizes(bv->t,bv->n,bv->N);
128: VecSetFromOptions(bv->t);
129: VecGetSize(bv->t,&bv->N);
130: VecGetLocalSize(bv->t,&bv->n);
131: if (bv->matrix) { /* check compatible dimensions of user-provided matrix */
132: MatGetLocalSize(bv->matrix,&ma,NULL);
133: if (bv->n!=ma) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Local dimension %D does not match that of matrix given at BVSetMatrix %D",bv->n,ma);
134: }
135: }
136: if (bv->ops->create) {
137: PetscLogEventBegin(BV_Create,bv,0,0,0);
138: (*bv->ops->create)(bv);
139: PetscLogEventEnd(BV_Create,bv,0,0,0);
140: bv->ops->create = 0;
141: bv->defersfo = PETSC_FALSE;
142: }
143: return(0);
144: }
146: /*@
147: BVSetSizesFromVec - Sets the local and global sizes, and the number of columns.
148: Local and global sizes are specified indirectly by passing a template vector.
150: Collective on BV152: Input Parameters:
153: + bv - the basis vectors
154: . t - the template vector
155: - m - the number of columns
157: Level: beginner
159: .seealso: BVSetSizes(), BVGetSizes(), BVResize()
160: @*/
161: PetscErrorCode BVSetSizesFromVec(BV bv,Vec t,PetscInt m)162: {
164: PetscInt ma;
171: if (m <= 0) SETERRQ1(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Number of columns %D must be positive",m);
172: if (bv->t) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Template vector was already set by a previous call to BVSetSizes/FromVec");
173: VecGetSize(t,&bv->N);
174: VecGetLocalSize(t,&bv->n);
175: if (bv->matrix) { /* check compatible dimensions of user-provided matrix */
176: MatGetLocalSize(bv->matrix,&ma,NULL);
177: if (bv->n!=ma) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Local dimension %D does not match that of matrix given at BVSetMatrix %D",bv->n,ma);
178: }
179: bv->m = m;
180: bv->k = m;
181: bv->t = t;
182: PetscObjectReference((PetscObject)t);
183: if (bv->ops->create) {
184: (*bv->ops->create)(bv);
185: bv->ops->create = 0;
186: bv->defersfo = PETSC_FALSE;
187: }
188: return(0);
189: }
191: /*@
192: BVGetSizes - Returns the local and global sizes, and the number of columns.
194: Not Collective
196: Input Parameter:
197: . bv - the basis vectors
199: Output Parameters:
200: + n - the local size
201: . N - the global size
202: - m - the number of columns
204: Note:
205: Normal usage requires that bv has already been given its sizes, otherwise
206: the call fails. However, this function can also be used to determine if
207: a BV object has been initialized completely (sizes and type). For this,
208: call with n=NULL and N=NULL, then a return value of m=0 indicates that
209: the BV object is not ready for use yet.
211: Level: beginner
213: .seealso: BVSetSizes(), BVSetSizesFromVec()
214: @*/
215: PetscErrorCode BVGetSizes(BV bv,PetscInt *n,PetscInt *N,PetscInt *m)216: {
218: if (!bv) {
219: if (m && !n && !N) *m = 0;
220: return(0);
221: }
223: if (n || N) BVCheckSizes(bv,1);
224: if (n) *n = bv->n;
225: if (N) *N = bv->N;
226: if (m) *m = bv->m;
227: if (m && !n && !N && !((PetscObject)bv)->type_name) *m = 0;
228: return(0);
229: }
231: /*@
232: BVSetNumConstraints - Set the number of constraints.
234: Logically Collective on BV236: Input Parameters:
237: + V - basis vectors
238: - nc - number of constraints
240: Notes:
241: This function sets the number of constraints to nc and marks all remaining
242: columns as regular. Normal user would call BVInsertConstraints() instead.
244: If nc is smaller than the previously set value, then some of the constraints
245: are discarded. In particular, using nc=0 removes all constraints preserving
246: the content of regular columns.
248: Level: developer
250: .seealso: BVInsertConstraints()
251: @*/
252: PetscErrorCode BVSetNumConstraints(BV V,PetscInt nc)253: {
255: PetscInt total,diff,i;
256: Vec x,y;
261: if (nc<0) SETERRQ1(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Number of constraints (given %D) cannot be negative",nc);
263: BVCheckSizes(V,1);
264: 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");
266: diff = nc-V->nc;
267: if (!diff) return(0);
268: total = V->nc+V->m;
269: if (total-nc<=0) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Not enough columns for the given nc value");
270: if (diff<0) { /* lessen constraints, shift contents of BV */
271: for (i=0;i<V->m;i++) {
272: BVGetColumn(V,i,&x);
273: BVGetColumn(V,i+diff,&y);
274: VecCopy(x,y);
275: BVRestoreColumn(V,i,&x);
276: BVRestoreColumn(V,i+diff,&y);
277: }
278: }
279: V->nc = nc;
280: V->ci[0] = -V->nc-1;
281: V->ci[1] = -V->nc-1;
282: V->l = 0;
283: V->m = total-nc;
284: V->k = V->m;
285: PetscObjectStateIncrease((PetscObject)V);
286: return(0);
287: }
289: /*@
290: BVGetNumConstraints - Returns the number of constraints.
292: Not Collective
294: Input Parameter:
295: . bv - the basis vectors
297: Output Parameters:
298: . nc - the number of constraints
300: Level: advanced
302: .seealso: BVGetSizes(), BVInsertConstraints()
303: @*/
304: PetscErrorCode BVGetNumConstraints(BV bv,PetscInt *nc)305: {
309: *nc = bv->nc;
310: return(0);
311: }
313: /*@
314: BVResize - Change the number of columns.
316: Collective on BV318: Input Parameters:
319: + bv - the basis vectors
320: . m - the new number of columns
321: - copy - a flag indicating whether current values should be kept
323: Note:
324: Internal storage is reallocated. If the copy flag is set to true, then
325: the contents are copied to the leading part of the new space.
327: Level: advanced
329: .seealso: BVSetSizes(), BVSetSizesFromVec()
330: @*/
331: PetscErrorCode BVResize(BV bv,PetscInt m,PetscBool copy)332: {
333: PetscErrorCode ierr;
334: PetscScalar *array;
335: const PetscScalar *omega;
336: Vec v;
343: if (m <= 0) SETERRQ1(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Number of columns %D must be positive",m);
344: if (bv->nc) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"Cannot resize a BV with constraints");
345: if (bv->m == m) return(0);
347: PetscLogEventBegin(BV_Create,bv,0,0,0);
348: (*bv->ops->resize)(bv,m,copy);
349: VecDestroy(&bv->buffer);
350: BVDestroy(&bv->cached);
351: PetscFree2(bv->h,bv->c);
352: if (bv->omega) {
353: if (bv->cuda) {
354: #if defined(PETSC_HAVE_VECCUDA)
355: VecCreateSeqCUDA(PETSC_COMM_SELF,m,&v);
356: #else
357: SETERRQ(PetscObjectComm((PetscObject)bv),1,"Something wrong happened");
358: #endif
359: } else {
360: VecCreateSeq(PETSC_COMM_SELF,m,&v);
361: }
362: PetscLogObjectParent((PetscObject)bv,(PetscObject)v);
363: if (copy) {
364: VecGetArray(v,&array);
365: VecGetArrayRead(bv->omega,&omega);
366: PetscMemcpy(array,omega,PetscMin(m,bv->m)*sizeof(PetscScalar));
367: VecRestoreArrayRead(bv->omega,&omega);
368: VecRestoreArray(v,&array);
369: } else {
370: VecSet(v,1.0);
371: }
372: VecDestroy(&bv->omega);
373: bv->omega = v;
374: }
375: bv->m = m;
376: bv->k = PetscMin(bv->k,m);
377: bv->l = PetscMin(bv->l,m);
378: PetscLogEventEnd(BV_Create,bv,0,0,0);
379: PetscObjectStateIncrease((PetscObject)bv);
380: return(0);
381: }
383: /*@
384: BVSetActiveColumns - Specify the columns that will be involved in operations.
386: Logically Collective on BV388: Input Parameters:
389: + bv - the basis vectors context
390: . l - number of leading columns
391: - k - number of active columns
393: Notes:
394: In operations such as BVMult() or BVDot(), only the first k columns are
395: considered. This is useful when the BV is filled from left to right, so
396: the last m-k columns do not have relevant information.
398: Also in operations such as BVMult() or BVDot(), the first l columns are
399: normally not included in the computation. See the manpage of each
400: operation.
402: In orthogonalization operations, the first l columns are treated
403: differently: they participate in the orthogonalization but the computed
404: coefficients are not stored.
406: Level: intermediate
408: .seealso: BVGetActiveColumns(), BVSetSizes()
409: @*/
410: PetscErrorCode BVSetActiveColumns(BV bv,PetscInt l,PetscInt k)411: {
416: BVCheckSizes(bv,1);
417: if (k==PETSC_DECIDE || k==PETSC_DEFAULT) {
418: bv->k = bv->m;
419: } else {
420: if (k<0 || k>bv->m) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of k. Must be between 0 and m");
421: bv->k = k;
422: }
423: if (l==PETSC_DECIDE || l==PETSC_DEFAULT) {
424: bv->l = 0;
425: } else {
426: if (l<0 || l>bv->k) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of l. Must be between 0 and k");
427: bv->l = l;
428: }
429: return(0);
430: }
432: /*@
433: BVGetActiveColumns - Returns the current active dimensions.
435: Not Collective
437: Input Parameter:
438: . bv - the basis vectors context
440: Output Parameter:
441: + l - number of leading columns
442: - k - number of active columns
444: Level: intermediate
446: .seealso: BVSetActiveColumns()
447: @*/
448: PetscErrorCode BVGetActiveColumns(BV bv,PetscInt *l,PetscInt *k)449: {
452: if (l) *l = bv->l;
453: if (k) *k = bv->k;
454: return(0);
455: }
457: /*@
458: BVSetMatrix - Specifies the inner product to be used in orthogonalization.
460: Collective on BV462: Input Parameters:
463: + bv - the basis vectors context
464: . B - a symmetric matrix (may be NULL)
465: - indef - a flag indicating if the matrix is indefinite
467: Notes:
468: This is used to specify a non-standard inner product, whose matrix
469: representation is given by B. Then, all inner products required during
470: orthogonalization are computed as (x,y)_B=y^H*B*x rather than the
471: standard form (x,y)=y^H*x.
473: Matrix B must be real symmetric (or complex Hermitian). A genuine inner
474: product requires that B is also positive (semi-)definite. However, we
475: also allow for an indefinite B (setting indef=PETSC_TRUE), in which
476: case the orthogonalization uses an indefinite inner product.
478: This affects operations BVDot(), BVNorm(), BVOrthogonalize(), and variants.
480: Setting B=NULL has the same effect as if the identity matrix was passed.
482: Level: advanced
484: .seealso: BVGetMatrix(), BVDot(), BVNorm(), BVOrthogonalize()
485: @*/
486: PetscErrorCode BVSetMatrix(BV bv,Mat B,PetscBool indef)487: {
489: PetscInt m,n;
494: if (B) {
496: MatGetLocalSize(B,&m,&n);
497: if (m!=n) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Matrix must be square");
498: if (bv->m && bv->n!=n) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension BV %D, Mat %D",bv->n,n);
499: }
500: MatDestroy(&bv->matrix);
501: if (B) PetscObjectReference((PetscObject)B);
502: bv->matrix = B;
503: bv->indef = indef;
504: if (B && !bv->Bx) {
505: MatCreateVecs(B,&bv->Bx,NULL);
506: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->Bx);
507: }
508: PetscObjectStateIncrease((PetscObject)bv);
509: return(0);
510: }
512: /*@
513: BVGetMatrix - Retrieves the matrix representation of the inner product.
515: Not collective, though a parallel Mat may be returned
517: Input Parameter:
518: . bv - the basis vectors context
520: Output Parameter:
521: + B - the matrix of the inner product (may be NULL)
522: - indef - the flag indicating if the matrix is indefinite
524: Level: advanced
526: .seealso: BVSetMatrix()
527: @*/
528: PetscErrorCode BVGetMatrix(BV bv,Mat *B,PetscBool *indef)529: {
532: if (B) *B = bv->matrix;
533: if (indef) *indef = bv->indef;
534: return(0);
535: }
537: /*@
538: BVApplyMatrix - Multiplies a vector by the matrix representation of the
539: inner product.
541: Neighbor-wise Collective on BV and Vec
543: Input Parameter:
544: + bv - the basis vectors context
545: - x - the vector
547: Output Parameter:
548: . y - the result
550: Note:
551: If no matrix was specified this function copies the vector.
553: Level: advanced
555: .seealso: BVSetMatrix(), BVApplyMatrixBV()
556: @*/
557: PetscErrorCode BVApplyMatrix(BV bv,Vec x,Vec y)558: {
565: if (bv->matrix) {
566: BV_IPMatMult(bv,x);
567: VecCopy(bv->Bx,y);
568: } else {
569: VecCopy(x,y);
570: }
571: return(0);
572: }
574: /*@
575: BVApplyMatrixBV - Multiplies the BV vectors by the matrix representation
576: of the inner product.
578: Neighbor-wise Collective on BV580: Input Parameter:
581: . X - the basis vectors context
583: Output Parameter:
584: . Y - the basis vectors to store the result (optional)
586: Note:
587: This function computes Y = B*X, where B is the matrix given with
588: BVSetMatrix(). This operation is computed as in BVMatMult().
589: If no matrix was specified, then it just copies Y = X.
591: If no Y is given, the result is stored internally in the cached BV.
593: Level: developer
595: .seealso: BVSetMatrix(), BVApplyMatrix(), BVMatMult(), BVGetCachedBV()
596: @*/
597: PetscErrorCode BVApplyMatrixBV(BV X,BV Y)598: {
603: if (Y) {
605: if (X->matrix) {
606: BVMatMult(X,X->matrix,Y);
607: } else {
608: BVCopy(X,Y);
609: }
610: } else {
611: BV_IPMatMultBV(X);
612: }
613: return(0);
614: }
616: /*@
617: BVGetCachedBV - Returns a BV object stored internally that holds the
618: result of B*X after a call to BVApplyMatrixBV().
620: Not collective
622: Input Parameter:
623: . bv - the basis vectors context
625: Output Parameter:
626: . cached - the cached BV628: Note:
629: The cached BV is created if not available previously.
631: Level: developer
633: .seealso: BVApplyMatrixBV()
634: @*/
635: PetscErrorCode BVGetCachedBV(BV bv,BV *cached)636: {
642: BVCheckSizes(bv,1);
643: if (!bv->cached) {
644: BVCreate(PetscObjectComm((PetscObject)bv),&bv->cached);
645: BVSetSizesFromVec(bv->cached,bv->t,bv->m);
646: BVSetType(bv->cached,((PetscObject)bv)->type_name);
647: BVSetOrthogonalization(bv->cached,bv->orthog_type,bv->orthog_ref,bv->orthog_eta,bv->orthog_block);
648: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->cached);
649: }
650: *cached = bv->cached;
651: return(0);
652: }
654: /*@
655: BVSetSignature - Sets the signature matrix to be used in orthogonalization.
657: Logically Collective on BV659: Input Parameter:
660: + bv - the basis vectors context
661: - omega - a vector representing the diagonal of the signature matrix
663: Note:
664: The signature matrix Omega = V'*B*V is relevant only for an indefinite B.
666: Level: developer
668: .seealso: BVSetMatrix(), BVGetSignature()
669: @*/
670: PetscErrorCode BVSetSignature(BV bv,Vec omega)671: {
672: PetscErrorCode ierr;
673: PetscInt i,n;
674: const PetscScalar *pomega;
675: PetscScalar *intern;
679: BVCheckSizes(bv,1);
683: VecGetSize(omega,&n);
684: if (n!=bv->k) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Vec argument has %D elements, should be %D",n,bv->k);
685: BV_AllocateSignature(bv);
686: if (bv->indef) {
687: VecGetArrayRead(omega,&pomega);
688: VecGetArray(bv->omega,&intern);
689: for (i=0;i<n;i++) intern[bv->nc+i] = pomega[i];
690: VecRestoreArray(bv->omega,&intern);
691: VecRestoreArrayRead(omega,&pomega);
692: } else {
693: PetscInfo(bv,"Ignoring signature because BV is not indefinite\n");
694: }
695: PetscObjectStateIncrease((PetscObject)bv);
696: return(0);
697: }
699: /*@
700: BVGetSignature - Retrieves the signature matrix from last orthogonalization.
702: Not collective
704: Input Parameter:
705: . bv - the basis vectors context
707: Output Parameter:
708: . omega - a vector representing the diagonal of the signature matrix
710: Note:
711: The signature matrix Omega = V'*B*V is relevant only for an indefinite B.
713: Level: developer
715: .seealso: BVSetMatrix(), BVSetSignature()
716: @*/
717: PetscErrorCode BVGetSignature(BV bv,Vec omega)718: {
719: PetscErrorCode ierr;
720: PetscInt i,n;
721: PetscScalar *pomega;
722: const PetscScalar *intern;
726: BVCheckSizes(bv,1);
730: VecGetSize(omega,&n);
731: if (n!=bv->k) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Vec argument has %D elements, should be %D",n,bv->k);
732: if (bv->indef && bv->omega) {
733: VecGetArray(omega,&pomega);
734: VecGetArrayRead(bv->omega,&intern);
735: for (i=0;i<n;i++) pomega[i] = intern[bv->nc+i];
736: VecRestoreArrayRead(bv->omega,&intern);
737: VecRestoreArray(omega,&pomega);
738: } else {
739: VecSet(omega,1.0);
740: }
741: return(0);
742: }
744: /*@
745: BVSetBufferVec - Attach a vector object to be used as buffer space for
746: several operations.
748: Collective on BV750: Input Parameters:
751: + bv - the basis vectors context)
752: - buffer - the vector
754: Notes:
755: Use BVGetBufferVec() to retrieve the vector (for example, to free it
756: at the end of the computations).
758: The vector must be sequential of length (nc+m)*m, where m is the number
759: of columns of bv and nc is the number of constraints.
761: Level: developer
763: .seealso: BVGetBufferVec(), BVSetSizes(), BVGetNumConstraints()
764: @*/
765: PetscErrorCode BVSetBufferVec(BV bv,Vec buffer)766: {
768: PetscInt ld,n;
769: PetscMPIInt size;
774: BVCheckSizes(bv,1);
775: VecGetSize(buffer,&n);
776: ld = bv->m+bv->nc;
777: if (n != ld*bv->m) SETERRQ1(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Buffer size must be %d",ld*bv->m);
778: MPI_Comm_size(PetscObjectComm((PetscObject)buffer),&size);
779: if (size>1) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Buffer must be a sequential vector");
780: 781: PetscObjectReference((PetscObject)buffer);
782: VecDestroy(&bv->buffer);
783: bv->buffer = buffer;
784: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->buffer);
785: return(0);
786: }
788: /*@
789: BVGetBufferVec - Obtain the buffer vector associated with the BV object.
791: Not Collective
793: Input Parameters:
794: . bv - the basis vectors context
796: Output Parameter:
797: . buffer - vector
799: Notes:
800: The vector is created if not available previously. It is a sequential vector
801: of length (nc+m)*m, where m is the number of columns of bv and nc is the number
802: of constraints.
804: Developer Notes:
805: The buffer vector is viewed as a column-major matrix with leading dimension
806: ld=nc+m, and m columns at most. In the most common usage, it has the structure
807: .vb
808: | | C |
809: |s|---|
810: | | H |
811: .ve
812: where H is an upper Hessenberg matrix of order m x (m-1), C contains coefficients
813: related to orthogonalization against constraints (first nc rows), and s is the
814: first column that contains scratch values computed during Gram-Schmidt
815: orthogonalization. In particular, BVDotColumn() and BVMultColumn() use s to
816: store the coefficients.
818: Level: developer
820: .seealso: BVSetBufferVec(), BVSetSizes(), BVGetNumConstraints(), BVDotColumn(), BVMultColumn()
821: @*/
822: PetscErrorCode BVGetBufferVec(BV bv,Vec *buffer)823: {
825: PetscInt ld;
830: BVCheckSizes(bv,1);
831: if (!bv->buffer) {
832: ld = bv->m+bv->nc;
833: VecCreate(PETSC_COMM_SELF,&bv->buffer);
834: VecSetSizes(bv->buffer,PETSC_DECIDE,ld*bv->m);
835: VecSetType(bv->buffer,((PetscObject)bv->t)->type_name);
836: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->buffer);
837: }
838: *buffer = bv->buffer;
839: return(0);
840: }
842: /*@
843: BVSetRandomContext - Sets the PetscRandom object associated with the BV,
844: to be used in operations that need random numbers.
846: Collective on BV848: Input Parameters:
849: + bv - the basis vectors context
850: - rand - the random number generator context
852: Level: advanced
854: .seealso: BVGetRandomContext(), BVSetRandom(), BVSetRandomColumn()
855: @*/
856: PetscErrorCode BVSetRandomContext(BV bv,PetscRandom rand)857: {
864: PetscObjectReference((PetscObject)rand);
865: PetscRandomDestroy(&bv->rand);
866: bv->rand = rand;
867: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->rand);
868: return(0);
869: }
871: /*@
872: BVGetRandomContext - Gets the PetscRandom object associated with the BV.
874: Not Collective
876: Input Parameter:
877: . bv - the basis vectors context
879: Output Parameter:
880: . rand - the random number generator context
882: Level: advanced
884: .seealso: BVSetRandomContext(), BVSetRandom(), BVSetRandomColumn()
885: @*/
886: PetscErrorCode BVGetRandomContext(BV bv,PetscRandom* rand)887: {
893: if (!bv->rand) {
894: PetscRandomCreate(PetscObjectComm((PetscObject)bv),&bv->rand);
895: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->rand);
896: if (bv->rrandom) {
897: PetscRandomSetSeed(bv->rand,0x12345678);
898: PetscRandomSeed(bv->rand);
899: }
900: }
901: *rand = bv->rand;
902: return(0);
903: }
905: /*@
906: BVSetFromOptions - Sets BV options from the options database.
908: Collective on BV910: Input Parameter:
911: . bv - the basis vectors context
913: Level: beginner
914: @*/
915: PetscErrorCode BVSetFromOptions(BV bv)916: {
917: PetscErrorCode ierr;
918: char type[256];
919: PetscBool flg1,flg2,flg3,flg4;
920: PetscReal r;
921: BVOrthogType otype;
922: BVOrthogRefineType orefine;
923: BVOrthogBlockType oblock;
927: BVRegisterAll();
928: PetscObjectOptionsBegin((PetscObject)bv);
929: PetscOptionsFList("-bv_type","Basis Vectors type","BVSetType",BVList,(char*)(((PetscObject)bv)->type_name?((PetscObject)bv)->type_name:BVSVEC),type,256,&flg1);
930: if (flg1) {
931: BVSetType(bv,type);
932: } else if (!((PetscObject)bv)->type_name) {
933: BVSetType(bv,BVSVEC);
934: }
936: otype = bv->orthog_type;
937: PetscOptionsEnum("-bv_orthog_type","Orthogonalization method","BVSetOrthogonalization",BVOrthogTypes,(PetscEnum)otype,(PetscEnum*)&otype,&flg1);
938: orefine = bv->orthog_ref;
939: PetscOptionsEnum("-bv_orthog_refine","Iterative refinement mode during orthogonalization","BVSetOrthogonalization",BVOrthogRefineTypes,(PetscEnum)orefine,(PetscEnum*)&orefine,&flg2);
940: r = bv->orthog_eta;
941: PetscOptionsReal("-bv_orthog_eta","Parameter of iterative refinement during orthogonalization","BVSetOrthogonalization",r,&r,&flg3);
942: oblock = bv->orthog_block;
943: PetscOptionsEnum("-bv_orthog_block","Block orthogonalization method","BVSetOrthogonalization",BVOrthogBlockTypes,(PetscEnum)oblock,(PetscEnum*)&oblock,&flg4);
944: if (flg1 || flg2 || flg3 || flg4) { BVSetOrthogonalization(bv,otype,orefine,r,oblock); }
946: PetscOptionsEnum("-bv_matmult","Method for BVMatMult","BVSetMatMultMethod",BVMatMultTypes,(PetscEnum)bv->vmm,(PetscEnum*)&bv->vmm,NULL);
948: /* undocumented option to generate random vectors that are independent of the number of processes */
949: PetscOptionsGetBool(NULL,NULL,"-bv_reproducible_random",&bv->rrandom,NULL);
951: if (bv->ops->create) bv->defersfo = PETSC_TRUE; /* defer call to setfromoptions */
952: else if (bv->ops->setfromoptions) {
953: (*bv->ops->setfromoptions)(PetscOptionsObject,bv);
954: }
955: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)bv);
956: PetscOptionsEnd();
958: if (!bv->rand) { BVGetRandomContext(bv,&bv->rand); }
959: PetscRandomSetFromOptions(bv->rand);
960: return(0);
961: }
963: /*@
964: BVSetOrthogonalization - Specifies the method used for the orthogonalization of
965: vectors (classical or modified Gram-Schmidt with or without refinement), and
966: for the block-orthogonalization (simultaneous orthogonalization of a set of
967: vectors).
969: Logically Collective on BV971: Input Parameters:
972: + bv - the basis vectors context
973: . type - the method of vector orthogonalization
974: . refine - type of refinement
975: . eta - parameter for selective refinement
976: - block - the method of block orthogonalization
978: Options Database Keys:
979: + -bv_orthog_type <type> - Where <type> is cgs for Classical Gram-Schmidt orthogonalization
980: (default) or mgs for Modified Gram-Schmidt orthogonalization
981: . -bv_orthog_refine <ref> - Where <ref> is one of never, ifneeded (default) or always
982: . -bv_orthog_eta <eta> - For setting the value of eta
983: - -bv_orthog_block <block> - Where <block> is the block-orthogonalization method
985: Notes:
986: The default settings work well for most problems.
988: The parameter eta should be a real value between 0 and 1 (or PETSC_DEFAULT).
989: The value of eta is used only when the refinement type is "ifneeded".
991: When using several processors, MGS is likely to result in bad scalability.
993: If the method set for block orthogonalization is GS, then the computation
994: is done column by column with the vector orthogonalization.
996: Level: advanced
998: .seealso: BVOrthogonalizeColumn(), BVGetOrthogonalization(), BVOrthogType, BVOrthogRefineType, BVOrthogBlockType999: @*/
1000: PetscErrorCode BVSetOrthogonalization(BV bv,BVOrthogType type,BVOrthogRefineType refine,PetscReal eta,BVOrthogBlockType block)1001: {
1008: switch (type) {
1009: case BV_ORTHOG_CGS:
1010: case BV_ORTHOG_MGS:
1011: bv->orthog_type = type;
1012: break;
1013: default:1014: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
1015: }
1016: switch (refine) {
1017: case BV_ORTHOG_REFINE_NEVER:
1018: case BV_ORTHOG_REFINE_IFNEEDED:
1019: case BV_ORTHOG_REFINE_ALWAYS:
1020: bv->orthog_ref = refine;
1021: break;
1022: default:1023: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown refinement type");
1024: }
1025: if (eta == PETSC_DEFAULT) {
1026: bv->orthog_eta = 0.7071;
1027: } else {
1028: if (eta <= 0.0 || eta > 1.0) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Invalid eta value");
1029: bv->orthog_eta = eta;
1030: }
1031: switch (block) {
1032: case BV_ORTHOG_BLOCK_GS:
1033: case BV_ORTHOG_BLOCK_CHOL:
1034: case BV_ORTHOG_BLOCK_TSQR:
1035: bv->orthog_block = block;
1036: break;
1037: default:1038: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown block orthogonalization type");
1039: }
1040: return(0);
1041: }
1043: /*@
1044: BVGetOrthogonalization - Gets the orthogonalization settings from the BV object.
1046: Not Collective
1048: Input Parameter:
1049: . bv - basis vectors context
1051: Output Parameter:
1052: + type - the method of vector orthogonalization
1053: . refine - type of refinement
1054: . eta - parameter for selective refinement
1055: - block - the method of block orthogonalization
1057: Level: advanced
1059: .seealso: BVOrthogonalizeColumn(), BVSetOrthogonalization(), BVOrthogType, BVOrthogRefineType, BVOrthogBlockType1060: @*/
1061: PetscErrorCode BVGetOrthogonalization(BV bv,BVOrthogType *type,BVOrthogRefineType *refine,PetscReal *eta,BVOrthogBlockType *block)1062: {
1065: if (type) *type = bv->orthog_type;
1066: if (refine) *refine = bv->orthog_ref;
1067: if (eta) *eta = bv->orthog_eta;
1068: if (block) *block = bv->orthog_block;
1069: return(0);
1070: }
1072: /*@
1073: BVSetMatMultMethod - Specifies the method used for the BVMatMult() operation.
1075: Logically Collective on BV1077: Input Parameters:
1078: + bv - the basis vectors context
1079: - method - the method for the BVMatMult() operation
1081: Options Database Keys:
1082: . -bv_matmult <meth> - choose one of the methods: vecs, mat, mat_save
1084: Note:
1085: Allowed values are:
1086: + BV_MATMULT_VECS - perform a matrix-vector multiply per each column
1087: . BV_MATMULT_MAT - carry out a MatMatMult() product with a dense matrix
1088: - BV_MATMULT_MAT_SAVE - call MatMatMult() and keep auxiliary matrices
1089: The default is BV_MATMULT_MAT.
1091: Level: advanced
1093: .seealso: BVMatMult(), BVGetMatMultMethod(), BVMatMultType1094: @*/
1095: PetscErrorCode BVSetMatMultMethod(BV bv,BVMatMultType method)1096: {
1100: switch (method) {
1101: case BV_MATMULT_VECS:
1102: case BV_MATMULT_MAT:
1103: case BV_MATMULT_MAT_SAVE:
1104: bv->vmm = method;
1105: break;
1106: default:1107: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown matmult method");
1108: }
1109: return(0);
1110: }
1112: /*@
1113: BVGetMatMultMethod - Gets the method used for the BVMatMult() operation.
1115: Not Collective
1117: Input Parameter:
1118: . bv - basis vectors context
1120: Output Parameter:
1121: . method - the method for the BVMatMult() operation
1123: Level: advanced
1125: .seealso: BVMatMult(), BVSetMatMultMethod(), BVMatMultType1126: @*/
1127: PetscErrorCode BVGetMatMultMethod(BV bv,BVMatMultType *method)1128: {
1132: *method = bv->vmm;
1133: return(0);
1134: }
1136: /*@
1137: BVGetColumn - Returns a Vec object that contains the entries of the
1138: requested column of the basis vectors object.
1140: Logically Collective on BV1142: Input Parameters:
1143: + bv - the basis vectors context
1144: - j - the index of the requested column
1146: Output Parameter:
1147: . v - vector containing the jth column
1149: Notes:
1150: The returned Vec must be seen as a reference (not a copy) of the BV1151: column, that is, modifying the Vec will change the BV entries as well.
1153: The returned Vec must not be destroyed. BVRestoreColumn() must be
1154: called when it is no longer needed. At most, two columns can be fetched,
1155: that is, this function can only be called twice before the corresponding
1156: BVRestoreColumn() is invoked.
1158: A negative index j selects the i-th constraint, where i=-j. Constraints
1159: should not be modified.
1161: Level: beginner
1163: .seealso: BVRestoreColumn(), BVInsertConstraints()
1164: @*/
1165: PetscErrorCode BVGetColumn(BV bv,PetscInt j,Vec *v)1166: {
1168: PetscInt l;
1173: BVCheckSizes(bv,1);
1175: 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);
1176: if (j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"You requested column %D but only %D are available",j,bv->m);
1177: 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);
1178: l = BVAvailableVec;
1179: 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");
1180: (*bv->ops->getcolumn)(bv,j,v);
1181: bv->ci[l] = j;
1182: PetscObjectStateGet((PetscObject)bv->cv[l],&bv->st[l]);
1183: PetscObjectGetId((PetscObject)bv->cv[l],&bv->id[l]);
1184: *v = bv->cv[l];
1185: return(0);
1186: }
1188: /*@
1189: BVRestoreColumn - Restore a column obtained with BVGetColumn().
1191: Logically Collective on BV1193: Input Parameters:
1194: + bv - the basis vectors context
1195: . j - the index of the column
1196: - v - vector obtained with BVGetColumn()
1198: Note:
1199: The arguments must match the corresponding call to BVGetColumn().
1201: Level: beginner
1203: .seealso: BVGetColumn()
1204: @*/
1205: PetscErrorCode BVRestoreColumn(BV bv,PetscInt j,Vec *v)1206: {
1207: PetscErrorCode ierr;
1208: PetscObjectId id;
1209: PetscObjectState st;
1210: PetscInt l;
1215: BVCheckSizes(bv,1);
1219: 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);
1220: if (j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"You requested column %D but only %D are available",j,bv->m);
1221: 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);
1222: l = (j==bv->ci[0])? 0: 1;
1223: PetscObjectGetId((PetscObject)*v,&id);
1224: 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");
1225: PetscObjectStateGet((PetscObject)*v,&st);
1226: if (st!=bv->st[l]) {
1227: PetscObjectStateIncrease((PetscObject)bv);
1228: }
1229: if (bv->ops->restorecolumn) {
1230: (*bv->ops->restorecolumn)(bv,j,v);
1231: } else bv->cv[l] = NULL;
1232: bv->ci[l] = -bv->nc-1;
1233: bv->st[l] = -1;
1234: bv->id[l] = 0;
1235: *v = NULL;
1236: return(0);
1237: }
1239: /*@C
1240: BVGetArray - Returns a pointer to a contiguous array that contains this
1241: processor's portion of the BV data.
1243: Logically Collective on BV1245: Input Parameters:
1246: . bv - the basis vectors context
1248: Output Parameter:
1249: . a - location to put pointer to the array
1251: Notes:
1252: BVRestoreArray() must be called when access to the array is no longer needed.
1253: This operation may imply a data copy, for BV types that do not store
1254: 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 BVGetArray(BV bv,PetscScalar **a)1264: {
1270: BVCheckSizes(bv,1);
1271: (*bv->ops->getarray)(bv,a);
1272: return(0);
1273: }
1275: /*@C
1276: BVRestoreArray - Restore the BV object after BVGetArray() has been called.
1278: Logically Collective on BV1280: Input Parameters:
1281: + bv - the basis vectors context
1282: - a - location of pointer to array obtained from BVGetArray()
1284: Note:
1285: This operation may imply a data copy, for BV types that do not store
1286: data contiguously in memory.
1288: Level: advanced
1290: .seealso: BVGetColumn()
1291: @*/
1292: PetscErrorCode BVRestoreArray(BV bv,PetscScalar **a)1293: {
1299: BVCheckSizes(bv,1);
1300: if (bv->ops->restorearray) {
1301: (*bv->ops->restorearray)(bv,a);
1302: }
1303: if (a) *a = NULL;
1304: PetscObjectStateIncrease((PetscObject)bv);
1305: return(0);
1306: }
1308: /*@C
1309: BVGetArrayRead - Returns a read-only pointer to a contiguous array that
1310: contains this processor's portion of the BV data.
1312: Not Collective
1314: Input Parameters:
1315: . bv - the basis vectors context
1317: Output Parameter:
1318: . a - location to put pointer to the array
1320: Notes:
1321: BVRestoreArrayRead() must be called when access to the array is no
1322: longer needed. This operation may imply a data copy, for BV types that
1323: do not store data contiguously in memory.
1325: The pointer will normally point to the first entry of the first column,
1326: but if the BV has constraints then these go before the regular columns.
1328: Level: advanced
1330: .seealso: BVRestoreArray(), BVInsertConstraints()
1331: @*/
1332: PetscErrorCode BVGetArrayRead(BV bv,const PetscScalar **a)1333: {
1339: BVCheckSizes(bv,1);
1340: (*bv->ops->getarrayread)(bv,a);
1341: return(0);
1342: }
1344: /*@C
1345: BVRestoreArrayRead - Restore the BV object after BVGetArrayRead() has
1346: been called.
1348: Not Collective
1350: Input Parameters:
1351: + bv - the basis vectors context
1352: - a - location of pointer to array obtained from BVGetArrayRead()
1354: Level: advanced
1356: .seealso: BVGetColumn()
1357: @*/
1358: PetscErrorCode BVRestoreArrayRead(BV bv,const PetscScalar **a)1359: {
1365: BVCheckSizes(bv,1);
1366: if (bv->ops->restorearrayread) {
1367: (*bv->ops->restorearrayread)(bv,a);
1368: }
1369: if (a) *a = NULL;
1370: return(0);
1371: }
1373: /*@
1374: BVCreateVec - Creates a new Vec object with the same type and dimensions
1375: as the columns of the basis vectors object.
1377: Collective on BV1379: Input Parameter:
1380: . bv - the basis vectors context
1382: Output Parameter:
1383: . v - the new vector
1385: Note:
1386: The user is responsible of destroying the returned vector.
1388: Level: beginner
1390: .seealso: BVCreateMat()
1391: @*/
1392: PetscErrorCode BVCreateVec(BV bv,Vec *v)1393: {
1398: BVCheckSizes(bv,1);
1400: VecDuplicate(bv->t,v);
1401: return(0);
1402: }
1404: /*@
1405: BVCreateMat - Creates a new Mat object of dense type and copies the contents
1406: of the BV object.
1408: Collective on BV1410: Input Parameter:
1411: . bv - the basis vectors context
1413: Output Parameter:
1414: . A - the new matrix
1416: Note:
1417: The user is responsible of destroying the returned matrix.
1419: Level: intermediate
1421: .seealso: BVCreateFromMat(), BVCreateVec()
1422: @*/
1423: PetscErrorCode BVCreateMat(BV bv,Mat *A)1424: {
1425: PetscErrorCode ierr;
1426: PetscScalar *aa;
1427: const PetscScalar *vv;
1431: BVCheckSizes(bv,1);
1434: MatCreateDense(PetscObjectComm((PetscObject)bv->t),bv->n,PETSC_DECIDE,bv->N,bv->m,NULL,A);
1435: MatDenseGetArray(*A,&aa);
1436: BVGetArrayRead(bv,&vv);
1437: PetscMemcpy(aa,vv,bv->m*bv->n*sizeof(PetscScalar));
1438: BVRestoreArrayRead(bv,&vv);
1439: MatDenseRestoreArray(*A,&aa);
1440: MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
1441: MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
1442: return(0);
1443: }
1445: PETSC_STATIC_INLINE PetscErrorCode BVDuplicate_Private(BV V,PetscInt m,BV *W)1446: {
1450: BVCreate(PetscObjectComm((PetscObject)V),W);
1451: BVSetSizesFromVec(*W,V->t,m);
1452: BVSetType(*W,((PetscObject)V)->type_name);
1453: BVSetMatrix(*W,V->matrix,V->indef);
1454: BVSetOrthogonalization(*W,V->orthog_type,V->orthog_ref,V->orthog_eta,V->orthog_block);
1455: (*W)->vmm = V->vmm;
1456: (*W)->rrandom = V->rrandom;
1457: if (V->ops->duplicate) { (*V->ops->duplicate)(V,W); }
1458: PetscObjectStateIncrease((PetscObject)*W);
1459: return(0);
1460: }
1462: /*@
1463: BVDuplicate - Creates a new basis vector object of the same type and
1464: dimensions as an existing one.
1466: Collective on BV1468: Input Parameter:
1469: . V - basis vectors context
1471: Output Parameter:
1472: . W - location to put the new BV1474: Notes:
1475: The new BV has the same type and dimensions as V, and it shares the same
1476: template vector. Also, the inner product matrix and orthogonalization
1477: options are copied.
1479: BVDuplicate() DOES NOT COPY the entries, but rather allocates storage
1480: for the new basis vectors. Use BVCopy() to copy the contents.
1482: Level: intermediate
1484: .seealso: BVDuplicateResize(), BVCreate(), BVSetSizesFromVec(), BVCopy()
1485: @*/
1486: PetscErrorCode BVDuplicate(BV V,BV *W)1487: {
1493: BVCheckSizes(V,1);
1495: BVDuplicate_Private(V,V->m,W);
1496: return(0);
1497: }
1499: /*@
1500: BVDuplicateResize - Creates a new basis vector object of the same type and
1501: dimensions as an existing one, but with possibly different number of columns.
1503: Collective on BV1505: Input Parameter:
1506: + V - basis vectors context
1507: - m - the new number of columns
1509: Output Parameter:
1510: . W - location to put the new BV1512: Note:
1513: This is equivalent of a call to BVDuplicate() followed by BVResize(). The
1514: contents of V are not copied to W.
1516: Level: intermediate
1518: .seealso: BVDuplicate(), BVResize()
1519: @*/
1520: PetscErrorCode BVDuplicateResize(BV V,PetscInt m,BV *W)1521: {
1527: BVCheckSizes(V,1);
1530: BVDuplicate_Private(V,m,W);
1531: return(0);
1532: }
1534: /*@
1535: BVCopy - Copies a basis vector object into another one, W <- V.
1537: Logically Collective on BV1539: Input Parameter:
1540: . V - basis vectors context
1542: Output Parameter:
1543: . W - the copy
1545: Note:
1546: Both V and W must be distributed in the same manner; local copies are
1547: done. Only active columns (excluding the leading ones) are copied.
1548: In the destination W, columns are overwritten starting from the leading ones.
1549: Constraints are not copied.
1551: Level: beginner
1553: .seealso: BVCopyVec(), BVCopyColumn(), BVDuplicate(), BVSetActiveColumns()
1554: @*/
1555: PetscErrorCode BVCopy(BV V,BV W)1556: {
1557: PetscErrorCode ierr;
1558: PetscScalar *womega;
1559: const PetscScalar *vomega;
1564: BVCheckSizes(V,1);
1567: BVCheckSizes(W,2);
1569: if (V->n!=W->n) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension V %D, W %D",V->n,W->n);
1570: if (V->k-V->l>W->m-W->l) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"W has %D non-leading columns, not enough to store %D columns",W->m-W->l,V->k-V->l);
1571: if (!V->n) return(0);
1573: PetscLogEventBegin(BV_Copy,V,W,0,0);
1574: if (V->indef && V->matrix && V->indef==W->indef && V->matrix==W->matrix) {
1575: /* copy signature */
1576: BV_AllocateSignature(W);
1577: VecGetArrayRead(V->omega,&vomega);
1578: VecGetArray(W->omega,&womega);
1579: PetscMemcpy(womega+W->nc+W->l,vomega+V->nc+V->l,(V->k-V->l)*sizeof(PetscScalar));
1580: VecRestoreArray(W->omega,&womega);
1581: VecRestoreArrayRead(V->omega,&vomega);
1582: }
1583: (*V->ops->copy)(V,W);
1584: PetscLogEventEnd(BV_Copy,V,W,0,0);
1585: PetscObjectStateIncrease((PetscObject)W);
1586: return(0);
1587: }
1589: /*@
1590: BVCopyVec - Copies one of the columns of a basis vectors object into a Vec.
1592: Logically Collective on BV1594: Input Parameter:
1595: + V - basis vectors context
1596: - j - the column number to be copied
1598: Output Parameter:
1599: . w - the copied column
1601: Note:
1602: Both V and w must be distributed in the same manner; local copies are done.
1604: Level: beginner
1606: .seealso: BVCopy(), BVCopyColumn(), BVInsertVec()
1607: @*/
1608: PetscErrorCode BVCopyVec(BV V,PetscInt j,Vec w)1609: {
1611: PetscInt n,N;
1612: Vec z;
1617: BVCheckSizes(V,1);
1622: VecGetSize(w,&N);
1623: VecGetLocalSize(w,&n);
1624: 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);
1626: PetscLogEventBegin(BV_Copy,V,w,0,0);
1627: BVGetColumn(V,j,&z);
1628: VecCopy(z,w);
1629: BVRestoreColumn(V,j,&z);
1630: PetscLogEventEnd(BV_Copy,V,w,0,0);
1631: return(0);
1632: }
1634: /*@
1635: BVCopyColumn - Copies the values from one of the columns to another one.
1637: Logically Collective on BV1639: Input Parameter:
1640: + V - basis vectors context
1641: . j - the number of the source column
1642: - i - the number of the destination column
1644: Level: beginner
1646: .seealso: BVCopy(), BVCopyVec()
1647: @*/
1648: PetscErrorCode BVCopyColumn(BV V,PetscInt j,PetscInt i)1649: {
1651: Vec z,w;
1652: PetscScalar *omega;
1657: BVCheckSizes(V,1);
1660: if (j==i) return(0);
1662: PetscLogEventBegin(BV_Copy,V,0,0,0);
1663: if (V->omega) {
1664: VecGetArray(V->omega,&omega);
1665: omega[i] = omega[j];
1666: VecRestoreArray(V->omega,&omega);
1667: }
1668: BVGetColumn(V,j,&z);
1669: BVGetColumn(V,i,&w);
1670: VecCopy(z,w);
1671: BVRestoreColumn(V,j,&z);
1672: BVRestoreColumn(V,i,&w);
1673: PetscLogEventEnd(BV_Copy,V,0,0,0);
1674: PetscObjectStateIncrease((PetscObject)V);
1675: return(0);
1676: }