Actual source code: cgeig.c
petsc-3.7.4 2016-10-02
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 = cgP->ned;
21: if (nmax < n) SETERRQ(PetscObjectComm((PetscObject)ksp),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: return(0);
27: }
28: d = cgP->d; e = cgP->e; ee = cgP->ee;
30: /* copy tridiagonal matrix to work space */
31: for (j=0; j<n; j++) {
32: r[j] = PetscRealPart(d[j]);
33: ee[j] = PetscRealPart(e[j]);
34: }
36: LINPACKcgtql1(&n,r,ee,&j);
37: if (j != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error from tql1(); eispack eigenvalue routine");
38: PetscSortReal(n,r);
39: return(0);
40: }
44: PetscErrorCode KSPComputeExtremeSingularValues_CG(KSP ksp,PetscReal *emax,PetscReal *emin)
45: {
46: KSP_CG *cgP = (KSP_CG*)ksp->data;
47: PetscScalar *d,*e;
48: PetscReal *dd,*ee;
49: PetscInt j,n = cgP->ned;
52: if (!n) {
53: *emax = *emin = 1.0;
54: return(0);
55: }
56: d = cgP->d; e = cgP->e; dd = cgP->dd; ee = cgP->ee;
58: /* copy tridiagonal matrix to work space */
59: for (j=0; j<n; j++) {
60: dd[j] = PetscRealPart(d[j]);
61: ee[j] = PetscRealPart(e[j]);
62: }
64: LINPACKcgtql1(&n,dd,ee,&j);
65: if (j != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error from tql1(); eispack eigenvalue routine");
66: *emin = dd[0]; *emax = dd[n-1];
67: return(0);
68: }
70: /* tql1.f -- translated by f2c (version of 25 March 1992 12:58:56).
71: By Barry Smith on March 27, 1994.
72: Eispack routine to determine eigenvalues of symmetric
73: tridiagonal matrix
75: Note that this routine always uses real numbers (not complex) even if the underlying
76: matrix is Hermitian. This is because the Lanczos process applied to Hermitian matrices
77: always produces a real, symmetric tridiagonal matrix.
78: */
80: static PetscReal LINPACKcgpthy(PetscReal*,PetscReal*);
84: static PetscErrorCode LINPACKcgtql1(PetscInt *n,PetscReal *d,PetscReal *e,PetscInt *ierr)
85: {
86: /* System generated locals */
87: PetscInt i__1,i__2;
88: PetscReal d__1,d__2,c_b10 = 1.0;
90: /* Local variables */
91: PetscReal c,f,g,h;
92: PetscInt i,j,l,m;
93: PetscReal p,r,s,c2,c3 = 0.0;
94: PetscInt l1,l2;
95: PetscReal s2 = 0.0;
96: PetscInt ii;
97: PetscReal dl1,el1;
98: PetscInt mml;
99: PetscReal tst1,tst2;
101: /* THIS SUBROUTINE IS A TRANSLATION OF THE ALGOL PROCEDURE TQL1, */
102: /* NUM. MATH. 11, 293-306(1968) BY BOWDLER, MARTIN, REINSCH, AND */
103: /* WILKINSON. */
104: /* HANDBOOK FOR AUTO. COMP., VOL.II-LINEAR ALGEBRA, 227-240(1971). */
106: /* THIS SUBROUTINE FINDS THE EIGENVALUES OF A SYMMETRIC */
107: /* TRIDIAGONAL MATRIX BY THE QL METHOD. */
109: /* ON INPUT */
111: /* N IS THE ORDER OF THE MATRIX. */
113: /* D CONTAINS THE DIAGONAL ELEMENTS OF THE INPUT MATRIX. */
115: /* E CONTAINS THE SUBDIAGONAL ELEMENTS OF THE INPUT MATRIX */
116: /* IN ITS LAST N-1 POSITIONS. E(1) IS ARBITRARY. */
118: /* ON OUTPUT */
120: /* D CONTAINS THE EIGENVALUES IN ASCENDING ORDER. IF AN */
121: /* ERROR EXIT IS MADE, THE EIGENVALUES ARE CORRECT AND */
122: /* ORDERED FOR INDICES 1,2,...IERR-1, BUT MAY NOT BE */
123: /* THE SMALLEST EIGENVALUES. */
125: /* E HAS BEEN DESTROYED. */
127: /* IERR IS SET TO */
128: /* ZERO FOR NORMAL RETURN, */
129: /* J IF THE J-TH EIGENVALUE HAS NOT BEEN */
130: /* DETERMINED AFTER 30 ITERATIONS. */
132: /* CALLS CGPTHY FOR DSQRT(A*A + B*B) . */
134: /* QUESTIONS AND COMMENTS SHOULD BE DIRECTED TO BURTON S. GARBOW, */
135: /* MATHEMATICS AND COMPUTER SCIENCE DIV, ARGONNE NATIONAL LABORATORY
136: */
138: /* THIS VERSION DATED AUGUST 1983. */
140: /* ------------------------------------------------------------------
141: */
142: PetscReal ds;
145: --e;
146: --d;
148: *0;
149: if (*n == 1) goto L1001;
152: i__1 = *n;
153: for (i = 2; i <= i__1; ++i) e[i - 1] = e[i];
155: f = 0.;
156: tst1 = 0.;
157: e[*n] = 0.;
159: i__1 = *n;
160: for (l = 1; l <= i__1; ++l) {
161: j = 0;
162: h = (d__1 = d[l],PetscAbsReal(d__1)) + (d__2 = e[l],PetscAbsReal(d__2));
163: if (tst1 < h) tst1 = h;
164: /* .......... LOOK FOR SMALL SUB-DIAGONAL ELEMENT .......... */
165: i__2 = *n;
166: for (m = l; m <= i__2; ++m) {
167: tst2 = tst1 + (d__1 = e[m],PetscAbsReal(d__1));
168: if (tst2 == tst1) goto L120;
169: /* .......... E(N) IS ALWAYS ZERO,SO THERE IS NO EXIT */
170: /* THROUGH THE BOTTOM OF THE LOOP .......... */
171: }
172: L120:
173: if (m == l) goto L210;
174: L130:
175: if (j == 30) goto L1000;
176: ++j;
177: /* .......... FORM SHIFT .......... */
178: l1 = l + 1;
179: l2 = l1 + 1;
180: g = d[l];
181: p = (d[l1] - g) / (e[l] * 2.);
182: r = LINPACKcgpthy(&p,&c_b10);
183: ds = 1.0; if (p < 0.0) ds = -1.0;
184: d[l] = e[l] / (p + ds*r);
185: d[l1] = e[l] * (p + ds*r);
186: dl1 = d[l1];
187: h = g - d[l];
188: if (l2 > *n) goto L145;
190: i__2 = *n;
191: for (i = l2; i <= i__2; ++i) d[i] -= h;
193: L145:
194: f += h;
195: /* .......... QL TRANSFORMATION .......... */
196: p = d[m];
197: c = 1.;
198: c2 = c;
199: el1 = e[l1];
200: s = 0.;
201: mml = m - l;
202: /* .......... FOR I=M-1 STEP -1 UNTIL L DO -- .......... */
203: i__2 = mml;
204: for (ii = 1; ii <= i__2; ++ii) {
205: c3 = c2;
206: c2 = c;
207: s2 = s;
208: i = m - ii;
209: g = c * e[i];
210: h = c * p;
211: r = LINPACKcgpthy(&p,&e[i]);
212: e[i + 1] = s * r;
213: s = e[i] / r;
214: c = p / r;
215: p = c * d[i] - s * g;
216: d[i + 1] = h + s * (c * g + s * d[i]);
217: }
219: p = -s * s2 * c3 * el1 * e[l] / dl1;
220: e[l] = s * p;
221: d[l] = c * p;
222: tst2 = tst1 + (d__1 = e[l],PetscAbsReal(d__1));
223: if (tst2 > tst1) goto L130;
224: L210:
225: p = d[l] + f;
226: /* .......... ORDER EIGENVALUES .......... */
227: if (l == 1) goto L250;
228: /* .......... FOR I=L STEP -1 UNTIL 2 DO -- .......... */
229: i__2 = l;
230: for (ii = 2; ii <= i__2; ++ii) {
231: i = l + 2 - ii;
232: if (p >= d[i - 1]) goto L270;
233: d[i] = d[i - 1];
234: }
236: L250:
237: i = 1;
238: L270:
239: d[i] = p;
240: }
242: goto L1001;
243: /* .......... SET ERROR -- NO CONVERGENCE TO AN */
244: /* EIGENVALUE AFTER 30 ITERATIONS .......... */
245: L1000:
246: *l;
247: L1001:
248: return(0);
249: } /* cgtql1_ */
253: static PetscReal LINPACKcgpthy(PetscReal *a,PetscReal *b)
254: {
255: /* System generated locals */
256: PetscReal ret_val,d__1,d__2,d__3;
258: /* Local variables */
259: PetscReal p,r,s,t,u;
262: /* FINDS DSQRT(A**2+B**2) WITHOUT OVERFLOW OR DESTRUCTIVE UNDERFLOW */
265: /* Computing MAX */
266: d__1 = PetscAbsReal(*a),d__2 = PetscAbsReal(*b);
267: p = PetscMax(d__1,d__2);
268: if (!p) goto L20;
269: /* Computing MIN */
270: d__2 = PetscAbsReal(*a),d__3 = PetscAbsReal(*b);
271: /* Computing 2nd power */
272: d__1 = PetscMin(d__2,d__3) / p;
273: r = d__1 * d__1;
274: L10:
275: t = r + 4.;
276: if (t == 4.) goto L20;
277: s = r / t;
278: u = s * 2. + 1.;
279: p = u * p;
280: /* Computing 2nd power */
281: d__1 = s / u;
282: r = d__1 * d__1 * r;
283: goto L10;
284: L20:
285: ret_val = p;
286: PetscFunctionReturn(ret_val);
287: } /* cgpthy_ */