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-2011, Universitat 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> /*I "slepcip.h" I*/
26: #include <slepcblaslapack.h>
28: /*
29: IPOrthogonalizeMGS1 - Compute one step of Modified Gram-Schmidt
30: */
33: static PetscErrorCode IPOrthogonalizeMGS1(IP ip,PetscInt n,PetscBool *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,PetscBool *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 = PetscSqrtReal(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 = PetscSqrtReal(*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,PetscBool *which,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
132: {
134: PetscInt i,k;
135: PetscReal onrm,nrm;
138: if (H) {
139: for (i=0;i<n;i++)
140: H[i] = 0;
141: }
142:
143: switch (ip->orthog_ref) {
144:
145: case IP_ORTHOG_REFINE_NEVER:
146: IPOrthogonalizeMGS1(ip,nds,PETSC_NULL,DS,v,PETSC_NULL);
147: IPOrthogonalizeMGS1(ip,n,which,V,v,H);
148: /* compute |v| */
149: if (norm) { IPNorm(ip,v,norm); }
150: /* linear dependence check does not work without refinement */
151: if (lindep) *lindep = PETSC_FALSE;
152: break;
153:
154: case IP_ORTHOG_REFINE_ALWAYS:
155: /* first step */
156: IPOrthogonalizeMGS1(ip,nds,PETSC_NULL,DS,v,PETSC_NULL);
157: IPOrthogonalizeMGS1(ip,n,which,V,v,H);
158: if (lindep) { IPNorm(ip,v,&onrm); }
159: /* second step */
160: IPOrthogonalizeMGS1(ip,nds,PETSC_NULL,DS,v,PETSC_NULL);
161: IPOrthogonalizeMGS1(ip,n,which,V,v,H);
162: if (lindep || norm) { IPNorm(ip,v,&nrm); }
163: if (lindep) {
164: if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
165: else *lindep = PETSC_FALSE;
166: }
167: if (norm) *norm = nrm;
168: break;
169:
170: case IP_ORTHOG_REFINE_IFNEEDED:
171: /* first step */
172: IPNorm(ip,v,&onrm);
173: IPOrthogonalizeMGS1(ip,nds,PETSC_NULL,DS,v,PETSC_NULL);
174: IPOrthogonalizeMGS1(ip,n,which,V,v,H);
175: IPNorm(ip,v,&nrm);
176: /* ||q|| < eta ||h|| */
177: k = 1;
178: while (k<3 && nrm < ip->orthog_eta * onrm) {
179: k++;
180: onrm = nrm;
181: IPOrthogonalizeMGS1(ip,nds,PETSC_NULL,DS,v,PETSC_NULL);
182: IPOrthogonalizeMGS1(ip,n,which,V,v,H);
183: IPNorm(ip,v,&nrm);
184: }
185: if (lindep) {
186: if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
187: else *lindep = PETSC_FALSE;
188: }
189: if (norm) *norm = nrm;
190: break;
191:
192: default:
193: SETERRQ(((PetscObject)ip)->comm,PETSC_ERR_ARG_WRONG,"Unknown orthogonalization refinement");
194: }
195: return(0);
196: }
198: /*
199: IPOrthogonalizeCGS - Orthogonalize with classical Gram-Schmidt
200: */
203: static PetscErrorCode IPOrthogonalizeCGS(IP ip,PetscInt nds,Vec *DS,PetscInt n,PetscBool *which,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
204: {
206: PetscScalar lh[100],*h,lc[100],*c;
207: PetscBool allocatedh = PETSC_FALSE,allocatedc = PETSC_FALSE;
208: PetscReal onrm,nrm;
209: PetscInt j,k;
212: /* allocate h and c if needed */
213: if (!H || nds>0) {
214: if (nds+n<=100) h = lh;
215: else {
216: PetscMalloc((nds+n)*sizeof(PetscScalar),&h);
217: allocatedh = PETSC_TRUE;
218: }
219: } else h = H;
220: if (ip->orthog_ref != IP_ORTHOG_REFINE_NEVER) {
221: if (nds+n<=100) c = lc;
222: else {
223: PetscMalloc((nds+n)*sizeof(PetscScalar),&c);
224: allocatedc = PETSC_TRUE;
225: }
226: }
228: /* orthogonalize and compute onorm */
229: switch (ip->orthog_ref) {
230:
231: case IP_ORTHOG_REFINE_NEVER:
232: IPOrthogonalizeCGS1(ip,nds,DS,n,which,V,v,h,PETSC_NULL,PETSC_NULL);
233: /* compute |v| */
234: if (norm) { IPNorm(ip,v,norm); }
235: /* linear dependence check does not work without refinement */
236: if (lindep) *lindep = PETSC_FALSE;
237: break;
238:
239: case IP_ORTHOG_REFINE_ALWAYS:
240: IPOrthogonalizeCGS1(ip,nds,DS,n,which,V,v,h,PETSC_NULL,PETSC_NULL);
241: if (lindep) {
242: IPOrthogonalizeCGS1(ip,nds,DS,n,which,V,v,c,&onrm,&nrm);
243: if (norm) *norm = nrm;
244: if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
245: else *lindep = PETSC_FALSE;
246: } else {
247: IPOrthogonalizeCGS1(ip,nds,DS,n,which,V,v,c,PETSC_NULL,norm);
248: }
249: for (j=0;j<n;j++)
250: if (!which || which[j]) h[nds+j] += c[nds+j];
251: break;
252:
253: case IP_ORTHOG_REFINE_IFNEEDED:
254: IPOrthogonalizeCGS1(ip,nds,DS,n,which,V,v,h,&onrm,&nrm);
255: /* ||q|| < eta ||h|| */
256: k = 1;
257: while (k<3 && nrm < ip->orthog_eta * onrm) {
258: k++;
259: if (!ip->matrix) {
260: IPOrthogonalizeCGS1(ip,nds,DS,n,which,V,v,c,&onrm,&nrm);
261: } else {
262: onrm = nrm;
263: IPOrthogonalizeCGS1(ip,nds,DS,n,which,V,v,c,PETSC_NULL,&nrm);
264: }
265: for (j=0;j<n;j++)
266: if (!which || which[j]) h[nds+j] += c[nds+j];
267: }
268: if (norm) *norm = nrm;
269: if (lindep) {
270: if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
271: else *lindep = PETSC_FALSE;
272: }
273: break;
275: default:
276: SETERRQ(((PetscObject)ip)->comm,PETSC_ERR_ARG_WRONG,"Unknown orthogonalization refinement");
277: }
279: /* recover H from workspace */
280: if (H && nds>0) {
281: for (j=0;j<n;j++)
282: if (!which || which[j]) H[j] = h[nds+j];
283: }
285: /* free work space */
286: if (allocatedc) { PetscFree(c); }
287: if (allocatedh) { PetscFree(h); }
288: return(0);
289: }
293: /*@
294: IPOrthogonalize - Orthogonalize a vector with respect to two set of vectors.
296: Collective on IP and Vec
298: Input Parameters:
299: + ip - the inner product (IP) context
300: . nds - number of columns of DS
301: . DS - first set of vectors
302: . n - number of columns of V
303: . which - logical array indicating columns of V to be used
304: - V - second set of vectors
306: Input/Output Parameter:
307: . v - (input) vector to be orthogonalized and (output) result of
308: orthogonalization
310: Output Parameter:
311: + H - coefficients computed during orthogonalization with V
312: . norm - norm of the vector after being orthogonalized
313: - lindep - flag indicating that refinement did not improve the quality
314: of orthogonalization
316: Notes:
317: This function applies an orthogonal projector to project vector v onto the
318: orthogonal complement of the span of the columns of DS and V.
320: On exit, v0 = [V v]*H, where v0 is the original vector v.
322: This routine does not normalize the resulting vector.
324: Level: developer
326: .seealso: IPSetOrthogonalization(), IPBiOrthogonalize()
327: @*/
328: PetscErrorCode IPOrthogonalize(IP ip,PetscInt nds,Vec *DS,PetscInt n,PetscBool *which,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
329: {
336: PetscLogEventBegin(IP_Orthogonalize,ip,0,0,0);
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_ORTHOG_CGS:
343: IPOrthogonalizeCGS(ip,nds,DS,n,which,V,v,H,norm,lindep);
344: break;
345: case IP_ORTHOG_MGS:
346: IPOrthogonalizeMGS(ip,nds,DS,n,which,V,v,H,norm,lindep);
347: break;
348: default:
349: SETERRQ(((PetscObject)ip)->comm,PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
350: }
351: }
352: PetscLogEventEnd(IP_Orthogonalize,ip,0,0,0);
353: return(0);
354: }
358: /*@
359: IPQRDecomposition - Compute the QR factorization of a set of vectors.
361: Collective on IP and Vec
363: Input Parameters:
364: + ip - the eigenproblem solver context
365: . V - set of vectors
366: . m - starting column
367: . n - ending column
368: - ldr - leading dimension of R
370: Output Parameter:
371: . R - triangular matrix of QR factorization
373: Notes:
374: This routine orthonormalizes the columns of V so that V'*V=I up to
375: working precision. It assumes that the first m columns (from 0 to m-1)
376: are already orthonormal. The coefficients of the orthogonalization are
377: stored in R so that V*R is equal to the original V.
379: In some cases, this routine makes V B-orthonormal, that is, V'*B*V=I.
381: If one of the vectors is linearly dependent wrt the rest (at working
382: precision) then it is replaced by a random vector.
384: Level: developer
386: .seealso: IPOrthogonalize(), IPNorm(), IPInnerProduct().
387: @*/
388: PetscErrorCode IPQRDecomposition(IP ip,Vec *V,PetscInt m,PetscInt n,PetscScalar *R,PetscInt ldr)
389: {
391: PetscInt k;
392: PetscReal norm;
393: PetscBool lindep;
394: PetscRandom rctx=PETSC_NULL;
395:
397: for (k=m; k<n; k++) {
399: /* orthogonalize v_k with respect to v_0, ..., v_{k-1} */
400: if (R) { IPOrthogonalize(ip,0,PETSC_NULL,k,PETSC_NULL,V,V[k],&R[0+ldr*k],&norm,&lindep); }
401: else { IPOrthogonalize(ip,0,PETSC_NULL,k,PETSC_NULL,V,V[k],PETSC_NULL,&norm,&lindep); }
403: /* normalize v_k: r_{k,k} = ||v_k||_2; v_k = v_k/r_{k,k} */
404: if (norm==0.0 || lindep) {
405: PetscInfo(ip,"Linearly dependent vector found, generating a new random vector\n");
406: if (!rctx) {
407: PetscRandomCreate(((PetscObject)ip)->comm,&rctx);
408: PetscRandomSetSeed(rctx,0x12345678);
409: PetscRandomSetFromOptions(rctx);
410: }
411: SlepcVecSetRandom(V[k],rctx);
412: IPNorm(ip,V[k],&norm);
413: }
414: VecScale(V[k],1.0/norm);
415: if (R) R[k+ldr*k] = norm;
417: }
418: PetscRandomDestroy(&rctx);
419: return(0);
420: }
422: /*
423: Biorthogonalization routine using classical Gram-Schmidt with refinement.
424: */
427: static PetscErrorCode IPCGSBiOrthogonalization(IP ip,PetscInt n_,Vec *V,Vec *W,Vec v,PetscScalar *H,PetscReal *hnorm,PetscReal *norm)
428: {
429: #if defined(SLEPC_MISSING_LAPACK_GELQF) || defined(SLEPC_MISSING_LAPACK_ORMLQ)
431: SETERRQ(((PetscObject)ip)->comm,PETSC_ERR_SUP,"xGELQF - Lapack routine is unavailable.");
432: #else
434: PetscBLASInt j,ione=1,lwork,info,n=n_;
435: PetscScalar shh[100],*lhh,*vw,*tau,one=1.0,*work;
438: /* Don't allocate small arrays */
439: if (n<=100) lhh = shh;
440: else { PetscMalloc(n*sizeof(PetscScalar),&lhh); }
441: PetscMalloc(n*n*sizeof(PetscScalar),&vw);
442:
443: for (j=0;j<n;j++) {
444: IPMInnerProduct(ip,V[j],n,W,vw+j*n);
445: }
446: lwork = n;
447: PetscMalloc(n*sizeof(PetscScalar),&tau);
448: PetscMalloc(lwork*sizeof(PetscScalar),&work);
449: LAPACKgelqf_(&n,&n,vw,&n,tau,work,&lwork,&info);
450: if (info) SETERRQ1(((PetscObject)ip)->comm,PETSC_ERR_LIB,"Error in Lapack xGELQF %d",info);
451:
452: /*** First orthogonalization ***/
454: /* h = W^* v */
455: /* q = v - V h */
456: IPMInnerProduct(ip,v,n,W,H);
457: BLAStrsm_("L","L","N","N",&n,&ione,&one,vw,&n,H,&n);
458: LAPACKormlq_("L","N",&n,&ione,&n,vw,&n,tau,H,&n,work,&lwork,&info);
459: if (info) SETERRQ1(((PetscObject)ip)->comm,PETSC_ERR_LIB,"Error in Lapack xORMLQ %d",info);
460: SlepcVecMAXPBY(v,1.0,-1.0,n,H,V);
462: /* compute norm of v */
463: if (norm) { IPNorm(ip,v,norm); }
464:
465: if (n>100) { PetscFree(lhh); }
466: PetscFree(vw);
467: PetscFree(tau);
468: PetscFree(work);
469: return(0);
470: #endif
471: }
475: /*@
476: IPBiOrthogonalize - Bi-orthogonalize a vector with respect to a set of vectors.
478: Collective on IP and Vec
480: Input Parameters:
481: + ip - the inner product context
482: . n - number of columns of V
483: . V - set of vectors
484: - W - set of vectors
486: Input/Output Parameter:
487: . v - vector to be orthogonalized
489: Output Parameter:
490: + H - coefficients computed during orthogonalization
491: - norm - norm of the vector after being orthogonalized
493: Notes:
494: This function applies an oblique projector to project vector v onto the
495: span of the columns of V along the orthogonal complement of the column
496: space of W.
498: On exit, v0 = [V v]*H, where v0 is the original vector v.
500: This routine does not normalize the resulting vector.
502: Level: developer
504: .seealso: IPSetOrthogonalization(), IPOrthogonalize()
505: @*/
506: PetscErrorCode IPBiOrthogonalize(IP ip,PetscInt n,Vec *V,Vec *W,Vec v,PetscScalar *H,PetscReal *norm)
507: {
509: PetscScalar lh[100],*h;
510: PetscBool allocated = PETSC_FALSE;
511: 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: /* allocate H if needed */
521: if (!H) {
522: if (n<=100) h = lh;
523: else {
524: PetscMalloc(n*sizeof(PetscScalar),&h);
525: allocated = PETSC_TRUE;
526: }
527: } else h = H;
528:
529: /* retrieve hnrm and nrm for linear dependence check or conditional refinement */
530: if (ip->orthog_ref == IP_ORTHOG_REFINE_IFNEEDED) {
531: hnrm = &lhnrm;
532: if (norm) nrm = norm;
533: else nrm = &lnrm;
534: } else {
535: hnrm = PETSC_NULL;
536: nrm = norm;
537: }
538:
539: switch (ip->orthog_type) {
540: case IP_ORTHOG_CGS:
541: IPCGSBiOrthogonalization(ip,n,V,W,v,h,hnrm,nrm);
542: break;
543: default:
544: SETERRQ(((PetscObject)ip)->comm,PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
545: }
546:
547: if (allocated) { PetscFree(h); }
548: PetscLogEventEnd(IP_Orthogonalize,ip,0,0,0);
549: }
550: return(0);
551: }