1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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: BV operations involving global communication
12: */
14: #include <slepc/private/bvimpl.h> 16: /*
17: BVDot for the particular case of non-standard inner product with
18: matrix B, which is assumed to be symmetric (or complex Hermitian)
19: */
20: PETSC_STATIC_INLINE PetscErrorCode BVDot_Private(BV X,BV Y,Mat M) 21: {
23: PetscObjectId idx,idy;
24: PetscInt i,j,jend,m;
25: PetscScalar *marray;
26: PetscBool symm=PETSC_FALSE;
27: Mat B;
28: Vec z;
31: BVCheckOp(Y,1,dotvec);
32: MatGetSize(M,&m,NULL);
33: MatDenseGetArray(M,&marray);
34: PetscObjectGetId((PetscObject)X,&idx);
35: PetscObjectGetId((PetscObject)Y,&idy);
36: B = Y->matrix;
37: Y->matrix = NULL;
38: if (idx==idy) symm=PETSC_TRUE; /* M=X'BX is symmetric */
39: jend = X->k;
40: for (j=X->l;j<jend;j++) {
41: if (symm) Y->k = j+1;
42: BVGetColumn(X->cached,j,&z);
43: (*Y->ops->dotvec)(Y,z,marray+j*m+Y->l);
44: BVRestoreColumn(X->cached,j,&z);
45: if (symm) {
46: for (i=X->l;i<j;i++)
47: marray[j+i*m] = PetscConj(marray[i+j*m]);
48: }
49: }
50: MatDenseRestoreArray(M,&marray);
51: Y->matrix = B;
52: return(0);
53: }
55: /*@
56: BVDot - Computes the 'block-dot' product of two basis vectors objects.
58: Collective on X
60: Input Parameters:
61: + X - first basis vectors
62: - Y - second basis vectors
64: Output Parameter:
65: . M - the resulting matrix
67: Notes:
68: This is the generalization of VecDot() for a collection of vectors, M = Y^H*X.
69: The result is a matrix M whose entry m_ij is equal to y_i^H x_j (where y_i^H
70: denotes the conjugate transpose of y_i).
72: If a non-standard inner product has been specified with BVSetMatrix(),
73: then the result is M = Y^H*B*X. In this case, both X and Y must have
74: the same associated matrix.
76: On entry, M must be a sequential dense Mat with dimensions m,n at least, where
77: m is the number of active columns of Y and n is the number of active columns of X.
78: Only rows (resp. columns) of M starting from ly (resp. lx) are computed,
79: where ly (resp. lx) is the number of leading columns of Y (resp. X).
81: X and Y need not be different objects.
83: Level: intermediate
85: .seealso: BVDotVec(), BVDotColumn(), BVSetActiveColumns(), BVSetMatrix()
86: @*/
87: PetscErrorCode BVDot(BV X,BV Y,Mat M) 88: {
90: PetscInt m,n;
97: BVCheckSizes(X,1);
99: BVCheckSizes(Y,2);
104: MatGetSize(M,&m,&n);
105: if (m<Y->k) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_SIZ,"Mat argument has %D rows, should have at least %D",m,Y->k);
106: if (n<X->k) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_SIZ,"Mat argument has %D columns, should have at least %D",n,X->k);
107: if (X->n!=Y->n) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension X %D, Y %D",X->n,Y->n);
108: if (X->matrix!=Y->matrix) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONGSTATE,"X and Y must have the same inner product matrix");
109: if (X->l==X->k || Y->l==Y->k) return(0);
111: PetscLogEventBegin(BV_Dot,X,Y,0,0);
112: if (X->matrix) { /* non-standard inner product */
113: /* compute BX first */
114: BV_IPMatMultBV(X);
115: if (X->vmm==BV_MATMULT_VECS) {
116: /* perform computation column by column */
117: BVDot_Private(X,Y,M);
118: } else {
119: (*X->ops->dot)(X->cached,Y,M);
120: }
121: } else {
122: (*X->ops->dot)(X,Y,M);
123: }
124: PetscLogEventEnd(BV_Dot,X,Y,0,0);
125: return(0);
126: }
128: /*@
129: BVDotVec - Computes multiple dot products of a vector against all the
130: column vectors of a BV.
132: Collective on X
134: Input Parameters:
135: + X - basis vectors
136: - y - a vector
138: Output Parameter:
139: . m - an array where the result must be placed
141: Notes:
142: This is analogue to VecMDot(), but using BV to represent a collection
143: of vectors. The result is m = X^H*y, so m_i is equal to x_j^H y. Note
144: that here X is transposed as opposed to BVDot().
146: If a non-standard inner product has been specified with BVSetMatrix(),
147: then the result is m = X^H*B*y.
149: The length of array m must be equal to the number of active columns of X
150: minus the number of leading columns, i.e. the first entry of m is the
151: product of the first non-leading column with y.
153: Level: intermediate
155: .seealso: BVDot(), BVDotColumn(), BVSetActiveColumns(), BVSetMatrix()
156: @*/
157: PetscErrorCode BVDotVec(BV X,Vec y,PetscScalar m[])158: {
160: PetscInt n;
166: BVCheckSizes(X,1);
167: BVCheckOp(X,1,dotvec);
171: VecGetLocalSize(y,&n);
172: if (X->n!=n) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension X %D, y %D",X->n,n);
174: PetscLogEventBegin(BV_DotVec,X,y,0,0);
175: (*X->ops->dotvec)(X,y,m);
176: PetscLogEventEnd(BV_DotVec,X,y,0,0);
177: return(0);
178: }
180: /*@
181: BVDotVecBegin - Starts a split phase dot product computation.
183: Input Parameters:
184: + X - basis vectors
185: . y - a vector
186: - m - an array where the result will go (can be NULL)
188: Note:
189: Each call to BVDotVecBegin() should be paired with a call to BVDotVecEnd().
191: Level: advanced
193: .seealso: BVDotVecEnd(), BVDotVec()
194: @*/
195: PetscErrorCode BVDotVecBegin(BV X,Vec y,PetscScalar *m)196: {
197: PetscErrorCode ierr;
198: PetscInt i,n,nv;
199: PetscSplitReduction *sr;
200: MPI_Comm comm;
206: BVCheckSizes(X,1);
210: VecGetLocalSize(y,&n);
211: if (X->n!=n) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension X %D, y %D",X->n,n);
213: if (X->ops->dotvec_begin) {
214: (*X->ops->dotvec_begin)(X,y,m);
215: } else {
216: BVCheckOp(X,1,dotvec_local);
217: nv = X->k-X->l;
218: PetscObjectGetComm((PetscObject)X,&comm);
219: PetscSplitReductionGet(comm,&sr);
220: if (sr->state != STATE_BEGIN) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ORDER,"Called before all VecxxxEnd() called");
221: for (i=0;i<nv;i++) {
222: if (sr->numopsbegin+i >= sr->maxops) {
223: PetscSplitReductionExtend(sr);
224: }
225: sr->reducetype[sr->numopsbegin+i] = PETSC_SR_REDUCE_SUM;
226: sr->invecs[sr->numopsbegin+i] = (void*)X;
227: }
228: PetscLogEventBegin(BV_DotVec,X,y,0,0);
229: (*X->ops->dotvec_local)(X,y,sr->lvalues+sr->numopsbegin);
230: sr->numopsbegin += nv;
231: PetscLogEventEnd(BV_DotVec,X,y,0,0);
232: }
233: return(0);
234: }
236: /*@
237: BVDotVecEnd - Ends a split phase dot product computation.
239: Input Parameters:
240: + X - basis vectors
241: . y - a vector
242: - m - an array where the result will go
244: Note:
245: Each call to BVDotVecBegin() should be paired with a call to BVDotVecEnd().
247: Level: advanced
249: .seealso: BVDotVecBegin(), BVDotVec()
250: @*/
251: PetscErrorCode BVDotVecEnd(BV X,Vec y,PetscScalar *m)252: {
253: PetscErrorCode ierr;
254: PetscInt i,nv;
255: PetscSplitReduction *sr;
256: MPI_Comm comm;
261: BVCheckSizes(X,1);
263: if (X->ops->dotvec_end) {
264: (*X->ops->dotvec_end)(X,y,m);
265: } else {
266: nv = X->k-X->l;
267: PetscObjectGetComm((PetscObject)X,&comm);
268: PetscSplitReductionGet(comm,&sr);
269: PetscSplitReductionEnd(sr);
271: if (sr->numopsend >= sr->numopsbegin) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONGSTATE,"Called VecxxxEnd() more times than VecxxxBegin()");
272: if ((void*)X != sr->invecs[sr->numopsend]) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONGSTATE,"Called BVxxxEnd() in a different order or with a different BV than BVxxxBegin()");
273: if (sr->reducetype[sr->numopsend] != PETSC_SR_REDUCE_SUM) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONGSTATE,"Wrong type of reduction");
274: for (i=0;i<nv;i++) m[i] = sr->gvalues[sr->numopsend++];
276: /* Finished getting all the results so reset to no outstanding requests */
277: if (sr->numopsend == sr->numopsbegin) {
278: sr->state = STATE_BEGIN;
279: sr->numopsend = 0;
280: sr->numopsbegin = 0;
281: }
282: }
283: return(0);
284: }
286: /*@
287: BVDotColumn - Computes multiple dot products of a column against all the
288: previous columns of a BV.
290: Collective on X
292: Input Parameters:
293: + X - basis vectors
294: - j - the column index
296: Output Parameter:
297: . q - an array where the result must be placed
299: Notes:
300: This operation is equivalent to BVDotVec() but it uses column j of X
301: rather than taking a Vec as an argument. The number of active columns of
302: X is set to j before the computation, and restored afterwards.
303: If X has leading columns specified, then these columns do not participate
304: in the computation. Therefore, the length of array q must be equal to j
305: minus the number of leading columns.
307: Developer Notes:
308: If q is NULL, then the result is written in position nc+l of the internal
309: buffer vector, see BVGetBufferVec().
311: Level: advanced
313: .seealso: BVDot(), BVDotVec(), BVSetActiveColumns(), BVSetMatrix()
314: @*/
315: PetscErrorCode BVDotColumn(BV X,PetscInt j,PetscScalar *q)316: {
318: PetscInt ksave;
319: Vec y;
325: BVCheckSizes(X,1);
326: BVCheckOp(X,1,dotvec);
328: if (j<0) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
329: if (j>=X->m) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_OUTOFRANGE,"Index j=%D but BV only has %D columns",j,X->m);
331: PetscLogEventBegin(BV_DotVec,X,0,0,0);
332: ksave = X->k;
333: X->k = j;
334: if (!q && !X->buffer) { BVGetBufferVec(X,&X->buffer); }
335: BVGetColumn(X,j,&y);
336: (*X->ops->dotvec)(X,y,q);
337: BVRestoreColumn(X,j,&y);
338: X->k = ksave;
339: PetscLogEventEnd(BV_DotVec,X,0,0,0);
340: return(0);
341: }
343: /*@
344: BVDotColumnBegin - Starts a split phase dot product computation.
346: Input Parameters:
347: + X - basis vectors
348: - j - the column index
349: - m - an array where the result will go (can be NULL)
351: Note:
352: Each call to BVDotColumnBegin() should be paired with a call to BVDotColumnEnd().
354: Level: advanced
356: .seealso: BVDotColumnEnd(), BVDotColumn()
357: @*/
358: PetscErrorCode BVDotColumnBegin(BV X,PetscInt j,PetscScalar *m)359: {
360: PetscErrorCode ierr;
361: PetscInt i,nv,ksave;
362: PetscSplitReduction *sr;
363: MPI_Comm comm;
364: Vec y;
370: BVCheckSizes(X,1);
372: if (j<0) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
373: if (j>=X->m) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_OUTOFRANGE,"Index j=%D but BV only has %D columns",j,X->m);
374: ksave = X->k;
375: X->k = j;
376: BVGetColumn(X,j,&y);
378: if (X->ops->dotvec_begin) {
379: (*X->ops->dotvec_begin)(X,y,m);
380: } else {
381: BVCheckOp(X,1,dotvec_local);
382: nv = X->k-X->l;
383: PetscObjectGetComm((PetscObject)X,&comm);
384: PetscSplitReductionGet(comm,&sr);
385: if (sr->state != STATE_BEGIN) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ORDER,"Called before all VecxxxEnd() called");
386: for (i=0;i<nv;i++) {
387: if (sr->numopsbegin+i >= sr->maxops) {
388: PetscSplitReductionExtend(sr);
389: }
390: sr->reducetype[sr->numopsbegin+i] = PETSC_SR_REDUCE_SUM;
391: sr->invecs[sr->numopsbegin+i] = (void*)X;
392: }
393: PetscLogEventBegin(BV_DotVec,X,0,0,0);
394: (*X->ops->dotvec_local)(X,y,sr->lvalues+sr->numopsbegin);
395: sr->numopsbegin += nv;
396: PetscLogEventEnd(BV_DotVec,X,0,0,0);
397: }
398: BVRestoreColumn(X,j,&y);
399: X->k = ksave;
400: return(0);
401: }
403: /*@
404: BVDotColumnEnd - Ends a split phase dot product computation.
406: Input Parameters:
407: + X - basis vectors
408: . j - the column index
409: - m - an array where the result will go
411: Notes:
412: Each call to BVDotColumnBegin() should be paired with a call to BVDotColumnEnd().
414: Level: advanced
416: .seealso: BVDotColumnBegin(), BVDotColumn()
417: @*/
418: PetscErrorCode BVDotColumnEnd(BV X,PetscInt j,PetscScalar *m)419: {
420: PetscErrorCode ierr;
421: PetscInt i,nv,ksave;
422: PetscSplitReduction *sr;
423: MPI_Comm comm;
424: Vec y;
430: BVCheckSizes(X,1);
432: if (j<0) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
433: if (j>=X->m) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_OUTOFRANGE,"Index j=%D but BV only has %D columns",j,X->m);
434: ksave = X->k;
435: X->k = j;
437: if (X->ops->dotvec_end) {
438: BVGetColumn(X,j,&y);
439: (*X->ops->dotvec_end)(X,y,m);
440: BVRestoreColumn(X,j,&y);
441: } else {
442: nv = X->k-X->l;
443: PetscObjectGetComm((PetscObject)X,&comm);
444: PetscSplitReductionGet(comm,&sr);
445: PetscSplitReductionEnd(sr);
447: if (sr->numopsend >= sr->numopsbegin) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONGSTATE,"Called VecxxxEnd() more times than VecxxxBegin()");
448: if ((void*)X != sr->invecs[sr->numopsend]) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONGSTATE,"Called BVxxxEnd() in a different order or with a different BV than BVxxxBegin()");
449: if (sr->reducetype[sr->numopsend] != PETSC_SR_REDUCE_SUM) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONGSTATE,"Wrong type of reduction");
450: for (i=0;i<nv;i++) m[i] = sr->gvalues[sr->numopsend++];
452: /* Finished getting all the results so reset to no outstanding requests */
453: if (sr->numopsend == sr->numopsbegin) {
454: sr->state = STATE_BEGIN;
455: sr->numopsend = 0;
456: sr->numopsbegin = 0;
457: }
458: }
459: X->k = ksave;
460: return(0);
461: }
463: PETSC_STATIC_INLINE PetscErrorCode BVNorm_Private(BV bv,Vec z,NormType type,PetscReal *val)464: {
466: PetscScalar p;
469: BV_IPMatMult(bv,z);
470: VecDot(bv->Bx,z,&p);
471: BV_SafeSqrt(bv,p,val);
472: return(0);
473: }
475: PETSC_STATIC_INLINE PetscErrorCode BVNorm_Begin_Private(BV bv,Vec z,NormType type,PetscReal *val)476: {
478: PetscScalar p;
481: BV_IPMatMult(bv,z);
482: VecDotBegin(bv->Bx,z,&p);
483: return(0);
484: }
486: PETSC_STATIC_INLINE PetscErrorCode BVNorm_End_Private(BV bv,Vec z,NormType type,PetscReal *val)487: {
489: PetscScalar p;
492: VecDotEnd(bv->Bx,z,&p);
493: BV_SafeSqrt(bv,p,val);
494: return(0);
495: }
497: /*@
498: BVNorm - Computes the matrix norm of the BV.
500: Collective on bv
502: Input Parameters:
503: + bv - basis vectors
504: - type - the norm type
506: Output Parameter:
507: . val - the norm
509: Notes:
510: All active columns (except the leading ones) are considered as a matrix.
511: The allowed norms are NORM_1, NORM_FROBENIUS, and NORM_INFINITY.
513: This operation fails if a non-standard inner product has been
514: specified with BVSetMatrix().
516: Level: intermediate
518: .seealso: BVNormVec(), BVNormColumn(), BVNormalize(), BVSetActiveColumns(), BVSetMatrix()
519: @*/
520: PetscErrorCode BVNorm(BV bv,NormType type,PetscReal *val)521: {
529: BVCheckSizes(bv,1);
531: if (type==NORM_2 || type==NORM_1_AND_2) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not available");
532: if (bv->matrix) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Matrix norm not available for non-standard inner product");
534: PetscLogEventBegin(BV_Norm,bv,0,0,0);
535: (*bv->ops->norm)(bv,-1,type,val);
536: PetscLogEventEnd(BV_Norm,bv,0,0,0);
537: return(0);
538: }
540: /*@
541: BVNormVec - Computes the norm of a given vector.
543: Collective on bv
545: Input Parameters:
546: + bv - basis vectors
547: . v - the vector
548: - type - the norm type
550: Output Parameter:
551: . val - the norm
553: Notes:
554: This is the analogue of BVNormColumn() but for a vector that is not in the BV.
555: If a non-standard inner product has been specified with BVSetMatrix(),
556: then the returned value is sqrt(v'*B*v), where B is the inner product
557: matrix (argument 'type' is ignored). Otherwise, VecNorm() is called.
559: Level: developer
561: .seealso: BVNorm(), BVNormColumn(), BVSetMatrix()
562: @*/
563: PetscErrorCode BVNormVec(BV bv,Vec v,NormType type,PetscReal *val)564: {
566: PetscInt n;
574: BVCheckSizes(bv,1);
578: if (type==NORM_1_AND_2 && !bv->matrix) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not available");
580: PetscLogEventBegin(BV_NormVec,bv,0,0,0);
581: if (bv->matrix) { /* non-standard inner product */
582: VecGetLocalSize(v,&n);
583: if (bv->n!=n) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension bv %D, v %D",bv->n,n);
584: BVNorm_Private(bv,v,type,val);
585: } else {
586: VecNorm(v,type,val);
587: }
588: PetscLogEventEnd(BV_NormVec,bv,0,0,0);
589: return(0);
590: }
592: /*@
593: BVNormVecBegin - Starts a split phase norm computation.
595: Input Parameters:
596: + bv - basis vectors
597: . v - the vector
598: . type - the norm type
599: - val - the norm
601: Note:
602: Each call to BVNormVecBegin() should be paired with a call to BVNormVecEnd().
604: Level: advanced
606: .seealso: BVNormVecEnd(), BVNormVec()
607: @*/
608: PetscErrorCode BVNormVecBegin(BV bv,Vec v,NormType type,PetscReal *val)609: {
611: PetscInt n;
619: BVCheckSizes(bv,1);
623: if (type==NORM_1_AND_2 && !bv->matrix) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not available");
625: PetscLogEventBegin(BV_NormVec,bv,0,0,0);
626: if (bv->matrix) { /* non-standard inner product */
627: VecGetLocalSize(v,&n);
628: if (bv->n!=n) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension bv %D, v %D",bv->n,n);
629: BVNorm_Begin_Private(bv,v,type,val);
630: } else {
631: VecNormBegin(v,type,val);
632: }
633: PetscLogEventEnd(BV_NormVec,bv,0,0,0);
634: return(0);
635: }
637: /*@
638: BVNormVecEnd - Ends a split phase norm computation.
640: Input Parameters:
641: + bv - basis vectors
642: . v - the vector
643: . type - the norm type
644: - val - the norm
646: Note:
647: Each call to BVNormVecBegin() should be paired with a call to BVNormVecEnd().
649: Level: advanced
651: .seealso: BVNormVecBegin(), BVNormVec()
652: @*/
653: PetscErrorCode BVNormVecEnd(BV bv,Vec v,NormType type,PetscReal *val)654: {
662: BVCheckSizes(bv,1);
664: if (type==NORM_1_AND_2) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not available");
666: if (bv->matrix) { /* non-standard inner product */
667: BVNorm_End_Private(bv,v,type,val);
668: } else {
669: VecNormEnd(v,type,val);
670: }
671: return(0);
672: }
674: /*@
675: BVNormColumn - Computes the vector norm of a selected column.
677: Collective on bv
679: Input Parameters:
680: + bv - basis vectors
681: . j - column number to be used
682: - type - the norm type
684: Output Parameter:
685: . val - the norm
687: Notes:
688: The norm of V[j] is computed (NORM_1, NORM_2, or NORM_INFINITY).
689: If a non-standard inner product has been specified with BVSetMatrix(),
690: then the returned value is sqrt(V[j]'*B*V[j]),
691: where B is the inner product matrix (argument 'type' is ignored).
693: Level: intermediate
695: .seealso: BVNorm(), BVNormVec(), BVNormalize(), BVSetActiveColumns(), BVSetMatrix()
696: @*/
697: PetscErrorCode BVNormColumn(BV bv,PetscInt j,NormType type,PetscReal *val)698: {
700: Vec z;
708: BVCheckSizes(bv,1);
710: if (j<0 || j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Argument j has wrong value %D, the number of columns is %D",j,bv->m);
711: if (type==NORM_1_AND_2 && !bv->matrix) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not available");
713: PetscLogEventBegin(BV_NormVec,bv,0,0,0);
714: if (bv->matrix) { /* non-standard inner product */
715: BVGetColumn(bv,j,&z);
716: BVNorm_Private(bv,z,type,val);
717: BVRestoreColumn(bv,j,&z);
718: } else {
719: (*bv->ops->norm)(bv,j,type,val);
720: }
721: PetscLogEventEnd(BV_NormVec,bv,0,0,0);
722: return(0);
723: }
725: /*@
726: BVNormColumnBegin - Starts a split phase norm computation.
728: Input Parameters:
729: + bv - basis vectors
730: . j - column number to be used
731: . type - the norm type
732: - val - the norm
734: Note:
735: Each call to BVNormColumnBegin() should be paired with a call to BVNormColumnEnd().
737: Level: advanced
739: .seealso: BVNormColumnEnd(), BVNormColumn()
740: @*/
741: PetscErrorCode BVNormColumnBegin(BV bv,PetscInt j,NormType type,PetscReal *val)742: {
743: PetscErrorCode ierr;
744: PetscSplitReduction *sr;
745: PetscReal lresult;
746: MPI_Comm comm;
747: Vec z;
755: BVCheckSizes(bv,1);
757: if (j<0 || j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Argument j has wrong value %D, the number of columns is %D",j,bv->m);
758: if (type==NORM_1_AND_2 && !bv->matrix) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not available");
760: PetscLogEventBegin(BV_NormVec,bv,0,0,0);
761: BVGetColumn(bv,j,&z);
762: if (bv->matrix) { /* non-standard inner product */
763: BVNorm_Begin_Private(bv,z,type,val);
764: } else if (bv->ops->norm_begin) {
765: (*bv->ops->norm_begin)(bv,j,type,val);
766: } else {
767: BVCheckOp(bv,1,norm_local);
768: PetscObjectGetComm((PetscObject)z,&comm);
769: PetscSplitReductionGet(comm,&sr);
770: if (sr->state != STATE_BEGIN) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ORDER,"Called before all VecxxxEnd() called");
771: if (sr->numopsbegin >= sr->maxops) {
772: PetscSplitReductionExtend(sr);
773: }
774: sr->invecs[sr->numopsbegin] = (void*)bv;
775: (*bv->ops->norm_local)(bv,j,type,&lresult);
776: if (type == NORM_2) lresult = lresult*lresult;
777: if (type == NORM_MAX) sr->reducetype[sr->numopsbegin] = PETSC_SR_REDUCE_MAX;
778: else sr->reducetype[sr->numopsbegin] = PETSC_SR_REDUCE_SUM;
779: sr->lvalues[sr->numopsbegin++] = lresult;
780: }
781: BVRestoreColumn(bv,j,&z);
782: PetscLogEventEnd(BV_NormVec,bv,0,0,0);
783: return(0);
784: }
786: /*@
787: BVNormColumnEnd - Ends a split phase norm computation.
789: Input Parameters:
790: + bv - basis vectors
791: . j - column number to be used
792: . type - the norm type
793: - val - the norm
795: Note:
796: Each call to BVNormColumnBegin() should be paired with a call to BVNormColumnEnd().
798: Level: advanced
800: .seealso: BVNormColumnBegin(), BVNormColumn()
801: @*/
802: PetscErrorCode BVNormColumnEnd(BV bv,PetscInt j,NormType type,PetscReal *val)803: {
804: PetscErrorCode ierr;
805: PetscSplitReduction *sr;
806: MPI_Comm comm;
807: Vec z;
815: BVCheckSizes(bv,1);
817: if (type==NORM_1_AND_2) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not available");
819: BVGetColumn(bv,j,&z);
820: if (bv->matrix) { /* non-standard inner product */
821: BVNorm_End_Private(bv,z,type,val);
822: } else if (bv->ops->norm_end) {
823: (*bv->ops->norm_end)(bv,j,type,val);
824: } else {
825: PetscObjectGetComm((PetscObject)z,&comm);
826: PetscSplitReductionGet(comm,&sr);
827: PetscSplitReductionEnd(sr);
829: if (sr->numopsend >= sr->numopsbegin) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"Called VecxxxEnd() more times then VecxxxBegin()");
830: if ((void*)bv != sr->invecs[sr->numopsend]) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"Called VecxxxEnd() in a different order or with a different vector than VecxxxBegin()");
831: if (sr->reducetype[sr->numopsend] != PETSC_SR_REDUCE_MAX && type == NORM_MAX) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"Called BVNormEnd(,NORM_MAX,) on a reduction started with VecDotBegin() or NORM_1 or NORM_2");
832: *val = PetscRealPart(sr->gvalues[sr->numopsend++]);
833: if (type == NORM_2) *val = PetscSqrtReal(*val);
834: if (sr->numopsend == sr->numopsbegin) {
835: sr->state = STATE_BEGIN;
836: sr->numopsend = 0;
837: sr->numopsbegin = 0;
838: }
839: }
840: BVRestoreColumn(bv,j,&z);
841: return(0);
842: }
844: /*@
845: BVNormalize - Normalize all columns (starting from the leading ones).
847: Collective on bv
849: Input Parameters:
850: + bv - basis vectors
851: - eigi - (optional) imaginary parts of eigenvalues
853: Notes:
854: On output, all columns will have unit norm. The normalization is done with
855: respect to the 2-norm, or to the B-norm if a non-standard inner product has
856: been specified with BVSetMatrix(), see BVNormColumn().
858: If the optional argument eigi is passed (taken into account only in real
859: scalars) it is interpreted as the imaginary parts of the eigenvalues and
860: the BV is supposed to contain the corresponding eigenvectors. Suppose the
861: first three values are eigi = { 0, alpha, -alpha }, then the first column
862: is normalized as usual, but the second and third ones are normalized assuming
863: that they contain the real and imaginary parts of a complex conjugate pair of
864: eigenvectors.
866: If eigi is passed, the inner-product matrix is ignored.
868: If there are leading columns, they are not modified (are assumed to be already
869: normalized).
871: Level: intermediate
873: .seealso: BVNormColumn()
874: @*/
875: PetscErrorCode BVNormalize(BV bv,PetscScalar *eigi)876: {
878: PetscReal norm;
879: PetscInt i;
880: Vec v;
881: #if !defined(PETSC_USE_COMPLEX)
882: Vec v1;
883: #endif
888: BVCheckSizes(bv,1);
890: PetscLogEventBegin(BV_Normalize,bv,0,0,0);
891: if (bv->matrix && !eigi) {
892: for (i=bv->l;i<bv->k;i++) {
893: BVNormColumn(bv,i,NORM_2,&norm);
894: BVScaleColumn(bv,i,1.0/norm);
895: }
896: } else if (bv->ops->normalize) {
897: (*bv->ops->normalize)(bv,eigi);
898: } else {
899: for (i=bv->l;i<bv->k;i++) {
900: #if !defined(PETSC_USE_COMPLEX)
901: if (eigi && eigi[i] != 0.0) {
902: BVGetColumn(bv,i,&v);
903: BVGetColumn(bv,i+1,&v1);
904: VecNormalizeComplex(v,v1,PETSC_TRUE,NULL);
905: BVRestoreColumn(bv,i,&v);
906: BVRestoreColumn(bv,i+1,&v1);
907: i++;
908: } else
909: #endif
910: {
911: BVGetColumn(bv,i,&v);
912: VecNormalize(v,NULL);
913: BVRestoreColumn(bv,i,&v);
914: }
915: }
916: }
917: PetscLogEventEnd(BV_Normalize,bv,0,0,0);
918: PetscObjectStateIncrease((PetscObject)bv);
919: return(0);
920: }
922: /*
923: Compute Y^H*A*X: right part column by column (with MatMult) and bottom
924: part row by row (with MatMultHermitianTranspose); result placed in marray[*,ldm]
925: */
926: PETSC_STATIC_INLINE PetscErrorCode BVMatProject_Vec(BV X,Mat A,BV Y,PetscScalar *marray,PetscInt ldm,PetscBool symm)927: {
929: PetscInt i,j,lx,ly,kx,ky,ulim;
930: Vec z,f;
933: lx = X->l; kx = X->k;
934: ly = Y->l; ky = Y->k;
935: BVCreateVec(X,&f);
936: BVCheckOp(Y,3,dotvec);
937: for (j=lx;j<kx;j++) {
938: BVGetColumn(X,j,&z);
939: MatMult(A,z,f);
940: BVRestoreColumn(X,j,&z);
941: ulim = PetscMin(ly+(j-lx)+1,ky);
942: Y->l = 0; Y->k = ulim;
943: (*Y->ops->dotvec)(Y,f,marray+j*ldm);
944: if (symm) {
945: for (i=0;i<j;i++) marray[j+i*ldm] = PetscConj(marray[i+j*ldm]);
946: }
947: }
948: if (!symm) {
949: BVCheckOp(X,1,dotvec);
950: BV_AllocateCoeffs(Y);
951: for (j=ly;j<ky;j++) {
952: BVGetColumn(Y,j,&z);
953: MatMultHermitianTranspose(A,z,f);
954: BVRestoreColumn(Y,j,&z);
955: ulim = PetscMin(lx+(j-ly),kx);
956: X->l = 0; X->k = ulim;
957: (*X->ops->dotvec)(X,f,Y->h);
958: for (i=0;i<ulim;i++) marray[j+i*ldm] = PetscConj(Y->h[i]);
959: }
960: }
961: VecDestroy(&f);
962: X->l = lx; X->k = kx;
963: Y->l = ly; Y->k = ky;
964: return(0);
965: }
967: /*
968: Compute Y^H*A*X= [ -- | Y0'*W1 ]
969: [ Y1'*W0 | Y1'*W1 ]
970: Allocates auxiliary BV to store the result of A*X, then one BVDot971: call for top-right part and another one for bottom part;
972: result placed in marray[*,ldm]
973: */
974: PETSC_STATIC_INLINE PetscErrorCode BVMatProject_MatMult(BV X,Mat A,BV Y,PetscScalar *marray,PetscInt ldm)975: {
976: PetscErrorCode ierr;
977: PetscInt j,lx,ly,kx,ky;
978: const PetscScalar *harray;
979: Mat H;
980: BV W;
983: lx = X->l; kx = X->k;
984: ly = Y->l; ky = Y->k;
985: BVDuplicate(X,&W);
986: X->l = 0; X->k = kx;
987: W->l = 0; W->k = kx;
988: BVMatMult(X,A,W);
990: /* top-right part, Y0'*AX1 */
991: if (ly>0 && lx<kx) {
992: MatCreateSeqDense(PETSC_COMM_SELF,ly,kx,NULL,&H);
993: W->l = lx; W->k = kx;
994: Y->l = 0; Y->k = ly;
995: BVDot(W,Y,H);
996: MatDenseGetArrayRead(H,&harray);
997: for (j=lx;j<kx;j++) {
998: PetscArraycpy(marray+j*ldm,harray+j*ly,ly);
999: }
1000: MatDenseRestoreArrayRead(H,&harray);
1001: MatDestroy(&H);
1002: }
1004: /* bottom part, Y1'*AX */
1005: if (kx>0 && ly<ky) {
1006: MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&H);
1007: W->l = 0; W->k = kx;
1008: Y->l = ly; Y->k = ky;
1009: BVDot(W,Y,H);
1010: MatDenseGetArrayRead(H,&harray);
1011: for (j=0;j<kx;j++) {
1012: PetscArraycpy(marray+j*ldm+ly,harray+j*ky+ly,ky-ly);
1013: }
1014: MatDenseRestoreArrayRead(H,&harray);
1015: MatDestroy(&H);
1016: }
1017: BVDestroy(&W);
1018: X->l = lx; X->k = kx;
1019: Y->l = ly; Y->k = ky;
1020: return(0);
1021: }
1023: /*
1024: Compute Y^H*A*X= [ -- | Y0'*W1 ]
1025: [ Y1'*W0 | Y1'*W1 ]
1026: First stage: allocate auxiliary BV to store A*X1, one BVDot for right part;
1027: Second stage: resize BV to accommodate A'*Y1, then call BVDot for transpose of
1028: bottom-left part; result placed in marray[*,ldm]
1029: */
1030: PETSC_STATIC_INLINE PetscErrorCode BVMatProject_MatMult_2(BV X,Mat A,BV Y,PetscScalar *marray,PetscInt ldm,PetscBool symm)1031: {
1032: PetscErrorCode ierr;
1033: PetscInt i,j,lx,ly,kx,ky;
1034: const PetscScalar *harray;
1035: Mat H;
1036: BV W;
1039: lx = X->l; kx = X->k;
1040: ly = Y->l; ky = Y->k;
1042: /* right part, Y'*AX1 */
1043: BVDuplicateResize(X,kx-lx,&W);
1044: if (ky>0 && lx<kx) {
1045: BVMatMult(X,A,W);
1046: MatCreateSeqDense(PETSC_COMM_SELF,ky,kx-lx,NULL,&H);
1047: Y->l = 0; Y->k = ky;
1048: BVDot(W,Y,H);
1049: MatDenseGetArrayRead(H,&harray);
1050: for (j=lx;j<kx;j++) {
1051: PetscArraycpy(marray+j*ldm,harray+(j-lx)*ky,ky);
1052: }
1053: MatDenseRestoreArrayRead(H,&harray);
1054: MatDestroy(&H);
1055: }
1057: /* bottom-left part, Y1'*AX0 */
1058: if (lx>0 && ly<ky) {
1059: if (symm) {
1060: /* do not compute, just copy symmetric elements */
1061: for (i=ly;i<ky;i++) {
1062: for (j=0;j<lx;j++) marray[i+j*ldm] = PetscConj(marray[j+i*ldm]);
1063: }
1064: } else {
1065: BVResize(W,ky-ly,PETSC_FALSE);
1066: Y->l = ly; Y->k = ky;
1067: BVMatMultHermitianTranspose(Y,A,W);
1068: MatCreateSeqDense(PETSC_COMM_SELF,lx,ky-ly,NULL,&H);
1069: X->l = 0; X->k = lx;
1070: BVDot(W,X,H);
1071: MatDenseGetArrayRead(H,&harray);
1072: for (i=0;i<ky-ly;i++) {
1073: for (j=0;j<lx;j++) {
1074: marray[i+j*ldm+ly] = PetscConj(harray[j+i*(ky-ly)]);
1075: }
1076: }
1077: MatDenseRestoreArrayRead(H,&harray);
1078: MatDestroy(&H);
1079: }
1080: }
1081: BVDestroy(&W);
1082: X->l = lx; X->k = kx;
1083: Y->l = ly; Y->k = ky;
1084: return(0);
1085: }
1087: /*
1088: Compute Y^H*X = [ -- | Y0'*X1 ] (X contains A*X):
1089: [ Y1'*X0 | Y1'*X1 ]
1090: one BVDot call for top-right part and another one for bottom part;
1091: result placed in marray[*,ldm]
1092: */
1093: PETSC_STATIC_INLINE PetscErrorCode BVMatProject_Dot(BV X,BV Y,PetscScalar *marray,PetscInt ldm)1094: {
1095: PetscErrorCode ierr;
1096: PetscInt j,lx,ly,kx,ky;
1097: const PetscScalar *harray;
1098: Mat H;
1101: lx = X->l; kx = X->k;
1102: ly = Y->l; ky = Y->k;
1104: /* top-right part, Y0'*X1 */
1105: if (ly>0 && lx<kx) {
1106: MatCreateSeqDense(PETSC_COMM_SELF,ly,kx,NULL,&H);
1107: X->l = lx; X->k = kx;
1108: Y->l = 0; Y->k = ly;
1109: BVDot(X,Y,H);
1110: MatDenseGetArrayRead(H,&harray);
1111: for (j=lx;j<kx;j++) {
1112: PetscArraycpy(marray+j*ldm,harray+j*ly,ly);
1113: }
1114: MatDenseRestoreArrayRead(H,&harray);
1115: MatDestroy(&H);
1116: }
1118: /* bottom part, Y1'*X */
1119: if (kx>0 && ly<ky) {
1120: MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&H);
1121: X->l = 0; X->k = kx;
1122: Y->l = ly; Y->k = ky;
1123: BVDot(X,Y,H);
1124: MatDenseGetArrayRead(H,&harray);
1125: for (j=0;j<kx;j++) {
1126: PetscArraycpy(marray+j*ldm+ly,harray+j*ky+ly,ky-ly);
1127: }
1128: MatDenseRestoreArrayRead(H,&harray);
1129: MatDestroy(&H);
1130: }
1131: X->l = lx; X->k = kx;
1132: Y->l = ly; Y->k = ky;
1133: return(0);
1134: }
1136: /*@
1137: BVMatProject - Computes the projection of a matrix onto a subspace.
1139: Collective on X
1141: Input Parameters:
1142: + X - basis vectors
1143: . A - (optional) matrix to be projected
1144: . Y - left basis vectors, can be equal to X
1145: - M - Mat object where the result must be placed
1147: Output Parameter:
1148: . M - the resulting matrix
1150: Notes:
1151: If A=NULL, then it is assumed that X already contains A*X.
1153: This operation is similar to BVDot(), with important differences.
1154: The goal is to compute the matrix resulting from the orthogonal projection
1155: of A onto the subspace spanned by the columns of X, M = X^H*A*X, or the
1156: oblique projection onto X along Y, M = Y^H*A*X.
1158: A difference with respect to BVDot() is that the standard inner product
1159: is always used, regardless of a non-standard inner product being specified
1160: with BVSetMatrix().
1162: On entry, M must be a sequential dense Mat with dimensions ky,kx at least,
1163: where ky (resp. kx) is the number of active columns of Y (resp. X).
1164: Another difference with respect to BVDot() is that all entries of M are
1165: computed except the leading ly,lx part, where ly (resp. lx) is the
1166: number of leading columns of Y (resp. X). Hence, the leading columns of
1167: X and Y participate in the computation, as opposed to BVDot().
1168: The leading part of M is assumed to be already available from previous
1169: computations.
1171: In the orthogonal projection case, Y=X, some computation can be saved if
1172: A is real symmetric (or complex Hermitian). In order to exploit this
1173: property, the symmetry flag of A must be set with MatSetOption().
1175: Level: intermediate
1177: .seealso: BVDot(), BVSetActiveColumns(), BVSetMatrix()
1178: @*/
1179: PetscErrorCode BVMatProject(BV X,Mat A,BV Y,Mat M)1180: {
1182: PetscBool set,flg,symm=PETSC_FALSE;
1183: PetscInt m,n;
1184: PetscScalar *marray;
1185: Mat Xmatrix,Ymatrix;
1186: PetscObjectId idx,idy;
1194: BVCheckSizes(X,1);
1195: if (A) {
1198: }
1200: BVCheckSizes(Y,3);
1205: MatGetSize(M,&m,&n);
1206: if (m<Y->k) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_SIZ,"Matrix M has %D rows, should have at least %D",m,Y->k);
1207: if (n<X->k) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_SIZ,"Matrix M has %D columns, should have at least %D",n,X->k);
1208: if (X->n!=Y->n) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension X %D, Y %D",X->n,Y->n);
1210: PetscLogEventBegin(BV_MatProject,X,A,Y,0);
1211: /* temporarily set standard inner product */
1212: Xmatrix = X->matrix;
1213: Ymatrix = Y->matrix;
1214: X->matrix = Y->matrix = NULL;
1216: PetscObjectGetId((PetscObject)X,&idx);
1217: PetscObjectGetId((PetscObject)Y,&idy);
1218: if (A && idx==idy) { /* check symmetry of M=X'AX */
1219: MatIsHermitianKnown(A,&set,&flg);
1220: symm = set? flg: PETSC_FALSE;
1221: }
1223: MatDenseGetArray(M,&marray);
1225: if (A) {
1226: if (X->vmm==BV_MATMULT_VECS) {
1227: /* perform computation column by column */
1228: BVMatProject_Vec(X,A,Y,marray,m,symm);
1229: } else {
1230: /* use BVMatMult, then BVDot */
1231: MatHasOperation(A,MATOP_MULT_TRANSPOSE,&flg);
1232: if (symm || (flg && X->l>=X->k/2 && Y->l>=Y->k/2)) {
1233: BVMatProject_MatMult_2(X,A,Y,marray,m,symm);
1234: } else {
1235: BVMatProject_MatMult(X,A,Y,marray,m);
1236: }
1237: }
1238: } else {
1239: /* use BVDot on subblocks */
1240: BVMatProject_Dot(X,Y,marray,m);
1241: }
1243: MatDenseRestoreArray(M,&marray);
1244: PetscLogEventEnd(BV_MatProject,X,A,Y,0);
1245: /* restore non-standard inner product */
1246: X->matrix = Xmatrix;
1247: Y->matrix = Ymatrix;
1248: return(0);
1249: }