Actual source code: cgeig.c

  1: /*$Id: cgeig.c,v 1.55 2001/08/07 21:30:43 bsmith Exp $*/
  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/cgctx.h
  8: static int LINPACKcgtql1(int *,PetscReal *,PetscReal *,int *);

 12: int KSPComputeEigenvalues_CG(KSP ksp,int nmax,PetscReal *r,PetscReal *c,int *neig)
 13: {
 14:   KSP_CG      *cgP = (KSP_CG*)ksp->data;
 15:   PetscScalar *d,*e;
 16:   PetscReal   *ee;
 17:   int          j,n = ksp->its,ierr;

 20:   if (nmax < n) SETERRQ(PETSC_ERR_ARG_SIZ,"Not enough room in work space r and c for eigenvalues");
 21:   *neig = n;

 23:   PetscMemzero(c,nmax*sizeof(PetscReal));
 24:   if (!n) {
 25:     *r = 0.0;
 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_ERR_LIB,"Error from tql1(); eispack eigenvalue routine");
 38:   PetscSortReal(n,r);
 39:   return(0);
 40: }

 44: int 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:   int         j,n = ksp->its;

 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_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 int LINPACKcgtql1(int *n,PetscReal *d,PetscReal *e,int *ierr)
 85: {
 86:     /* System generated locals */
 87:     int    i__1,i__2;
 88:     PetscReal d__1,d__2,c_b10 = 1.0;

 90:     /* Local variables */
 91:      PetscReal c,f,g,h;
 92:      int    i,j,l,m;
 93:      PetscReal p,r,s,c2,c3 = 0.0;
 94:      int    l1,l2;
 95:      PetscReal s2 = 0.0;
 96:      int    ii;
 97:      PetscReal dl1,el1;
 98:      int    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) {
150:         goto L1001;
151:     }

153:     i__1 = *n;
154:     for (i = 2; i <= i__1; ++i) {
155:         e[i - 1] = e[i];
156:     }

158:     f = 0.;
159:     tst1 = 0.;
160:     e[*n] = 0.;

162:     i__1 = *n;
163:     for (l = 1; l <= i__1; ++l) {
164:         j = 0;
165:         h = (d__1 = d[l],PetscAbsReal(d__1)) + (d__2 = e[l],PetscAbsReal(d__2));
166:         if (tst1 < h) {
167:             tst1 = h;
168:         }
169: /*     .......... LOOK FOR SMALL SUB-DIAGONAL ELEMENT .......... */
170:         i__2 = *n;
171:         for (m = l; m <= i__2; ++m) {
172:             tst2 = tst1 + (d__1 = e[m],PetscAbsReal(d__1));
173:             if (tst2 == tst1) {
174:                 goto L120;
175:             }
176: /*     .......... E(N) IS ALWAYS ZERO,SO THERE IS NO EXIT */
177: /*                THROUGH THE BOTTOM OF THE LOOP .......... */
178:         }
179: L120:
180:         if (m == l) {
181:             goto L210;
182:         }
183: L130:
184:         if (j == 30) {
185:             goto L1000;
186:         }
187:         ++j;
188: /*     .......... FORM SHIFT .......... */
189:         l1 = l + 1;
190:         l2 = l1 + 1;
191:         g = d[l];
192:         p = (d[l1] - g) / (e[l] * 2.);
193:         r = LINPACKcgpthy(&p,&c_b10);
194:         ds = 1.0; if (p < 0.0) ds = -1.0;
195:         d[l] = e[l] / (p + ds*r);
196:         d[l1] = e[l] * (p + ds*r);
197:         dl1 = d[l1];
198:         h = g - d[l];
199:         if (l2 > *n) {
200:             goto L145;
201:         }

203:         i__2 = *n;
204:         for (i = l2; i <= i__2; ++i) {
205:             d[i] -= h;
206:         }

208: L145:
209:         f += h;
210: /*     .......... QL TRANSFORMATION .......... */
211:         p = d[m];
212:         c = 1.;
213:         c2 = c;
214:         el1 = e[l1];
215:         s = 0.;
216:         mml = m - l;
217: /*     .......... FOR I=M-1 STEP -1 UNTIL L DO -- .......... */
218:         i__2 = mml;
219:         for (ii = 1; ii <= i__2; ++ii) {
220:             c3 = c2;
221:             c2 = c;
222:             s2 = s;
223:             i = m - ii;
224:             g = c * e[i];
225:             h = c * p;
226:             r = LINPACKcgpthy(&p,&e[i]);
227:             e[i + 1] = s * r;
228:             s = e[i] / r;
229:             c = p / r;
230:             p = c * d[i] - s * g;
231:             d[i + 1] = h + s * (c * g + s * d[i]);
232:         }
233: 
234:         p = -s * s2 * c3 * el1 * e[l] / dl1;
235:         e[l] = s * p;
236:         d[l] = c * p;
237:         tst2 = tst1 + (d__1 = e[l],PetscAbsReal(d__1));
238:         if (tst2 > tst1) {
239:             goto L130;
240:         }
241: L210:
242:         p = d[l] + f;
243: /*     .......... ORDER EIGENVALUES .......... */
244:         if (l == 1) {
245:             goto L250;
246:         }
247: /*     .......... FOR I=L STEP -1 UNTIL 2 DO -- .......... */
248:         i__2 = l;
249:         for (ii = 2; ii <= i__2; ++ii) {
250:             i = l + 2 - ii;
251:             if (p >= d[i - 1]) {
252:                 goto L270;
253:             }
254:             d[i] = d[i - 1];
255:         }

257: L250:
258:         i = 1;
259: L270:
260:         d[i] = p;
261:     }

263:     goto L1001;
264: /*     .......... SET ERROR -- NO CONVERGENCE TO AN */
265: /*                EIGENVALUE AFTER 30 ITERATIONS .......... */
266: L1000:
267:     *l;
268: L1001:
269:     return(0);
270: } /* cgtql1_ */

274: static PetscReal LINPACKcgpthy(PetscReal *a,PetscReal *b)
275: {
276:     /* System generated locals */
277:     PetscReal ret_val,d__1,d__2,d__3;
278: 
279:     /* Local variables */
280:     PetscReal p,r,s,t,u;

283: /*     FINDS DSQRT(A**2+B**2) WITHOUT OVERFLOW OR DESTRUCTIVE UNDERFLOW */


286: /* Computing MAX */
287:     d__1 = PetscAbsReal(*a),d__2 = PetscAbsReal(*b);
288:     p = PetscMax(d__1,d__2);
289:     if (!p) {
290:         goto L20;
291:     }
292: /* Computing MIN */
293:     d__2 = PetscAbsReal(*a),d__3 = PetscAbsReal(*b);
294: /* Computing 2nd power */
295:     d__1 = PetscMin(d__2,d__3) / p;
296:     r = d__1 * d__1;
297: L10:
298:     t = r + 4.;
299:     if (t == 4.) {
300:         goto L20;
301:     }
302:     s = r / t;
303:     u = s * 2. + 1.;
304:     p = u * p;
305: /* Computing 2nd power */
306:     d__1 = s / u;
307:     r = d__1 * d__1 * r;
308:     goto L10;
309: L20:
310:     ret_val = p;
311:     PetscFunctionReturn(ret_val);
312: } /* cgpthy_ */