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-2013, Universitat Politecnica de Valencia, Spain
9: This file is part of SLEPc.
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 <slepc-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;
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 *defl,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 && !which && !onorm && !norm) {
66: /* use simpler function */
67: IPMInnerProduct(ip,v,n,V,H);
68: } else {
69: /* merge comunications */
70: IPMInnerProductBegin(ip,v,nds,defl,H);
71: if (which) { /* use select array */
72: for (j=0;j<n;j++)
73: if (which[j]) {
74: IPInnerProductBegin(ip,v,V[j],H+nds+j);
75: }
76: } else {
77: IPMInnerProductBegin(ip,v,n,V,H+nds);
78: }
79: if (onorm || (norm && !ip->matrix)) {
80: IPInnerProductBegin(ip,v,v,&alpha);
81: }
83: IPMInnerProductEnd(ip,v,nds,defl,H);
84: if (which) { /* use select array */
85: for (j=0;j<n;j++)
86: if (which[j]) {
87: IPInnerProductEnd(ip,v,V[j],H+nds+j);
88: }
89: } else {
90: IPMInnerProductEnd(ip,v,n,V,H+nds);
91: }
92: if (onorm || (norm && !ip->matrix)) {
93: IPInnerProductEnd(ip,v,v,&alpha);
94: }
95: }
97: /* q = v - V h */
98: SlepcVecMAXPBY(v,1.0,-1.0,nds,H,defl);
99: if (which) {
100: for (j=0;j<n;j++)
101: if (which[j]) {
102: VecAXPBY(v,-H[nds+j],1.0,V[j]);
103: }
104: } else {
105: SlepcVecMAXPBY(v,1.0,-1.0,n,H+nds,V);
106: }
108: /* compute |v| */
109: if (onorm) *onorm = PetscSqrtReal(PetscRealPart(alpha));
111: if (norm) {
112: if (!ip->matrix) {
113: /* estimate |v'| from |v| */
114: sum = 0.0;
115: for (j=0; j<nds; j++)
116: sum += PetscRealPart(H[j] * PetscConj(H[j]));
117: for (j=0; j<n; j++)
118: if (!which || which[j])
119: sum += PetscRealPart(H[nds+j] * PetscConj(H[nds+j]));
120: *norm = PetscRealPart(alpha)-sum;
121: if (*norm <= 0.0) {
122: IPNorm(ip,v,norm);
123: } else *norm = PetscSqrtReal(*norm);
124: } else {
125: /* compute |v'| */
126: IPNorm(ip,v,norm);
127: }
128: }
129: return(0);
130: }
132: /*
133: IPOrthogonalizeMGS - Orthogonalize with modified Gram-Schmidt
134: */
137: static PetscErrorCode IPOrthogonalizeMGS(IP ip,PetscInt nds,Vec *defl,PetscInt n,PetscBool *which,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
138: {
140: PetscInt i,k;
141: PetscReal onrm,nrm;
144: if (H) {
145: for (i=0;i<n;i++) H[i] = 0;
146: }
148: switch (ip->orthog_ref) {
150: case IP_ORTHOG_REFINE_NEVER:
151: IPOrthogonalizeMGS1(ip,nds,NULL,defl,v,NULL);
152: IPOrthogonalizeMGS1(ip,n,which,V,v,H);
153: /* compute |v| */
154: if (norm) { IPNorm(ip,v,norm); }
155: /* linear dependence check does not work without refinement */
156: if (lindep) *lindep = PETSC_FALSE;
157: break;
159: case IP_ORTHOG_REFINE_ALWAYS:
160: /* first step */
161: IPOrthogonalizeMGS1(ip,nds,NULL,defl,v,NULL);
162: IPOrthogonalizeMGS1(ip,n,which,V,v,H);
163: if (lindep) { IPNorm(ip,v,&onrm); }
164: /* second step */
165: IPOrthogonalizeMGS1(ip,nds,NULL,defl,v,NULL);
166: IPOrthogonalizeMGS1(ip,n,which,V,v,H);
167: if (lindep || norm) { IPNorm(ip,v,&nrm); }
168: if (lindep) {
169: if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
170: else *lindep = PETSC_FALSE;
171: }
172: if (norm) *norm = nrm;
173: break;
175: case IP_ORTHOG_REFINE_IFNEEDED:
176: /* first step */
177: IPNorm(ip,v,&onrm);
178: IPOrthogonalizeMGS1(ip,nds,NULL,defl,v,NULL);
179: IPOrthogonalizeMGS1(ip,n,which,V,v,H);
180: IPNorm(ip,v,&nrm);
181: /* ||q|| < eta ||h|| */
182: k = 1;
183: while (k<3 && nrm < ip->orthog_eta * onrm) {
184: k++;
185: onrm = nrm;
186: IPOrthogonalizeMGS1(ip,nds,NULL,defl,v,NULL);
187: IPOrthogonalizeMGS1(ip,n,which,V,v,H);
188: IPNorm(ip,v,&nrm);
189: }
190: if (lindep) {
191: if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
192: else *lindep = PETSC_FALSE;
193: }
194: if (norm) *norm = nrm;
195: break;
197: default:
198: SETERRQ(PetscObjectComm((PetscObject)ip),PETSC_ERR_ARG_WRONG,"Unknown orthogonalization refinement");
199: }
200: return(0);
201: }
203: /*
204: IPOrthogonalizeCGS - Orthogonalize with classical Gram-Schmidt
205: */
208: static PetscErrorCode IPOrthogonalizeCGS(IP ip,PetscInt nds,Vec *defl,PetscInt n,PetscBool *which,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
209: {
211: PetscScalar *h,*c;
212: PetscReal onrm,nrm;
213: PetscInt sz=0,sz1,j,k;
216: /* allocate h and c if needed */
217: if (!H || nds>0) sz = nds+n;
218: sz1 = sz;
219: if (ip->orthog_ref != IP_ORTHOG_REFINE_NEVER) sz += nds+n;
220: if (sz>ip->lwork) {
221: PetscFree(ip->work);
222: PetscMalloc(sz*sizeof(PetscScalar),&ip->work);
223: PetscLogObjectMemory(ip,(sz-ip->lwork)*sizeof(PetscScalar));
224: ip->lwork = sz;
225: }
226: if (!H || nds>0) h = ip->work;
227: else h = H;
228: if (ip->orthog_ref != IP_ORTHOG_REFINE_NEVER) c = ip->work + sz1;
230: /* orthogonalize and compute onorm */
231: switch (ip->orthog_ref) {
233: case IP_ORTHOG_REFINE_NEVER:
234: IPOrthogonalizeCGS1(ip,nds,defl,n,which,V,v,h,NULL,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;
241: case IP_ORTHOG_REFINE_ALWAYS:
242: IPOrthogonalizeCGS1(ip,nds,defl,n,which,V,v,h,NULL,NULL);
243: if (lindep) {
244: IPOrthogonalizeCGS1(ip,nds,defl,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,defl,n,which,V,v,c,NULL,norm);
250: }
251: for (j=0;j<n;j++)
252: if (!which || which[j]) h[nds+j] += c[nds+j];
253: break;
255: case IP_ORTHOG_REFINE_IFNEEDED:
256: IPOrthogonalizeCGS1(ip,nds,defl,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,defl,n,which,V,v,c,&onrm,&nrm);
263: } else {
264: onrm = nrm;
265: IPOrthogonalizeCGS1(ip,nds,defl,n,which,V,v,c,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(PetscObjectComm((PetscObject)ip),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: }
286: return(0);
287: }
291: /*@
292: IPOrthogonalize - Orthogonalize a vector with respect to a set of vectors.
294: Collective on IP and Vec
296: Input Parameters:
297: + ip - the inner product (IP) context
298: . nds - number of columns of defl
299: . defl - first set of vectors
300: . n - number of columns of V
301: . which - logical array indicating columns of V to be used
302: - V - second set of vectors
304: Input/Output Parameter:
305: . v - (input) vector to be orthogonalized and (output) result of
306: orthogonalization
308: Output Parameter:
309: + H - coefficients computed during orthogonalization with V
310: . norm - norm of the vector after being orthogonalized
311: - lindep - flag indicating that refinement did not improve the quality
312: of orthogonalization
314: Notes:
315: This function applies an orthogonal projector to project vector v onto the
316: orthogonal complement of the span of the columns of defl and V. The columns
317: of defl and V are assumed to be mutually orthonormal.
319: On exit, v = v0 - V*H, where v0 is the original vector v.
321: This routine does not normalize the resulting vector.
323: Level: developer
325: .seealso: IPSetOrthogonalization(), IPBiOrthogonalize()
326: @*/
327: PetscErrorCode IPOrthogonalize(IP ip,PetscInt nds,Vec *defl,PetscInt n,PetscBool *which,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
328: {
335: PetscLogEventBegin(IP_Orthogonalize,ip,0,0,0);
336: if (!nds && !n) {
337: if (norm) { IPNorm(ip,v,norm); }
338: if (lindep) *lindep = PETSC_FALSE;
339: } else {
340: switch (ip->orthog_type) {
341: case IP_ORTHOG_CGS:
342: IPOrthogonalizeCGS(ip,nds,defl,n,which,V,v,H,norm,lindep);
343: break;
344: case IP_ORTHOG_MGS:
345: IPOrthogonalizeMGS(ip,nds,defl,n,which,V,v,H,norm,lindep);
346: break;
347: default:
348: SETERRQ(PetscObjectComm((PetscObject)ip),PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
349: }
350: }
351: PetscLogEventEnd(IP_Orthogonalize,ip,0,0,0);
352: return(0);
353: }
357: /*@
358: IPQRDecomposition - Compute the QR factorization of a set of vectors.
360: Collective on IP and Vec
362: Input Parameters:
363: + ip - the eigenproblem solver context
364: . V - set of vectors
365: . m - starting column
366: . n - ending column
367: - ldr - leading dimension of R
369: Output Parameter:
370: . R - triangular matrix of QR factorization
372: Notes:
373: This routine orthonormalizes the columns of V so that V'*V=I up to
374: working precision. It assumes that the first m columns (from 0 to m-1)
375: are already orthonormal. The coefficients of the orthogonalization are
376: stored in R so that V*R is equal to the original V.
378: In some cases, this routine makes V B-orthonormal, that is, V'*B*V=I.
380: If one of the vectors is linearly dependent wrt the rest (at working
381: precision) then it is replaced by a random vector.
383: Level: developer
385: .seealso: IPOrthogonalize(), IPNorm(), IPInnerProduct().
386: @*/
387: PetscErrorCode IPQRDecomposition(IP ip,Vec *V,PetscInt m,PetscInt n,PetscScalar *R,PetscInt ldr)
388: {
390: PetscInt k;
391: PetscReal norm;
392: PetscBool lindep;
393: PetscRandom rctx=NULL;
396: for (k=m; k<n; k++) {
398: /* orthogonalize v_k with respect to v_0, ..., v_{k-1} */
399: if (R) {
400: IPOrthogonalize(ip,0,NULL,k,NULL,V,V[k],&R[0+ldr*k],&norm,&lindep);
401: } else {
402: IPOrthogonalize(ip,0,NULL,k,NULL,V,V[k],NULL,&norm,&lindep);
403: }
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(PetscObjectComm((PetscObject)ip),&rctx);
410: PetscRandomSetSeed(rctx,0x12345678);
411: PetscRandomSetFromOptions(rctx);
412: }
413: SlepcVecSetRandom(V[k],rctx);
414: IPNorm(ip,V[k],&norm);
415: }
416: VecScale(V[k],1.0/norm);
417: if (R) R[k+ldr*k] = norm;
419: }
420: PetscRandomDestroy(&rctx);
421: return(0);
422: }
426: /*
427: Given m vectors in W, this function transfers them to V while
428: orthonormalizing them with respect to ip.
429: The original vectors in W are destroyed.
430: */
431: PetscErrorCode IPOrthonormalizeBasis_Private(IP ip,PetscInt *m,Vec **W,Vec *V)
432: {
434: PetscInt i,k;
435: PetscBool lindep;
436: PetscReal norm;
439: k = 0;
440: for (i=0;i<*m;i++) {
441: VecCopy((*W)[i],V[k]);
442: VecDestroy(&(*W)[i]);
443: IPOrthogonalize(ip,0,NULL,k,NULL,V,V[k],NULL,&norm,&lindep);
444: if (norm==0.0 || lindep) { PetscInfo(ip,"Linearly dependent vector found, removing...\n"); }
445: else {
446: VecScale(V[k],1.0/norm);
447: k++;
448: }
449: }
450: *m = k;
451: PetscFree(*W);
452: return(0);
453: }