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/sles/ksp/impls/cg/cgctx.h
  8: static int LINPACKcgtql1(int *,PetscReal *,PetscReal *,int *);

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

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

 21:   PetscMemzero(c,nmax*sizeof(PetscReal));
 22:   if (!n) {
 23:     *r = 0.0;
 24:     return(0);
 25:   }
 26:   d = cgP->d; e = cgP->e; ee = cgP->ee;

 28:   /* copy tridiagonal matrix to work space */
 29:   for (j=0; j<n ; j++) {
 30:     r[j]  = PetscRealPart(d[j]);
 31:     ee[j] = PetscRealPart(e[j]);
 32:   }

 34:   LINPACKcgtql1(&n,r,ee,&j);
 35:   if (j != 0) SETERRQ(PETSC_ERR_LIB,"Error from tql1(); eispack eigenvalue routine");
 36:   PetscSortReal(n,r);
 37:   return(0);
 38: }

 40: int KSPComputeExtremeSingularValues_CG(KSP ksp,PetscReal *emax,PetscReal *emin)
 41: {
 42:   KSP_CG      *cgP = (KSP_CG*)ksp->data;
 43:   PetscScalar *d,*e;
 44:   PetscReal   *dd,*ee;
 45:   int         j,n = ksp->its;

 48:   if (!n) {
 49:     *emax = *emin = 1.0;
 50:     return(0);
 51:   }
 52:   d = cgP->d; e = cgP->e; dd = cgP->dd; ee = cgP->ee;

 54:   /* copy tridiagonal matrix to work space */
 55:   for (j=0; j<n ; j++) {
 56:     dd[j] = PetscRealPart(d[j]);
 57:     ee[j] = PetscRealPart(e[j]);
 58:   }

 60:   LINPACKcgtql1(&n,dd,ee,&j);
 61:   if (j != 0) SETERRQ(PETSC_ERR_LIB,"Error from tql1(); eispack eigenvalue routine");
 62:   *emin = dd[0]; *emax = dd[n-1];
 63:   return(0);
 64: }

 66: /* tql1.f -- translated by f2c (version of 25 March 1992  12:58:56).
 67:    By Barry Smith on March 27, 1994. 
 68:    Eispack routine to determine eigenvalues of symmetric 
 69:    tridiagonal matrix 

 71:   Note that this routine always uses real numbers (not complex) even if the underlying 
 72:   matrix is Hermitian. This is because the Lanczos process applied to Hermitian matrices
 73:   always produces a real, symmetric tridiagonal matrix.
 74: */

 76: static PetscReal LINPACKcgpthy(PetscReal*,PetscReal*);

 78: static int LINPACKcgtql1(int *n,PetscReal *d,PetscReal *e,int *ierr)
 79: {
 80:     /* System generated locals */
 81:     int    i__1,i__2;
 82:     PetscReal d__1,d__2,c_b10 = 1.0;

 84:     /* Local variables */
 85:      PetscReal c,f,g,h;
 86:      int    i,j,l,m;
 87:      PetscReal p,r,s,c2,c3 = 0.0;
 88:      int    l1,l2;
 89:      PetscReal s2 = 0.0;
 90:      int    ii;
 91:      PetscReal dl1,el1;
 92:      int    mml;
 93:      PetscReal tst1,tst2;

 95: /*     THIS SUBROUTINE IS A TRANSLATION OF THE ALGOL PROCEDURE TQL1, */
 96: /*     NUM. MATH. 11, 293-306(1968) BY BOWDLER, MARTIN, REINSCH, AND */
 97: /*     WILKINSON. */
 98: /*     HANDBOOK FOR AUTO. COMP., VOL.II-LINEAR ALGEBRA, 227-240(1971). */

100: /*     THIS SUBROUTINE FINDS THE EIGENVALUES OF A SYMMETRIC */
101: /*     TRIDIAGONAL MATRIX BY THE QL METHOD. */

103: /*     ON INPUT */

105: /*        N IS THE ORDER OF THE MATRIX. */

107: /*        D CONTAINS THE DIAGONAL ELEMENTS OF THE INPUT MATRIX. */

109: /*        E CONTAINS THE SUBDIAGONAL ELEMENTS OF THE INPUT MATRIX */
110: /*          IN ITS LAST N-1 POSITIONS.  E(1) IS ARBITRARY. */

112: /*      ON OUTPUT */

114: /*        D CONTAINS THE EIGENVALUES IN ASCENDING ORDER.  IF AN */
115: /*          ERROR EXIT IS MADE, THE EIGENVALUES ARE CORRECT AND */
116: /*          ORDERED FOR INDICES 1,2,...IERR-1, BUT MAY NOT BE */
117: /*          THE SMALLEST EIGENVALUES. */

119: /*        E HAS BEEN DESTROYED. */

121: /*        IERR IS SET TO */
122: /*          ZERO       FOR NORMAL RETURN, */
123: /*          J          IF THE J-TH EIGENVALUE HAS NOT BEEN */
124: /*                     DETERMINED AFTER 30 ITERATIONS. */

126: /*     CALLS CGPTHY FOR  DSQRT(A*A + B*B) . */

128: /*     QUESTIONS AND COMMENTS SHOULD BE DIRECTED TO BURTON S. GARBOW, */
129: /*     MATHEMATICS AND COMPUTER SCIENCE DIV, ARGONNE NATIONAL LABORATORY 
130: */

132: /*     THIS VERSION DATED AUGUST 1983. */

134: /*     ------------------------------------------------------------------ 
135: */
136:     PetscReal ds;

139:     --e;
140:     --d;

142:     *0;
143:     if (*n == 1) {
144:         goto L1001;
145:     }

147:     i__1 = *n;
148:     for (i = 2; i <= i__1; ++i) {
149:         e[i - 1] = e[i];
150:     }

152:     f = 0.;
153:     tst1 = 0.;
154:     e[*n] = 0.;

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

197:         i__2 = *n;
198:         for (i = l2; i <= i__2; ++i) {
199:             d[i] -= h;
200:         }

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

251: L250:
252:         i = 1;
253: L270:
254:         d[i] = p;
255:     }

257:     goto L1001;
258: /*     .......... SET ERROR -- NO CONVERGENCE TO AN */
259: /*                EIGENVALUE AFTER 30 ITERATIONS .......... */
260: L1000:
261:     *l;
262: L1001:
263:     return(0);
264: } /* cgtql1_ */

266: static PetscReal LINPACKcgpthy(PetscReal *a,PetscReal *b)
267: {
268:     /* System generated locals */
269:     PetscReal ret_val,d__1,d__2,d__3;
270: 
271:     /* Local variables */
272:     PetscReal p,r,s,t,u;

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


278: /* Computing MAX */
279:     d__1 = PetscAbsReal(*a),d__2 = PetscAbsReal(*b);
280:     p = PetscMax(d__1,d__2);
281:     if (!p) {
282:         goto L20;
283:     }
284: /* Computing MIN */
285:     d__2 = PetscAbsReal(*a),d__3 = PetscAbsReal(*b);
286: /* Computing 2nd power */
287:     d__1 = PetscMin(d__2,d__3) / p;
288:     r = d__1 * d__1;
289: L10:
290:     t = r + 4.;
291:     if (t == 4.) {
292:         goto L20;
293:     }
294:     s = r / t;
295:     u = s * 2. + 1.;
296:     p = u * p;
297: /* Computing 2nd power */
298:     d__1 = s / u;
299:     r = d__1 * d__1 * r;
300:     goto L10;
301: L20:
302:     ret_val = p;
303:     PetscFunctionReturn(ret_val);
304: } /* cgpthy_ */