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-2010, Universidad Politecnica de Valencia, Spain
9: This file is part of SLEPc.
10:
11: SLEPc is free software: you can redistribute it and/or modify it under the
12: terms of version 3 of the GNU Lesser General Public License as published by
13: the Free Software Foundation.
15: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
16: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
17: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
18: more details.
20: You should have received a copy of the GNU Lesser General Public License
21: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
22: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23: */
25: #include private/ipimpl.h
26: #include slepcblaslapack.h
28: /*
29: IPOrthogonalizeMGS1 - Compute one step of Modified Gram-Schmidt
30: */
33: static PetscErrorCode IPOrthogonalizeMGS1(IP ip,PetscInt n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H)
34: {
36: PetscInt j;
37: PetscScalar dot;
38:
40: for (j=0; j<n; j++)
41: if (!which || which[j]) {
42: /* h_j = ( v, v_j ) */
43: IPInnerProduct(ip,v,V[j],&dot);
44: /* v <- v - h_j v_j */
45: VecAXPY(v,-dot,V[j]);
46: if (H) H[j] += dot;
47: }
48: return(0);
49: }
51: /*
52: IPOrthogonalizeCGS1 - Compute |v'| (estimated), |v| and one step of CGS with only one global synchronization
53: */
56: PetscErrorCode IPOrthogonalizeCGS1(IP ip,PetscInt nds,Vec *DS,PetscInt n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H,PetscReal *onorm,PetscReal *norm)
57: {
59: PetscInt j;
60: PetscScalar alpha;
61: PetscReal sum;
64: /* h = W^* v ; alpha = (v , v) */
65: if (nds==0 && !which && !onorm && !norm) {
66: /* use simpler function */
67: IPMInnerProduct(ip,v,n,V,H);
68: } else {
69: /* merge comunications */
70: IPMInnerProductBegin(ip,v,nds,DS,H);
71: if (which) { /* use select array */
72: for (j=0; j<n; j++)
73: if (which[j]) { IPInnerProductBegin(ip,v,V[j],H+nds+j); }
74: } else {
75: IPMInnerProductBegin(ip,v,n,V,H+nds);
76: }
77: if (onorm || (norm && !ip->matrix)) {
78: IPInnerProductBegin(ip,v,v,&alpha);
79: }
81: IPMInnerProductEnd(ip,v,nds,DS,H);
82: if (which) { /* use select array */
83: for (j=0; j<n; j++)
84: if (which[j]) { IPInnerProductEnd(ip,v,V[j],H+nds+j); }
85: } else {
86: IPMInnerProductEnd(ip,v,n,V,H+nds);
87: }
88: if (onorm || (norm && !ip->matrix)) {
89: IPInnerProductEnd(ip,v,v,&alpha);
90: }
91: }
93: /* q = v - V h */
94: SlepcVecMAXPBY(v,1.0,-1.0,nds,H,DS);
95: if (which) {
96: for (j=0; j<n; j++)
97: if (which[j]) { VecAXPBY(v,-H[nds+j],1.0,V[j]); }
98: } else {
99: SlepcVecMAXPBY(v,1.0,-1.0,n,H+nds,V);
100: }
101:
102: /* compute |v| */
103: if (onorm) *onorm = sqrt(PetscRealPart(alpha));
105: if (norm) {
106: if (!ip->matrix) {
107: /* estimate |v'| from |v| */
108: sum = 0.0;
109: for (j=0; j<nds; j++)
110: sum += PetscRealPart(H[j] * PetscConj(H[j]));
111: for (j=0; j<n; j++)
112: if (!which || which[j])
113: sum += PetscRealPart(H[nds+j] * PetscConj(H[nds+j]));
114: *norm = PetscRealPart(alpha)-sum;
115: if (*norm <= 0.0) {
116: IPNorm(ip,v,norm);
117: } else *norm = sqrt(*norm);
118: } else {
119: /* compute |v'| */
120: IPNorm(ip,v,norm);
121: }
122: }
123: return(0);
124: }
126: /*
127: IPOrthogonalizeMGS - Orthogonalize with modified Gram-Schmidt
128: */
131: static PetscErrorCode IPOrthogonalizeMGS(IP ip,PetscInt nds,Vec *DS,PetscInt n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscTruth *lindep)
132: {
134: PetscInt i,k;
135: PetscReal onrm,nrm;
138:
139: if (H) {
140: for (i=0;i<n;i++)
141: H[i] = 0;
142: }
143:
144: switch (ip->orthog_ref) {
145:
146: case IP_ORTH_REFINE_NEVER:
147: IPOrthogonalizeMGS1(ip,nds,PETSC_NULL,DS,v,PETSC_NULL);
148: IPOrthogonalizeMGS1(ip,n,which,V,v,H);
149: /* compute |v| */
150: if (norm) { IPNorm(ip,v,norm); }
151: /* linear dependence check does not work without refinement */
152: if (lindep) *lindep = PETSC_FALSE;
153: break;
154:
155: case IP_ORTH_REFINE_ALWAYS:
156: /* first step */
157: IPOrthogonalizeMGS1(ip,nds,PETSC_NULL,DS,v,PETSC_NULL);
158: IPOrthogonalizeMGS1(ip,n,which,V,v,H);
159: if (lindep) { IPNorm(ip,v,&onrm); }
160: /* second step */
161: IPOrthogonalizeMGS1(ip,nds,PETSC_NULL,DS,v,PETSC_NULL);
162: IPOrthogonalizeMGS1(ip,n,which,V,v,H);
163: if (lindep || norm) { IPNorm(ip,v,&nrm); }
164: if (lindep) {
165: if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
166: else *lindep = PETSC_FALSE;
167: }
168: if (norm) *norm = nrm;
169: break;
170:
171: case IP_ORTH_REFINE_IFNEEDED:
172: /* first step */
173: IPNorm(ip,v,&onrm);
174: IPOrthogonalizeMGS1(ip,nds,PETSC_NULL,DS,v,PETSC_NULL);
175: IPOrthogonalizeMGS1(ip,n,which,V,v,H);
176: IPNorm(ip,v,&nrm);
177: /* ||q|| < eta ||h|| */
178: k = 1;
179: while (k<3 && nrm < ip->orthog_eta * onrm) {
180: k++;
181: onrm = nrm;
182: IPOrthogonalizeMGS1(ip,nds,PETSC_NULL,DS,v,PETSC_NULL);
183: IPOrthogonalizeMGS1(ip,n,which,V,v,H);
184: IPNorm(ip,v,&nrm);
185: }
186: if (lindep) {
187: if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
188: else *lindep = PETSC_FALSE;
189: }
190: if (norm) *norm = nrm;
191: break;
192:
193: default:
194: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization refinement");
195: }
197: return(0);
198: }
200: /*
201: IPOrthogonalizeCGS - Orthogonalize with classical Gram-Schmidt
202: */
205: static PetscErrorCode IPOrthogonalizeCGS(IP ip,PetscInt nds,Vec *DS,PetscInt n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscTruth *lindep)
206: {
208: PetscScalar lh[100],*h,lc[100],*c;
209: PetscTruth allocatedh = PETSC_FALSE,allocatedc = PETSC_FALSE;
210: PetscReal onrm,nrm;
211: PetscInt j,k;
214: /* allocate h and c if needed */
215: if (!H || nds>0) {
216: if (nds+n<=100) h = lh;
217: else {
218: PetscMalloc((nds+n)*sizeof(PetscScalar),&h);
219: allocatedh = PETSC_TRUE;
220: }
221: } else h = H;
222: if (ip->orthog_ref != IP_ORTH_REFINE_NEVER) {
223: if (nds+n<=100) c = lc;
224: else {
225: PetscMalloc((nds+n)*sizeof(PetscScalar),&c);
226: allocatedc = PETSC_TRUE;
227: }
228: }
230: /* orthogonalize and compute onorm */
231: switch (ip->orthog_ref) {
232:
233: case IP_ORTH_REFINE_NEVER:
234: IPOrthogonalizeCGS1(ip,nds,DS,n,which,V,v,h,PETSC_NULL,PETSC_NULL);
235: /* compute |v| */
236: if (norm) { IPNorm(ip,v,norm); }
237: /* linear dependence check does not work without refinement */
238: if (lindep) *lindep = PETSC_FALSE;
239: break;
240:
241: case IP_ORTH_REFINE_ALWAYS:
242: IPOrthogonalizeCGS1(ip,nds,DS,n,which,V,v,h,PETSC_NULL,PETSC_NULL);
243: if (lindep) {
244: IPOrthogonalizeCGS1(ip,nds,DS,n,which,V,v,c,&onrm,&nrm);
245: if (norm) *norm = nrm;
246: if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
247: else *lindep = PETSC_FALSE;
248: } else {
249: IPOrthogonalizeCGS1(ip,nds,DS,n,which,V,v,c,PETSC_NULL,norm);
250: }
251: for (j=0;j<n;j++)
252: if (!which || which[j]) h[nds+j] += c[nds+j];
253: break;
254:
255: case IP_ORTH_REFINE_IFNEEDED:
256: IPOrthogonalizeCGS1(ip,nds,DS,n,which,V,v,h,&onrm,&nrm);
257: /* ||q|| < eta ||h|| */
258: k = 1;
259: while (k<3 && nrm < ip->orthog_eta * onrm) {
260: k++;
261: if (!ip->matrix) {
262: IPOrthogonalizeCGS1(ip,nds,DS,n,which,V,v,c,&onrm,&nrm);
263: } else {
264: onrm = nrm;
265: IPOrthogonalizeCGS1(ip,nds,DS,n,which,V,v,c,PETSC_NULL,&nrm);
266: }
267: for (j=0;j<n;j++)
268: if (!which || which[j]) h[nds+j] += c[nds+j];
269: }
270: if (norm) *norm = nrm;
271: if (lindep) {
272: if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
273: else *lindep = PETSC_FALSE;
274: }
275: break;
277: default:
278: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization refinement");
279: }
281: /* recover H from workspace */
282: if (H && nds>0) {
283: for (j=0;j<n;j++)
284: if (!which || which[j]) H[j] = h[nds+j];
285: }
287: /* free work space */
288: if (allocatedc) { PetscFree(c); }
289: if (allocatedh) { PetscFree(h); }
290: return(0);
291: }
295: /*@
296: IPOrthogonalize - Orthogonalize a vector with respect to two set of vectors.
298: Collective on IP
300: Input Parameters:
301: + ip - the inner product (IP) context
302: . nds - number of columns of DS
303: . DS - first set of vectors
304: . n - number of columns of V
305: . which - logical array indicating columns of V to be used
306: - V - second set of vectors
308: Input/Output Parameter:
309: . v - (input) vector to be orthogonalized and (output) result of
310: orthogonalization
312: Output Parameter:
313: + H - coefficients computed during orthogonalization with V
314: . norm - norm of the vector after being orthogonalized
315: - lindep - flag indicating that refinement did not improve the quality
316: of orthogonalization
318: Notes:
319: This function applies an orthogonal projector to project vector v onto the
320: orthogonal complement of the span of the columns of DS and V.
322: On exit, v0 = [V v]*H, where v0 is the original vector v.
324: This routine does not normalize the resulting vector.
326: Level: developer
328: .seealso: IPSetOrthogonalization(), IPBiOrthogonalize()
329: @*/
330: PetscErrorCode IPOrthogonalize(IP ip,PetscInt nds,Vec *DS,PetscInt n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscTruth *lindep)
331: {
335: PetscLogEventBegin(IP_Orthogonalize,ip,0,0,0);
336:
337: if (nds==0 && n==0) {
338: if (norm) { IPNorm(ip,v,norm); }
339: if (lindep) *lindep = PETSC_FALSE;
340: } else {
341: switch (ip->orthog_type) {
342: case IP_ORTH_CGS:
343: IPOrthogonalizeCGS(ip,nds,DS,n,which,V,v,H,norm,lindep);
344: break;
345: case IP_ORTH_MGS:
346: IPOrthogonalizeMGS(ip,nds,DS,n,which,V,v,H,norm,lindep);
347: break;
348: default:
349: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
350: }
351: }
352:
353: PetscLogEventEnd(IP_Orthogonalize,ip,0,0,0);
354: return(0);
355: }
359: /*@
360: IPQRDecomposition - Compute the QR factorization of a set of vectors.
362: Collective on IP
364: Input Parameters:
365: + ip - the eigenproblem solver context
366: . V - set of vectors
367: . m - starting column
368: . n - ending column
369: - ldr - leading dimension of R
371: Output Parameter:
372: . R - triangular matrix of QR factorization
374: Notes:
375: This routine orthonormalizes the columns of V so that V'*V=I up to
376: working precision. It assumes that the first m columns (from 0 to m-1)
377: are already orthonormal. The coefficients of the orthogonalization are
378: stored in R so that V*R is equal to the original V.
380: In some cases, this routine makes V B-orthonormal, that is, V'*B*V=I.
382: If one of the vectors is linearly dependent wrt the rest (at working
383: precision) then it is replaced by a random vector.
385: Level: developer
387: .seealso: IPOrthogonalize(), IPNorm(), IPInnerProduct().
388: @*/
389: PetscErrorCode IPQRDecomposition(IP ip,Vec *V,PetscInt m,PetscInt n,PetscScalar *R,PetscInt ldr)
390: {
392: PetscInt k;
393: PetscReal norm;
394: PetscTruth lindep;
395: PetscRandom rctx=PETSC_NULL;
396:
399: for (k=m; k<n; k++) {
401: /* orthogonalize v_k with respect to v_0, ..., v_{k-1} */
402: if (R) { IPOrthogonalize(ip,0,PETSC_NULL,k,PETSC_NULL,V,V[k],&R[0+ldr*k],&norm,&lindep); }
403: else { IPOrthogonalize(ip,0,PETSC_NULL,k,PETSC_NULL,V,V[k],PETSC_NULL,&norm,&lindep); }
405: /* normalize v_k: r_{k,k} = ||v_k||_2; v_k = v_k/r_{k,k} */
406: if (norm==0.0 || lindep) {
407: PetscInfo(ip,"Linearly dependent vector found, generating a new random vector\n");
408: if (!rctx) {
409: PetscRandomCreate(((PetscObject)ip)->comm,&rctx);
410: PetscRandomSetFromOptions(rctx);
411: }
412: SlepcVecSetRandom(V[k],rctx);
413: IPNorm(ip,V[k],&norm);
414: }
415: VecScale(V[k],1.0/norm);
416: if (R) R[k+ldr*k] = norm;
418: }
419: if (rctx) { PetscRandomDestroy(rctx); }
421: return(0);
422: }
424: /*
425: Biorthogonalization routine using classical Gram-Schmidt with refinement.
426: */
429: static PetscErrorCode IPCGSBiOrthogonalization(IP ip,PetscInt n_,Vec *V,Vec *W,Vec v,PetscScalar *H,PetscReal *hnorm,PetscReal *norm)
430: {
431: #if defined(SLEPC_MISSING_LAPACK_GELQF) || defined(SLEPC_MISSING_LAPACK_ORMLQ)
433: SETERRQ(PETSC_ERR_SUP,"xGELQF - Lapack routine is unavailable.");
434: #else
436: PetscBLASInt j,ione=1,lwork,info,n=n_;
437: PetscScalar shh[100],*lhh,*vw,*tau,one=1.0,*work;
441: /* Don't allocate small arrays */
442: if (n<=100) lhh = shh;
443: else { PetscMalloc(n*sizeof(PetscScalar),&lhh); }
444: PetscMalloc(n*n*sizeof(PetscScalar),&vw);
445:
446: for (j=0;j<n;j++) {
447: IPMInnerProduct(ip,V[j],n,W,vw+j*n);
448: }
449: lwork = n;
450: PetscMalloc(n*sizeof(PetscScalar),&tau);
451: PetscMalloc(lwork*sizeof(PetscScalar),&work);
452: LAPACKgelqf_(&n,&n,vw,&n,tau,work,&lwork,&info);
453: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGELQF %i",info);
454:
455: /*** First orthogonalization ***/
457: /* h = W^* v */
458: /* q = v - V h */
459: IPMInnerProduct(ip,v,n,W,H);
460: BLAStrsm_("L","L","N","N",&n,&ione,&one,vw,&n,H,&n);
461: LAPACKormlq_("L","N",&n,&ione,&n,vw,&n,tau,H,&n,work,&lwork,&info);
462: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xORMLQ %i",info);
463: SlepcVecMAXPBY(v,1.0,-1.0,n,H,V);
465: /* compute norm of v */
466: if (norm) { IPNorm(ip,v,norm); }
467:
468: if (n>100) { PetscFree(lhh); }
469: PetscFree(vw);
470: PetscFree(tau);
471: PetscFree(work);
472: return(0);
473: #endif
474: }
478: /*@
479: IPBiOrthogonalize - Bi-orthogonalize a vector with respect to a set of vectors.
481: Collective on IP
483: Input Parameters:
484: + ip - the inner product context
485: . n - number of columns of V
486: . V - set of vectors
487: - W - set of vectors
489: Input/Output Parameter:
490: . v - vector to be orthogonalized
492: Output Parameter:
493: + H - coefficients computed during orthogonalization
494: - norm - norm of the vector after being orthogonalized
496: Notes:
497: This function applies an oblique projector to project vector v onto the
498: span of the columns of V along the orthogonal complement of the column
499: space of W.
501: On exit, v0 = [V v]*H, where v0 is the original vector v.
503: This routine does not normalize the resulting vector.
505: Level: developer
507: .seealso: IPSetOrthogonalization(), IPOrthogonalize()
508: @*/
509: PetscErrorCode IPBiOrthogonalize(IP ip,PetscInt n,Vec *V,Vec *W,Vec v,PetscScalar *H,PetscReal *norm)
510: {
512: PetscScalar lh[100],*h;
513: PetscTruth allocated = PETSC_FALSE;
514: PetscReal lhnrm,*hnrm,lnrm,*nrm;
516: if (n==0) {
517: if (norm) { IPNorm(ip,v,norm); }
518: } else {
519: PetscLogEventBegin(IP_Orthogonalize,ip,0,0,0);
520:
521: /* allocate H if needed */
522: if (!H) {
523: if (n<=100) h = lh;
524: else {
525: PetscMalloc(n*sizeof(PetscScalar),&h);
526: allocated = PETSC_TRUE;
527: }
528: } else h = H;
529:
530: /* retrieve hnrm and nrm for linear dependence check or conditional refinement */
531: if (ip->orthog_ref == IP_ORTH_REFINE_IFNEEDED) {
532: hnrm = &lhnrm;
533: if (norm) nrm = norm;
534: else nrm = &lnrm;
535: } else {
536: hnrm = PETSC_NULL;
537: nrm = norm;
538: }
539:
540: switch (ip->orthog_type) {
541: case IP_ORTH_CGS:
542: IPCGSBiOrthogonalization(ip,n,V,W,v,h,hnrm,nrm);
543: break;
544: default:
545: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
546: }
547:
548: if (allocated) { PetscFree(h); }
549:
550: PetscLogEventEnd(IP_Orthogonalize,ip,0,0,0);
551: }
552: return(0);
553: }