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: return(0);
380: }
382: /*@
383: BVSetActiveColumns - Specify the columns that will be involved in operations.
385: Logically Collective on BV387: Input Parameters:
388: + bv - the basis vectors context
389: . l - number of leading columns
390: - k - number of active columns
392: Notes:
393: In operations such as BVMult() or BVDot(), only the first k columns are
394: considered. This is useful when the BV is filled from left to right, so
395: the last m-k columns do not have relevant information.
397: Also in operations such as BVMult() or BVDot(), the first l columns are
398: normally not included in the computation. See the manpage of each
399: operation.
401: In orthogonalization operations, the first l columns are treated
402: differently: they participate in the orthogonalization but the computed
403: coefficients are not stored.
405: Level: intermediate
407: .seealso: BVGetActiveColumns(), BVSetSizes()
408: @*/
409: PetscErrorCode BVSetActiveColumns(BV bv,PetscInt l,PetscInt k)410: {
415: BVCheckSizes(bv,1);
416: if (k==PETSC_DECIDE || k==PETSC_DEFAULT) {
417: bv->k = bv->m;
418: } else {
419: if (k<0 || k>bv->m) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of k. Must be between 0 and m");
420: bv->k = k;
421: }
422: if (l==PETSC_DECIDE || l==PETSC_DEFAULT) {
423: bv->l = 0;
424: } else {
425: if (l<0 || l>bv->k) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of l. Must be between 0 and k");
426: bv->l = l;
427: }
428: return(0);
429: }
431: /*@
432: BVGetActiveColumns - Returns the current active dimensions.
434: Not Collective
436: Input Parameter:
437: . bv - the basis vectors context
439: Output Parameter:
440: + l - number of leading columns
441: - k - number of active columns
443: Level: intermediate
445: .seealso: BVSetActiveColumns()
446: @*/
447: PetscErrorCode BVGetActiveColumns(BV bv,PetscInt *l,PetscInt *k)448: {
451: if (l) *l = bv->l;
452: if (k) *k = bv->k;
453: return(0);
454: }
456: /*@
457: BVSetMatrix - Specifies the inner product to be used in orthogonalization.
459: Collective on BV461: Input Parameters:
462: + bv - the basis vectors context
463: . B - a symmetric matrix (may be NULL)
464: - indef - a flag indicating if the matrix is indefinite
466: Notes:
467: This is used to specify a non-standard inner product, whose matrix
468: representation is given by B. Then, all inner products required during
469: orthogonalization are computed as (x,y)_B=y^H*B*x rather than the
470: standard form (x,y)=y^H*x.
472: Matrix B must be real symmetric (or complex Hermitian). A genuine inner
473: product requires that B is also positive (semi-)definite. However, we
474: also allow for an indefinite B (setting indef=PETSC_TRUE), in which
475: case the orthogonalization uses an indefinite inner product.
477: This affects operations BVDot(), BVNorm(), BVOrthogonalize(), and variants.
479: Setting B=NULL has the same effect as if the identity matrix was passed.
481: Level: advanced
483: .seealso: BVGetMatrix(), BVDot(), BVNorm(), BVOrthogonalize()
484: @*/
485: PetscErrorCode BVSetMatrix(BV bv,Mat B,PetscBool indef)486: {
488: PetscInt m,n;
493: if (B) {
495: MatGetLocalSize(B,&m,&n);
496: if (m!=n) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Matrix must be square");
497: if (bv->m && bv->n!=n) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension BV %D, Mat %D",bv->n,n);
498: }
499: MatDestroy(&bv->matrix);
500: if (B) PetscObjectReference((PetscObject)B);
501: bv->matrix = B;
502: bv->indef = indef;
503: if (B && !bv->Bx) {
504: MatCreateVecs(B,&bv->Bx,NULL);
505: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->Bx);
506: }
507: return(0);
508: }
510: /*@
511: BVGetMatrix - Retrieves the matrix representation of the inner product.
513: Not collective, though a parallel Mat may be returned
515: Input Parameter:
516: . bv - the basis vectors context
518: Output Parameter:
519: + B - the matrix of the inner product (may be NULL)
520: - indef - the flag indicating if the matrix is indefinite
522: Level: advanced
524: .seealso: BVSetMatrix()
525: @*/
526: PetscErrorCode BVGetMatrix(BV bv,Mat *B,PetscBool *indef)527: {
530: if (B) *B = bv->matrix;
531: if (indef) *indef = bv->indef;
532: return(0);
533: }
535: /*@
536: BVApplyMatrix - Multiplies a vector by the matrix representation of the
537: inner product.
539: Neighbor-wise Collective on BV and Vec
541: Input Parameter:
542: + bv - the basis vectors context
543: - x - the vector
545: Output Parameter:
546: . y - the result
548: Note:
549: If no matrix was specified this function copies the vector.
551: Level: advanced
553: .seealso: BVSetMatrix(), BVApplyMatrixBV()
554: @*/
555: PetscErrorCode BVApplyMatrix(BV bv,Vec x,Vec y)556: {
563: if (bv->matrix) {
564: BV_IPMatMult(bv,x);
565: VecCopy(bv->Bx,y);
566: } else {
567: VecCopy(x,y);
568: }
569: return(0);
570: }
572: /*@
573: BVApplyMatrixBV - Multiplies the BV vectors by the matrix representation
574: of the inner product.
576: Neighbor-wise Collective on BV578: Input Parameter:
579: . X - the basis vectors context
581: Output Parameter:
582: . Y - the basis vectors to store the result (optional)
584: Note:
585: This function computes Y = B*X, where B is the matrix given with
586: BVSetMatrix(). This operation is computed as in BVMatMult().
587: If no matrix was specified, then it just copies Y = X.
589: If no Y is given, the result is stored internally in the cached BV.
591: Level: developer
593: .seealso: BVSetMatrix(), BVApplyMatrix(), BVMatMult(), BVGetCachedBV()
594: @*/
595: PetscErrorCode BVApplyMatrixBV(BV X,BV Y)596: {
601: if (Y) {
603: if (X->matrix) {
604: BVMatMult(X,X->matrix,Y);
605: } else {
606: BVCopy(X,Y);
607: }
608: } else {
609: BV_IPMatMultBV(X);
610: }
611: return(0);
612: }
614: /*@
615: BVGetCachedBV - Returns a BV object stored internally that holds the
616: result of B*X after a call to BVApplyMatrixBV().
618: Not collective
620: Input Parameter:
621: . bv - the basis vectors context
623: Output Parameter:
624: . cached - the cached BV626: Note:
627: The cached BV is created if not available previously.
629: Level: developer
631: .seealso: BVApplyMatrixBV()
632: @*/
633: PetscErrorCode BVGetCachedBV(BV bv,BV *cached)634: {
640: BVCheckSizes(bv,1);
641: if (!bv->cached) {
642: BVCreate(PetscObjectComm((PetscObject)bv),&bv->cached);
643: BVSetSizesFromVec(bv->cached,bv->t,bv->m);
644: BVSetType(bv->cached,((PetscObject)bv)->type_name);
645: BVSetOrthogonalization(bv->cached,bv->orthog_type,bv->orthog_ref,bv->orthog_eta,bv->orthog_block);
646: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->cached);
647: }
648: *cached = bv->cached;
649: return(0);
650: }
652: /*@
653: BVSetSignature - Sets the signature matrix to be used in orthogonalization.
655: Logically Collective on BV657: Input Parameter:
658: + bv - the basis vectors context
659: - omega - a vector representing the diagonal of the signature matrix
661: Note:
662: The signature matrix Omega = V'*B*V is relevant only for an indefinite B.
664: Level: developer
666: .seealso: BVSetMatrix(), BVGetSignature()
667: @*/
668: PetscErrorCode BVSetSignature(BV bv,Vec omega)669: {
670: PetscErrorCode ierr;
671: PetscInt i,n;
672: const PetscScalar *pomega;
673: PetscScalar *intern;
677: BVCheckSizes(bv,1);
681: VecGetSize(omega,&n);
682: if (n!=bv->k) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Vec argument has %D elements, should be %D",n,bv->k);
683: BV_AllocateSignature(bv);
684: if (bv->indef) {
685: VecGetArrayRead(omega,&pomega);
686: VecGetArray(bv->omega,&intern);
687: for (i=0;i<n;i++) intern[bv->nc+i] = pomega[i];
688: VecRestoreArray(bv->omega,&intern);
689: VecRestoreArrayRead(omega,&pomega);
690: } else {
691: PetscInfo(bv,"Ignoring signature because BV is not indefinite\n");
692: }
693: return(0);
694: }
696: /*@
697: BVGetSignature - Retrieves the signature matrix from last orthogonalization.
699: Not collective
701: Input Parameter:
702: . bv - the basis vectors context
704: Output Parameter:
705: . omega - a vector representing the diagonal of the signature matrix
707: Note:
708: The signature matrix Omega = V'*B*V is relevant only for an indefinite B.
710: Level: developer
712: .seealso: BVSetMatrix(), BVSetSignature()
713: @*/
714: PetscErrorCode BVGetSignature(BV bv,Vec omega)715: {
716: PetscErrorCode ierr;
717: PetscInt i,n;
718: PetscScalar *pomega;
719: const PetscScalar *intern;
723: BVCheckSizes(bv,1);
727: VecGetSize(omega,&n);
728: if (n!=bv->k) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Vec argument has %D elements, should be %D",n,bv->k);
729: if (bv->indef && bv->omega) {
730: VecGetArray(omega,&pomega);
731: VecGetArrayRead(bv->omega,&intern);
732: for (i=0;i<n;i++) pomega[i] = intern[bv->nc+i];
733: VecRestoreArrayRead(bv->omega,&intern);
734: VecRestoreArray(omega,&pomega);
735: } else {
736: VecSet(omega,1.0);
737: }
738: return(0);
739: }
741: /*@
742: BVSetBufferVec - Attach a vector object to be used as buffer space for
743: several operations.
745: Collective on BV747: Input Parameters:
748: + bv - the basis vectors context)
749: - buffer - the vector
751: Notes:
752: Use BVGetBufferVec() to retrieve the vector (for example, to free it
753: at the end of the computations).
755: The vector must be sequential of length (nc+m)*m, where m is the number
756: of columns of bv and nc is the number of constraints.
758: Level: developer
760: .seealso: BVGetBufferVec(), BVSetSizes(), BVGetNumConstraints()
761: @*/
762: PetscErrorCode BVSetBufferVec(BV bv,Vec buffer)763: {
765: PetscInt ld,n;
766: PetscMPIInt size;
771: BVCheckSizes(bv,1);
772: VecGetSize(buffer,&n);
773: ld = bv->m+bv->nc;
774: if (n != ld*bv->m) SETERRQ1(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_SIZ,"Buffer size must be %d",ld*bv->m);
775: MPI_Comm_size(PetscObjectComm((PetscObject)buffer),&size);
776: if (size>1) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Buffer must be a sequential vector");
777: 778: PetscObjectReference((PetscObject)buffer);
779: VecDestroy(&bv->buffer);
780: bv->buffer = buffer;
781: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->buffer);
782: return(0);
783: }
785: /*@
786: BVGetBufferVec - Obtain the buffer vector associated with the BV object.
788: Not Collective
790: Input Parameters:
791: . bv - the basis vectors context
793: Output Parameter:
794: . buffer - vector
796: Notes:
797: The vector is created if not available previously. It is a sequential vector
798: of length (nc+m)*m, where m is the number of columns of bv and nc is the number
799: of constraints.
801: Developer Notes:
802: The buffer vector is viewed as a column-major matrix with leading dimension
803: ld=nc+m, and m columns at most. In the most common usage, it has the structure
804: .vb
805: | | C |
806: |s|---|
807: | | H |
808: .ve
809: where H is an upper Hessenberg matrix of order m x (m-1), C contains coefficients
810: related to orthogonalization against constraints (first nc rows), and s is the
811: first column that contains scratch values computed during Gram-Schmidt
812: orthogonalization. In particular, BVDotColumn() and BVMultColumn() use s to
813: store the coefficients.
815: Level: developer
817: .seealso: BVSetBufferVec(), BVSetSizes(), BVGetNumConstraints(), BVDotColumn(), BVMultColumn()
818: @*/
819: PetscErrorCode BVGetBufferVec(BV bv,Vec *buffer)820: {
822: PetscInt ld;
827: BVCheckSizes(bv,1);
828: if (!bv->buffer) {
829: ld = bv->m+bv->nc;
830: VecCreate(PETSC_COMM_SELF,&bv->buffer);
831: VecSetSizes(bv->buffer,PETSC_DECIDE,ld*bv->m);
832: VecSetType(bv->buffer,((PetscObject)bv->t)->type_name);
833: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->buffer);
834: }
835: *buffer = bv->buffer;
836: return(0);
837: }
839: /*@
840: BVSetRandomContext - Sets the PetscRandom object associated with the BV,
841: to be used in operations that need random numbers.
843: Collective on BV845: Input Parameters:
846: + bv - the basis vectors context
847: - rand - the random number generator context
849: Level: advanced
851: .seealso: BVGetRandomContext(), BVSetRandom(), BVSetRandomColumn()
852: @*/
853: PetscErrorCode BVSetRandomContext(BV bv,PetscRandom rand)854: {
861: PetscObjectReference((PetscObject)rand);
862: PetscRandomDestroy(&bv->rand);
863: bv->rand = rand;
864: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->rand);
865: return(0);
866: }
868: /*@
869: BVGetRandomContext - Gets the PetscRandom object associated with the BV.
871: Not Collective
873: Input Parameter:
874: . bv - the basis vectors context
876: Output Parameter:
877: . rand - the random number generator context
879: Level: advanced
881: .seealso: BVSetRandomContext(), BVSetRandom(), BVSetRandomColumn()
882: @*/
883: PetscErrorCode BVGetRandomContext(BV bv,PetscRandom* rand)884: {
890: if (!bv->rand) {
891: PetscRandomCreate(PetscObjectComm((PetscObject)bv),&bv->rand);
892: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->rand);
893: if (bv->rrandom) {
894: PetscRandomSetSeed(bv->rand,0x12345678);
895: PetscRandomSeed(bv->rand);
896: }
897: }
898: *rand = bv->rand;
899: return(0);
900: }
902: /*@
903: BVSetFromOptions - Sets BV options from the options database.
905: Collective on BV907: Input Parameter:
908: . bv - the basis vectors context
910: Level: beginner
911: @*/
912: PetscErrorCode BVSetFromOptions(BV bv)913: {
914: PetscErrorCode ierr;
915: char type[256];
916: PetscBool flg1,flg2,flg3,flg4;
917: PetscReal r;
918: BVOrthogType otype;
919: BVOrthogRefineType orefine;
920: BVOrthogBlockType oblock;
924: BVRegisterAll();
925: PetscObjectOptionsBegin((PetscObject)bv);
926: PetscOptionsFList("-bv_type","Basis Vectors type","BVSetType",BVList,(char*)(((PetscObject)bv)->type_name?((PetscObject)bv)->type_name:BVSVEC),type,256,&flg1);
927: if (flg1) {
928: BVSetType(bv,type);
929: } else if (!((PetscObject)bv)->type_name) {
930: BVSetType(bv,BVSVEC);
931: }
933: otype = bv->orthog_type;
934: PetscOptionsEnum("-bv_orthog_type","Orthogonalization method","BVSetOrthogonalization",BVOrthogTypes,(PetscEnum)otype,(PetscEnum*)&otype,&flg1);
935: orefine = bv->orthog_ref;
936: PetscOptionsEnum("-bv_orthog_refine","Iterative refinement mode during orthogonalization","BVSetOrthogonalization",BVOrthogRefineTypes,(PetscEnum)orefine,(PetscEnum*)&orefine,&flg2);
937: r = bv->orthog_eta;
938: PetscOptionsReal("-bv_orthog_eta","Parameter of iterative refinement during orthogonalization","BVSetOrthogonalization",r,&r,&flg3);
939: oblock = bv->orthog_block;
940: PetscOptionsEnum("-bv_orthog_block","Block orthogonalization method","BVSetOrthogonalization",BVOrthogBlockTypes,(PetscEnum)oblock,(PetscEnum*)&oblock,&flg4);
941: if (flg1 || flg2 || flg3 || flg4) { BVSetOrthogonalization(bv,otype,orefine,r,oblock); }
943: PetscOptionsEnum("-bv_matmult","Method for BVMatMult","BVSetMatMultMethod",BVMatMultTypes,(PetscEnum)bv->vmm,(PetscEnum*)&bv->vmm,NULL);
945: /* undocumented option to generate random vectors that are independent of the number of processes */
946: PetscOptionsGetBool(NULL,NULL,"-bv_reproducible_random",&bv->rrandom,NULL);
948: if (bv->ops->create) bv->defersfo = PETSC_TRUE; /* defer call to setfromoptions */
949: else if (bv->ops->setfromoptions) {
950: (*bv->ops->setfromoptions)(PetscOptionsObject,bv);
951: }
952: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)bv);
953: PetscOptionsEnd();
955: if (!bv->rand) { BVGetRandomContext(bv,&bv->rand); }
956: PetscRandomSetFromOptions(bv->rand);
957: return(0);
958: }
960: /*@
961: BVSetOrthogonalization - Specifies the method used for the orthogonalization of
962: vectors (classical or modified Gram-Schmidt with or without refinement), and
963: for the block-orthogonalization (simultaneous orthogonalization of a set of
964: vectors).
966: Logically Collective on BV968: Input Parameters:
969: + bv - the basis vectors context
970: . type - the method of vector orthogonalization
971: . refine - type of refinement
972: . eta - parameter for selective refinement
973: - block - the method of block orthogonalization
975: Options Database Keys:
976: + -bv_orthog_type <type> - Where <type> is cgs for Classical Gram-Schmidt orthogonalization
977: (default) or mgs for Modified Gram-Schmidt orthogonalization
978: . -bv_orthog_refine <ref> - Where <ref> is one of never, ifneeded (default) or always
979: . -bv_orthog_eta <eta> - For setting the value of eta
980: - -bv_orthog_block <block> - Where <block> is the block-orthogonalization method
982: Notes:
983: The default settings work well for most problems.
985: The parameter eta should be a real value between 0 and 1 (or PETSC_DEFAULT).
986: The value of eta is used only when the refinement type is "ifneeded".
988: When using several processors, MGS is likely to result in bad scalability.
990: If the method set for block orthogonalization is GS, then the computation
991: is done column by column with the vector orthogonalization.
993: Level: advanced
995: .seealso: BVOrthogonalizeColumn(), BVGetOrthogonalization(), BVOrthogType, BVOrthogRefineType, BVOrthogBlockType996: @*/
997: PetscErrorCode BVSetOrthogonalization(BV bv,BVOrthogType type,BVOrthogRefineType refine,PetscReal eta,BVOrthogBlockType block)998: {
1005: switch (type) {
1006: case BV_ORTHOG_CGS:
1007: case BV_ORTHOG_MGS:
1008: bv->orthog_type = type;
1009: break;
1010: default:1011: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
1012: }
1013: switch (refine) {
1014: case BV_ORTHOG_REFINE_NEVER:
1015: case BV_ORTHOG_REFINE_IFNEEDED:
1016: case BV_ORTHOG_REFINE_ALWAYS:
1017: bv->orthog_ref = refine;
1018: break;
1019: default:1020: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown refinement type");
1021: }
1022: if (eta == PETSC_DEFAULT) {
1023: bv->orthog_eta = 0.7071;
1024: } else {
1025: if (eta <= 0.0 || eta > 1.0) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Invalid eta value");
1026: bv->orthog_eta = eta;
1027: }
1028: switch (block) {
1029: case BV_ORTHOG_BLOCK_GS:
1030: case BV_ORTHOG_BLOCK_CHOL:
1031: case BV_ORTHOG_BLOCK_TSQR:
1032: bv->orthog_block = block;
1033: break;
1034: default:1035: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown block orthogonalization type");
1036: }
1037: return(0);
1038: }
1040: /*@
1041: BVGetOrthogonalization - Gets the orthogonalization settings from the BV object.
1043: Not Collective
1045: Input Parameter:
1046: . bv - basis vectors context
1048: Output Parameter:
1049: + type - the method of vector orthogonalization
1050: . refine - type of refinement
1051: . eta - parameter for selective refinement
1052: - block - the method of block orthogonalization
1054: Level: advanced
1056: .seealso: BVOrthogonalizeColumn(), BVSetOrthogonalization(), BVOrthogType, BVOrthogRefineType, BVOrthogBlockType1057: @*/
1058: PetscErrorCode BVGetOrthogonalization(BV bv,BVOrthogType *type,BVOrthogRefineType *refine,PetscReal *eta,BVOrthogBlockType *block)1059: {
1062: if (type) *type = bv->orthog_type;
1063: if (refine) *refine = bv->orthog_ref;
1064: if (eta) *eta = bv->orthog_eta;
1065: if (block) *block = bv->orthog_block;
1066: return(0);
1067: }
1069: /*@
1070: BVSetMatMultMethod - Specifies the method used for the BVMatMult() operation.
1072: Logically Collective on BV1074: Input Parameters:
1075: + bv - the basis vectors context
1076: - method - the method for the BVMatMult() operation
1078: Options Database Keys:
1079: . -bv_matmult <meth> - choose one of the methods: vecs, mat, mat_save
1081: Note:
1082: Allowed values are:
1083: + BV_MATMULT_VECS - perform a matrix-vector multiply per each column
1084: . BV_MATMULT_MAT - carry out a MatMatMult() product with a dense matrix
1085: - BV_MATMULT_MAT_SAVE - call MatMatMult() and keep auxiliary matrices
1086: The default is BV_MATMULT_MAT.
1088: Level: advanced
1090: .seealso: BVMatMult(), BVGetMatMultMethod(), BVMatMultType1091: @*/
1092: PetscErrorCode BVSetMatMultMethod(BV bv,BVMatMultType method)1093: {
1097: switch (method) {
1098: case BV_MATMULT_VECS:
1099: case BV_MATMULT_MAT:
1100: case BV_MATMULT_MAT_SAVE:
1101: bv->vmm = method;
1102: break;
1103: default:1104: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONG,"Unknown matmult method");
1105: }
1106: return(0);
1107: }
1109: /*@
1110: BVGetMatMultMethod - Gets the method used for the BVMatMult() operation.
1112: Not Collective
1114: Input Parameter:
1115: . bv - basis vectors context
1117: Output Parameter:
1118: . method - the method for the BVMatMult() operation
1120: Level: advanced
1122: .seealso: BVMatMult(), BVSetMatMultMethod(), BVMatMultType1123: @*/
1124: PetscErrorCode BVGetMatMultMethod(BV bv,BVMatMultType *method)1125: {
1129: *method = bv->vmm;
1130: return(0);
1131: }
1133: /*@
1134: BVGetColumn - Returns a Vec object that contains the entries of the
1135: requested column of the basis vectors object.
1137: Logically Collective on BV1139: Input Parameters:
1140: + bv - the basis vectors context
1141: - j - the index of the requested column
1143: Output Parameter:
1144: . v - vector containing the jth column
1146: Notes:
1147: The returned Vec must be seen as a reference (not a copy) of the BV1148: column, that is, modifying the Vec will change the BV entries as well.
1150: The returned Vec must not be destroyed. BVRestoreColumn() must be
1151: called when it is no longer needed. At most, two columns can be fetched,
1152: that is, this function can only be called twice before the corresponding
1153: BVRestoreColumn() is invoked.
1155: A negative index j selects the i-th constraint, where i=-j. Constraints
1156: should not be modified.
1158: Level: beginner
1160: .seealso: BVRestoreColumn(), BVInsertConstraints()
1161: @*/
1162: PetscErrorCode BVGetColumn(BV bv,PetscInt j,Vec *v)1163: {
1165: PetscInt l;
1170: BVCheckSizes(bv,1);
1172: 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);
1173: if (j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"You requested column %D but only %D are available",j,bv->m);
1174: 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);
1175: l = BVAvailableVec;
1176: 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");
1177: (*bv->ops->getcolumn)(bv,j,v);
1178: bv->ci[l] = j;
1179: PetscObjectStateGet((PetscObject)bv->cv[l],&bv->st[l]);
1180: PetscObjectGetId((PetscObject)bv->cv[l],&bv->id[l]);
1181: *v = bv->cv[l];
1182: return(0);
1183: }
1185: /*@
1186: BVRestoreColumn - Restore a column obtained with BVGetColumn().
1188: Logically Collective on BV1190: Input Parameters:
1191: + bv - the basis vectors context
1192: . j - the index of the column
1193: - v - vector obtained with BVGetColumn()
1195: Note:
1196: The arguments must match the corresponding call to BVGetColumn().
1198: Level: beginner
1200: .seealso: BVGetColumn()
1201: @*/
1202: PetscErrorCode BVRestoreColumn(BV bv,PetscInt j,Vec *v)1203: {
1204: PetscErrorCode ierr;
1205: PetscObjectId id;
1206: PetscObjectState st;
1207: PetscInt l;
1212: BVCheckSizes(bv,1);
1216: 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);
1217: if (j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"You requested column %D but only %D are available",j,bv->m);
1218: 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);
1219: l = (j==bv->ci[0])? 0: 1;
1220: PetscObjectGetId((PetscObject)*v,&id);
1221: 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");
1222: PetscObjectStateGet((PetscObject)*v,&st);
1223: if (st!=bv->st[l]) {
1224: PetscObjectStateIncrease((PetscObject)bv);
1225: }
1226: if (bv->ops->restorecolumn) {
1227: (*bv->ops->restorecolumn)(bv,j,v);
1228: } else bv->cv[l] = NULL;
1229: bv->ci[l] = -bv->nc-1;
1230: bv->st[l] = -1;
1231: bv->id[l] = 0;
1232: *v = NULL;
1233: return(0);
1234: }
1236: /*@C
1237: BVGetArray - Returns a pointer to a contiguous array that contains this
1238: processor's portion of the BV data.
1240: Logically Collective on BV1242: Input Parameters:
1243: . bv - the basis vectors context
1245: Output Parameter:
1246: . a - location to put pointer to the array
1248: Notes:
1249: BVRestoreArray() must be called when access to the array is no longer needed.
1250: This operation may imply a data copy, for BV types that do not store
1251: data contiguously in memory.
1253: The pointer will normally point to the first entry of the first column,
1254: but if the BV has constraints then these go before the regular columns.
1256: Level: advanced
1258: .seealso: BVRestoreArray(), BVInsertConstraints()
1259: @*/
1260: PetscErrorCode BVGetArray(BV bv,PetscScalar **a)1261: {
1267: BVCheckSizes(bv,1);
1268: (*bv->ops->getarray)(bv,a);
1269: return(0);
1270: }
1272: /*@C
1273: BVRestoreArray - Restore the BV object after BVGetArray() has been called.
1275: Logically Collective on BV1277: Input Parameters:
1278: + bv - the basis vectors context
1279: - a - location of pointer to array obtained from BVGetArray()
1281: Note:
1282: This operation may imply a data copy, for BV types that do not store
1283: data contiguously in memory.
1285: Level: advanced
1287: .seealso: BVGetColumn()
1288: @*/
1289: PetscErrorCode BVRestoreArray(BV bv,PetscScalar **a)1290: {
1296: BVCheckSizes(bv,1);
1297: if (bv->ops->restorearray) {
1298: (*bv->ops->restorearray)(bv,a);
1299: }
1300: if (a) *a = NULL;
1301: PetscObjectStateIncrease((PetscObject)bv);
1302: return(0);
1303: }
1305: /*@C
1306: BVGetArrayRead - Returns a read-only pointer to a contiguous array that
1307: contains this processor's portion of the BV data.
1309: Not Collective
1311: Input Parameters:
1312: . bv - the basis vectors context
1314: Output Parameter:
1315: . a - location to put pointer to the array
1317: Notes:
1318: BVRestoreArrayRead() must be called when access to the array is no
1319: longer needed. This operation may imply a data copy, for BV types that
1320: do not store data contiguously in memory.
1322: The pointer will normally point to the first entry of the first column,
1323: but if the BV has constraints then these go before the regular columns.
1325: Level: advanced
1327: .seealso: BVRestoreArray(), BVInsertConstraints()
1328: @*/
1329: PetscErrorCode BVGetArrayRead(BV bv,const PetscScalar **a)1330: {
1336: BVCheckSizes(bv,1);
1337: (*bv->ops->getarrayread)(bv,a);
1338: return(0);
1339: }
1341: /*@C
1342: BVRestoreArrayRead - Restore the BV object after BVGetArrayRead() has
1343: been called.
1345: Not Collective
1347: Input Parameters:
1348: + bv - the basis vectors context
1349: - a - location of pointer to array obtained from BVGetArrayRead()
1351: Level: advanced
1353: .seealso: BVGetColumn()
1354: @*/
1355: PetscErrorCode BVRestoreArrayRead(BV bv,const PetscScalar **a)1356: {
1362: BVCheckSizes(bv,1);
1363: if (bv->ops->restorearrayread) {
1364: (*bv->ops->restorearrayread)(bv,a);
1365: }
1366: if (a) *a = NULL;
1367: return(0);
1368: }
1370: /*@
1371: BVCreateVec - Creates a new Vec object with the same type and dimensions
1372: as the columns of the basis vectors object.
1374: Collective on BV1376: Input Parameter:
1377: . bv - the basis vectors context
1379: Output Parameter:
1380: . v - the new vector
1382: Note:
1383: The user is responsible of destroying the returned vector.
1385: Level: beginner
1387: .seealso: BVCreateMat()
1388: @*/
1389: PetscErrorCode BVCreateVec(BV bv,Vec *v)1390: {
1395: BVCheckSizes(bv,1);
1397: VecDuplicate(bv->t,v);
1398: return(0);
1399: }
1401: /*@
1402: BVCreateMat - Creates a new Mat object of dense type and copies the contents
1403: of the BV object.
1405: Collective on BV1407: Input Parameter:
1408: . bv - the basis vectors context
1410: Output Parameter:
1411: . A - the new matrix
1413: Note:
1414: The user is responsible of destroying the returned matrix.
1416: Level: intermediate
1418: .seealso: BVCreateFromMat(), BVCreateVec()
1419: @*/
1420: PetscErrorCode BVCreateMat(BV bv,Mat *A)1421: {
1422: PetscErrorCode ierr;
1423: PetscScalar *aa;
1424: const PetscScalar *vv;
1428: BVCheckSizes(bv,1);
1431: MatCreateDense(PetscObjectComm((PetscObject)bv->t),bv->n,PETSC_DECIDE,bv->N,bv->m,NULL,A);
1432: MatDenseGetArray(*A,&aa);
1433: BVGetArrayRead(bv,&vv);
1434: PetscMemcpy(aa,vv,bv->m*bv->n*sizeof(PetscScalar));
1435: BVRestoreArrayRead(bv,&vv);
1436: MatDenseRestoreArray(*A,&aa);
1437: MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
1438: MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
1439: return(0);
1440: }
1442: PETSC_STATIC_INLINE PetscErrorCode BVDuplicate_Private(BV V,PetscInt m,BV *W)1443: {
1447: BVCreate(PetscObjectComm((PetscObject)V),W);
1448: BVSetSizesFromVec(*W,V->t,m);
1449: BVSetType(*W,((PetscObject)V)->type_name);
1450: BVSetMatrix(*W,V->matrix,V->indef);
1451: BVSetOrthogonalization(*W,V->orthog_type,V->orthog_ref,V->orthog_eta,V->orthog_block);
1452: (*W)->vmm = V->vmm;
1453: (*W)->rrandom = V->rrandom;
1454: if (V->ops->duplicate) { (*V->ops->duplicate)(V,W); }
1455: PetscObjectStateIncrease((PetscObject)*W);
1456: return(0);
1457: }
1459: /*@
1460: BVDuplicate - Creates a new basis vector object of the same type and
1461: dimensions as an existing one.
1463: Collective on BV1465: Input Parameter:
1466: . V - basis vectors context
1468: Output Parameter:
1469: . W - location to put the new BV1471: Notes:
1472: The new BV has the same type and dimensions as V, and it shares the same
1473: template vector. Also, the inner product matrix and orthogonalization
1474: options are copied.
1476: BVDuplicate() DOES NOT COPY the entries, but rather allocates storage
1477: for the new basis vectors. Use BVCopy() to copy the contents.
1479: Level: intermediate
1481: .seealso: BVDuplicateResize(), BVCreate(), BVSetSizesFromVec(), BVCopy()
1482: @*/
1483: PetscErrorCode BVDuplicate(BV V,BV *W)1484: {
1490: BVCheckSizes(V,1);
1492: BVDuplicate_Private(V,V->m,W);
1493: return(0);
1494: }
1496: /*@
1497: BVDuplicateResize - Creates a new basis vector object of the same type and
1498: dimensions as an existing one, but with possibly different number of columns.
1500: Collective on BV1502: Input Parameter:
1503: + V - basis vectors context
1504: - m - the new number of columns
1506: Output Parameter:
1507: . W - location to put the new BV1509: Note:
1510: This is equivalent of a call to BVDuplicate() followed by BVResize(). The
1511: contents of V are not copied to W.
1513: Level: intermediate
1515: .seealso: BVDuplicate(), BVResize()
1516: @*/
1517: PetscErrorCode BVDuplicateResize(BV V,PetscInt m,BV *W)1518: {
1524: BVCheckSizes(V,1);
1527: BVDuplicate_Private(V,m,W);
1528: return(0);
1529: }
1531: /*@
1532: BVCopy - Copies a basis vector object into another one, W <- V.
1534: Logically Collective on BV1536: Input Parameter:
1537: . V - basis vectors context
1539: Output Parameter:
1540: . W - the copy
1542: Note:
1543: Both V and W must be distributed in the same manner; local copies are
1544: done. Only active columns (excluding the leading ones) are copied.
1545: In the destination W, columns are overwritten starting from the leading ones.
1546: Constraints are not copied.
1548: Level: beginner
1550: .seealso: BVCopyVec(), BVCopyColumn(), BVDuplicate(), BVSetActiveColumns()
1551: @*/
1552: PetscErrorCode BVCopy(BV V,BV W)1553: {
1554: PetscErrorCode ierr;
1555: PetscScalar *womega;
1556: const PetscScalar *vomega;
1561: BVCheckSizes(V,1);
1564: BVCheckSizes(W,2);
1566: if (V->n!=W->n) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension V %D, W %D",V->n,W->n);
1567: 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);
1568: if (!V->n) return(0);
1570: PetscLogEventBegin(BV_Copy,V,W,0,0);
1571: if (V->indef && V->matrix && V->indef==W->indef && V->matrix==W->matrix) {
1572: /* copy signature */
1573: BV_AllocateSignature(W);
1574: VecGetArrayRead(V->omega,&vomega);
1575: VecGetArray(W->omega,&womega);
1576: PetscMemcpy(womega+W->nc+W->l,vomega+V->nc+V->l,(V->k-V->l)*sizeof(PetscScalar));
1577: VecRestoreArray(W->omega,&womega);
1578: VecRestoreArrayRead(V->omega,&vomega);
1579: }
1580: (*V->ops->copy)(V,W);
1581: PetscLogEventEnd(BV_Copy,V,W,0,0);
1582: return(0);
1583: }
1585: /*@
1586: BVCopyVec - Copies one of the columns of a basis vectors object into a Vec.
1588: Logically Collective on BV1590: Input Parameter:
1591: + V - basis vectors context
1592: - j - the column number to be copied
1594: Output Parameter:
1595: . w - the copied column
1597: Note:
1598: Both V and w must be distributed in the same manner; local copies are done.
1600: Level: beginner
1602: .seealso: BVCopy(), BVCopyColumn(), BVInsertVec()
1603: @*/
1604: PetscErrorCode BVCopyVec(BV V,PetscInt j,Vec w)1605: {
1607: PetscInt n,N;
1608: Vec z;
1613: BVCheckSizes(V,1);
1618: VecGetSize(w,&N);
1619: VecGetLocalSize(w,&n);
1620: 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);
1622: PetscLogEventBegin(BV_Copy,V,w,0,0);
1623: BVGetColumn(V,j,&z);
1624: VecCopy(z,w);
1625: BVRestoreColumn(V,j,&z);
1626: PetscLogEventEnd(BV_Copy,V,w,0,0);
1627: return(0);
1628: }
1630: /*@
1631: BVCopyColumn - Copies the values from one of the columns to another one.
1633: Logically Collective on BV1635: Input Parameter:
1636: + V - basis vectors context
1637: . j - the number of the source column
1638: - i - the number of the destination column
1640: Level: beginner
1642: .seealso: BVCopy(), BVCopyVec()
1643: @*/
1644: PetscErrorCode BVCopyColumn(BV V,PetscInt j,PetscInt i)1645: {
1647: Vec z,w;
1648: PetscScalar *omega;
1653: BVCheckSizes(V,1);
1656: if (j==i) return(0);
1658: PetscLogEventBegin(BV_Copy,V,0,0,0);
1659: if (V->omega) {
1660: VecGetArray(V->omega,&omega);
1661: omega[i] = omega[j];
1662: VecRestoreArray(V->omega,&omega);
1663: }
1664: BVGetColumn(V,j,&z);
1665: BVGetColumn(V,i,&w);
1666: VecCopy(z,w);
1667: BVRestoreColumn(V,j,&z);
1668: BVRestoreColumn(V,i,&w);
1669: PetscLogEventEnd(BV_Copy,V,0,0,0);
1670: return(0);
1671: }