Actual source code: bvbiorthog.c
slepc-3.17.2 2022-08-09
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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 bi-orthogonalization routines
12: */
14: #include <slepc/private/bvimpl.h>
16: /*
17: BVBiorthogonalizeMGS1 - Compute one step of Modified Gram-Schmidt bi-orthogonalization
18: */
19: static PetscErrorCode BVBiorthogonalizeMGS1(BV V,BV W,Vec v,PetscScalar *h,PetscScalar *c)
20: {
21: PetscInt i;
22: PetscScalar dot;
23: Vec vi,wi;
25: for (i=-V->nc;i<V->k;i++) {
26: BVGetColumn(W,i,&wi);
27: /* h_i = (v, w_i) */
28: VecDot(v,wi,&dot);
29: BVRestoreColumn(W,i,&wi);
30: /* v <- v - h_i v_i */
31: BV_SetValue(V,i,0,c,dot);
32: BVGetColumn(V,i,&vi);
33: VecAXPY(v,-dot,vi);
34: BVRestoreColumn(V,i,&vi);
35: }
36: BV_AddCoefficients(V,V->k,h,c);
37: PetscFunctionReturn(0);
38: }
40: /*
41: BVBiorthogonalizeCGS1 - Compute one step of CGS bi-orthogonalization: v = (I-V*W')*v
42: */
43: static PetscErrorCode BVBiorthogonalizeCGS1(BV V,BV W,Vec v,PetscScalar *h,PetscScalar *c)
44: {
45: /* h = W'*v */
46: BVDotVec(W,v,c);
48: /* v = v - V h */
49: BVMultVec(V,-1.0,1.0,v,c);
51: BV_AddCoefficients(V,V->k,h,c);
52: PetscFunctionReturn(0);
53: }
55: #define BVBiorthogonalizeGS1(a,b,c,d,e) ((V->orthog_type==BV_ORTHOG_MGS)?BVBiorthogonalizeMGS1:BVBiorthogonalizeCGS1)(a,b,c,d,e)
57: /*
58: BVBiorthogonalizeGS - Orthogonalize with (classical or modified) Gram-Schmidt
60: V, W - the two basis vectors objects
61: v - the vector to bi-orthogonalize
62: */
63: static PetscErrorCode BVBiorthogonalizeGS(BV V,BV W,Vec v)
64: {
65: PetscScalar *h,*c;
67: h = V->h;
68: c = V->c;
69: BV_CleanCoefficients(V,V->k,h);
70: BVBiorthogonalizeGS1(V,W,v,h,c);
71: if (V->orthog_ref!=BV_ORTHOG_REFINE_NEVER) BVBiorthogonalizeGS1(V,W,v,h,c);
72: PetscFunctionReturn(0);
73: }
75: /*@
76: BVBiorthogonalizeColumn - Bi-orthogonalize a column of two BV objects.
78: Collective on V
80: Input Parameters:
81: + V - first basis vectors context
82: . W - second basis vectors context
83: - j - index of column to be bi-orthonormalized
85: Notes:
86: This function bi-orthogonalizes vectors V[j],W[j] against W[0..j-1],
87: and V[0..j-1], respectively, so that W[0..j]'*V[0..j] = diagonal.
89: Level: advanced
91: .seealso: BVOrthogonalizeColumn(), BVBiorthonormalizeColumn()
92: @*/
93: PetscErrorCode BVBiorthogonalizeColumn(BV V,BV W,PetscInt j)
94: {
95: PetscInt ksavev,lsavev,ksavew,lsavew;
96: Vec y,z;
102: BVCheckSizes(V,1);
104: BVCheckSizes(W,2);
114: /* bi-orthogonalize */
115: PetscLogEventBegin(BV_OrthogonalizeVec,V,0,0,0);
116: ksavev = V->k;
117: lsavev = V->l;
118: ksavew = W->k;
119: lsavew = W->l;
120: V->k = j;
121: V->l = -V->nc; /* must also bi-orthogonalize against constraints and leading columns */
122: W->k = j;
123: W->l = -W->nc;
124: BV_AllocateCoeffs(V);
125: BV_AllocateCoeffs(W);
126: BVGetColumn(V,j,&y);
127: BVBiorthogonalizeGS(V,W,y);
128: BVRestoreColumn(V,j,&y);
129: BVGetColumn(W,j,&z);
130: BVBiorthogonalizeGS(W,V,z);
131: BVRestoreColumn(W,j,&z);
132: V->k = ksavev;
133: V->l = lsavev;
134: W->k = ksavew;
135: W->l = lsavew;
136: PetscLogEventEnd(BV_OrthogonalizeVec,V,0,0,0);
137: PetscObjectStateIncrease((PetscObject)V);
138: PetscObjectStateIncrease((PetscObject)W);
139: PetscFunctionReturn(0);
140: }
142: /*@
143: BVBiorthonormalizeColumn - Bi-orthonormalize a column of two BV objects.
145: Collective on V
147: Input Parameters:
148: + V - first basis vectors context
149: . W - second basis vectors context
150: - j - index of column to be bi-orthonormalized
152: Output Parameters:
153: . delta - (optional) value used for normalization
155: Notes:
156: This function first bi-orthogonalizes vectors V[j],W[j] against W[0..j-1],
157: and V[0..j-1], respectively. Then, it scales the vectors with 1/delta, so
158: that the resulting vectors satisfy W[j]'*V[j] = 1.
160: Level: advanced
162: .seealso: BVOrthonormalizeColumn(), BVBiorthogonalizeColumn()
163: @*/
164: PetscErrorCode BVBiorthonormalizeColumn(BV V,BV W,PetscInt j,PetscReal *delta)
165: {
166: PetscScalar alpha;
167: PetscReal deltat;
168: PetscInt ksavev,lsavev,ksavew,lsavew;
169: Vec y,z;
175: BVCheckSizes(V,1);
177: BVCheckSizes(W,2);
187: /* bi-orthogonalize */
188: PetscLogEventBegin(BV_OrthogonalizeVec,V,0,0,0);
189: ksavev = V->k;
190: lsavev = V->l;
191: ksavew = W->k;
192: lsavew = W->l;
193: V->k = j;
194: V->l = -V->nc; /* must also bi-orthogonalize against constraints and leading columns */
195: W->k = j;
196: W->l = -W->nc;
197: BV_AllocateCoeffs(V);
198: BV_AllocateCoeffs(W);
199: BVGetColumn(V,j,&y);
200: BVBiorthogonalizeGS(V,W,y);
201: BVRestoreColumn(V,j,&y);
202: BVGetColumn(W,j,&z);
203: BVBiorthogonalizeGS(W,V,z);
204: BVRestoreColumn(W,j,&z);
205: V->k = ksavev;
206: V->l = lsavev;
207: W->k = ksavew;
208: W->l = lsavew;
209: PetscLogEventEnd(BV_OrthogonalizeVec,V,0,0,0);
211: /* scale */
212: PetscLogEventBegin(BV_Scale,V,0,0,0);
213: BVGetColumn(V,j,&y);
214: BVGetColumn(W,j,&z);
215: VecDot(z,y,&alpha);
216: BVRestoreColumn(V,j,&y);
217: BVRestoreColumn(W,j,&z);
218: deltat = PetscSqrtReal(PetscAbsScalar(alpha));
219: if (V->n) (*V->ops->scale)(V,j,1.0/PetscConj(alpha/deltat));
220: if (W->n) (*W->ops->scale)(W,j,1.0/deltat);
221: PetscLogEventEnd(BV_Scale,V,0,0,0);
222: if (delta) *delta = deltat;
223: PetscObjectStateIncrease((PetscObject)V);
224: PetscObjectStateIncrease((PetscObject)W);
225: PetscFunctionReturn(0);
226: }