Actual source code: orthog.c
1: /*
2: EPS routines related to orthogonalization.
3: See the SLEPc users manual for a detailed explanation.
4: */
5: #include src/eps/epsimpl.h
9: /*@
10: EPSQRDecomposition - Compute the QR factorization of a set of vectors.
12: Collective on EPS
14: Input Parameters:
15: + eps - the eigenproblem solver context
16: . V - set of vectors
17: . m - starting column
18: . n - ending column
19: - ldr - leading dimension of R
21: Output Parameter:
22: . R - triangular matrix of QR factorization
24: Notes:
25: This routine orthonormalizes the columns of V so that V'*V=I up to
26: working precision. It assumes that the first m columns (from 0 to m-1)
27: are already orthonormal. The coefficients of the orthogonalization are
28: stored in R so that V*R is equal to the original V.
30: In some cases, this routine makes V B-orthonormal, that is, V'*B*V=I.
32: If one of the vectors is linearly dependent wrt the rest (at working
33: precision) then it is replaced by a random vector.
35: Level: developer
37: .seealso: EPSOrthogonalize(), STNorm(), STInnerProduct().
38: @*/
39: PetscErrorCode EPSQRDecomposition(EPS eps,Vec *V,int m,int n,PetscScalar *R,int ldr)
40: {
41: PetscErrorCode ierr;
42: int k;
43: PetscScalar alpha;
44: PetscReal norm;
45: PetscTruth breakdown;
46:
49: for (k=m; k<n; k++) {
51: /* orthogonalize v_k with respect to v_0, ..., v_{k-1} */
52: if (R) { EPSOrthogonalize(eps,k,V,V[k],&R[0+ldr*k],&norm,&breakdown); }
53: else { EPSOrthogonalize(eps,k,V,V[k],PETSC_NULL,&norm,&breakdown); }
55: /* normalize v_k: r_{k,k} = ||v_k||_2; v_k = v_k/r_{k,k} */
56: if (norm==0.0 || breakdown) {
57: PetscLogInfo(eps,"EPSQRDecomposition: Linearly dependent vector found, generating a new random vector\n" );
58: SlepcVecSetRandom(V[k]);
59: STNorm(eps->OP,V[k],&norm);
60: }
61: alpha = 1.0/norm;
62: VecScale(&alpha,V[k]);
63: if (R) R[k+ldr*k] = norm;
65: }
67: return(0);
68: }
70: /*
71: Orthogonalization routine using classical Gram-Schmidt with refinement.
72: */
75: static PetscErrorCode EPSClassicalGramSchmidtOrthogonalization(EPS eps,int n,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscTruth *breakdown)
76: {
77: PetscErrorCode ierr;
78: int j;
79: PetscScalar shh[100],shh2[100],*lhh;
80: PetscTruth alloc = PETSC_FALSE;
81: PetscReal hnorm = 0,lnorm;
85: if (!H) {
86: if (n<=100) H = shh2; /* Don't allocate small arrays */
87: else {
88: PetscMalloc(n*sizeof(PetscScalar),&H);
89: alloc = PETSC_TRUE;
90: }
91: }
92: /* Don't allocate small arrays */
93: if (n<=100) lhh = shh;
94: else { PetscMalloc(n*sizeof(PetscScalar),&lhh); }
95: if (!norm) { norm = &lnorm; }
96:
97: /*** First orthogonalization ***/
99: /* h = V^* v */
100: /* q = v - V h */
101: for (j=0;j<n;j++) {
102: STInnerProduct(eps->OP,v,V[j],&H[j]);
103: lhh[j] = -H[j];
104: }
105: VecMAXPY(n,lhh,v,V);
106: STNorm(eps->OP,v,norm);
108: /* compute hnorm */
109: if (eps->orthog_ref == EPS_ORTH_REFINE_IFNEEDED || breakdown) {
110: hnorm = 0.0;
111: for (j=0; j<n; j++) {
112: hnorm += PetscRealPart(H[j] * PetscConj(H[j]));
113: }
114: hnorm = sqrt(hnorm);
115: }
116:
117: /*** Second orthogonalization if necessary ***/
118:
119: /* if ||q|| < eta ||h|| */
120: if (eps->orthog_ref == EPS_ORTH_REFINE_ALWAYS ||
121: (eps->orthog_ref == EPS_ORTH_REFINE_IFNEEDED &&
122: *norm < eps->orthog_eta * hnorm)) {
123: PetscLogInfo(eps,"EPSClassicalGramSchmidtOrthogonalization:Performing iterative refinement wnorm %g hnorm %g\n",*norm,hnorm);
125: /* s = V^* q */
126: /* q = q - V s ; h = h + s */
127: for (j=0;j<n;j++) {
128: STInnerProduct(eps->OP,v,V[j],&lhh[j]);
129: H[j] += lhh[j];
130: lhh[j] = -lhh[j];
131: }
132: VecMAXPY(n,lhh,v,V);
134: hnorm = *norm;
135: STNorm(eps->OP,v,norm);
136: }
137:
138: /* check breakdown */
139: if (breakdown) {
140: if (*norm < eps->orthog_eta * hnorm) *breakdown = PETSC_TRUE;
141: else *breakdown = PETSC_FALSE;
142: }
144: if (n>100) { PetscFree(lhh); }
145: if (alloc) { PetscFree(H); }
146: return(0);
147: }
149: /*
150: Orthogonalization routine using modified Gram-Schmidt with refinement.
151: */
154: PetscErrorCode EPSModifiedGramSchmidtOrthogonalization(EPS eps,int n,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscTruth *breakdown)
155: {
156: PetscErrorCode ierr;
157: int j;
158: PetscScalar alpha;
159: PetscTruth allocated;
160: PetscScalar lh[100],*h;
161: PetscReal hnorm = 0,lnorm;
162:
164:
165: allocated = PETSC_FALSE;
166: if (!H) {
167: if (n<=100) h = lh;
168: else {
169: PetscMalloc(n*sizeof(PetscScalar),&h);
170: allocated = PETSC_TRUE;
171: }
172: } else h = H;
173:
174: if (!norm) { norm = &lnorm; }
176: for (j=0; j<n; j++) {
177: /* alpha = ( v, v_j ) */
178: STInnerProduct(eps->OP,v,V[j],&alpha);
179: /* store coefficients if requested */
180: h[j] = alpha;
181: /* v <- v - alpha v_j */
182: alpha = -alpha;
183: VecAXPY(&alpha,V[j],v);
184: }
185: STNorm(eps->OP,v,norm);
187: if (eps->orthog_ref == EPS_ORTH_REFINE_IFNEEDED || breakdown) {
188: hnorm = 0.0;
189: for (j=0; j<n; j++) {
190: hnorm += PetscRealPart(h[j] * PetscConj(h[j]));
191: }
192: hnorm = sqrt(hnorm);
193: }
194:
195: if (eps->orthog_ref == EPS_ORTH_REFINE_ALWAYS ||
196: (eps->orthog_ref == EPS_ORTH_REFINE_IFNEEDED &&
197: *norm < eps->orthog_eta * hnorm)) {
198: PetscLogInfo(eps,"EPSModifiedGramSchmidtOrthogonalization:Performing iterative refinement wnorm %g hnorm %g\n",*norm,hnorm);
199: for (j=0; j<n; j++) {
200: /* alpha = ( v, v_j ) */
201: STInnerProduct(eps->OP,v,V[j],&alpha);
202: /* store coefficients if requested */
203: h[j] += alpha;
204: /* v <- v - alpha v_j */
205: alpha = -alpha;
206: VecAXPY(&alpha,V[j],v);
207: }
208: hnorm = *norm;
209: STNorm(eps->OP,v,norm);
210: }
212: /* check breakdown */
213: if (breakdown) {
214: if (*norm < eps->orthog_eta * hnorm) *breakdown = PETSC_TRUE;
215: else *breakdown = PETSC_FALSE;
216: }
217:
218: if (allocated) {
219: PetscFree(h);
220: }
221: return(0);
222: }
226: /*@
227: EPSPurge - Purge a vector of all converged vectors.
229: Collective on EPS
231: Input Parameters:
232: . eps - the eigenproblem solver context
234: Input/Output Parameter:
235: . v - vector to be purged
237: Notes:
238: On exit, v is orthogonal to all the basis vectors of the currently
239: converged invariant subspace as well as all the deflation vectors
240: provided by the user.
242: This routine does not normalize the resulting vector.
244: Level: developer
246: .seealso: EPSOrthogonalize(), EPSAttachDeflationSpace()
247: @*/
248: PetscErrorCode EPSPurge(EPS eps,Vec v)
249: {
250: PetscErrorCode ierr;
252: EPSOrthogonalize(eps,eps->nds+eps->nconv,eps->DSV,v,PETSC_NULL,PETSC_NULL,PETSC_NULL);
253: return(0);
254: }
258: /*@
259: EPSOrthogonalize - Orthogonalize a vector with respect to a set of vectors.
261: Collective on EPS
263: Input Parameters:
264: + eps - the eigenproblem solver context
265: . n - number of columns of V
266: - V - set of vectors
268: Input/Output Parameter:
269: . v - vector to be orthogonalized
271: Output Parameter:
272: + H - coefficients computed during orthogonalization
273: . norm - norm of the vector ofter being orthogonalized
274: - breakdown - flag indicating that refinement did not improve the quality
275: of orthogonalization
277: Notes:
278: On exit, v0 = [V v]*H, where v0 is the original vector v.
280: This routine does not normalize the resulting vector.
282: Level: developer
284: .seealso: EPSSetOrthogonalization()
285: @*/
286: PetscErrorCode EPSOrthogonalize(EPS eps,int n,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscTruth *breakdown)
287: {
288: PetscErrorCode ierr;
290: if (n==0) {
291: if (norm) { STNorm(eps->OP,v,norm); }
292: if (breakdown) *breakdown = PETSC_FALSE;
293: } else {
294: PetscLogEventBegin(EPS_Orthogonalization,eps,0,0,0);
295: switch (eps->orthog_type) {
296: case EPS_CGS_ORTH:
297: EPSClassicalGramSchmidtOrthogonalization(eps,n,V,v,H,norm,breakdown);
298: break;
299: case EPS_MGS_ORTH:
300: EPSModifiedGramSchmidtOrthogonalization(eps,n,V,v,H,norm,breakdown);
301: break;
302: default:
303: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
304: }
305: PetscLogEventEnd(EPS_Orthogonalization,eps,0,0,0);
306: }
307: return(0);
308: }