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: Oct 2006

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

 24:       This file is part of SLEPc. See the README file for conditions of use
 25:       and additional information.
 26:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 27: */

 29:  #include src/eps/epsimpl.h
 30:  #include slepcblaslapack.h

 32: typedef struct {
 33:   EPSLanczosReorthogType reorthog;
 34: } EPS_LANCZOS;

 38: PetscErrorCode EPSSetUp_LANCZOS(EPS eps)
 39: {
 41:   PetscInt       N;

 44:   VecGetSize(eps->vec_initial,&N);
 45:   if (eps->ncv) {
 46:     if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
 47:   }
 48:   else eps->ncv = PetscMin(N,PetscMax(2*eps->nev,eps->nev+15));
 49:   if (!eps->max_it) eps->max_it = PetscMax(100,2*N/eps->ncv);

 51:   if (eps->solverclass==EPS_ONE_SIDE) {
 52:     if (eps->which == EPS_LARGEST_IMAGINARY || eps->which == EPS_SMALLEST_IMAGINARY)
 53:       SETERRQ(1,"Wrong value of eps->which");
 54:     if (!eps->ishermitian)
 55:       SETERRQ(PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
 56:   } else {
 57:     if (eps->which != EPS_LARGEST_MAGNITUDE)
 58:       SETERRQ(1,"Wrong value of eps->which");
 59:   }
 60:   EPSAllocateSolution(eps);
 61:   PetscFree(eps->T);
 62:   PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);
 63:   if (eps->solverclass==EPS_TWO_SIDE) {
 64:     PetscFree(eps->Tl);
 65:     PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->Tl);
 66:   }
 67:   EPSDefaultGetWork(eps,2);
 68:   return(0);
 69: }

 73: /*
 74:    EPSLocalLanczos - Local reorthogonalization.

 76:    This is the simplest variant. At each Lanczos step, the corresponding Lanczos vector 
 77:    is orthogonalized with respect to the two previous Lanczos vectors, according to
 78:    the three term Lanczos recurrence. WARNING: This variant does not track the loss of 
 79:    orthogonality that occurs in finite-precision arithmetic and, therefore, the 
 80:    generated vectors are not guaranteed to be (semi-)orthogonal.
 81: */
 82: static PetscErrorCode EPSLocalLanczos(EPS eps,PetscScalar *T,Vec *V,int k,int *M,Vec f,PetscReal *beta,PetscTruth *breakdown)
 83: {
 85:   int            i,j,m = *M;
 86:   PetscReal      norm;
 87:   PetscTruth     *which,lwhich[100];
 88: 
 90:   if (m>100) {
 91:     PetscMalloc(sizeof(PetscTruth)*m,&which);
 92:   } else which = lwhich;
 93:   for (i=0;i<k;i++)
 94:     which[i] = PETSC_TRUE;

 96:   for (j=k;j<m;j++) {
 97:     STApply(eps->OP,V[j],f);
 98:     IPOrthogonalize(eps->ip,eps->nds,PETSC_NULL,eps->DS,f,PETSC_NULL,PETSC_NULL,PETSC_NULL,eps->work[0]);
 99:     which[j] = PETSC_TRUE;
100:     if (j-2>=k) which[j-2] = PETSC_FALSE;
101:     IPOrthogonalize(eps->ip,j+1,which,V,f,T+m*j,&norm,breakdown,eps->work[0]);
102:     if (*breakdown) {
103:       *M = j+1;
104:       break;
105:     }
106:     if (j<m-1) {
107:       T[m*j+j+1] = norm;
108:       VecScale(f,1.0/norm);
109:       VecCopy(f,V[j+1]);
110:     }
111:   }
112:   *beta = norm;

114:   if (m>100) { PetscFree(which); }
115:   return(0);
116: }

120: /*
121:    EPSSelectiveLanczos - Selective reorthogonalization.
122: */
123: static PetscErrorCode EPSSelectiveLanczos(EPS eps,PetscScalar *T,Vec *V,int k,int *M,Vec f,PetscReal *beta,PetscTruth *breakdown,PetscReal anorm)
124: {
126:   int            i,j,m = *M,n,nritz=0,nritzo;
127:   PetscReal      *ritz,norm;
128:   PetscScalar    *Y;
129:   PetscTruth     *which,lwhich[100];

132:   PetscMalloc(m*sizeof(PetscReal),&ritz);
133:   PetscMalloc(m*m*sizeof(PetscScalar),&Y);

135:   if (m>100) {
136:     PetscMalloc(sizeof(PetscTruth)*m,&which);
137:   } else which = lwhich;
138:   for (i=0;i<k;i++)
139:     which[i] = PETSC_TRUE;

141:   for (j=k;j<m;j++) {
142:     /* Lanczos step */
143:     STApply(eps->OP,V[j],f);
144:     IPOrthogonalize(eps->ip,eps->nds,PETSC_NULL,eps->DS,f,PETSC_NULL,PETSC_NULL,PETSC_NULL,eps->work[0]);
145:     which[j] = PETSC_TRUE;
146:     if (j-2>=k) which[j-2] = PETSC_FALSE;
147:     IPOrthogonalize(eps->ip,j+1,which,V,f,T+m*j,&norm,breakdown,eps->work[0]);
148:     if (*breakdown) {
149:       *M = j+1;
150:       break;
151:     }

153:     /* Compute eigenvalues and eigenvectors Y of the tridiagonal block */
154:     n = j-k+1;
155:     EPSDenseTridiagonal(n,T+k*(m+1),m,ritz,Y);
156: 
157:     /* Estimate ||A|| */
158:     for (i=0;i<n;i++)
159:       if (PetscAbsReal(ritz[i]) > anorm) anorm = PetscAbsReal(ritz[i]);

161:     /* Compute nearly converged Ritz vectors */
162:     nritzo = 0;
163:     for (i=0;i<n;i++)
164:       if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm)
165:         nritzo++;

167:     if (nritzo>nritz) {
168:       nritz = 0;
169:       for (i=0;i<n;i++) {
170:         if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm) {
171:           VecSet(eps->AV[nritz],0.0);
172:           VecMAXPY(eps->AV[nritz],n,Y+i*n,V+k);
173:           nritz++;
174:         }
175:       }
176:     }

178:     if (nritz > 0) {
179:       IPOrthogonalize(eps->ip,nritz,PETSC_NULL,eps->AV,f,PETSC_NULL,&norm,breakdown,eps->work[0]);
180:       if (*breakdown) {
181:         *M = j+1;
182:         break;
183:       }
184:     }
185: 
186:     if (j<m-1) {
187:       T[m*j+j+1] = norm;
188:       VecScale(f,1.0 / norm);
189:       VecCopy(f,V[j+1]);
190:     }
191:   }
192:   *beta = norm;
193: 
194:   PetscFree(ritz);
195:   PetscFree(Y);
196:   if (m>100) { PetscFree(which); }
197:   return(0);
198: }

202: static void update_omega(PetscReal *omega,PetscReal *omega_old,int j,PetscReal *alpha,PetscReal *beta,PetscReal eps1,PetscReal anorm)
203: {
204:   int            k;
205:   PetscReal      T,binv,temp;

208:   /* Estimate of contribution to roundoff errors from A*v 
209:        fl(A*v) = A*v + f, 
210:      where ||f|| \approx eps1*||A||.
211:      For a full matrix A, a rule-of-thumb estimate is eps1 = sqrt(n)*eps. */
212:   T = eps1*anorm;
213:   binv = 1.0/beta[j+1];

215:   /* Update omega(1) using omega(0)==0. */
216:   omega_old[0]= beta[1]*omega[1] + (alpha[0]-alpha[j])*omega[0] -
217:                 beta[j]*omega_old[0];
218:   if (omega_old[0] > 0)
219:     omega_old[0] = binv*(omega_old[0] + T);
220:   else
221:     omega_old[0] = binv*(omega_old[0] - T);
222: 
223:   /* Update remaining components. */
224:   for (k=1;k<j-1;k++) {
225:     omega_old[k] = beta[k+1]*omega[k+1] + (alpha[k]-alpha[j])*omega[k] +
226:                    beta[k]*omega[k-1] - beta[j]*omega_old[k];
227:     if (omega_old[k] > 0)
228:       omega_old[k] = binv*(omega_old[k] + T);
229:     else
230:       omega_old[k] = binv*(omega_old[k] - T);
231:   }
232:   omega_old[j-1] = binv*T;
233: 
234:   /* Swap omega and omega_old. */
235:   for (k=0;k<j;k++) {
236:     temp = omega[k];
237:     omega[k] = omega_old[k];
238:     omega_old[k] = omega[k];
239:   }
240:   omega[j] = eps1;
241:   PetscFunctionReturnVoid();
242: }

246: static void compute_int(PetscTruth *which,PetscReal *mu,int j,PetscReal delta,PetscReal eta)
247: {
248:   int        i,k,maxpos;
249:   PetscReal  max;
250:   PetscTruth found;
251: 
253:   /* initialize which */
254:   found = PETSC_FALSE;
255:   maxpos = 0;
256:   max = 0.0;
257:   for (i=0;i<j;i++) {
258:     if (PetscAbsReal(mu[i]) >= delta) {
259:       which[i] = PETSC_TRUE;
260:       found = PETSC_TRUE;
261:     } else which[i] = PETSC_FALSE;
262:     if (PetscAbsReal(mu[i]) > max) {
263:       maxpos = i;
264:       max = PetscAbsReal(mu[i]);
265:     }
266:   }
267:   if (!found) which[maxpos] = PETSC_TRUE;
268: 
269:   for (i=0;i<j;i++)
270:     if (which[i]) {
271:       /* find left interval */
272:       for (k=i;k>=0;k--) {
273:         if (PetscAbsReal(mu[k])<eta || which[k]) break;
274:         else which[k] = PETSC_TRUE;
275:       }
276:       /* find right interval */
277:       for (k=i+1;k<j;k++) {
278:         if (PetscAbsReal(mu[k])<eta || which[k]) break;
279:         else which[k] = PETSC_TRUE;
280:       }
281:     }
282:   PetscFunctionReturnVoid();
283: }

287: /*
288:    EPSPartialLanczos - Partial reorthogonalization.
289: */
290: static PetscErrorCode EPSPartialLanczos(EPS eps,PetscScalar *T,Vec *V,int k,int *M,Vec f,PetscReal *beta, PetscTruth *breakdown,PetscReal anorm)
291: {
292:   EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
294:   Mat            A;
295:   int            i,j,m = *M;
296:   PetscInt       n;
297:   PetscReal      norm,*omega,lomega[100],*omega_old,lomega_old[100],eps1,delta,eta,*b,lb[101],*a,la[100];
298:   PetscTruth     *which,lwhich[100],*which2,lwhich2[100],
299:                  reorth = PETSC_FALSE,force_reorth = PETSC_FALSE,fro = PETSC_FALSE,estimate_anorm = PETSC_FALSE;

302:   if (m>100) {
303:     PetscMalloc(m*sizeof(PetscReal),&a);
304:     PetscMalloc((m+1)*sizeof(PetscReal),&b);
305:     PetscMalloc(m*sizeof(PetscReal),&omega);
306:     PetscMalloc(m*sizeof(PetscReal),&omega_old);
307:     PetscMalloc(sizeof(PetscTruth)*m,&which);
308:     PetscMalloc(sizeof(PetscTruth)*m,&which2);
309:   } else {
310:     a = la;
311:     b = lb;
312:     omega = lomega;
313:     omega_old = lomega_old;
314:     which = lwhich;
315:     which2 = lwhich2;
316:   }
317: 
318:   STGetOperators(eps->OP,&A,PETSC_NULL);
319:   MatGetSize(A,&n,PETSC_NULL);
320:   eps1 = sqrt((PetscReal)n)*PETSC_MACHINE_EPSILON/2;
321:   delta = PETSC_SQRT_MACHINE_EPSILON/sqrt((PetscReal)eps->ncv);
322:   eta = pow(PETSC_MACHINE_EPSILON,3.0/4.0)/sqrt((PetscReal)eps->ncv);
323:   if (anorm < 0.0) {
324:     anorm = 1.0;
325:     estimate_anorm = PETSC_TRUE;
326:   }
327:   for (i=0;i<m-k;i++)
328:     omega[i] = omega_old[i] = 0.0;
329:   for (i=0;i<k;i++)
330:     which[i] = PETSC_TRUE;
331: 
332:   for (j=k;j<m;j++) {
333:     STApply(eps->OP,V[j],f);
334:     IPOrthogonalize(eps->ip,eps->nds,PETSC_NULL,eps->DS,f,PETSC_NULL,PETSC_NULL,PETSC_NULL,eps->work[0]);
335:     if (fro) {
336:       /* Lanczos step with full reorthogonalization */
337:       IPOrthogonalize(eps->ip,j+1,PETSC_NULL,V,f,T+m*j,&norm,breakdown,eps->work[0]);
338:     } else {
339:       /* Lanczos step */
340:       which[j] = PETSC_TRUE;
341:       if (j-2>=k) which[j-2] = PETSC_FALSE;
342:       IPOrthogonalize(eps->ip,j+1,which,V,f,T+m*j,&norm,breakdown,eps->work[0]);
343:       a[j-k] = PetscRealPart(T[m*j+j]);
344:       b[j-k+1] = norm;
345: 
346:       /* Estimate ||A|| if needed */
347:       if (estimate_anorm) {
348:         if (j>k) anorm = PetscMax(anorm,PetscAbsReal(a[j-k])+norm+b[j-k]);
349:         else anorm = PetscMax(anorm,PetscAbsReal(a[j-k])+norm);
350:       }

352:       /* Check if reorthogonalization is needed */
353:       reorth = PETSC_FALSE;
354:       if (j>k) {
355:         update_omega(omega,omega_old,j-k,a,b,eps1,anorm);
356:         for (i=0;i<j-k;i++)
357:           if (PetscAbsScalar(omega[i]) > delta) reorth = PETSC_TRUE;
358:       }

360:       if (reorth || force_reorth) {
361:         if (lanczos->reorthog == EPSLANCZOS_REORTHOG_PERIODIC) {
362:           /* Periodic reorthogonalization */
363:           if (force_reorth) force_reorth = PETSC_FALSE;
364:           else force_reorth = PETSC_TRUE;
365:           IPOrthogonalize(eps->ip,j-k,PETSC_NULL,V+k,f,PETSC_NULL,&norm,breakdown,eps->work[0]);
366:           for (i=0;i<j-k;i++)
367:             omega[i] = eps1;
368:         } else {
369:           /* Partial reorthogonalization */
370:           if (force_reorth) force_reorth = PETSC_FALSE;
371:           else {
372:             force_reorth = PETSC_TRUE;
373:             compute_int(which2,omega,j-k,delta,eta);
374:             for (i=0;i<j-k;i++)
375:               if (which2[i]) omega[i] = eps1;
376:           }
377:           IPOrthogonalize(eps->ip,j-k,which2,V+k,f,PETSC_NULL,&norm,breakdown,eps->work[0]);
378:         }
379:       }
380:     }
381: 
382:     if (*breakdown || norm < n*anorm*PETSC_MACHINE_EPSILON) {
383:       *M = j+1;
384:       break;
385:     }
386:     if (!fro && norm*delta < anorm*eps1) {
387:       fro = PETSC_TRUE;
388:       PetscInfo1(eps,"Switching to full reorthogonalization at iteration %i\n",eps->its);
389:     }
390:     if (j<m-1) {
391:       T[m*j+j+1] = b[j-k+1] = norm;
392:       VecScale(f,1.0/norm);
393:       VecCopy(f,V[j+1]);
394:     }
395:   }
396:   *beta = norm;

398:   if (m>100) {
399:     PetscFree(a);
400:     PetscFree(b);
401:     PetscFree(omega);
402:     PetscFree(omega_old);
403:     PetscFree(which);
404:     PetscFree(which2);
405:   }
406:   return(0);
407: }

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

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

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

425:    This function simply calls another function which depends on the selected
426:    reorthogonalization strategy.
427: */
428: static PetscErrorCode EPSBasicLanczos(EPS eps,PetscScalar *T,Vec *V,int k,int *m,Vec f,PetscReal *beta,PetscTruth *breakdown,PetscReal anorm)
429: {
430:   EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
432:   IPOrthogonalizationRefinementType orthog_ref;

435:   switch (lanczos->reorthog) {
436:     case EPSLANCZOS_REORTHOG_LOCAL:
437:       EPSLocalLanczos(eps,T,V,k,m,f,beta,breakdown);
438:       break;
439:     case EPSLANCZOS_REORTHOG_SELECTIVE:
440:       EPSSelectiveLanczos(eps,T,V,k,m,f,beta,breakdown,anorm);
441:       break;
442:     case EPSLANCZOS_REORTHOG_PARTIAL:
443:     case EPSLANCZOS_REORTHOG_PERIODIC:
444:       EPSPartialLanczos(eps,T,V,k,m,f,beta,breakdown,anorm);
445:       break;
446:     case EPSLANCZOS_REORTHOG_FULL:
447:       EPSBasicArnoldi(eps,PETSC_FALSE,T,V,k,m,f,beta,breakdown);
448:       break;
449:     case EPSLANCZOS_REORTHOG_DELAYED:
450:       IPGetOrthogonalization(eps->ip,PETSC_NULL,&orthog_ref,PETSC_NULL);
451:       if (orthog_ref == IP_ORTH_REFINE_NEVER) {
452:         EPSDelayedArnoldi1(eps,T,V,k,m,f,beta,breakdown);
453:       } else {
454:         EPSDelayedArnoldi(eps,T,V,k,m,f,beta,breakdown);
455:       }
456:       break;
457:     default:
458:       SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
459:   }
460:   return(0);
461: }

465: PetscErrorCode EPSSolve_LANCZOS(EPS eps)
466: {
467:   EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
469:   int            nconv,i,j,k,n,m,*perm,restart,ncv=eps->ncv;
470:   Vec            f=eps->work[1];
471:   PetscScalar    *T=eps->T,*Y;
472:   PetscReal      *ritz,*bnd,anorm,beta,norm;
473:   PetscTruth     breakdown;
474:   char           *conv;

477:   PetscMalloc(ncv*sizeof(PetscReal),&ritz);
478:   PetscMalloc(ncv*ncv*sizeof(PetscScalar),&Y);
479:   PetscMalloc(ncv*sizeof(PetscReal),&bnd);
480:   PetscMalloc(ncv*sizeof(int),&perm);
481:   PetscMalloc(ncv*sizeof(char),&conv);

483:   /* The first Lanczos vector is the normalized initial vector */
484:   EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
485: 
486:   anorm = -1.0;
487:   nconv = 0;
488: 
489:   /* Restart loop */
490:   while (eps->reason == EPS_CONVERGED_ITERATING) {
491:     eps->its++;
492:     /* Compute an ncv-step Lanczos factorization */
493:     m = ncv;
494:     EPSBasicLanczos(eps,T,eps->V,nconv,&m,f,&beta,&breakdown,anorm);

496:     /* Compute eigenvalues and eigenvectors Y of the tridiagonal block */
497:     n = m - nconv;
498:     EPSDenseTridiagonal(n,T+nconv*(ncv+1),ncv,ritz,Y);

500:     /* Estimate ||A|| */
501:     for (i=0;i<n;i++)
502:       if (PetscAbsReal(ritz[i]) > anorm) anorm = PetscAbsReal(ritz[i]);
503: 
504:     /* Compute residual norm estimates as beta*abs(Y(m,:)) + eps*||A|| */
505:     for (i=0;i<n;i++)
506:       bnd[i] = beta*PetscAbsScalar(Y[i*n+n-1]) + PETSC_MACHINE_EPSILON*anorm;

508:     /* Sort eigenvalues according to eps->which */
509:     if (eps->which == EPS_SMALLEST_REAL) {
510:       /* LAPACK function has already ordered the eigenvalues and eigenvectors */
511:       for (i=0;i<n;i++)
512:         perm[i] = i;
513:     } else {
514: #ifdef PETSC_USE_COMPLEX
515:       for (i=0;i<n;i++)
516:         eps->eigr[i+nconv] = ritz[i];
517:       EPSSortEigenvalues(n,eps->eigr+nconv,eps->eigi,eps->which,n,perm);
518: #else
519:       EPSSortEigenvalues(n,ritz,eps->eigi,eps->which,n,perm);
520: #endif
521:     }

523:     /* Look for converged eigenpairs */
524:     k = nconv;
525:     for (i=0;i<n;i++) {
526:       eps->eigr[k] = ritz[perm[i]];
527:       eps->errest[k] = bnd[perm[i]] / PetscAbsScalar(eps->eigr[k]);
528:       if (eps->errest[k] < eps->tol) {
529: 
530:         if (lanczos->reorthog == EPSLANCZOS_REORTHOG_LOCAL) {
531:           if (i>0 && PetscAbsScalar((eps->eigr[k]-ritz[perm[i-1]])/eps->eigr[k]) < eps->tol) {
532:               /* Discard repeated eigenvalues */
533:             conv[i] = 'R';
534:             continue;
535:            }
536:         }
537: 
538:         VecSet(eps->AV[k],0.0);
539:         VecMAXPY(eps->AV[k],n,Y+perm[i]*n,eps->V+nconv);

541:         if (lanczos->reorthog == EPSLANCZOS_REORTHOG_LOCAL) {
542:           /* normalize locked vector and compute residual norm */
543:           VecNorm(eps->AV[k],NORM_2,&norm);
544:           VecScale(eps->AV[k],1.0/norm);
545:           STApply(eps->OP,eps->AV[k],f);
546:           VecAXPY(f,-eps->eigr[k],eps->AV[k]);
547:           VecNorm(f,NORM_2,&norm);
548:           eps->errest[k] = norm / PetscAbsScalar(eps->eigr[k]);
549:           if (eps->errest[k] >= eps->tol) {
550:             conv[i] = 'S';
551:             continue;
552:           }
553:         }
554: 
555:         conv[i] = 'C';
556:         k++;
557:       } else conv[i] = 'N';
558:     }

560:     /* Look for non-converged eigenpairs */
561:     j = k;
562:     restart = -1;
563:     for (i=0;i<n;i++) {
564:       if (conv[i] != 'C') {
565:         if (restart == -1 && conv[i] == 'N') restart = i;
566:         eps->eigr[j] = ritz[perm[i]];
567:         eps->errest[j] = bnd[perm[i]] / ritz[perm[i]];
568:         j++;
569:       }
570:     }

572:     if (breakdown) {
573:       restart = -1;
574:       PetscInfo2(eps,"Breakdown in Lanczos method (it=%i norm=%g)\n",eps->its,beta);
575:     }
576: 
577:     if (k<eps->nev) {
578:       if (restart != -1) {
579:         /* Use first non-converged vector for restarting */
580:         VecSet(eps->AV[k],0.0);
581:         VecMAXPY(eps->AV[k],n,Y+perm[restart]*n,eps->V+nconv);
582:         VecCopy(eps->AV[k],eps->V[k]);
583:       }
584:     }
585: 
586:     /* Copy converged vectors to V */
587:     for (i=nconv;i<k;i++) {
588:       VecCopy(eps->AV[i],eps->V[i]);
589:     }

591:     if (k<eps->nev) {
592:       if (restart == -1) {
593:         /* Use random vector for restarting */
594:         PetscInfo(eps,"Using random vector for restart\n");
595:         EPSGetStartVector(eps,k,eps->V[k],&breakdown);
596:       } else if (lanczos->reorthog == EPSLANCZOS_REORTHOG_LOCAL) {
597:         /* Reorthonormalize restart vector */
598:         IPOrthogonalize(eps->ip,eps->nds+k,PETSC_NULL,eps->DSV,eps->V[k],PETSC_NULL,&norm,&breakdown,eps->work[0]);
599:         VecScale(eps->V[k],1.0/norm);
600:       } else breakdown = PETSC_FALSE;
601:       if (breakdown) {
602:         eps->reason = EPS_DIVERGED_BREAKDOWN;
603:         PetscInfo(eps,"Unable to generate more start vectors\n");
604:       }
605:     }

607:     EPSMonitor(eps,eps->its,k,eps->eigr,eps->eigi,eps->errest,nconv+n);
608:     nconv = k;
609:     if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
610:     if (nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
611:   }
612: 
613:   eps->nconv = nconv;

615:   PetscFree(ritz);
616:   PetscFree(Y);
617:   PetscFree(bnd);
618:   PetscFree(perm);
619:   PetscFree(conv);
620:   return(0);
621: }

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

627: PetscErrorCode EPSSetFromOptions_LANCZOS(EPS eps)
628: {
630:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS *)eps->data;
631:   PetscTruth     flg;
632:   PetscInt       i;

635:   PetscOptionsHead("LANCZOS options");
636:   PetscOptionsEList("-eps_lanczos_reorthog","Lanczos reorthogonalization","EPSLanczosSetReorthog",lanczoslist,6,lanczoslist[lanczos->reorthog],&i,&flg);
637:   if (flg) lanczos->reorthog = (EPSLanczosReorthogType)i;
638:   PetscOptionsTail();
639:   return(0);
640: }

645: PetscErrorCode EPSLanczosSetReorthog_LANCZOS(EPS eps,EPSLanczosReorthogType reorthog)
646: {
647:   EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;

650:   switch (reorthog) {
651:     case EPSLANCZOS_REORTHOG_LOCAL:
652:     case EPSLANCZOS_REORTHOG_FULL:
653:     case EPSLANCZOS_REORTHOG_DELAYED:
654:     case EPSLANCZOS_REORTHOG_SELECTIVE:
655:     case EPSLANCZOS_REORTHOG_PERIODIC:
656:     case EPSLANCZOS_REORTHOG_PARTIAL:
657:       lanczos->reorthog = reorthog;
658:       break;
659:     default:
660:       SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
661:   }
662:   return(0);
663: }

668: /*@
669:    EPSLanczosSetReorthog - Sets the type of reorthogonalization used during the Lanczos
670:    iteration. 

672:    Collective on EPS

674:    Input Parameters:
675: +  eps - the eigenproblem solver context
676: -  reorthog - the type of reorthogonalization

678:    Options Database Key:
679: .  -eps_lanczos_orthog - Sets the reorthogonalization type (either 'local', 'selective',
680:                          'periodic', 'partial' or 'full')
681:    
682:    Level: advanced

684: .seealso: EPSLanczosGetReorthog(), EPSLanczosReorthogType
685: @*/
686: PetscErrorCode EPSLanczosSetReorthog(EPS eps,EPSLanczosReorthogType reorthog)
687: {
688:   PetscErrorCode ierr, (*f)(EPS,EPSLanczosReorthogType);

692:   PetscObjectQueryFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",(void (**)())&f);
693:   if (f) {
694:     (*f)(eps,reorthog);
695:   }
696:   return(0);
697: }

702: PetscErrorCode EPSLanczosGetReorthog_LANCZOS(EPS eps,EPSLanczosReorthogType *reorthog)
703: {
704:   EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
706:   *reorthog = lanczos->reorthog;
707:   return(0);
708: }

713: /*@C
714:    EPSLanczosGetReorthog - Gets the type of reorthogonalization used during the Lanczos
715:    iteration. 

717:    Collective on EPS

719:    Input Parameter:
720: .  eps - the eigenproblem solver context

722:    Input Parameter:
723: .  reorthog - the type of reorthogonalization

725:    Level: advanced

727: .seealso: EPSLanczosSetReorthog(), EPSLanczosReorthogType
728: @*/
729: PetscErrorCode EPSLanczosGetReorthog(EPS eps,EPSLanczosReorthogType *reorthog)
730: {
731:   PetscErrorCode ierr, (*f)(EPS,EPSLanczosReorthogType*);

735:   PetscObjectQueryFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",(void (**)())&f);
736:   if (f) {
737:     (*f)(eps,reorthog);
738:   }
739:   return(0);
740: }

744: PetscErrorCode EPSView_LANCZOS(EPS eps,PetscViewer viewer)
745: {
747:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS *)eps->data;
748:   PetscTruth     isascii;

751:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
752:   if (!isascii) {
753:     SETERRQ1(1,"Viewer type %s not supported for EPSLANCZOS",((PetscObject)viewer)->type_name);
754:   }
755:   PetscViewerASCIIPrintf(viewer,"reorthogonalization: %s\n",lanczoslist[lanczos->reorthog]);
756:   return(0);
757: }

759: /*
760: EXTERN PetscErrorCode EPSSolve_TS_LANCZOS(EPS);
761: */

766: PetscErrorCode EPSCreate_LANCZOS(EPS eps)
767: {
769:   EPS_LANCZOS    *lanczos;

772:   PetscNew(EPS_LANCZOS,&lanczos);
773:   PetscLogObjectMemory(eps,sizeof(EPS_LANCZOS));
774:   eps->data                      = (void *) lanczos;
775:   eps->ops->solve                = EPSSolve_LANCZOS;
776: /*  eps->ops->solvets              = EPSSolve_TS_LANCZOS;*/
777:   eps->ops->setup                = EPSSetUp_LANCZOS;
778:   eps->ops->setfromoptions       = EPSSetFromOptions_LANCZOS;
779:   eps->ops->destroy              = EPSDestroy_Default;
780:   eps->ops->view                 = EPSView_LANCZOS;
781:   eps->ops->backtransform        = EPSBackTransform_Default;
782:   /*if (eps->solverclass==EPS_TWO_SIDE)
783:        eps->ops->computevectors       = EPSComputeVectors_Schur;
784:   else*/ eps->ops->computevectors       = EPSComputeVectors_Hermitian;
785:   lanczos->reorthog              = EPSLANCZOS_REORTHOG_LOCAL;
786:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosSetReorthog_C","EPSLanczosSetReorthog_LANCZOS",EPSLanczosSetReorthog_LANCZOS);
787:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosGetReorthog_C","EPSLanczosGetReorthog_LANCZOS",EPSLanczosGetReorthog_LANCZOS);
788:   return(0);
789: }