Actual source code: orthog.c
1: /*
2: EPS routines related to orthogonalization.
3: See the SLEPc Technical Report STR-1 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,PETSC_NULL,V,V[k],&R[0+ldr*k],&norm,&lindep); }
53: else { EPSOrthogonalize(eps,k,PETSC_NULL,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: }
69: /*
70: EPSOrthogonalizeGS - Compute |v'|, |v| and one step of CGS or MGS
71: */
74: static PetscErrorCode EPSOrthogonalizeGS(EPS eps,int n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H,PetscReal *onorm,PetscReal *norm,Vec w)
75: {
77: int j;
78: PetscScalar alpha;
79: PetscReal sum;
80:
82: switch (eps->orthog_type) {
83:
84: case EPS_CGS_ORTH:
85: /* h = W^* v ; alpha = (v , v) */
86: if (which) { /* use select array */
87: for (j=0; j<n; j++)
88: if (which[j]) {
89: STInnerProductBegin(eps->OP,v,V[j],&H[j]);
90: }
91: if (onorm || norm) {
92: STInnerProductBegin(eps->OP,v,v,&alpha);
93: }
94: for (j=0; j<n; j++)
95: if (which[j]) { STInnerProductEnd(eps->OP,v,V[j],&H[j]); }
96: if (onorm || norm) { STInnerProductEnd(eps->OP,v,v,&alpha); }
97: } else { /* merge comunications */
98: if (onorm || norm) {
99: STMInnerProductBegin(eps->OP,n,v,V,H);
100: STInnerProductBegin(eps->OP,v,v,&alpha);
101: STMInnerProductEnd(eps->OP,n,v,V,H);
102: STInnerProductEnd(eps->OP,v,v,&alpha);
103: } else { /* use simpler function */
104: STMInnerProduct(eps->OP,n,v,V,H);
105: }
106: }
108: /* q = v - V h */
109: VecSet(w,0.0);
110: if (which) {
111: for (j=0; j<n; j++)
112: if (which[j]) { VecAXPY(w,H[j],V[j]); }
113: } else {
114: VecMAXPY(w,n,H,V);
115: }
116: VecAXPY(v,-1.0,w);
117:
118: /* compute |v| and |v'| */
119: if (onorm) *onorm = sqrt(PetscRealPart(alpha));
120: if (norm) {
121: sum = 0.0;
122: for (j=0; j<n; j++)
123: if (!which || which[j])
124: sum += PetscRealPart(H[j] * PetscConj(H[j]));
125: *norm = PetscRealPart(alpha)-sum;
126: if (*norm < 0.0) {
127: STNorm(eps->OP,v,norm);
128: } else *norm = sqrt(*norm);
129: }
130: break;
131:
132: case EPS_MGS_ORTH:
133: /* compute |v| */
134: if (onorm) { STNorm(eps->OP,v,onorm); }
135: for (j=0; j<n; j++)
136: if (!which || which[j]) {
137: /* h_j = ( v, v_j ) */
138: STInnerProduct(eps->OP,v,V[j],&H[j]);
139: /* v <- v - h_j v_j */
140: VecAXPY(v,-H[j],V[j]);
141: }
142: /* compute |v'| */
143: if (norm) { STNorm(eps->OP,v,norm); }
144: break;
145:
146: default:
147: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
148: }
149: return(0);
150: }
154: /*@
155: EPSOrthogonalize - Orthogonalize a vector with respect to a set of vectors.
157: Collective on EPS
159: Input Parameters:
160: + eps - the eigenproblem solver context
161: . n - number of columns of V
162: . which - logical array indicating columns of V to be used
163: - V - set of vectors
165: Input/Output Parameter:
166: . v - (input) vector to be orthogonalized and (output) result of
167: orthogonalization
169: Output Parameter:
170: + H - coefficients computed during orthogonalization
171: . norm - norm of the vector after being orthogonalized
172: - lindep - flag indicating that refinement did not improve the quality
173: of orthogonalization
175: Notes:
176: This function applies an orthogonal projector to project vector v onto the
177: orthogonal complement of the span of the columns of V.
179: On exit, v0 = [V v]*H, where v0 is the original vector v.
181: This routine does not normalize the resulting vector.
183: Level: developer
185: .seealso: EPSSetOrthogonalization(), EPSBiOrthogonalize()
186: @*/
187: PetscErrorCode EPSOrthogonalize(EPS eps,int n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscTruth *lindep)
188: {
190: Vec w = PETSC_NULL;
191: PetscScalar lh[100],*h,lc[100],*c;
192: PetscTruth allocatedh = PETSC_FALSE,allocatedc = PETSC_FALSE;
193: PetscReal onrm,nrm;
194: int j,k;
196: if (n==0) {
197: if (norm) { STNorm(eps->OP,v,norm); }
198: if (lindep) *lindep = PETSC_FALSE;
199: return(0);
200: }
201: PetscLogEventBegin(EPS_Orthogonalize,eps,0,0,0);
203: /* allocate H, c and w if needed */
204: if (!H) {
205: if (n<=100) h = lh;
206: else {
207: PetscMalloc(n*sizeof(PetscScalar),&h);
208: allocatedh = PETSC_TRUE;
209: }
210: } else h = H;
211: if (eps->orthog_ref != EPS_ORTH_REFINE_NEVER) {
212: if (n<=100) c = lc;
213: else {
214: PetscMalloc(n*sizeof(PetscScalar),&c);
215: allocatedc = PETSC_TRUE;
216: }
217: }
218: if (eps->orthog_type != EPS_MGS_ORTH) {
219: VecDuplicate(v,&w);
220: }
222: /* orthogonalize and compute onorm */
223: switch (eps->orthog_ref) {
224:
225: case EPS_ORTH_REFINE_NEVER:
226: EPSOrthogonalizeGS(eps,n,which,V,v,h,PETSC_NULL,PETSC_NULL,w);
227: /* compute |v| */
228: if (norm) { STNorm(eps->OP,v,norm); }
229: /* linear dependence check does not work without refinement */
230: if (lindep) *lindep = PETSC_FALSE;
231: break;
232:
233: case EPS_ORTH_REFINE_ALWAYS:
234: EPSOrthogonalizeGS(eps,n,which,V,v,h,PETSC_NULL,PETSC_NULL,w);
235: if (lindep) {
236: EPSOrthogonalizeGS(eps,n,which,V,v,c,&onrm,&nrm,w);
237: if (norm) *norm = nrm;
238: if (nrm < eps->orthog_eta * onrm) *lindep = PETSC_TRUE;
239: else *lindep = PETSC_FALSE;
240: } else {
241: EPSOrthogonalizeGS(eps,n,which,V,v,c,PETSC_NULL,norm,w);
242: }
243: for (j=0;j<n;j++)
244: if (!which || which[j]) h[j] += c[j];
245: break;
246:
247: case EPS_ORTH_REFINE_IFNEEDED:
248: EPSOrthogonalizeGS(eps,n,which,V,v,h,&onrm,&nrm,w);
249: /* ||q|| < eta ||h|| */
250: k = 1;
251: while (k<3 && nrm < eps->orthog_eta * onrm) {
252: k++;
253: switch (eps->orthog_type) {
254: case EPS_CGS_ORTH:
255: EPSOrthogonalizeGS(eps,n,which,V,v,c,&onrm,&nrm,w);
256: break;
257: case EPS_MGS_ORTH:
258: onrm = nrm;
259: EPSOrthogonalizeGS(eps,n,which,V,v,c,PETSC_NULL,&nrm,w);
260: break;
261: default:
262: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
263: }
264: for (j=0;j<n;j++)
265: if (!which || which[j]) h[j] += c[j];
266: }
267: if (norm) *norm = nrm;
268: if (lindep) {
269: if (nrm < eps->orthog_eta * onrm) *lindep = PETSC_TRUE;
270: else *lindep = PETSC_FALSE;
271: }
272:
273: break;
274: default:
275: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization refinement");
276: }
278: /* free work space */
279: if (allocatedc) { PetscFree(c); }
280: if (allocatedh) { PetscFree(h); }
281: if (w) { VecDestroy(w); }
282:
283: PetscLogEventEnd(EPS_Orthogonalize,eps,0,0,0);
284: return(0);
285: }
287: /*
288: Biorthogonalization routine using classical Gram-Schmidt with refinement.
289: */
292: static PetscErrorCode EPSCGSBiOrthogonalization(EPS eps,int n,Vec *V,Vec *W,Vec v,PetscScalar *H,PetscReal *hnorm,PetscReal *norm)
293: {
294: #if defined(SLEPC_MISSING_LAPACK_GELQF) || defined(SLEPC_MISSING_LAPACK_ORMLQ)
296: SETERRQ(PETSC_ERR_SUP,"xGELQF - Lapack routine is unavailable.");
297: #else
299: int j,ione=1,lwork,info;
300: PetscScalar shh[100],*lhh,*vw,*tau,one=1.0,*work;
301: Vec w;
305: /* Don't allocate small arrays */
306: if (n<=100) lhh = shh;
307: else { PetscMalloc(n*sizeof(PetscScalar),&lhh); }
308: PetscMalloc(n*n*sizeof(PetscScalar),&vw);
309: VecDuplicate(v,&w);
310:
311: for (j=0;j<n;j++) {
312: STMInnerProduct(eps->OP,n,V[j],W,vw+j*n);
313: }
314: lwork = n;
315: PetscMalloc(n*sizeof(PetscScalar),&tau);
316: PetscMalloc(lwork*sizeof(PetscScalar),&work);
317: LAPACKgelqf_(&n,&n,vw,&n,tau,work,&lwork,&info);
318: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGELQF %i",info);
319:
320: /*** First orthogonalization ***/
322: /* h = W^* v */
323: /* q = v - V h */
324: STMInnerProduct(eps->OP,n,v,W,H);
325: BLAStrsm_("L","L","N","N",&n,&ione,&one,vw,&n,H,&n,1,1,1,1);
326: LAPACKormlq_("L","N",&n,&ione,&n,vw,&n,tau,H,&n,work,&lwork,&info,1,1);
327: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xORMLQ %i",info);
328: VecSet(w,0.0);
329: VecMAXPY(w,n,H,V);
330: VecAXPY(v,-1.0,w);
332: /* compute norm of v */
333: if (norm) { STNorm(eps->OP,v,norm); }
334:
335: if (n>100) { PetscFree(lhh); }
336: PetscFree(vw);
337: PetscFree(tau);
338: PetscFree(work);
339: VecDestroy(w);
340: return(0);
341: #endif
342: }
346: /*@
347: EPSBiOrthogonalize - Bi-orthogonalize a vector with respect to a set of vectors.
349: Collective on EPS
351: Input Parameters:
352: + eps - the eigenproblem solver context
353: . n - number of columns of V
354: . V - set of vectors
355: - W - set of vectors
357: Input/Output Parameter:
358: . v - vector to be orthogonalized
360: Output Parameter:
361: + H - coefficients computed during orthogonalization
362: - norm - norm of the vector after being orthogonalized
364: Notes:
365: This function applies an oblique projector to project vector v onto the
366: span of the columns of V along the orthogonal complement of the column
367: space of W.
369: On exit, v0 = [V v]*H, where v0 is the original vector v.
371: This routine does not normalize the resulting vector.
373: Level: developer
375: .seealso: EPSSetOrthogonalization(), EPSOrthogonalize()
376: @*/
377: PetscErrorCode EPSBiOrthogonalize(EPS eps,int n,Vec *V,Vec *W,Vec v,PetscScalar *H,PetscReal *norm)
378: {
380: PetscScalar lh[100],*h;
381: PetscTruth allocated = PETSC_FALSE;
382: PetscReal lhnrm,*hnrm,lnrm,*nrm;
384: if (n==0) {
385: if (norm) { STNorm(eps->OP,v,norm); }
386: } else {
387: PetscLogEventBegin(EPS_Orthogonalize,eps,0,0,0);
388:
389: /* allocate H if needed */
390: if (!H) {
391: if (n<=100) h = lh;
392: else {
393: PetscMalloc(n*sizeof(PetscScalar),&h);
394: allocated = PETSC_TRUE;
395: }
396: } else h = H;
397:
398: /* retrieve hnrm and nrm for linear dependence check or conditional refinement */
399: if (eps->orthog_ref == EPS_ORTH_REFINE_IFNEEDED) {
400: hnrm = &lhnrm;
401: if (norm) nrm = norm;
402: else nrm = &lnrm;
403: } else {
404: hnrm = PETSC_NULL;
405: nrm = norm;
406: }
407:
408: switch (eps->orthog_type) {
409: case EPS_CGS_ORTH:
410: EPSCGSBiOrthogonalization(eps,n,V,W,v,h,hnrm,nrm);
411: break;
412: default:
413: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
414: }
415:
416: if (allocated) { PetscFree(h); }
417:
418: PetscLogEventEnd(EPS_Orthogonalize,eps,0,0,0);
419: }
420: return(0);
421: }