Actual source code: lanczos.c

  1: /*                       

  3:    SLEPc eigensolver: "lanczos"

  5:    Method: Explicitly Restarted Symmetric/Hermitian Lanczos

  7:    Algorithm:

  9:        Lanczos method for symmetric (Hermitian) problems, with explicit 
 10:        restart and deflation. Several reorthogonalization strategies can
 11:        be selected.

 13:    References:

 15:        [1] "Lanczos Methods in SLEPc", SLEPc Technical Report STR-5, 
 16:            available at http://www.grycap.upv.es/slepc.

 18:    Last update: Feb 2009

 20:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 21:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 22:    Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain

 24:    This file is part of SLEPc.
 25:       
 26:    SLEPc is free software: you can redistribute it and/or modify it under  the
 27:    terms of version 3 of the GNU Lesser General Public License as published by
 28:    the Free Software Foundation.

 30:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 31:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 32:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 33:    more details.

 35:    You  should have received a copy of the GNU Lesser General  Public  License
 36:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 37:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 38: */

 40:  #include private/epsimpl.h

 42: PetscErrorCode EPSSolve_LANCZOS(EPS);

 44: typedef struct {
 45:   EPSLanczosReorthogType reorthog;
 46:   Vec *AV;
 47: } EPS_LANCZOS;

 51: PetscErrorCode EPSSetUp_LANCZOS(EPS eps)
 52: {
 53:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS *)eps->data;

 57:   if (eps->ncv) { /* ncv set */
 58:     if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
 59:   }
 60:   else if (eps->mpd) { /* mpd set */
 61:     eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd);
 62:   }
 63:   else { /* neither set: defaults depend on nev being small or large */
 64:     if (eps->nev<500) eps->ncv = PetscMin(eps->n,PetscMax(2*eps->nev,eps->nev+15));
 65:     else { eps->mpd = 500; eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd); }
 66:   }
 67:   if (!eps->mpd) eps->mpd = eps->ncv;
 68:   if (eps->ncv>eps->nev+eps->mpd) SETERRQ(1,"The value of ncv must not be larger than nev+mpd");
 69:   if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);

 71:   if (!eps->which) eps->which = EPS_LARGEST_MAGNITUDE;
 72:   switch (eps->which) {
 73:     case EPS_LARGEST_IMAGINARY:
 74:     case EPS_SMALLEST_IMAGINARY:
 75:     case EPS_TARGET_IMAGINARY:
 76:       SETERRQ(1,"Wrong value of eps->which");
 77:     default: ; /* default case to remove warning */
 78:   }
 79:   if (!eps->ishermitian)
 80:     SETERRQ(PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
 81:   if (!eps->extraction) {
 82:     EPSSetExtraction(eps,EPS_RITZ);
 83:   } else if (eps->extraction!=EPS_RITZ) {
 84:     SETERRQ(PETSC_ERR_SUP,"Unsupported extraction type\n");
 85:   }

 87:   EPSAllocateSolution(eps);
 88:   if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_SELECTIVE) {
 89:     VecDuplicateVecs(eps->V[0],eps->ncv,&lanczos->AV);
 90:   }
 91:   if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
 92:     EPSDefaultGetWork(eps,2);
 93:   } else {
 94:     EPSDefaultGetWork(eps,1);
 95:   }

 97:   /* dispatch solve method */
 98:   if (eps->leftvecs) SETERRQ(PETSC_ERR_SUP,"Left vectors not supported in this solver");
 99:   eps->ops->solve = EPSSolve_LANCZOS;
100:   return(0);
101: }

105: /*
106:    EPSLocalLanczos - Local reorthogonalization.

108:    This is the simplest variant. At each Lanczos step, the corresponding Lanczos vector 
109:    is orthogonalized with respect to the two previous Lanczos vectors, according to
110:    the three term Lanczos recurrence. WARNING: This variant does not track the loss of 
111:    orthogonality that occurs in finite-precision arithmetic and, therefore, the 
112:    generated vectors are not guaranteed to be (semi-)orthogonal.
113: */
114: static PetscErrorCode EPSLocalLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscTruth *breakdown)
115: {
117:   PetscInt       i,j,m = *M;
118:   PetscReal      norm;
119:   PetscTruth     *which,lwhich[100];
120:   PetscScalar    *hwork,lhwork[100];
121: 
123:   if (m > 100) {
124:     PetscMalloc(sizeof(PetscTruth)*m,&which);
125:     PetscMalloc(m*sizeof(PetscScalar),&hwork);
126:   } else {
127:     which = lwhich;
128:     hwork = lhwork;
129:   }
130:   for (i=0;i<k;i++)
131:     which[i] = PETSC_TRUE;

133:   for (j=k;j<m-1;j++) {
134:     STApply(eps->OP,V[j],V[j+1]);
135:     which[j] = PETSC_TRUE;
136:     if (j-2>=k) which[j-2] = PETSC_FALSE;
137:     IPOrthogonalize(eps->ip,eps->nds,eps->DS,j+1,which,V,V[j+1],hwork,&norm,breakdown);
138:     alpha[j-k] = PetscRealPart(hwork[j]);
139:     beta[j-k] = norm;
140:     if (*breakdown) {
141:       *M = j+1;
142:       if (m > 100) {
143:         PetscFree(which);
144:         PetscFree(hwork);
145:       }
146:       return(0);
147:     } else {
148:       VecScale(V[j+1],1.0/norm);
149:     }
150:   }
151:   STApply(eps->OP,V[m-1],f);
152:   IPOrthogonalize(eps->ip,eps->nds,eps->DS,m,PETSC_NULL,V,f,hwork,&norm,PETSC_NULL);
153:   alpha[m-1-k] = PetscRealPart(hwork[m-1]);
154:   beta[m-1-k] = norm;

156:   if (m > 100) {
157:     PetscFree(which);
158:     PetscFree(hwork);
159:   }
160:   return(0);
161: }

165: /*
166:    EPSSelectiveLanczos - Selective reorthogonalization.
167: */
168: static PetscErrorCode EPSSelectiveLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscTruth *breakdown,PetscReal anorm)
169: {
171:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS *)eps->data;
172:   PetscInt       i,j,m = *M,n,nritz=0,nritzo;
173:   PetscReal      *d,*e,*ritz,norm;
174:   PetscScalar    *Y,*hwork,lhwork[100];
175:   PetscTruth     *which,lwhich[100];

178:   PetscMalloc(m*sizeof(PetscReal),&d);
179:   PetscMalloc(m*sizeof(PetscReal),&e);
180:   PetscMalloc(m*sizeof(PetscReal),&ritz);
181:   PetscMalloc(m*m*sizeof(PetscScalar),&Y);
182:   if (m > 100) {
183:     PetscMalloc(sizeof(PetscTruth)*m,&which);
184:     PetscMalloc(m*sizeof(PetscScalar),&hwork);
185:   } else {
186:     which = lwhich;
187:     hwork = lhwork;
188:   }
189:   for (i=0;i<k;i++)
190:     which[i] = PETSC_TRUE;

192:   for (j=k;j<m;j++) {
193:     /* Lanczos step */
194:     STApply(eps->OP,V[j],f);
195:     which[j] = PETSC_TRUE;
196:     if (j-2>=k) which[j-2] = PETSC_FALSE;
197:     IPOrthogonalize(eps->ip,eps->nds,eps->DS,j+1,which,V,f,hwork,&norm,breakdown);
198:     alpha[j-k] = PetscRealPart(hwork[j]);
199:     beta[j-k] = norm;
200:     if (*breakdown) {
201:       *M = j+1;
202:       break;
203:     }

205:     /* Compute eigenvalues and eigenvectors Y of the tridiagonal block */
206:     n = j-k+1;
207:     for (i=0;i<n;i++) { d[i] = alpha[i]; e[i] = beta[i]; }
208:     EPSDenseTridiagonal(n,d,e,ritz,Y);
209: 
210:     /* Estimate ||A|| */
211:     for (i=0;i<n;i++)
212:       if (PetscAbsReal(ritz[i]) > anorm) anorm = PetscAbsReal(ritz[i]);

214:     /* Compute nearly converged Ritz vectors */
215:     nritzo = 0;
216:     for (i=0;i<n;i++)
217:       if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm)
218:         nritzo++;

220:     if (nritzo>nritz) {
221:       nritz = 0;
222:       for (i=0;i<n;i++) {
223:         if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm) {
224:           SlepcVecMAXPBY(lanczos->AV[nritz],0.0,1.0,n,Y+i*n,V+k);
225:           nritz++;
226:         }
227:       }
228:     }

230:     if (nritz > 0) {
231:       IPOrthogonalize(eps->ip,0,PETSC_NULL,nritz,PETSC_NULL,lanczos->AV,f,hwork,&norm,breakdown);
232:       if (*breakdown) {
233:         *M = j+1;
234:         break;
235:       }
236:     }
237: 
238:     if (j<m-1) {
239:       VecScale(f,1.0 / norm);
240:       VecCopy(f,V[j+1]);
241:     }
242:   }
243: 
244:   PetscFree(d);
245:   PetscFree(e);
246:   PetscFree(ritz);
247:   PetscFree(Y);
248:   if (m > 100) {
249:     PetscFree(which);
250:     PetscFree(hwork);
251:   }
252:   return(0);
253: }

257: static void update_omega(PetscReal *omega,PetscReal *omega_old,PetscInt j,PetscReal *alpha,PetscReal *beta,PetscReal eps1,PetscReal anorm)
258: {
259:   PetscInt       k;
260:   PetscReal      T,binv;

263:   /* Estimate of contribution to roundoff errors from A*v 
264:        fl(A*v) = A*v + f, 
265:      where ||f|| \approx eps1*||A||.
266:      For a full matrix A, a rule-of-thumb estimate is eps1 = sqrt(n)*eps. */
267:   T = eps1*anorm;
268:   binv = 1.0/beta[j+1];

270:   /* Update omega(1) using omega(0)==0. */
271:   omega_old[0]= beta[1]*omega[1] + (alpha[0]-alpha[j])*omega[0] -
272:                 beta[j]*omega_old[0];
273:   if (omega_old[0] > 0)
274:     omega_old[0] = binv*(omega_old[0] + T);
275:   else
276:     omega_old[0] = binv*(omega_old[0] - T);
277: 
278:   /* Update remaining components. */
279:   for (k=1;k<j-1;k++) {
280:     omega_old[k] = beta[k+1]*omega[k+1] + (alpha[k]-alpha[j])*omega[k] +
281:                    beta[k]*omega[k-1] - beta[j]*omega_old[k];
282:     if (omega_old[k] > 0)
283:       omega_old[k] = binv*(omega_old[k] + T);
284:     else
285:       omega_old[k] = binv*(omega_old[k] - T);
286:   }
287:   omega_old[j-1] = binv*T;
288: 
289:   /* Swap omega and omega_old. */
290:   for (k=0;k<j;k++) {
291:     omega[k] = omega_old[k];
292:     omega_old[k] = omega[k];
293:   }
294:   omega[j] = eps1;
295:   PetscFunctionReturnVoid();
296: }

300: static void compute_int(PetscTruth *which,PetscReal *mu,PetscInt j,PetscReal delta,PetscReal eta)
301: {
302:   PetscInt   i,k,maxpos;
303:   PetscReal  max;
304:   PetscTruth found;
305: 
307:   /* initialize which */
308:   found = PETSC_FALSE;
309:   maxpos = 0;
310:   max = 0.0;
311:   for (i=0;i<j;i++) {
312:     if (PetscAbsReal(mu[i]) >= delta) {
313:       which[i] = PETSC_TRUE;
314:       found = PETSC_TRUE;
315:     } else which[i] = PETSC_FALSE;
316:     if (PetscAbsReal(mu[i]) > max) {
317:       maxpos = i;
318:       max = PetscAbsReal(mu[i]);
319:     }
320:   }
321:   if (!found) which[maxpos] = PETSC_TRUE;
322: 
323:   for (i=0;i<j;i++)
324:     if (which[i]) {
325:       /* find left interval */
326:       for (k=i;k>=0;k--) {
327:         if (PetscAbsReal(mu[k])<eta || which[k]) break;
328:         else which[k] = PETSC_TRUE;
329:       }
330:       /* find right interval */
331:       for (k=i+1;k<j;k++) {
332:         if (PetscAbsReal(mu[k])<eta || which[k]) break;
333:         else which[k] = PETSC_TRUE;
334:       }
335:     }
336:   PetscFunctionReturnVoid();
337: }

341: /*
342:    EPSPartialLanczos - Partial reorthogonalization.
343: */
344: static PetscErrorCode EPSPartialLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscTruth *breakdown,PetscReal anorm)
345: {
346:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS *)eps->data;
348:   PetscInt       i,j,m = *M;
349:   PetscReal      norm,*omega,lomega[100],*omega_old,lomega_old[100],eps1,delta,eta;
350:   PetscTruth     *which,lwhich[100],*which2,lwhich2[100],
351:                  reorth = PETSC_FALSE,force_reorth = PETSC_FALSE,fro = PETSC_FALSE,estimate_anorm = PETSC_FALSE;
352:   PetscScalar    *hwork,lhwork[100];

355:   if (m>100) {
356:     PetscMalloc(m*sizeof(PetscReal),&omega);
357:     PetscMalloc(m*sizeof(PetscReal),&omega_old);
358:   } else {
359:     omega = lomega;
360:     omega_old = lomega_old;
361:   }
362:   if (m > 100) {
363:     PetscMalloc(sizeof(PetscTruth)*m,&which);
364:     PetscMalloc(sizeof(PetscTruth)*m,&which2);
365:     PetscMalloc(m*sizeof(PetscScalar),&hwork);
366:   } else {
367:     which = lwhich;
368:     which2 = lwhich2;
369:     hwork = lhwork;
370:   }

372:   eps1 = sqrt((PetscReal)eps->n)*PETSC_MACHINE_EPSILON/2;
373:   delta = PETSC_SQRT_MACHINE_EPSILON/sqrt((PetscReal)eps->ncv);
374:   eta = pow(PETSC_MACHINE_EPSILON,3.0/4.0)/sqrt((PetscReal)eps->ncv);
375:   if (anorm < 0.0) {
376:     anorm = 1.0;
377:     estimate_anorm = PETSC_TRUE;
378:   }
379:   for (i=0;i<m-k;i++)
380:     omega[i] = omega_old[i] = 0.0;
381:   for (i=0;i<k;i++)
382:     which[i] = PETSC_TRUE;
383: 
384:   for (j=k;j<m;j++) {
385:     STApply(eps->OP,V[j],f);
386:     if (fro) {
387:       /* Lanczos step with full reorthogonalization */
388:       IPOrthogonalize(eps->ip,eps->nds,eps->DS,j+1,PETSC_NULL,V,f,hwork,&norm,breakdown);
389:       alpha[j-k] = PetscRealPart(hwork[j]);
390:     } else {
391:       /* Lanczos step */
392:       which[j] = PETSC_TRUE;
393:       if (j-2>=k) which[j-2] = PETSC_FALSE;
394:       IPOrthogonalize(eps->ip,eps->nds,eps->DS,j+1,which,V,f,hwork,&norm,breakdown);
395:       alpha[j-k] = PetscRealPart(hwork[j]);
396:       beta[j-k] = norm;
397: 
398:       /* Estimate ||A|| if needed */
399:       if (estimate_anorm) {
400:         if (j>k) anorm = PetscMax(anorm,PetscAbsReal(alpha[j-k])+norm+beta[j-k-1]);
401:         else anorm = PetscMax(anorm,PetscAbsReal(alpha[j-k])+norm);
402:       }

404:       /* Check if reorthogonalization is needed */
405:       reorth = PETSC_FALSE;
406:       if (j>k) {
407:         update_omega(omega,omega_old,j-k,alpha,beta-1,eps1,anorm);
408:         for (i=0;i<j-k;i++)
409:           if (PetscAbsScalar(omega[i]) > delta) reorth = PETSC_TRUE;
410:       }

412:       if (reorth || force_reorth) {
413:         if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_PERIODIC) {
414:           /* Periodic reorthogonalization */
415:           if (force_reorth) force_reorth = PETSC_FALSE;
416:           else force_reorth = PETSC_TRUE;
417:           IPOrthogonalize(eps->ip,0,PETSC_NULL,j-k,PETSC_NULL,V+k,f,hwork,&norm,breakdown);
418:           for (i=0;i<j-k;i++)
419:             omega[i] = eps1;
420:         } else {
421:           /* Partial reorthogonalization */
422:           if (force_reorth) force_reorth = PETSC_FALSE;
423:           else {
424:             force_reorth = PETSC_TRUE;
425:             compute_int(which2,omega,j-k,delta,eta);
426:             for (i=0;i<j-k;i++)
427:               if (which2[i]) omega[i] = eps1;
428:           }
429:           IPOrthogonalize(eps->ip,0,PETSC_NULL,j-k,which2,V+k,f,hwork,&norm,breakdown);
430:         }
431:       }
432:     }
433: 
434:     if (*breakdown || norm < eps->n*anorm*PETSC_MACHINE_EPSILON) {
435:       *M = j+1;
436:       break;
437:     }
438:     if (!fro && norm*delta < anorm*eps1) {
439:       fro = PETSC_TRUE;
440:       PetscInfo1(eps,"Switching to full reorthogonalization at iteration %i\n",eps->its);
441:     }
442: 
443:     beta[j-k] = norm;
444:     if (j<m-1) {
445:       VecScale(f,1.0/norm);
446:       VecCopy(f,V[j+1]);
447:     }
448:   }

450:   if (m>100) {
451:     PetscFree(omega);
452:     PetscFree(omega_old);
453:   }
454:   if (m > 100) {
455:     PetscFree(which);
456:     PetscFree(which2);
457:     PetscFree(hwork);
458:   }
459:   return(0);
460: }

464: /*
465:    EPSBasicLanczos - Computes an m-step Lanczos factorization. The first k
466:    columns are assumed to be locked and therefore they are not modified. On
467:    exit, the following relation is satisfied:

469:                     OP * V - V * T = f * e_m^T

471:    where the columns of V are the Lanczos vectors, T is a tridiagonal matrix, 
472:    f is the residual vector and e_m is the m-th vector of the canonical basis. 
473:    The Lanczos vectors (together with vector f) are B-orthogonal (to working
474:    accuracy) if full reorthogonalization is being used, otherwise they are
475:    (B-)semi-orthogonal. On exit, beta contains the B-norm of f and the next 
476:    Lanczos vector can be computed as v_{m+1} = f / beta. 

478:    This function simply calls another function which depends on the selected
479:    reorthogonalization strategy.
480: */
481: static PetscErrorCode EPSBasicLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *m,Vec f,PetscTruth *breakdown,PetscReal anorm)
482: {
483:   EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
485:   IPOrthogonalizationRefinementType orthog_ref;
486:   PetscScalar *T;
487:   PetscInt i,n=*m;
488:   PetscReal betam;

491:   switch (lanczos->reorthog) {
492:     case EPS_LANCZOS_REORTHOG_LOCAL:
493:       EPSLocalLanczos(eps,alpha,beta,V,k,m,f,breakdown);
494:       break;
495:     case EPS_LANCZOS_REORTHOG_SELECTIVE:
496:       EPSSelectiveLanczos(eps,alpha,beta,V,k,m,f,breakdown,anorm);
497:       break;
498:     case EPS_LANCZOS_REORTHOG_FULL:
499:       EPSFullLanczos(eps,alpha,beta,V,k,m,f,breakdown);
500:       break;
501:     case EPS_LANCZOS_REORTHOG_PARTIAL:
502:     case EPS_LANCZOS_REORTHOG_PERIODIC:
503:       EPSPartialLanczos(eps,alpha,beta,V,k,m,f,breakdown,anorm);
504:       break;
505:     case EPS_LANCZOS_REORTHOG_DELAYED:
506:       PetscMalloc(n*n*sizeof(PetscScalar),&T);
507:       IPGetOrthogonalization(eps->ip,PETSC_NULL,&orthog_ref,PETSC_NULL);
508:       if (orthog_ref == IP_ORTH_REFINE_NEVER) {
509:         EPSDelayedArnoldi1(eps,T,n,V,k,m,f,&betam,breakdown);
510:       } else {
511:         EPSDelayedArnoldi(eps,T,n,V,k,m,f,&betam,breakdown);
512:       }
513:       for (i=k;i<n-1;i++) { alpha[i-k] = PetscRealPart(T[n*i+i]); beta[i-k] = PetscRealPart(T[n*i+i+1]); }
514:       alpha[n-1] = PetscRealPart(T[n*(n-1)+n-1]);
515:       beta[n-1] = betam;
516:       PetscFree(T);
517:       break;
518:     default:
519:       SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
520:   }
521:   return(0);
522: }

526: PetscErrorCode EPSSolve_LANCZOS(EPS eps)
527: {
528:   EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
530:   PetscInt       nconv,i,j,k,l,x,n,m,*perm,restart,ncv=eps->ncv,r;
531:   Vec            w=eps->work[1],f=eps->work[0];
532:   PetscScalar    *Y,stmp;
533:   PetscReal      *d,*e,*ritz,*bnd,anorm,beta,norm,rtmp,resnorm;
534:   PetscTruth     breakdown;
535:   char           *conv,ctmp;

538:   PetscMalloc(ncv*sizeof(PetscReal),&d);
539:   PetscMalloc(ncv*sizeof(PetscReal),&e);
540:   PetscMalloc(ncv*sizeof(PetscReal),&ritz);
541:   PetscMalloc(ncv*ncv*sizeof(PetscScalar),&Y);
542:   PetscMalloc(ncv*sizeof(PetscReal),&bnd);
543:   PetscMalloc(ncv*sizeof(PetscInt),&perm);
544:   PetscMalloc(ncv*sizeof(char),&conv);

546:   /* The first Lanczos vector is the normalized initial vector */
547:   EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
548: 
549:   anorm = -1.0;
550:   nconv = 0;
551: 
552:   /* Restart loop */
553:   while (eps->reason == EPS_CONVERGED_ITERATING) {
554:     eps->its++;
555:     /* Compute an ncv-step Lanczos factorization */
556:     m = PetscMin(nconv+eps->mpd,ncv);
557:     EPSBasicLanczos(eps,d,e,eps->V,nconv,&m,f,&breakdown,anorm);

559:     /* Compute eigenvalues and eigenvectors Y of the tridiagonal block */
560:     n = m - nconv;
561:     beta = e[n-1];
562:     EPSDenseTridiagonal(n,d,e,ritz,Y);
563: 
564:     /* Estimate ||A|| */
565:     for (i=0;i<n;i++)
566:       if (PetscAbsReal(ritz[i]) > anorm) anorm = PetscAbsReal(ritz[i]);
567: 
568:     /* Compute residual norm estimates as beta*abs(Y(m,:)) + eps*||A|| */
569:     for (i=0;i<n;i++) {
570:       resnorm = beta*PetscAbsScalar(Y[i*n+n-1]) + PETSC_MACHINE_EPSILON*anorm;
571:       (*eps->conv_func)(eps,ritz[i],eps->eigi[i],resnorm,&bnd[i],eps->conv_ctx);
572:       if (bnd[i]<eps->tol) {
573:         conv[i] = 'C';
574:       } else {
575:         conv[i] = 'N';
576:       }
577:     }

579:     /* purge repeated ritz values */
580:     if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL)
581:       for (i=1;i<n;i++)
582:         if (conv[i] == 'C')
583:           if (PetscAbsScalar((ritz[i]-ritz[i-1])/ritz[i]) < eps->tol)
584:             conv[i] = 'R';

586:     /* Compute restart vector */
587:     if (breakdown) {
588:       PetscInfo2(eps,"Breakdown in Lanczos method (it=%i norm=%g)\n",eps->its,beta);
589:     } else {
590:       restart = 0;
591:       while (restart<n && conv[restart] != 'N') restart++;
592:       if (restart >= n) {
593:         breakdown = PETSC_TRUE;
594:       } else {
595:         for (i=restart+1;i<n;i++)
596:           if (conv[i] == 'N') {
597:             EPSCompareEigenvalues(eps,ritz[restart],0.0,ritz[i],0.0,&r);
598:             if (r>0) restart = i;
599:           }
600:         SlepcVecMAXPBY(f,0.0,1.0,n,Y+restart*n,eps->V+nconv);
601:       }
602:     }

604:     /* Count and put converged eigenvalues first */
605:     for (i=0;i<n;i++) perm[i] = i;
606:     for (k=0;k<n;k++)
607:       if (conv[perm[k]] != 'C') {
608:         j = k + 1;
609:         while (j<n && conv[perm[j]] != 'C') j++;
610:         if (j>=n) break;
611:         l = perm[k]; perm[k] = perm[j]; perm[j] = l;
612:       }

614:     /* Sort eigenvectors according to permutation */
615:     for (i=0;i<k;i++) {
616:       x = perm[i];
617:       if (x != i) {
618:         j = i + 1;
619:         while (perm[j] != i) j++;
620:         /* swap eigenvalues i and j */
621:         rtmp = ritz[x]; ritz[x] = ritz[i]; ritz[i] = rtmp;
622:         rtmp = bnd[x]; bnd[x] = bnd[i]; bnd[i] = rtmp;
623:         ctmp = conv[x]; conv[x] = conv[i]; conv[i] = ctmp;
624:         perm[j] = x; perm[i] = i;
625:         /* swap eigenvectors i and j */
626:         for (l=0;l<n;l++) {
627:           stmp = Y[x*n+l]; Y[x*n+l] = Y[i*n+l]; Y[i*n+l] = stmp;
628:         }
629:       }
630:     }
631: 
632:     /* compute converged eigenvectors */
633:     SlepcUpdateVectors(n,eps->V+nconv,0,k,Y,n,PETSC_FALSE);
634: 
635:     /* purge spurious ritz values */
636:     if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
637:       for (i=0;i<k;i++) {
638:         VecNorm(eps->V[nconv+i],NORM_2,&norm);
639:         VecScale(eps->V[nconv+i],1.0/norm);
640:         STApply(eps->OP,eps->V[nconv+i],w);
641:         VecAXPY(w,-ritz[i],eps->V[nconv+i]);
642:         VecNorm(w,NORM_2,&norm);
643:         (*eps->conv_func)(eps,ritz[i],eps->eigi[i],norm,&bnd[i],eps->conv_ctx);
644:         if (bnd[i]>=eps->tol) conv[i] = 'S';
645:       }
646:       for (i=0;i<k;i++)
647:         if (conv[i] != 'C') {
648:           j = i + 1;
649:           while (j<k && conv[j] != 'C') j++;
650:           if (j>=k) break;
651:           /* swap eigenvalues i and j */
652:           rtmp = ritz[j]; ritz[j] = ritz[i]; ritz[i] = rtmp;
653:           rtmp = bnd[j]; bnd[j] = bnd[i]; bnd[i] = rtmp;
654:           ctmp = conv[j]; conv[j] = conv[i]; conv[i] = ctmp;
655:           /* swap eigenvectors i and j */
656:           VecSwap(eps->V[nconv+i],eps->V[nconv+j]);
657:         }
658:       k = i;
659:     }
660: 
661:     /* store ritz values and estimated errors */
662:     for (i=0;i<n;i++) {
663:       eps->eigr[nconv+i] = ritz[i];
664:       eps->errest[nconv+i] = bnd[i];
665:     }
666:     EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nconv+n);
667:     nconv = nconv+k;
668:     if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
669:     if (nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
670: 
671:      if (eps->reason == EPS_CONVERGED_ITERATING) { /* copy restart vector */
672:       if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL && !breakdown) {
673:         /* Reorthonormalize restart vector */
674:         IPOrthogonalize(eps->ip,eps->nds,eps->DS,nconv,PETSC_NULL,eps->V,f,PETSC_NULL,&norm,&breakdown);
675:         VecScale(f,1.0/norm);
676:       }
677:       if (breakdown) {
678:         /* Use random vector for restarting */
679:         PetscInfo(eps,"Using random vector for restart\n");
680:         EPSGetStartVector(eps,nconv,f,&breakdown);
681:       }
682:       if (breakdown) { /* give up */
683:         eps->reason = EPS_DIVERGED_BREAKDOWN;
684:         PetscInfo(eps,"Unable to generate more start vectors\n");
685:       } else {
686:         VecCopy(f,eps->V[nconv]);
687:       }
688:     }
689:   }
690: 
691:   eps->nconv = nconv;

693:   PetscFree(d);
694:   PetscFree(e);
695:   PetscFree(ritz);
696:   PetscFree(Y);
697:   PetscFree(bnd);
698:   PetscFree(perm);
699:   PetscFree(conv);
700:   return(0);
701: }

703: static const char *lanczoslist[6] = { "local", "full", "selective", "periodic", "partial" , "delayed" };

707: PetscErrorCode EPSSetFromOptions_LANCZOS(EPS eps)
708: {
710:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS *)eps->data;
711:   PetscTruth     flg;
712:   PetscInt       i;

715:   PetscOptionsBegin(((PetscObject)eps)->comm,((PetscObject)eps)->prefix,"LANCZOS Options","EPS");
716:   PetscOptionsEList("-eps_lanczos_reorthog","Lanczos reorthogonalization","EPSLanczosSetReorthog",lanczoslist,6,lanczoslist[lanczos->reorthog],&i,&flg);
717:   if (flg) lanczos->reorthog = (EPSLanczosReorthogType)i;
718:   PetscOptionsEnd();
719:   return(0);
720: }

725: PetscErrorCode EPSLanczosSetReorthog_LANCZOS(EPS eps,EPSLanczosReorthogType reorthog)
726: {
727:   EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;

730:   switch (reorthog) {
731:     case EPS_LANCZOS_REORTHOG_LOCAL:
732:     case EPS_LANCZOS_REORTHOG_FULL:
733:     case EPS_LANCZOS_REORTHOG_DELAYED:
734:     case EPS_LANCZOS_REORTHOG_SELECTIVE:
735:     case EPS_LANCZOS_REORTHOG_PERIODIC:
736:     case EPS_LANCZOS_REORTHOG_PARTIAL:
737:       lanczos->reorthog = reorthog;
738:       break;
739:     default:
740:       SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
741:   }
742:   return(0);
743: }

748: /*@
749:    EPSLanczosSetReorthog - Sets the type of reorthogonalization used during the Lanczos
750:    iteration. 

752:    Collective on EPS

754:    Input Parameters:
755: +  eps - the eigenproblem solver context
756: -  reorthog - the type of reorthogonalization

758:    Options Database Key:
759: .  -eps_lanczos_reorthog - Sets the reorthogonalization type (either 'local', 'selective',
760:                          'periodic', 'partial', 'full' or 'delayed')
761:    
762:    Level: advanced

764: .seealso: EPSLanczosGetReorthog(), EPSLanczosReorthogType
765: @*/
766: PetscErrorCode EPSLanczosSetReorthog(EPS eps,EPSLanczosReorthogType reorthog)
767: {
768:   PetscErrorCode ierr, (*f)(EPS,EPSLanczosReorthogType);

772:   PetscObjectQueryFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",(void (**)())&f);
773:   if (f) {
774:     (*f)(eps,reorthog);
775:   }
776:   return(0);
777: }

782: PetscErrorCode EPSLanczosGetReorthog_LANCZOS(EPS eps,EPSLanczosReorthogType *reorthog)
783: {
784:   EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
786:   *reorthog = lanczos->reorthog;
787:   return(0);
788: }

793: /*@C
794:    EPSLanczosGetReorthog - Gets the type of reorthogonalization used during the Lanczos
795:    iteration. 

797:    Collective on EPS

799:    Input Parameter:
800: .  eps - the eigenproblem solver context

802:    Input Parameter:
803: .  reorthog - the type of reorthogonalization

805:    Level: advanced

807: .seealso: EPSLanczosSetReorthog(), EPSLanczosReorthogType
808: @*/
809: PetscErrorCode EPSLanczosGetReorthog(EPS eps,EPSLanczosReorthogType *reorthog)
810: {
811:   PetscErrorCode ierr, (*f)(EPS,EPSLanczosReorthogType*);

815:   PetscObjectQueryFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",(void (**)())&f);
816:   if (f) {
817:     (*f)(eps,reorthog);
818:   }
819:   return(0);
820: }

824: PetscErrorCode EPSDestroy_LANCZOS(EPS eps)
825: {
827:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS *)eps->data;

831:   if (lanczos->AV) { VecDestroyVecs(lanczos->AV,eps->ncv); }
832:   EPSDestroy_Default(eps);
833:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosSetReorthog_C","",PETSC_NULL);
834:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosGetReorthog_C","",PETSC_NULL);
835:   return(0);
836: }

840: PetscErrorCode EPSView_LANCZOS(EPS eps,PetscViewer viewer)
841: {
843:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS *)eps->data;
844:   PetscTruth     isascii;

847:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
848:   if (!isascii) {
849:     SETERRQ1(1,"Viewer type %s not supported for EPSLANCZOS",((PetscObject)viewer)->type_name);
850:   }
851:   PetscViewerASCIIPrintf(viewer,"reorthogonalization: %s\n",lanczoslist[lanczos->reorthog]);
852:   return(0);
853: }

858: PetscErrorCode EPSCreate_LANCZOS(EPS eps)
859: {
861:   EPS_LANCZOS    *lanczos;

864:   PetscNew(EPS_LANCZOS,&lanczos);
865:   PetscLogObjectMemory(eps,sizeof(EPS_LANCZOS));
866:   eps->data                      = (void *) lanczos;
867:   eps->ops->setup                = EPSSetUp_LANCZOS;
868:   eps->ops->setfromoptions       = EPSSetFromOptions_LANCZOS;
869:   eps->ops->destroy              = EPSDestroy_LANCZOS;
870:   eps->ops->view                 = EPSView_LANCZOS;
871:   eps->ops->backtransform        = EPSBackTransform_Default;
872:   eps->ops->computevectors       = EPSComputeVectors_Hermitian;
873:   lanczos->reorthog              = EPS_LANCZOS_REORTHOG_LOCAL;
874:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosSetReorthog_C","EPSLanczosSetReorthog_LANCZOS",EPSLanczosSetReorthog_LANCZOS);
875:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosGetReorthog_C","EPSLanczosGetReorthog_LANCZOS",EPSLanczosGetReorthog_LANCZOS);
876:   return(0);
877: }