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
6: #include slepcblaslapack.h
10: /*@
11: EPSQRDecomposition - Compute the QR factorization of a set of vectors.
13: Collective on EPS
15: Input Parameters:
16: + eps - the eigenproblem solver context
17: . V - set of vectors
18: . m - starting column
19: . n - ending column
20: - ldr - leading dimension of R
22: Output Parameter:
23: . R - triangular matrix of QR factorization
25: Notes:
26: This routine orthonormalizes the columns of V so that V'*V=I up to
27: working precision. It assumes that the first m columns (from 0 to m-1)
28: are already orthonormal. The coefficients of the orthogonalization are
29: stored in R so that V*R is equal to the original V.
31: In some cases, this routine makes V B-orthonormal, that is, V'*B*V=I.
33: If one of the vectors is linearly dependent wrt the rest (at working
34: precision) then it is replaced by a random vector.
36: Level: developer
38: .seealso: EPSOrthogonalize(), STNorm(), STInnerProduct().
39: @*/
40: PetscErrorCode EPSQRDecomposition(EPS eps,Vec *V,int m,int n,PetscScalar *R,int ldr)
41: {
43: int k;
44: PetscReal norm;
45: PetscTruth lindep;
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,&lindep); }
53: else { EPSOrthogonalize(eps,k,V,V[k],PETSC_NULL,&norm,&lindep); }
55: /* normalize v_k: r_{k,k} = ||v_k||_2; v_k = v_k/r_{k,k} */
56: if (norm==0.0 || lindep) {
57: PetscInfo(eps,"Linearly dependent vector found, generating a new random vector\n");
58: SlepcVecSetRandom(V[k]);
59: STNorm(eps->OP,V[k],&norm);
60: }
61: VecScale(V[k],1.0/norm);
62: if (R) R[k+ldr*k] = norm;
64: }
66: return(0);
67: }
71: PetscErrorCode EPSOrthogonalizeGS(EPS eps,int n,Vec *V,Vec v,PetscScalar *H,PetscReal *onorm,PetscReal *norm,Vec w)
72: {
74: int j;
75: PetscScalar alpha;
76: PetscReal sum;
77:
79: /* orthogonalize */
80: switch (eps->orthog_type) {
81: case EPS_CGS_ORTH:
82: /* h = W^* v ; alpha = (v , v) */
83: STMInnerProductBegin(eps->OP,n,v,V,H);
84: if (onorm || norm) { STInnerProductBegin(eps->OP,v,v,&alpha); }
85: STMInnerProductEnd(eps->OP,n,v,V,H);
86: if (onorm || norm) { STInnerProductEnd(eps->OP,v,v,&alpha); }
87: /* q = v - V h */
88: VecSet(w,0.0);
89: VecMAXPY(w,n,H,V);
90: VecAXPY(v,-1.0,w);
91: break;
92: case EPS_MGS_ORTH:
93: if (onorm || norm) { STInnerProduct(eps->OP,v,v,&alpha); }
94: for (j=0; j<n; j++) {
95: /* h_j = ( v, v_j ) */
96: STInnerProduct(eps->OP,v,V[j],&H[j]);
97: /* v <- v - h_j v_j */
98: VecAXPY(v,-H[j],V[j]);
99: }
100: break;
101: default:
102: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
103: }
104: /* compute |v| and |v'| */
105: if (onorm) *onorm = sqrt(PetscRealPart(alpha));
106: if (norm) {
107: sum = 0.0;
108: for (j=0; j<n; j++)
109: sum += PetscRealPart(H[j] * PetscConj(H[j]));
110: *norm = PetscRealPart(alpha)-sum;
111: if (*norm < 0.0) {
112: STNorm(eps->OP,v,norm);
113: } else *norm = sqrt(*norm);
114: }
115: return(0);
116: }
120: /*@
121: EPSOrthogonalize - Orthogonalize a vector with respect to a set of vectors.
123: Collective on EPS
125: Input Parameters:
126: + eps - the eigenproblem solver context
127: . n - number of columns of V
128: - V - set of vectors
130: Input/Output Parameter:
131: . v - vector to be orthogonalized
133: Output Parameter:
134: + H - coefficients computed during orthogonalization
135: . norm - norm of the vector after being orthogonalized
136: - lindep - flag indicating that refinement did not improve the quality
137: of orthogonalization
139: Notes:
140: This function applies an orthogonal projector to project vector v onto the
141: orthogonal complement of the span of the columns of V.
143: On exit, v0 = [V v]*H, where v0 is the original vector v.
145: This routine does not normalize the resulting vector.
147: Level: developer
149: .seealso: EPSSetOrthogonalization(), EPSBiOrthogonalize()
150: @*/
151: PetscErrorCode EPSOrthogonalize(EPS eps,int n,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscTruth *lindep)
152: {
154: Vec w = PETSC_NULL;
155: PetscScalar lh[100],*h,lc[100],*c;
156: PetscTruth allocatedh = PETSC_FALSE,allocatedc = PETSC_FALSE;
157: PetscReal lonrm,*onrm,lnrm,*nrm;
158: int j,k;
160: if (n==0) {
161: if (norm) { STNorm(eps->OP,v,norm); }
162: if (lindep) *lindep = PETSC_FALSE;
163: return(0);
164: }
165: PetscLogEventBegin(EPS_Orthogonalize,eps,0,0,0);
167: /* allocate H, c and w if needed */
168: if (!H) {
169: if (n<=100) h = lh;
170: else {
171: PetscMalloc(n*sizeof(PetscScalar),&h);
172: allocatedh = PETSC_TRUE;
173: }
174: } else h = H;
175: if (eps->orthog_ref !=EPS_ORTH_REFINE_NEVER) {
176: if (n<=100) c = lc;
177: else {
178: PetscMalloc(n*sizeof(PetscScalar),&c);
179: allocatedc = PETSC_TRUE;
180: }
181: }
182: if (eps->orthog_type != EPS_MGS_ORTH) {
183: VecDuplicate(v,&w);
184: }
186: /* retrieve onrm and nrm for linear dependence check or conditional refinement */
187: if (lindep || eps->orthog_ref == EPS_ORTH_REFINE_IFNEEDED) {
188: onrm = &lonrm;
189: if (norm) nrm = norm;
190: else nrm = &lnrm;
191: } else {
192: onrm = PETSC_NULL;
193: nrm = norm;
194: }
196: /* orthogonalize and compute onorm */
197: switch (eps->orthog_ref) {
198: case EPS_ORTH_REFINE_NEVER:
199: EPSOrthogonalizeGS(eps,n,V,v,h,onrm,PETSC_NULL,w);
200: /* compute |v| */
201: if (nrm) { STNorm(eps->OP,v,nrm); }
202: break;
203: case EPS_ORTH_REFINE_ALWAYS:
204: EPSOrthogonalizeGS(eps,n,V,v,h,PETSC_NULL,PETSC_NULL,w);
205: eps->count_reorthog++;
206: EPSOrthogonalizeGS(eps,n,V,v,c,onrm,nrm,w);
207: for (j=0;j<n;j++)
208: h[j] += c[j];
209: break;
210: case EPS_ORTH_REFINE_IFNEEDED:
211: EPSOrthogonalizeGS(eps,n,V,v,h,onrm,nrm,w);
212: /* ||q|| < eta ||h|| */
213: k = 1;
214: while (k<3 && *nrm < eps->orthog_eta * *onrm) {
215: k++;
216: eps->count_reorthog++;
217: EPSOrthogonalizeGS(eps,n,V,v,c,onrm,nrm,w);
218: for (j=0;j<n;j++)
219: h[j] += c[j];
220: }
221: break;
222: default:
223: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization refinement");
224: }
226: /* check linear dependence */
227: if (lindep) {
228: if (*nrm < eps->orthog_eta * *onrm) *lindep = PETSC_TRUE;
229: else *lindep = PETSC_FALSE;
230: }
232: /* free work space */
233: if (allocatedc) { PetscFree(c); }
234: if (allocatedh) { PetscFree(h); }
235: if (w) { VecDestroy(w); }
236:
237: PetscLogEventEnd(EPS_Orthogonalize,eps,0,0,0);
238: return(0);
239: }
241: /*
242: Biorthogonalization routine using classical Gram-Schmidt with refinement.
243: */
246: static PetscErrorCode EPSCGSBiOrthogonalization(EPS eps,int n,Vec *V,Vec *W,Vec v,PetscScalar *H,PetscReal *hnorm,PetscReal *norm)
247: {
248: #if defined(SLEPC_MISSING_LAPACK_GELQF) || defined(SLEPC_MISSING_LAPACK_ORMLQ)
250: SETERRQ(PETSC_ERR_SUP,"xGELQF - Lapack routine is unavailable.");
251: #else
253: int j,ione=1,lwork,info;
254: PetscScalar shh[100],*lhh,*vw,*tau,one=1.0,*work;
255: Vec w;
259: /* Don't allocate small arrays */
260: if (n<=100) lhh = shh;
261: else { PetscMalloc(n*sizeof(PetscScalar),&lhh); }
262: PetscMalloc(n*n*sizeof(PetscScalar),&vw);
263: VecDuplicate(v,&w);
264:
265: for (j=0;j<n;j++) {
266: STMInnerProduct(eps->OP,n,V[j],W,vw+j*n);
267: }
268: lwork = n;
269: PetscMalloc(n*sizeof(PetscScalar),&tau);
270: PetscMalloc(lwork*sizeof(PetscScalar),&work);
271: LAPACKgelqf_(&n,&n,vw,&n,tau,work,&lwork,&info);
272: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGELQF %i",info);
273:
274: /*** First orthogonalization ***/
276: /* h = W^* v */
277: /* q = v - V h */
278: STMInnerProduct(eps->OP,n,v,W,H);
279: BLAStrsm_("L","L","N","N",&n,&ione,&one,vw,&n,H,&n,1,1,1,1);
280: LAPACKormlq_("L","N",&n,&ione,&n,vw,&n,tau,H,&n,work,&lwork,&info,1,1);
281: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xORMLQ %i",info);
282: VecSet(w,0.0);
283: VecMAXPY(w,n,H,V);
284: VecAXPY(v,-1.0,w);
286: /* compute norm of v */
287: if (norm) { STNorm(eps->OP,v,norm); }
288:
289: if (n>100) { PetscFree(lhh); }
290: PetscFree(vw);
291: PetscFree(tau);
292: PetscFree(work);
293: VecDestroy(w);
294: return(0);
295: #endif
296: }
300: /*@
301: EPSBiOrthogonalize - Bi-orthogonalize a vector with respect to a set of vectors.
303: Collective on EPS
305: Input Parameters:
306: + eps - the eigenproblem solver context
307: . n - number of columns of V
308: . V - set of vectors
309: - W - set of vectors
311: Input/Output Parameter:
312: . v - vector to be orthogonalized
314: Output Parameter:
315: + H - coefficients computed during orthogonalization
316: - norm - norm of the vector after being orthogonalized
318: Notes:
319: This function applies an oblique projector to project vector v onto the
320: span of the columns of V along the orthogonal complement of the column
321: space of W.
323: On exit, v0 = [V v]*H, where v0 is the original vector v.
325: This routine does not normalize the resulting vector.
327: Level: developer
329: .seealso: EPSSetOrthogonalization(), EPSOrthogonalize()
330: @*/
331: PetscErrorCode EPSBiOrthogonalize(EPS eps,int n,Vec *V,Vec *W,Vec v,PetscScalar *H,PetscReal *norm)
332: {
334: PetscScalar lh[100],*h;
335: PetscTruth allocated = PETSC_FALSE;
336: PetscReal lhnrm,*hnrm,lnrm,*nrm;
338: if (n==0) {
339: if (norm) { STNorm(eps->OP,v,norm); }
340: } else {
341: PetscLogEventBegin(EPS_Orthogonalize,eps,0,0,0);
342:
343: /* allocate H if needed */
344: if (!H) {
345: if (n<=100) h = lh;
346: else {
347: PetscMalloc(n*sizeof(PetscScalar),&h);
348: allocated = PETSC_TRUE;
349: }
350: } else h = H;
351:
352: /* retrieve hnrm and nrm for linear dependence check or conditional refinement */
353: if (eps->orthog_ref == EPS_ORTH_REFINE_IFNEEDED) {
354: hnrm = &lhnrm;
355: if (norm) nrm = norm;
356: else nrm = &lnrm;
357: } else {
358: hnrm = PETSC_NULL;
359: nrm = norm;
360: }
361:
362: switch (eps->orthog_type) {
363: case EPS_CGS_ORTH:
364: EPSCGSBiOrthogonalization(eps,n,V,W,v,h,hnrm,nrm);
365: break;
366: default:
367: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
368: }
369:
370: if (allocated) { PetscFree(h); }
371:
372: PetscLogEventEnd(EPS_Orthogonalize,eps,0,0,0);
373: }
374: return(0);
375: }