Actual source code: cgeig.c

petsc-3.4.2 2013-07-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 = ksp->its;

 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:     *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) goto L1001;


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

156:   f     = 0.;
157:   tst1  = 0.;
158:   e[*n] = 0.;

160:   i__1 = *n;
161:   for (l = 1; l <= i__1; ++l) {
162:     j = 0;
163:     h = (d__1 = d[l],PetscAbsReal(d__1)) + (d__2 = e[l],PetscAbsReal(d__2));
164:     if (tst1 < h) tst1 = h;
165: /*     .......... LOOK FOR SMALL SUB-DIAGONAL ELEMENT .......... */
166:     i__2 = *n;
167:     for (m = l; m <= i__2; ++m) {
168:       tst2 = tst1 + (d__1 = e[m],PetscAbsReal(d__1));
169:       if (tst2 == tst1) goto L120;
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) goto L210;
175: L130:
176:     if (j == 30) goto L1000;
177:     ++j;
178: /*     .......... FORM SHIFT .......... */
179:     l1    = l + 1;
180:     l2    = l1 + 1;
181:     g     = d[l];
182:     p     = (d[l1] - g) / (e[l] * 2.);
183:     r     = LINPACKcgpthy(&p,&c_b10);
184:     ds    = 1.0; if (p < 0.0) ds = -1.0;
185:     d[l]  = e[l] / (p + ds*r);
186:     d[l1] = e[l] * (p + ds*r);
187:     dl1   = d[l1];
188:     h     = g - d[l];
189:     if (l2 > *n) goto L145;

191:     i__2 = *n;
192:     for (i = l2; i <= i__2; ++i) d[i] -= h;

194: L145:
195:     f += h;
196: /*     .......... QL TRANSFORMATION .......... */
197:     p   = d[m];
198:     c   = 1.;
199:     c2  = c;
200:     el1 = e[l1];
201:     s   = 0.;
202:     mml = m - l;
203: /*     .......... FOR I=M-1 STEP -1 UNTIL L DO -- .......... */
204:     i__2 = mml;
205:     for (ii = 1; ii <= i__2; ++ii) {
206:       c3       = c2;
207:       c2       = c;
208:       s2       = s;
209:       i        = m - ii;
210:       g        = c * e[i];
211:       h        = c * p;
212:       r        = LINPACKcgpthy(&p,&e[i]);
213:       e[i + 1] = s * r;
214:       s        = e[i] / r;
215:       c        = p / r;
216:       p        = c * d[i] - s * g;
217:       d[i + 1] = h + s * (c * g + s * d[i]);
218:     }

220:     p    = -s * s2 * c3 * el1 * e[l] / dl1;
221:     e[l] = s * p;
222:     d[l] = c * p;
223:     tst2 = tst1 + (d__1 = e[l],PetscAbsReal(d__1));
224:     if (tst2 > tst1) goto L130;
225: L210:
226:     p = d[l] + f;
227: /*     .......... ORDER EIGENVALUES .......... */
228:     if (l == 1) goto L250;
229: /*     .......... FOR I=L STEP -1 UNTIL 2 DO -- .......... */
230:     i__2 = l;
231:     for (ii = 2; ii <= i__2; ++ii) {
232:       i = l + 2 - ii;
233:       if (p >= d[i - 1]) goto L270;
234:       d[i] = d[i - 1];
235:     }

237: L250:
238:     i = 1;
239: L270:
240:     d[i] = p;
241:   }

243:   goto L1001;
244: /*     .......... SET ERROR -- NO CONVERGENCE TO AN */
245: /*                EIGENVALUE AFTER 30 ITERATIONS .......... */
246: L1000:
247:   *l;
248: L1001:
249:   return(0);
250: } /* cgtql1_ */

254: static PetscReal LINPACKcgpthy(PetscReal *a,PetscReal *b)
255: {
256:   /* System generated locals */
257:   PetscReal ret_val,d__1,d__2,d__3;

259:   /* Local variables */
260:   PetscReal p,r,s,t,u;

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


266: /* Computing MAX */
267:   d__1 = PetscAbsReal(*a),d__2 = PetscAbsReal(*b);
268:   p    = PetscMax(d__1,d__2);
269:   if (!p) goto L20;
270: /* Computing MIN */
271:   d__2 = PetscAbsReal(*a),d__3 = PetscAbsReal(*b);
272: /* Computing 2nd power */
273:   d__1 = PetscMin(d__2,d__3) / p;
274:   r    = d__1 * d__1;
275: L10:
276:   t = r + 4.;
277:   if (t == 4.) goto L20;
278:   s = r / t;
279:   u = s * 2. + 1.;
280:   p = u * p;
281: /* Computing 2nd power */
282:   d__1 = s / u;
283:   r    = d__1 * d__1 * r;
284:   goto L10;
285: L20:
286:   ret_val = p;
287:   PetscFunctionReturn(ret_val);
288: } /* cgpthy_ */