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_MAGNITUDE:
 76:     case EPS_TARGET_REAL:
 77:     case EPS_TARGET_IMAGINARY:
 78:       SETERRQ(1,"Wrong value of eps->which");
 79:     default: ; /* default case to remove warning */
 80:   }
 81:   if (!eps->ishermitian)
 82:     SETERRQ(PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
 83:   if (!eps->extraction) {
 84:     EPSSetExtraction(eps,EPS_RITZ);
 85:   } else if (eps->extraction!=EPS_RITZ) {
 86:     SETERRQ(PETSC_ERR_SUP,"Unsupported extraction type\n");
 87:   }

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

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

107: /*
108:    EPSLocalLanczos - Local reorthogonalization.

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

588:     /* Compute restart vector */
589:     if (breakdown) {
590:       PetscInfo2(eps,"Breakdown in Lanczos method (it=%i norm=%g)\n",eps->its,beta);
591:     } else {
592:       restart = 0;
593:       while (restart<n && conv[restart] != 'N') restart++;
594:       if (restart >= n) {
595:         breakdown = PETSC_TRUE;
596:       } else {
597:         switch (eps->which) {
598:         case EPS_LARGEST_MAGNITUDE:
599:         case EPS_SMALLEST_MAGNITUDE:
600:           rtmp = PetscAbsReal(ritz[restart]);
601:           break;
602:         case EPS_LARGEST_REAL:
603:         case EPS_SMALLEST_REAL:
604:           rtmp = ritz[restart];
605:           break;
606:         default: SETERRQ(1,"Wrong value of which");
607:         }
608:         for (i=restart+1;i<n;i++)
609:           if (conv[i] == 'N') {
610:             switch (eps->which) {
611:             case EPS_LARGEST_MAGNITUDE:
612:               if (rtmp < PetscAbsReal(ritz[i])) { rtmp = PetscAbsReal(ritz[i]); restart = i; }
613:               break;
614:             case EPS_SMALLEST_MAGNITUDE:
615:               if (rtmp > PetscAbsReal(ritz[i])) { rtmp = PetscAbsReal(ritz[i]); restart = i; }
616:               break;
617:             case EPS_LARGEST_REAL:
618:               if (rtmp < ritz[i]) { rtmp = ritz[i]; restart = i; }
619:               break;
620:             case EPS_SMALLEST_REAL:
621:               if (rtmp > ritz[i]) { rtmp = ritz[i]; restart = i; }
622:               break;
623:             default: SETERRQ(1,"Wrong value of which");
624:             }
625:           }
626:         SlepcVecMAXPBY(f,0.0,1.0,n,Y+restart*n,eps->V+nconv);
627:       }
628:     }

630:     /* Count and put converged eigenvalues first */
631:     for (i=0;i<n;i++) perm[i] = i;
632:     for (k=0;k<n;k++)
633:       if (conv[perm[k]] != 'C') {
634:         j = k + 1;
635:         while (j<n && conv[perm[j]] != 'C') j++;
636:         if (j>=n) break;
637:         l = perm[k]; perm[k] = perm[j]; perm[j] = l;
638:       }

640:     /* Sort eigenvectors according to permutation */
641:     for (i=0;i<k;i++) {
642:       x = perm[i];
643:       if (x != i) {
644:         j = i + 1;
645:         while (perm[j] != i) j++;
646:         /* swap eigenvalues i and j */
647:         rtmp = ritz[x]; ritz[x] = ritz[i]; ritz[i] = rtmp;
648:         rtmp = bnd[x]; bnd[x] = bnd[i]; bnd[i] = rtmp;
649:         ctmp = conv[x]; conv[x] = conv[i]; conv[i] = ctmp;
650:         perm[j] = x; perm[i] = i;
651:         /* swap eigenvectors i and j */
652:         for (l=0;l<n;l++) {
653:           stmp = Y[x*n+l]; Y[x*n+l] = Y[i*n+l]; Y[i*n+l] = stmp;
654:         }
655:       }
656:     }
657: 
658:     /* compute converged eigenvectors */
659:     SlepcUpdateVectors(n,eps->V+nconv,0,k,Y,n,PETSC_FALSE);
660: 
661:     /* purge spurious ritz values */
662:     if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
663:       for (i=0;i<k;i++) {
664:         VecNorm(eps->V[nconv+i],NORM_2,&norm);
665:         VecScale(eps->V[nconv+i],1.0/norm);
666:         STApply(eps->OP,eps->V[nconv+i],w);
667:         VecAXPY(w,-ritz[i],eps->V[nconv+i]);
668:         VecNorm(w,NORM_2,&norm);
669:         (*eps->conv_func)(eps,ritz[i],eps->eigi[i],norm,&bnd[i],eps->conv_ctx);
670:         if (bnd[i]>=eps->tol) conv[i] = 'S';
671:       }
672:       for (i=0;i<k;i++)
673:         if (conv[i] != 'C') {
674:           j = i + 1;
675:           while (j<k && conv[j] != 'C') j++;
676:           if (j>=k) break;
677:           /* swap eigenvalues i and j */
678:           rtmp = ritz[j]; ritz[j] = ritz[i]; ritz[i] = rtmp;
679:           rtmp = bnd[j]; bnd[j] = bnd[i]; bnd[i] = rtmp;
680:           ctmp = conv[j]; conv[j] = conv[i]; conv[i] = ctmp;
681:           /* swap eigenvectors i and j */
682:           VecSwap(eps->V[nconv+i],eps->V[nconv+j]);
683:         }
684:       k = i;
685:     }
686: 
687:     /* store ritz values and estimated errors */
688:     for (i=0;i<n;i++) {
689:       eps->eigr[nconv+i] = ritz[i];
690:       eps->errest[nconv+i] = bnd[i];
691:     }
692:     EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nconv+n);
693:     nconv = nconv+k;
694:     if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
695:     if (nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
696: 
697:      if (eps->reason == EPS_CONVERGED_ITERATING) { /* copy restart vector */
698:       if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL && !breakdown) {
699:         /* Reorthonormalize restart vector */
700:         IPOrthogonalize(eps->ip,eps->nds,eps->DS,nconv,PETSC_NULL,eps->V,f,PETSC_NULL,&norm,&breakdown);
701:         VecScale(f,1.0/norm);
702:       }
703:       if (breakdown) {
704:         /* Use random vector for restarting */
705:         PetscInfo(eps,"Using random vector for restart\n");
706:         EPSGetStartVector(eps,nconv,f,&breakdown);
707:       }
708:       if (breakdown) { /* give up */
709:         eps->reason = EPS_DIVERGED_BREAKDOWN;
710:         PetscInfo(eps,"Unable to generate more start vectors\n");
711:       } else {
712:         VecCopy(f,eps->V[nconv]);
713:       }
714:     }
715:   }
716: 
717:   eps->nconv = nconv;

719:   PetscFree(d);
720:   PetscFree(e);
721:   PetscFree(ritz);
722:   PetscFree(Y);
723:   PetscFree(bnd);
724:   PetscFree(perm);
725:   PetscFree(conv);
726:   return(0);
727: }

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

733: PetscErrorCode EPSSetFromOptions_LANCZOS(EPS eps)
734: {
736:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS *)eps->data;
737:   PetscTruth     flg;
738:   PetscInt       i;

741:   PetscOptionsBegin(((PetscObject)eps)->comm,((PetscObject)eps)->prefix,"LANCZOS Options","EPS");
742:   PetscOptionsEList("-eps_lanczos_reorthog","Lanczos reorthogonalization","EPSLanczosSetReorthog",lanczoslist,6,lanczoslist[lanczos->reorthog],&i,&flg);
743:   if (flg) lanczos->reorthog = (EPSLanczosReorthogType)i;
744:   PetscOptionsEnd();
745:   return(0);
746: }

751: PetscErrorCode EPSLanczosSetReorthog_LANCZOS(EPS eps,EPSLanczosReorthogType reorthog)
752: {
753:   EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;

756:   switch (reorthog) {
757:     case EPS_LANCZOS_REORTHOG_LOCAL:
758:     case EPS_LANCZOS_REORTHOG_FULL:
759:     case EPS_LANCZOS_REORTHOG_DELAYED:
760:     case EPS_LANCZOS_REORTHOG_SELECTIVE:
761:     case EPS_LANCZOS_REORTHOG_PERIODIC:
762:     case EPS_LANCZOS_REORTHOG_PARTIAL:
763:       lanczos->reorthog = reorthog;
764:       break;
765:     default:
766:       SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
767:   }
768:   return(0);
769: }

774: /*@
775:    EPSLanczosSetReorthog - Sets the type of reorthogonalization used during the Lanczos
776:    iteration. 

778:    Collective on EPS

780:    Input Parameters:
781: +  eps - the eigenproblem solver context
782: -  reorthog - the type of reorthogonalization

784:    Options Database Key:
785: .  -eps_lanczos_reorthog - Sets the reorthogonalization type (either 'local', 'selective',
786:                          'periodic', 'partial', 'full' or 'delayed')
787:    
788:    Level: advanced

790: .seealso: EPSLanczosGetReorthog(), EPSLanczosReorthogType
791: @*/
792: PetscErrorCode EPSLanczosSetReorthog(EPS eps,EPSLanczosReorthogType reorthog)
793: {
794:   PetscErrorCode ierr, (*f)(EPS,EPSLanczosReorthogType);

798:   PetscObjectQueryFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",(void (**)())&f);
799:   if (f) {
800:     (*f)(eps,reorthog);
801:   }
802:   return(0);
803: }

808: PetscErrorCode EPSLanczosGetReorthog_LANCZOS(EPS eps,EPSLanczosReorthogType *reorthog)
809: {
810:   EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
812:   *reorthog = lanczos->reorthog;
813:   return(0);
814: }

819: /*@C
820:    EPSLanczosGetReorthog - Gets the type of reorthogonalization used during the Lanczos
821:    iteration. 

823:    Collective on EPS

825:    Input Parameter:
826: .  eps - the eigenproblem solver context

828:    Input Parameter:
829: .  reorthog - the type of reorthogonalization

831:    Level: advanced

833: .seealso: EPSLanczosSetReorthog(), EPSLanczosReorthogType
834: @*/
835: PetscErrorCode EPSLanczosGetReorthog(EPS eps,EPSLanczosReorthogType *reorthog)
836: {
837:   PetscErrorCode ierr, (*f)(EPS,EPSLanczosReorthogType*);

841:   PetscObjectQueryFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",(void (**)())&f);
842:   if (f) {
843:     (*f)(eps,reorthog);
844:   }
845:   return(0);
846: }

850: PetscErrorCode EPSDestroy_LANCZOS(EPS eps)
851: {
853:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS *)eps->data;

857:   if (lanczos->AV) { VecDestroyVecs(lanczos->AV,eps->ncv); }
858:   EPSDestroy_Default(eps);
859:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosSetReorthog_C","",PETSC_NULL);
860:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosGetReorthog_C","",PETSC_NULL);
861:   return(0);
862: }

866: PetscErrorCode EPSView_LANCZOS(EPS eps,PetscViewer viewer)
867: {
869:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS *)eps->data;
870:   PetscTruth     isascii;

873:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
874:   if (!isascii) {
875:     SETERRQ1(1,"Viewer type %s not supported for EPSLANCZOS",((PetscObject)viewer)->type_name);
876:   }
877:   PetscViewerASCIIPrintf(viewer,"reorthogonalization: %s\n",lanczoslist[lanczos->reorthog]);
878:   return(0);
879: }

884: PetscErrorCode EPSCreate_LANCZOS(EPS eps)
885: {
887:   EPS_LANCZOS    *lanczos;

890:   PetscNew(EPS_LANCZOS,&lanczos);
891:   PetscLogObjectMemory(eps,sizeof(EPS_LANCZOS));
892:   eps->data                      = (void *) lanczos;
893:   eps->ops->setup                = EPSSetUp_LANCZOS;
894:   eps->ops->setfromoptions       = EPSSetFromOptions_LANCZOS;
895:   eps->ops->destroy              = EPSDestroy_LANCZOS;
896:   eps->ops->view                 = EPSView_LANCZOS;
897:   eps->ops->backtransform        = EPSBackTransform_Default;
898:   eps->ops->computevectors       = EPSComputeVectors_Hermitian;
899:   lanczos->reorthog              = EPS_LANCZOS_REORTHOG_LOCAL;
900:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosSetReorthog_C","EPSLanczosSetReorthog_LANCZOS",EPSLanczosSetReorthog_LANCZOS);
901:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosGetReorthog_C","EPSLanczosGetReorthog_LANCZOS",EPSLanczosGetReorthog_LANCZOS);
902:   return(0);
903: }