Actual source code: cgeig.c
2: /*
3: Code for calculating extreme eigenvalues via the Lanczo method
4: running with CG. Note this only works for symmetric real and Hermitian
5: matrices (not complex matrices that are symmetric).
6: */
7: #include <../src/ksp/ksp/impls/cg/cgimpl.h>
8: static PetscErrorCode LINPACKcgtql1(PetscInt*,PetscReal *,PetscReal *,PetscInt *);
12: PetscErrorCode KSPComputeEigenvalues_CG(KSP ksp,PetscInt nmax,PetscReal *r,PetscReal *c,PetscInt *neig)
13: {
14: KSP_CG *cgP = (KSP_CG*)ksp->data;
15: PetscScalar *d,*e;
16: PetscReal *ee;
18: PetscInt j,n = ksp->its;
21: if (nmax < n) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_SIZ,"Not enough room in work space r and c for eigenvalues");
22: *neig = n;
24: PetscMemzero(c,nmax*sizeof(PetscReal));
25: if (!n) {
26: *r = 0.0;
27: return(0);
28: }
29: d = cgP->d; e = cgP->e; ee = cgP->ee;
31: /* copy tridiagonal matrix to work space */
32: for (j=0; j<n ; j++) {
33: r[j] = PetscRealPart(d[j]);
34: ee[j] = PetscRealPart(e[j]);
35: }
37: LINPACKcgtql1(&n,r,ee,&j);
38: if (j != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error from tql1(); eispack eigenvalue routine");
39: PetscSortReal(n,r);
40: return(0);
41: }
45: PetscErrorCode KSPComputeExtremeSingularValues_CG(KSP ksp,PetscReal *emax,PetscReal *emin)
46: {
47: KSP_CG *cgP = (KSP_CG*)ksp->data;
48: PetscScalar *d,*e;
49: PetscReal *dd,*ee;
50: PetscInt j,n = ksp->its;
53: if (!n) {
54: *emax = *emin = 1.0;
55: return(0);
56: }
57: d = cgP->d; e = cgP->e; dd = cgP->dd; ee = cgP->ee;
59: /* copy tridiagonal matrix to work space */
60: for (j=0; j<n ; j++) {
61: dd[j] = PetscRealPart(d[j]);
62: ee[j] = PetscRealPart(e[j]);
63: }
65: LINPACKcgtql1(&n,dd,ee,&j);
66: if (j != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error from tql1(); eispack eigenvalue routine");
67: *emin = dd[0]; *emax = dd[n-1];
68: return(0);
69: }
71: /* tql1.f -- translated by f2c (version of 25 March 1992 12:58:56).
72: By Barry Smith on March 27, 1994.
73: Eispack routine to determine eigenvalues of symmetric
74: tridiagonal matrix
76: Note that this routine always uses real numbers (not complex) even if the underlying
77: matrix is Hermitian. This is because the Lanczos process applied to Hermitian matrices
78: always produces a real, symmetric tridiagonal matrix.
79: */
81: static PetscReal LINPACKcgpthy(PetscReal*,PetscReal*);
85: static PetscErrorCode LINPACKcgtql1(PetscInt *n,PetscReal *d,PetscReal *e,PetscInt *ierr)
86: {
87: /* System generated locals */
88: PetscInt i__1,i__2;
89: PetscReal d__1,d__2,c_b10 = 1.0;
91: /* Local variables */
92: PetscReal c,f,g,h;
93: PetscInt i,j,l,m;
94: PetscReal p,r,s,c2,c3 = 0.0;
95: PetscInt l1,l2;
96: PetscReal s2 = 0.0;
97: PetscInt ii;
98: PetscReal dl1,el1;
99: PetscInt mml;
100: PetscReal tst1,tst2;
102: /* THIS SUBROUTINE IS A TRANSLATION OF THE ALGOL PROCEDURE TQL1, */
103: /* NUM. MATH. 11, 293-306(1968) BY BOWDLER, MARTIN, REINSCH, AND */
104: /* WILKINSON. */
105: /* HANDBOOK FOR AUTO. COMP., VOL.II-LINEAR ALGEBRA, 227-240(1971). */
107: /* THIS SUBROUTINE FINDS THE EIGENVALUES OF A SYMMETRIC */
108: /* TRIDIAGONAL MATRIX BY THE QL METHOD. */
110: /* ON INPUT */
112: /* N IS THE ORDER OF THE MATRIX. */
114: /* D CONTAINS THE DIAGONAL ELEMENTS OF THE INPUT MATRIX. */
116: /* E CONTAINS THE SUBDIAGONAL ELEMENTS OF THE INPUT MATRIX */
117: /* IN ITS LAST N-1 POSITIONS. E(1) IS ARBITRARY. */
119: /* ON OUTPUT */
121: /* D CONTAINS THE EIGENVALUES IN ASCENDING ORDER. IF AN */
122: /* ERROR EXIT IS MADE, THE EIGENVALUES ARE CORRECT AND */
123: /* ORDERED FOR INDICES 1,2,...IERR-1, BUT MAY NOT BE */
124: /* THE SMALLEST EIGENVALUES. */
126: /* E HAS BEEN DESTROYED. */
128: /* IERR IS SET TO */
129: /* ZERO FOR NORMAL RETURN, */
130: /* J IF THE J-TH EIGENVALUE HAS NOT BEEN */
131: /* DETERMINED AFTER 30 ITERATIONS. */
133: /* CALLS CGPTHY FOR DSQRT(A*A + B*B) . */
135: /* QUESTIONS AND COMMENTS SHOULD BE DIRECTED TO BURTON S. GARBOW, */
136: /* MATHEMATICS AND COMPUTER SCIENCE DIV, ARGONNE NATIONAL LABORATORY
137: */
139: /* THIS VERSION DATED AUGUST 1983. */
141: /* ------------------------------------------------------------------
142: */
143: PetscReal ds;
146: --e;
147: --d;
149: *0;
150: if (*n == 1) {
151: goto L1001;
152: }
154: i__1 = *n;
155: for (i = 2; i <= i__1; ++i) {
156: e[i - 1] = e[i];
157: }
159: f = 0.;
160: tst1 = 0.;
161: e[*n] = 0.;
163: i__1 = *n;
164: for (l = 1; l <= i__1; ++l) {
165: j = 0;
166: h = (d__1 = d[l],PetscAbsReal(d__1)) + (d__2 = e[l],PetscAbsReal(d__2));
167: if (tst1 < h) {
168: tst1 = h;
169: }
170: /* .......... LOOK FOR SMALL SUB-DIAGONAL ELEMENT .......... */
171: i__2 = *n;
172: for (m = l; m <= i__2; ++m) {
173: tst2 = tst1 + (d__1 = e[m],PetscAbsReal(d__1));
174: if (tst2 == tst1) {
175: goto L120;
176: }
177: /* .......... E(N) IS ALWAYS ZERO,SO THERE IS NO EXIT */
178: /* THROUGH THE BOTTOM OF THE LOOP .......... */
179: }
180: L120:
181: if (m == l) {
182: goto L210;
183: }
184: L130:
185: if (j == 30) {
186: goto L1000;
187: }
188: ++j;
189: /* .......... FORM SHIFT .......... */
190: l1 = l + 1;
191: l2 = l1 + 1;
192: g = d[l];
193: p = (d[l1] - g) / (e[l] * 2.);
194: r = LINPACKcgpthy(&p,&c_b10);
195: ds = 1.0; if (p < 0.0) ds = -1.0;
196: d[l] = e[l] / (p + ds*r);
197: d[l1] = e[l] * (p + ds*r);
198: dl1 = d[l1];
199: h = g - d[l];
200: if (l2 > *n) {
201: goto L145;
202: }
204: i__2 = *n;
205: for (i = l2; i <= i__2; ++i) {
206: d[i] -= h;
207: }
209: L145:
210: f += h;
211: /* .......... QL TRANSFORMATION .......... */
212: p = d[m];
213: c = 1.;
214: c2 = c;
215: el1 = e[l1];
216: s = 0.;
217: mml = m - l;
218: /* .......... FOR I=M-1 STEP -1 UNTIL L DO -- .......... */
219: i__2 = mml;
220: for (ii = 1; ii <= i__2; ++ii) {
221: c3 = c2;
222: c2 = c;
223: s2 = s;
224: i = m - ii;
225: g = c * e[i];
226: h = c * p;
227: r = LINPACKcgpthy(&p,&e[i]);
228: e[i + 1] = s * r;
229: s = e[i] / r;
230: c = p / r;
231: p = c * d[i] - s * g;
232: d[i + 1] = h + s * (c * g + s * d[i]);
233: }
234:
235: p = -s * s2 * c3 * el1 * e[l] / dl1;
236: e[l] = s * p;
237: d[l] = c * p;
238: tst2 = tst1 + (d__1 = e[l],PetscAbsReal(d__1));
239: if (tst2 > tst1) {
240: goto L130;
241: }
242: L210:
243: p = d[l] + f;
244: /* .......... ORDER EIGENVALUES .......... */
245: if (l == 1) {
246: goto L250;
247: }
248: /* .......... FOR I=L STEP -1 UNTIL 2 DO -- .......... */
249: i__2 = l;
250: for (ii = 2; ii <= i__2; ++ii) {
251: i = l + 2 - ii;
252: if (p >= d[i - 1]) {
253: goto L270;
254: }
255: d[i] = d[i - 1];
256: }
258: L250:
259: i = 1;
260: L270:
261: d[i] = p;
262: }
264: goto L1001;
265: /* .......... SET ERROR -- NO CONVERGENCE TO AN */
266: /* EIGENVALUE AFTER 30 ITERATIONS .......... */
267: L1000:
268: *l;
269: L1001:
270: return(0);
271: } /* cgtql1_ */
275: static PetscReal LINPACKcgpthy(PetscReal *a,PetscReal *b)
276: {
277: /* System generated locals */
278: PetscReal ret_val,d__1,d__2,d__3;
279:
280: /* Local variables */
281: PetscReal p,r,s,t,u;
284: /* FINDS DSQRT(A**2+B**2) WITHOUT OVERFLOW OR DESTRUCTIVE UNDERFLOW */
287: /* Computing MAX */
288: d__1 = PetscAbsReal(*a),d__2 = PetscAbsReal(*b);
289: p = PetscMax(d__1,d__2);
290: if (!p) {
291: goto L20;
292: }
293: /* Computing MIN */
294: d__2 = PetscAbsReal(*a),d__3 = PetscAbsReal(*b);
295: /* Computing 2nd power */
296: d__1 = PetscMin(d__2,d__3) / p;
297: r = d__1 * d__1;
298: L10:
299: t = r + 4.;
300: if (t == 4.) {
301: goto L20;
302: }
303: s = r / t;
304: u = s * 2. + 1.;
305: p = u * p;
306: /* Computing 2nd power */
307: d__1 = s / u;
308: r = d__1 * d__1 * r;
309: goto L10;
310: L20:
311: ret_val = p;
312: PetscFunctionReturn(ret_val);
313: } /* cgpthy_ */