Actual source code: iporthog.c
1: /*
2: Routines related to orthogonalization.
3: See the SLEPc Technical Report STR-1 for a detailed explanation.
5: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
6: SLEPc - Scalable Library for Eigenvalue Problem Computations
7: Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
9: This file is part of SLEPc. See the README file for conditions of use
10: and additional information.
11: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
12: */
14: #include "src/ip/ipimpl.h" /*I "slepcip.h" I*/
15: #include "slepcblaslapack.h"
17: /*
18: IPOrthogonalizeMGS - Compute one step of Modified Gram-Schmidt
19: */
22: static PetscErrorCode IPOrthogonalizeMGS(IP ip,int n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H)
23: {
25: int j;
26:
28: for (j=0; j<n; j++)
29: if (!which || which[j]) {
30: /* h_j = ( v, v_j ) */
31: IPInnerProduct(ip,v,V[j],&H[j]);
32: /* v <- v - h_j v_j */
33: VecAXPY(v,-H[j],V[j]);
34: }
35: return(0);
36: }
38: /*
39: IPOrthogonalizeCGS - Compute |v'| (estimated), |v| and one step of CGS with only one global synchronization
40: */
43: PetscErrorCode IPOrthogonalizeCGS(IP ip,int n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H,PetscReal *onorm,PetscReal *norm,Vec work)
44: {
46: int j;
47: PetscScalar alpha;
48: PetscReal sum;
49:
51: /* h = W^* v ; alpha = (v , v) */
52: if (which) {
53: /* use select array */
54: for (j=0; j<n; j++)
55: if (which[j]) { IPInnerProductBegin(ip,v,V[j],&H[j]); }
56: if (onorm || norm) { IPInnerProductBegin(ip,v,v,&alpha); }
57: for (j=0; j<n; j++)
58: if (which[j]) { IPInnerProductEnd(ip,v,V[j],&H[j]); }
59: if (onorm || norm) { IPInnerProductEnd(ip,v,v,&alpha); }
60: } else { /* merge comunications */
61: if (onorm || norm) {
62: IPMInnerProductBegin(ip,v,n,V,H);
63: IPInnerProductBegin(ip,v,v,&alpha);
64: IPMInnerProductEnd(ip,v,n,V,H);
65: IPInnerProductEnd(ip,v,v,&alpha);
66: } else { /* use simpler function */
67: IPMInnerProduct(ip,v,n,V,H);
68: }
69: }
71: /* q = v - V h */
72: VecSet(work,0.0);
73: if (which) {
74: for (j=0; j<n; j++)
75: if (which[j]) { VecAXPY(work,H[j],V[j]); }
76: } else {
77: VecMAXPY(work,n,H,V);
78: }
79: VecAXPY(v,-1.0,work);
80:
81: /* compute |v| */
82: if (onorm) *onorm = sqrt(PetscRealPart(alpha));
84: /* compute |v'| */
85: if (norm) {
86: sum = 0.0;
87: for (j=0; j<n; j++)
88: if (!which || which[j])
89: sum += PetscRealPart(H[j] * PetscConj(H[j]));
90: *norm = PetscRealPart(alpha)-sum;
91: if (*norm <= 0.0) {
92: IPNorm(ip,v,norm);
93: } else *norm = sqrt(*norm);
94: }
95: return(0);
96: }
98: /*
99: IPOrthogonalizeGS - Compute |v'|, |v| and one step of CGS or MGS
100: */
103: static PetscErrorCode IPOrthogonalizeGS(IP ip,int n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H,PetscReal *onorm,PetscReal *norm,Vec work)
104: {
106:
108: switch (ip->orthog_type) {
109: case IP_CGS_ORTH:
110: IPOrthogonalizeCGS(ip,n,which,V,v,H,onorm,norm,work);
111: break;
112: case IP_MGS_ORTH:
113: /* compute |v| */
114: if (onorm) { IPNorm(ip,v,onorm); }
115: IPOrthogonalizeMGS(ip,n,which,V,v,H);
116: /* compute |v'| */
117: if (norm) { IPNorm(ip,v,norm); }
118: break;
119: default:
120: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
121: }
122: return(0);
123: }
127: /*@
128: IPOrthogonalize - Orthogonalize a vector with respect to a set of vectors.
130: Collective on IP
132: Input Parameters:
133: + ip - the inner product (IP) context
134: . n - number of columns of V
135: . which - logical array indicating columns of V to be used
136: . V - set of vectors
137: - work - workspace vector
139: Input/Output Parameter:
140: . v - (input) vector to be orthogonalized and (output) result of
141: orthogonalization
143: Output Parameter:
144: + H - coefficients computed during orthogonalization
145: . norm - norm of the vector after being orthogonalized
146: - lindep - flag indicating that refinement did not improve the quality
147: of orthogonalization
149: Notes:
150: This function applies an orthogonal projector to project vector v onto the
151: orthogonal complement of the span of the columns of V.
153: On exit, v0 = [V v]*H, where v0 is the original vector v.
155: This routine does not normalize the resulting vector.
157: Level: developer
159: .seealso: IPSetOrthogonalization(), IPBiOrthogonalize()
160: @*/
161: PetscErrorCode IPOrthogonalize(IP ip,int n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscTruth *lindep,Vec work)
162: {
164: PetscScalar lh[100],*h,lc[100],*c;
165: PetscTruth allocatedh = PETSC_FALSE,allocatedc = PETSC_FALSE,allocatedw = PETSC_FALSE;
166: PetscReal onrm,nrm;
167: int j,k;
169: if (n==0) {
170: if (norm) { IPNorm(ip,v,norm); }
171: if (lindep) *lindep = PETSC_FALSE;
172: return(0);
173: }
174: PetscLogEventBegin(IP_Orthogonalize,ip,0,0,0);
176: /* allocate H, c and work if needed */
177: if (!H) {
178: if (n<=100) h = lh;
179: else {
180: PetscMalloc(n*sizeof(PetscScalar),&h);
181: allocatedh = PETSC_TRUE;
182: }
183: } else h = H;
184: if (ip->orthog_ref != IP_ORTH_REFINE_NEVER) {
185: if (n<=100) c = lc;
186: else {
187: PetscMalloc(n*sizeof(PetscScalar),&c);
188: allocatedc = PETSC_TRUE;
189: }
190: }
191: if (ip->orthog_type == IP_CGS_ORTH) {
192: VecDuplicate(v,&work);
193: allocatedw = PETSC_TRUE;
194: }
196: /* orthogonalize and compute onorm */
197: switch (ip->orthog_ref) {
198:
199: case IP_ORTH_REFINE_NEVER:
200: IPOrthogonalizeGS(ip,n,which,V,v,h,PETSC_NULL,PETSC_NULL,work);
201: /* compute |v| */
202: if (norm) { IPNorm(ip,v,norm); }
203: /* linear dependence check does not work without refinement */
204: if (lindep) *lindep = PETSC_FALSE;
205: break;
206:
207: case IP_ORTH_REFINE_ALWAYS:
208: IPOrthogonalizeGS(ip,n,which,V,v,h,PETSC_NULL,PETSC_NULL,work);
209: if (lindep) {
210: IPOrthogonalizeGS(ip,n,which,V,v,c,&onrm,&nrm,work);
211: if (norm) *norm = nrm;
212: if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
213: else *lindep = PETSC_FALSE;
214: } else {
215: IPOrthogonalizeGS(ip,n,which,V,v,c,PETSC_NULL,norm,work);
216: }
217: for (j=0;j<n;j++)
218: if (!which || which[j]) h[j] += c[j];
219: break;
220:
221: case IP_ORTH_REFINE_IFNEEDED:
222: IPOrthogonalizeGS(ip,n,which,V,v,h,&onrm,&nrm,work);
223: /* ||q|| < eta ||h|| */
224: k = 1;
225: while (k<3 && nrm < ip->orthog_eta * onrm) {
226: k++;
227: switch (ip->orthog_type) {
228: case IP_CGS_ORTH:
229: IPOrthogonalizeGS(ip,n,which,V,v,c,&onrm,&nrm,work);
230: break;
231: case IP_MGS_ORTH:
232: onrm = nrm;
233: IPOrthogonalizeGS(ip,n,which,V,v,c,PETSC_NULL,&nrm,work);
234: break;
235: default:
236: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
237: }
238: for (j=0;j<n;j++)
239: if (!which || which[j]) h[j] += c[j];
240: }
241: if (norm) *norm = nrm;
242: if (lindep) {
243: if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
244: else *lindep = PETSC_FALSE;
245: }
246:
247: break;
248: default:
249: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization refinement");
250: }
252: /* free work space */
253: if (allocatedc) { PetscFree(c); }
254: if (allocatedh) { PetscFree(h); }
255: if (allocatedw) { VecDestroy(work); }
256:
257: PetscLogEventEnd(IP_Orthogonalize,ip,0,0,0);
258: return(0);
259: }
263: /*@
264: IPQRDecomposition - Compute the QR factorization of a set of vectors.
266: Collective on IP
268: Input Parameters:
269: + ip - the eigenproblem solver context
270: . V - set of vectors
271: . m - starting column
272: . n - ending column
273: . ldr - leading dimension of R
274: - work - workspace vector
276: Output Parameter:
277: . R - triangular matrix of QR factorization
279: Notes:
280: This routine orthonormalizes the columns of V so that V'*V=I up to
281: working precision. It assumes that the first m columns (from 0 to m-1)
282: are already orthonormal. The coefficients of the orthogonalization are
283: stored in R so that V*R is equal to the original V.
285: In some cases, this routine makes V B-orthonormal, that is, V'*B*V=I.
287: If one of the vectors is linearly dependent wrt the rest (at working
288: precision) then it is replaced by a random vector.
290: Level: developer
292: .seealso: IPOrthogonalize(), IPNorm(), IPInnerProduct().
293: @*/
294: PetscErrorCode IPQRDecomposition(IP ip,Vec *V,int m,int n,PetscScalar *R,int ldr,Vec work)
295: {
297: int k;
298: PetscReal norm;
299: PetscTruth lindep;
300:
303: for (k=m; k<n; k++) {
305: /* orthogonalize v_k with respect to v_0, ..., v_{k-1} */
306: if (R) { IPOrthogonalize(ip,k,PETSC_NULL,V,V[k],&R[0+ldr*k],&norm,&lindep,work); }
307: else { IPOrthogonalize(ip,k,PETSC_NULL,V,V[k],PETSC_NULL,&norm,&lindep,work); }
309: /* normalize v_k: r_{k,k} = ||v_k||_2; v_k = v_k/r_{k,k} */
310: if (norm==0.0 || lindep) {
311: PetscInfo(ip,"Linearly dependent vector found, generating a new random vector\n");
312: SlepcVecSetRandom(V[k]);
313: IPNorm(ip,V[k],&norm);
314: }
315: VecScale(V[k],1.0/norm);
316: if (R) R[k+ldr*k] = norm;
318: }
320: return(0);
321: }
323: /*
324: Biorthogonalization routine using classical Gram-Schmidt with refinement.
325: */
328: static PetscErrorCode IPCGSBiOrthogonalization(IP ip,int n,Vec *V,Vec *W,Vec v,PetscScalar *H,PetscReal *hnorm,PetscReal *norm)
329: {
330: #if defined(SLEPC_MISSING_LAPACK_GELQF) || defined(SLEPC_MISSING_LAPACK_ORMLQ)
332: SETERRQ(PETSC_ERR_SUP,"xGELQF - Lapack routine is unavailable.");
333: #else
335: int j,ione=1,lwork,info;
336: PetscScalar shh[100],*lhh,*vw,*tau,one=1.0,*work;
337: Vec w;
341: /* Don't allocate small arrays */
342: if (n<=100) lhh = shh;
343: else { PetscMalloc(n*sizeof(PetscScalar),&lhh); }
344: PetscMalloc(n*n*sizeof(PetscScalar),&vw);
345: VecDuplicate(v,&w);
346:
347: for (j=0;j<n;j++) {
348: IPMInnerProduct(ip,V[j],n,W,vw+j*n);
349: }
350: lwork = n;
351: PetscMalloc(n*sizeof(PetscScalar),&tau);
352: PetscMalloc(lwork*sizeof(PetscScalar),&work);
353: LAPACKgelqf_(&n,&n,vw,&n,tau,work,&lwork,&info);
354: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGELQF %i",info);
355:
356: /*** First orthogonalization ***/
358: /* h = W^* v */
359: /* q = v - V h */
360: IPMInnerProduct(ip,v,n,W,H);
361: BLAStrsm_("L","L","N","N",&n,&ione,&one,vw,&n,H,&n,1,1,1,1);
362: LAPACKormlq_("L","N",&n,&ione,&n,vw,&n,tau,H,&n,work,&lwork,&info,1,1);
363: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xORMLQ %i",info);
364: VecSet(w,0.0);
365: VecMAXPY(w,n,H,V);
366: VecAXPY(v,-1.0,w);
368: /* compute norm of v */
369: if (norm) { IPNorm(ip,v,norm); }
370:
371: if (n>100) { PetscFree(lhh); }
372: PetscFree(vw);
373: PetscFree(tau);
374: PetscFree(work);
375: VecDestroy(w);
376: return(0);
377: #endif
378: }
382: /*@
383: IPBiOrthogonalize - Bi-orthogonalize a vector with respect to a set of vectors.
385: Collective on IP
387: Input Parameters:
388: + ip - the inner product context
389: . n - number of columns of V
390: . V - set of vectors
391: - W - set of vectors
393: Input/Output Parameter:
394: . v - vector to be orthogonalized
396: Output Parameter:
397: + H - coefficients computed during orthogonalization
398: - norm - norm of the vector after being orthogonalized
400: Notes:
401: This function applies an oblique projector to project vector v onto the
402: span of the columns of V along the orthogonal complement of the column
403: space of W.
405: On exit, v0 = [V v]*H, where v0 is the original vector v.
407: This routine does not normalize the resulting vector.
409: Level: developer
411: .seealso: IPSetOrthogonalization(), IPOrthogonalize()
412: @*/
413: PetscErrorCode IPBiOrthogonalize(IP ip,int n,Vec *V,Vec *W,Vec v,PetscScalar *H,PetscReal *norm)
414: {
416: PetscScalar lh[100],*h;
417: PetscTruth allocated = PETSC_FALSE;
418: PetscReal lhnrm,*hnrm,lnrm,*nrm;
420: if (n==0) {
421: if (norm) { IPNorm(ip,v,norm); }
422: } else {
423: PetscLogEventBegin(IP_Orthogonalize,ip,0,0,0);
424:
425: /* allocate H if needed */
426: if (!H) {
427: if (n<=100) h = lh;
428: else {
429: PetscMalloc(n*sizeof(PetscScalar),&h);
430: allocated = PETSC_TRUE;
431: }
432: } else h = H;
433:
434: /* retrieve hnrm and nrm for linear dependence check or conditional refinement */
435: if (ip->orthog_ref == IP_ORTH_REFINE_IFNEEDED) {
436: hnrm = &lhnrm;
437: if (norm) nrm = norm;
438: else nrm = &lnrm;
439: } else {
440: hnrm = PETSC_NULL;
441: nrm = norm;
442: }
443:
444: switch (ip->orthog_type) {
445: case IP_CGS_ORTH:
446: IPCGSBiOrthogonalization(ip,n,V,W,v,h,hnrm,nrm);
447: break;
448: default:
449: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
450: }
451:
452: if (allocated) { PetscFree(h); }
453:
454: PetscLogEventEnd(IP_Orthogonalize,ip,0,0,0);
455: }
456: return(0);
457: }