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-2013, Universitat Politecnica de Valencia, Spain

 24:    This file is part of SLEPc.

 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 <slepc-private/epsimpl.h>                /*I "slepceps.h" I*/
 41: #include <slepcblaslapack.h>

 43: PetscErrorCode EPSSolve_Lanczos(EPS);

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

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

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

 73:   if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
 74:   switch (eps->which) {
 75:     case EPS_LARGEST_IMAGINARY:
 76:     case EPS_SMALLEST_IMAGINARY:
 77:     case EPS_TARGET_IMAGINARY:
 78:       SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
 79:     default: ; /* default case to remove warning */
 80:   }
 81:   if (!eps->extraction) {
 82:     EPSSetExtraction(eps,EPS_RITZ);
 83:   } else if (eps->extraction!=EPS_RITZ) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
 84:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");

 86:   EPSAllocateSolution(eps);
 87:   if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_SELECTIVE) {
 88:     VecDuplicateVecs(eps->t,eps->ncv,&lanczos->AV);
 89:     PetscLogObjectParents(eps,eps->ncv,lanczos->AV);
 90:   }
 91:   DSSetType(eps->ds,DSHEP);
 92:   DSSetCompact(eps->ds,PETSC_TRUE);
 93:   DSAllocate(eps->ds,eps->ncv);
 94:   if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
 95:     EPSSetWorkVecs(eps,2);
 96:   } else {
 97:     EPSSetWorkVecs(eps,1);
 98:   }

100:   /* dispatch solve method */
101:   if (eps->leftvecs) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Left vectors not supported in this solver");
102:   if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
103:   if (eps->isgeneralized && eps->ishermitian && !eps->ispositive) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method does not work for indefinite problems");
104:   eps->ops->solve = EPSSolve_Lanczos;
105:   return(0);
106: }

110: /*
111:    EPSLocalLanczos - Local reorthogonalization.

113:    This is the simplest variant. At each Lanczos step, the corresponding Lanczos vector
114:    is orthogonalized with respect to the two previous Lanczos vectors, according to
115:    the three term Lanczos recurrence. WARNING: This variant does not track the loss of
116:    orthogonality that occurs in finite-precision arithmetic and, therefore, the
117:    generated vectors are not guaranteed to be (semi-)orthogonal.
118: */
119: static PetscErrorCode EPSLocalLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscBool *breakdown)
120: {
122:   PetscInt       i,j,m = *M;
123:   PetscReal      norm;
124:   PetscBool      *which,lwhich[100];
125:   PetscScalar    *hwork,lhwork[100];

128:   if (m > 100) {
129:     PetscMalloc(sizeof(PetscBool)*m,&which);
130:     PetscMalloc(m*sizeof(PetscScalar),&hwork);
131:   } else {
132:     which = lwhich;
133:     hwork = lhwork;
134:   }
135:   for (i=0;i<k;i++)
136:     which[i] = PETSC_TRUE;

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

161:   if (m > 100) {
162:     PetscFree(which);
163:     PetscFree(hwork);
164:   }
165:   return(0);
166: }

170: /*
171:    DenseTridiagonal - Solves a real tridiagonal Hermitian Eigenvalue Problem.

173:    Input Parameters:
174: +  n   - dimension of the eigenproblem
175: .  D   - pointer to the array containing the diagonal elements
176: -  E   - pointer to the array containing the off-diagonal elements

178:    Output Parameters:
179: +  w  - pointer to the array to store the computed eigenvalues
180: -  V  - pointer to the array to store the eigenvectors

182:    Notes:
183:    If V is NULL then the eigenvectors are not computed.

185:    This routine use LAPACK routines xSTEVR.

187: */
188: static PetscErrorCode DenseTridiagonal(PetscInt n_,PetscReal *D,PetscReal *E,PetscReal *w,PetscScalar *V)
189: {
190: #if defined(SLEPC_MISSING_LAPACK_STEVR)
192:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"STEVR - Lapack routine is unavailable");
193: #else
195:   PetscReal      abstol = 0.0,vl,vu,*work;
196:   PetscBLASInt   il,iu,m,*isuppz,n,lwork,*iwork,liwork,info;
197:   const char     *jobz;
198: #if defined(PETSC_USE_COMPLEX)
199:   PetscInt       i,j;
200:   PetscReal      *VV;
201: #endif

204:   PetscBLASIntCast(n_,&n);
205:   PetscBLASIntCast(20*n_,&lwork);
206:   PetscBLASIntCast(10*n_,&liwork);
207:   if (V) {
208:     jobz = "V";
209: #if defined(PETSC_USE_COMPLEX)
210:     PetscMalloc(n*n*sizeof(PetscReal),&VV);
211: #endif
212:   } else jobz = "N";
213:   PetscMalloc(2*n*sizeof(PetscBLASInt),&isuppz);
214:   PetscMalloc(lwork*sizeof(PetscReal),&work);
215:   PetscMalloc(liwork*sizeof(PetscBLASInt),&iwork);
216:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
217: #if defined(PETSC_USE_COMPLEX)
218:   PetscStackCallBLAS("LAPACKstevr",LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,VV,&n,isuppz,work,&lwork,iwork,&liwork,&info));
219: #else
220:   PetscStackCallBLAS("LAPACKstevr",LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info));
221: #endif
222:   PetscFPTrapPop();
223:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack DSTEVR %d",info);
224: #if defined(PETSC_USE_COMPLEX)
225:   if (V) {
226:     for (i=0;i<n;i++)
227:       for (j=0;j<n;j++)
228:         V[i*n+j] = VV[i*n+j];
229:     PetscFree(VV);
230:   }
231: #endif
232:   PetscFree(isuppz);
233:   PetscFree(work);
234:   PetscFree(iwork);
235:   return(0);
236: #endif
237: }

241: /*
242:    EPSSelectiveLanczos - Selective reorthogonalization.
243: */
244: static PetscErrorCode EPSSelectiveLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscBool *breakdown,PetscReal anorm)
245: {
247:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS*)eps->data;
248:   PetscInt       i,j,m = *M,n,nritz=0,nritzo;
249:   PetscReal      *d,*e,*ritz,norm;
250:   PetscScalar    *Y,*hwork,lhwork[100];
251:   PetscBool      *which,lwhich[100];

254:   PetscMalloc(m*sizeof(PetscReal),&d);
255:   PetscMalloc(m*sizeof(PetscReal),&e);
256:   PetscMalloc(m*sizeof(PetscReal),&ritz);
257:   PetscMalloc(m*m*sizeof(PetscScalar),&Y);
258:   if (m > 100) {
259:     PetscMalloc(sizeof(PetscBool)*m,&which);
260:     PetscMalloc(m*sizeof(PetscScalar),&hwork);
261:   } else {
262:     which = lwhich;
263:     hwork = lhwork;
264:   }
265:   for (i=0;i<k;i++)
266:     which[i] = PETSC_TRUE;

268:   for (j=k;j<m;j++) {
269:     /* Lanczos step */
270:     STApply(eps->st,V[j],f);
271:     which[j] = PETSC_TRUE;
272:     if (j-2>=k) which[j-2] = PETSC_FALSE;
273:     IPOrthogonalize(eps->ip,eps->nds,eps->defl,j+1,which,V,f,hwork,&norm,breakdown);
274:     alpha[j] = PetscRealPart(hwork[j]);
275:     beta[j] = norm;
276:     if (*breakdown) {
277:       *M = j+1;
278:       break;
279:     }

281:     /* Compute eigenvalues and eigenvectors Y of the tridiagonal block */
282:     n = j-k+1;
283:     for (i=0;i<n;i++) {
284:       d[i] = alpha[i+k];
285:       e[i] = beta[i+k];
286:     }
287:     DenseTridiagonal(n,d,e,ritz,Y);

289:     /* Estimate ||A|| */
290:     for (i=0;i<n;i++)
291:       if (PetscAbsReal(ritz[i]) > anorm) anorm = PetscAbsReal(ritz[i]);

293:     /* Compute nearly converged Ritz vectors */
294:     nritzo = 0;
295:     for (i=0;i<n;i++)
296:       if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm)
297:         nritzo++;

299:     if (nritzo>nritz) {
300:       nritz = 0;
301:       for (i=0;i<n;i++) {
302:         if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm) {
303:           SlepcVecMAXPBY(lanczos->AV[nritz],0.0,1.0,n,Y+i*n,V+k);
304:           nritz++;
305:         }
306:       }
307:     }

309:     if (nritz > 0) {
310:       IPOrthogonalize(eps->ip,0,NULL,nritz,NULL,lanczos->AV,f,hwork,&norm,breakdown);
311:       if (*breakdown) {
312:         *M = j+1;
313:         break;
314:       }
315:     }

317:     if (j<m-1) {
318:       VecScale(f,1.0 / norm);
319:       VecCopy(f,V[j+1]);
320:     }
321:   }

323:   PetscFree(d);
324:   PetscFree(e);
325:   PetscFree(ritz);
326:   PetscFree(Y);
327:   if (m > 100) {
328:     PetscFree(which);
329:     PetscFree(hwork);
330:   }
331:   return(0);
332: }

336: static void update_omega(PetscReal *omega,PetscReal *omega_old,PetscInt j,PetscReal *alpha,PetscReal *beta,PetscReal eps1,PetscReal anorm)
337: {
338:   PetscInt       k;
339:   PetscReal      T,binv;

342:   /* Estimate of contribution to roundoff errors from A*v
343:        fl(A*v) = A*v + f,
344:      where ||f|| \approx eps1*||A||.
345:      For a full matrix A, a rule-of-thumb estimate is eps1 = sqrt(n)*eps. */
346:   T = eps1*anorm;
347:   binv = 1.0/beta[j+1];

349:   /* Update omega(1) using omega(0)==0. */
350:   omega_old[0]= beta[1]*omega[1] + (alpha[0]-alpha[j])*omega[0] -
351:                 beta[j]*omega_old[0];
352:   if (omega_old[0] > 0) omega_old[0] = binv*(omega_old[0] + T);
353:   else omega_old[0] = binv*(omega_old[0] - T);

355:   /* Update remaining components. */
356:   for (k=1;k<j-1;k++) {
357:     omega_old[k] = beta[k+1]*omega[k+1] + (alpha[k]-alpha[j])*omega[k] + beta[k]*omega[k-1] - beta[j]*omega_old[k];
358:     if (omega_old[k] > 0) omega_old[k] = binv*(omega_old[k] + T);
359:     else omega_old[k] = binv*(omega_old[k] - T);
360:   }
361:   omega_old[j-1] = binv*T;

363:   /* Swap omega and omega_old. */
364:   for (k=0;k<j;k++) {
365:     omega[k] = omega_old[k];
366:     omega_old[k] = omega[k];
367:   }
368:   omega[j] = eps1;
369:   PetscFunctionReturnVoid();
370: }

374: static void compute_int(PetscBool *which,PetscReal *mu,PetscInt j,PetscReal delta,PetscReal eta)
375: {
376:   PetscInt  i,k,maxpos;
377:   PetscReal max;
378:   PetscBool found;

381:   /* initialize which */
382:   found = PETSC_FALSE;
383:   maxpos = 0;
384:   max = 0.0;
385:   for (i=0;i<j;i++) {
386:     if (PetscAbsReal(mu[i]) >= delta) {
387:       which[i] = PETSC_TRUE;
388:       found = PETSC_TRUE;
389:     } else which[i] = PETSC_FALSE;
390:     if (PetscAbsReal(mu[i]) > max) {
391:       maxpos = i;
392:       max = PetscAbsReal(mu[i]);
393:     }
394:   }
395:   if (!found) which[maxpos] = PETSC_TRUE;

397:   for (i=0;i<j;i++)
398:     if (which[i]) {
399:       /* find left interval */
400:       for (k=i;k>=0;k--) {
401:         if (PetscAbsReal(mu[k])<eta || which[k]) break;
402:         else which[k] = PETSC_TRUE;
403:       }
404:       /* find right interval */
405:       for (k=i+1;k<j;k++) {
406:         if (PetscAbsReal(mu[k])<eta || which[k]) break;
407:         else which[k] = PETSC_TRUE;
408:       }
409:     }
410:   PetscFunctionReturnVoid();
411: }

415: /*
416:    EPSPartialLanczos - Partial reorthogonalization.
417: */
418: static PetscErrorCode EPSPartialLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscBool *breakdown,PetscReal anorm)
419: {
420:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS*)eps->data;
422:   PetscInt       i,j,m = *M;
423:   PetscReal      norm,*omega,lomega[100],*omega_old,lomega_old[100],eps1,delta,eta;
424:   PetscBool      *which,lwhich[100],*which2,lwhich2[100];
425:   PetscBool      reorth = PETSC_FALSE,force_reorth = PETSC_FALSE;
426:   PetscBool      fro = PETSC_FALSE,estimate_anorm = PETSC_FALSE;
427:   PetscScalar    *hwork,lhwork[100];

430:   if (m>100) {
431:     PetscMalloc(m*sizeof(PetscReal),&omega);
432:     PetscMalloc(m*sizeof(PetscReal),&omega_old);
433:   } else {
434:     omega = lomega;
435:     omega_old = lomega_old;
436:   }
437:   if (m > 100) {
438:     PetscMalloc(sizeof(PetscBool)*m,&which);
439:     PetscMalloc(sizeof(PetscBool)*m,&which2);
440:     PetscMalloc(m*sizeof(PetscScalar),&hwork);
441:   } else {
442:     which = lwhich;
443:     which2 = lwhich2;
444:     hwork = lhwork;
445:   }

447:   eps1 = PetscSqrtReal((PetscReal)eps->n)*PETSC_MACHINE_EPSILON/2;
448:   delta = PETSC_SQRT_MACHINE_EPSILON/PetscSqrtReal((PetscReal)eps->ncv);
449:   eta = PetscPowReal(PETSC_MACHINE_EPSILON,3.0/4.0)/PetscSqrtReal((PetscReal)eps->ncv);
450:   if (anorm < 0.0) {
451:     anorm = 1.0;
452:     estimate_anorm = PETSC_TRUE;
453:   }
454:   for (i=0;i<m-k;i++)
455:     omega[i] = omega_old[i] = 0.0;
456:   for (i=0;i<k;i++)
457:     which[i] = PETSC_TRUE;

459:   for (j=k;j<m;j++) {
460:     STApply(eps->st,V[j],f);
461:     if (fro) {
462:       /* Lanczos step with full reorthogonalization */
463:       IPOrthogonalize(eps->ip,eps->nds,eps->defl,j+1,NULL,V,f,hwork,&norm,breakdown);
464:       alpha[j] = PetscRealPart(hwork[j]);
465:     } else {
466:       /* Lanczos step */
467:       which[j] = PETSC_TRUE;
468:       if (j-2>=k) which[j-2] = PETSC_FALSE;
469:       IPOrthogonalize(eps->ip,eps->nds,eps->defl,j+1,which,V,f,hwork,&norm,breakdown);
470:       alpha[j] = PetscRealPart(hwork[j]);
471:       beta[j] = norm;

473:       /* Estimate ||A|| if needed */
474:       if (estimate_anorm) {
475:         if (j>k) anorm = PetscMax(anorm,PetscAbsReal(alpha[j])+norm+beta[j-1]);
476:         else anorm = PetscMax(anorm,PetscAbsReal(alpha[j])+norm);
477:       }

479:       /* Check if reorthogonalization is needed */
480:       reorth = PETSC_FALSE;
481:       if (j>k) {
482:         update_omega(omega,omega_old,j,alpha,beta-1,eps1,anorm);
483:         for (i=0;i<j-k;i++)
484:           if (PetscAbsScalar(omega[i]) > delta) reorth = PETSC_TRUE;
485:       }

487:       if (reorth || force_reorth) {
488:         if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_PERIODIC) {
489:           /* Periodic reorthogonalization */
490:           if (force_reorth) force_reorth = PETSC_FALSE;
491:           else force_reorth = PETSC_TRUE;
492:           IPOrthogonalize(eps->ip,0,NULL,j-k,NULL,V+k,f,hwork,&norm,breakdown);
493:           for (i=0;i<j-k;i++)
494:             omega[i] = eps1;
495:         } else {
496:           /* Partial reorthogonalization */
497:           if (force_reorth) force_reorth = PETSC_FALSE;
498:           else {
499:             force_reorth = PETSC_TRUE;
500:             compute_int(which2,omega,j-k,delta,eta);
501:             for (i=0;i<j-k;i++)
502:               if (which2[i]) omega[i] = eps1;
503:           }
504:           IPOrthogonalize(eps->ip,0,NULL,j-k,which2,V+k,f,hwork,&norm,breakdown);
505:         }
506:       }
507:     }

509:     if (*breakdown || norm < eps->n*anorm*PETSC_MACHINE_EPSILON) {
510:       *M = j+1;
511:       break;
512:     }
513:     if (!fro && norm*delta < anorm*eps1) {
514:       fro = PETSC_TRUE;
515:       PetscInfo1(eps,"Switching to full reorthogonalization at iteration %D\n",eps->its);
516:     }

518:     beta[j] = norm;
519:     if (j<m-1) {
520:       VecScale(f,1.0/norm);
521:       VecCopy(f,V[j+1]);
522:     }
523:   }

525:   if (m>100) {
526:     PetscFree(omega);
527:     PetscFree(omega_old);
528:     PetscFree(which);
529:     PetscFree(which2);
530:     PetscFree(hwork);
531:   }
532:   return(0);
533: }

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

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

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

551:    This function simply calls another function which depends on the selected
552:    reorthogonalization strategy.
553: */
554: static PetscErrorCode EPSBasicLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *m,Vec f,PetscBool *breakdown,PetscReal anorm)
555: {
556:   EPS_LANCZOS        *lanczos = (EPS_LANCZOS*)eps->data;
557:   PetscScalar        *T;
558:   PetscInt           i,n=*m;
559:   PetscReal          betam;
560:   PetscErrorCode     ierr;
561:   IPOrthogRefineType orthog_ref;

564:   switch (lanczos->reorthog) {
565:     case EPS_LANCZOS_REORTHOG_LOCAL:
566:       EPSLocalLanczos(eps,alpha,beta,V,k,m,f,breakdown);
567:       break;
568:     case EPS_LANCZOS_REORTHOG_SELECTIVE:
569:       EPSSelectiveLanczos(eps,alpha,beta,V,k,m,f,breakdown,anorm);
570:       break;
571:     case EPS_LANCZOS_REORTHOG_FULL:
572:       EPSFullLanczos(eps,alpha,beta,V,k,m,f,breakdown);
573:       break;
574:     case EPS_LANCZOS_REORTHOG_PARTIAL:
575:     case EPS_LANCZOS_REORTHOG_PERIODIC:
576:       EPSPartialLanczos(eps,alpha,beta,V,k,m,f,breakdown,anorm);
577:       break;
578:     case EPS_LANCZOS_REORTHOG_DELAYED:
579:       PetscMalloc(n*n*sizeof(PetscScalar),&T);
580:       IPGetOrthogonalization(eps->ip,NULL,&orthog_ref,NULL);
581:       if (orthog_ref == IP_ORTHOG_REFINE_NEVER) {
582:         EPSDelayedArnoldi1(eps,T,n,V,k,m,f,&betam,breakdown);
583:       } else {
584:         EPSDelayedArnoldi(eps,T,n,V,k,m,f,&betam,breakdown);
585:       }
586:       for (i=k;i<n-1;i++) {
587:         alpha[i] = PetscRealPart(T[n*i+i]);
588:         beta[i] = PetscRealPart(T[n*i+i+1]);
589:       }
590:       alpha[n-1] = PetscRealPart(T[n*(n-1)+n-1]);
591:       beta[n-1] = betam;
592:       PetscFree(T);
593:       break;
594:     default:
595:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
596:   }
597:   return(0);
598: }

602: PetscErrorCode EPSSolve_Lanczos(EPS eps)
603: {
604:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS*)eps->data;
606:   PetscInt       nconv,i,j,k,l,x,n,*perm,restart,ncv=eps->ncv,r,ld;
607:   Vec            w=eps->work[1],f=eps->work[0];
608:   PetscScalar    *Y,*ritz,stmp;
609:   PetscReal      *d,*e,*bnd,anorm,beta,norm,rtmp,resnorm;
610:   PetscBool      breakdown;
611:   char           *conv,ctmp;

614:   DSGetLeadingDimension(eps->ds,&ld);
615:   PetscMalloc(ncv*sizeof(PetscScalar),&ritz);
616:   PetscMalloc(ncv*sizeof(PetscReal),&bnd);
617:   PetscMalloc(ncv*sizeof(PetscInt),&perm);
618:   PetscMalloc(ncv*sizeof(char),&conv);

620:   /* The first Lanczos vector is the normalized initial vector */
621:   EPSGetStartVector(eps,0,eps->V[0],NULL);

623:   anorm = -1.0;
624:   nconv = 0;

626:   /* Restart loop */
627:   while (eps->reason == EPS_CONVERGED_ITERATING) {
628:     eps->its++;

630:     /* Compute an ncv-step Lanczos factorization */
631:     n = PetscMin(nconv+eps->mpd,ncv);
632:     DSGetArrayReal(eps->ds,DS_MAT_T,&d);
633:     e = d + ld;
634:     EPSBasicLanczos(eps,d,e,eps->V,nconv,&n,f,&breakdown,anorm);
635:     beta = e[n-1];
636:     DSRestoreArrayReal(eps->ds,DS_MAT_T,&d);
637:     DSSetDimensions(eps->ds,n,0,nconv,0);
638:     DSSetState(eps->ds,DS_STATE_INTERMEDIATE);

640:     /* Solve projected problem */
641:     DSSolve(eps->ds,ritz,NULL);
642:     DSSort(eps->ds,ritz,NULL,NULL,NULL,NULL);

644:     /* Estimate ||A|| */
645:     for (i=nconv;i<n;i++)
646:       anorm = PetscMax(anorm,PetscAbsReal(PetscRealPart(ritz[i])));

648:     /* Compute residual norm estimates as beta*abs(Y(m,:)) + eps*||A|| */
649:     DSGetArray(eps->ds,DS_MAT_Q,&Y);
650:     for (i=nconv;i<n;i++) {
651:       resnorm = beta*PetscAbsScalar(Y[n-1+i*ld]) + PETSC_MACHINE_EPSILON*anorm;
652:       (*eps->converged)(eps,ritz[i],eps->eigi[i],resnorm,&bnd[i],eps->convergedctx);
653:       if (bnd[i]<eps->tol) conv[i] = 'C';
654:       else conv[i] = 'N';
655:     }
656:     DSRestoreArray(eps->ds,DS_MAT_Q,&Y);

658:     /* purge repeated ritz values */
659:     if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL)
660:       for (i=nconv+1;i<n;i++)
661:         if (conv[i] == 'C')
662:           if (PetscAbsScalar((ritz[i]-ritz[i-1])/ritz[i]) < eps->tol)
663:             conv[i] = 'R';

665:     /* Compute restart vector */
666:     if (breakdown) {
667:       PetscInfo2(eps,"Breakdown in Lanczos method (it=%D norm=%G)\n",eps->its,beta);
668:     } else {
669:       restart = nconv;
670:       while (restart<n && conv[restart] != 'N') restart++;
671:       if (restart >= n) {
672:         breakdown = PETSC_TRUE;
673:       } else {
674:         for (i=restart+1;i<n;i++)
675:           if (conv[i] == 'N') {
676:             (*eps->comparison)(ritz[restart],0.0,ritz[i],0.0,&r,eps->comparisonctx);
677:             if (r>0) restart = i;
678:           }
679:         DSGetArray(eps->ds,DS_MAT_Q,&Y);
680:         SlepcVecMAXPBY(f,0.0,1.0,n-nconv,Y+restart*ld+nconv,eps->V+nconv);
681:         DSRestoreArray(eps->ds,DS_MAT_Q,&Y);
682:       }
683:     }

685:     /* Count and put converged eigenvalues first */
686:     for (i=nconv;i<n;i++) perm[i] = i;
687:     for (k=nconv;k<n;k++)
688:       if (conv[perm[k]] != 'C') {
689:         j = k + 1;
690:         while (j<n && conv[perm[j]] != 'C') j++;
691:         if (j>=n) break;
692:         l = perm[k]; perm[k] = perm[j]; perm[j] = l;
693:       }

695:     /* Sort eigenvectors according to permutation */
696:     DSGetArray(eps->ds,DS_MAT_Q,&Y);
697:     for (i=nconv;i<k;i++) {
698:       x = perm[i];
699:       if (x != i) {
700:         j = i + 1;
701:         while (perm[j] != i) j++;
702:         /* swap eigenvalues i and j */
703:         stmp = ritz[x]; ritz[x] = ritz[i]; ritz[i] = stmp;
704:         rtmp = bnd[x]; bnd[x] = bnd[i]; bnd[i] = rtmp;
705:         ctmp = conv[x]; conv[x] = conv[i]; conv[i] = ctmp;
706:         perm[j] = x; perm[i] = i;
707:         /* swap eigenvectors i and j */
708:         for (l=0;l<n;l++) {
709:           stmp = Y[l+x*ld]; Y[l+x*ld] = Y[l+i*ld]; Y[l+i*ld] = stmp;
710:         }
711:       }
712:     }

714:     /* compute converged eigenvectors */
715:     SlepcUpdateVectors(n,eps->V,nconv,k,Y,ld,PETSC_FALSE);
716:     DSRestoreArray(eps->ds,DS_MAT_Q,&Y);

718:     /* purge spurious ritz values */
719:     if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
720:       for (i=nconv;i<k;i++) {
721:         VecNorm(eps->V[i],NORM_2,&norm);
722:         VecScale(eps->V[i],1.0/norm);
723:         STApply(eps->st,eps->V[i],w);
724:         VecAXPY(w,-ritz[i],eps->V[i]);
725:         VecNorm(w,NORM_2,&norm);
726:         (*eps->converged)(eps,ritz[i],eps->eigi[i],norm,&bnd[i],eps->convergedctx);
727:         if (bnd[i]>=eps->tol) conv[i] = 'S';
728:       }
729:       for (i=nconv;i<k;i++)
730:         if (conv[i] != 'C') {
731:           j = i + 1;
732:           while (j<k && conv[j] != 'C') j++;
733:           if (j>=k) break;
734:           /* swap eigenvalues i and j */
735:           stmp = ritz[j]; ritz[j] = ritz[i]; ritz[i] = stmp;
736:           rtmp = bnd[j]; bnd[j] = bnd[i]; bnd[i] = rtmp;
737:           ctmp = conv[j]; conv[j] = conv[i]; conv[i] = ctmp;
738:           /* swap eigenvectors i and j */
739:           VecSwap(eps->V[i],eps->V[j]);
740:         }
741:       k = i;
742:     }

744:     /* store ritz values and estimated errors */
745:     for (i=nconv;i<n;i++) {
746:       eps->eigr[i] = ritz[i];
747:       eps->errest[i] = bnd[i];
748:     }
749:     EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,n);
750:     nconv = k;
751:     if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
752:     if (nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;

754:     if (eps->reason == EPS_CONVERGED_ITERATING) { /* copy restart vector */
755:       if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL && !breakdown) {
756:         /* Reorthonormalize restart vector */
757:         IPOrthogonalize(eps->ip,eps->nds,eps->defl,nconv,NULL,eps->V,f,NULL,&norm,&breakdown);
758:         VecScale(f,1.0/norm);
759:       }
760:       if (breakdown) {
761:         /* Use random vector for restarting */
762:         PetscInfo(eps,"Using random vector for restart\n");
763:         EPSGetStartVector(eps,nconv,f,&breakdown);
764:       }
765:       if (breakdown) { /* give up */
766:         eps->reason = EPS_DIVERGED_BREAKDOWN;
767:         PetscInfo(eps,"Unable to generate more start vectors\n");
768:       } else {
769:         VecCopy(f,eps->V[nconv]);
770:       }
771:     }
772:   }

774:   eps->nconv = nconv;

776:   PetscFree(ritz);
777:   PetscFree(bnd);
778:   PetscFree(perm);
779:   PetscFree(conv);
780:   return(0);
781: }

785: PetscErrorCode EPSSetFromOptions_Lanczos(EPS eps)
786: {
787:   PetscErrorCode         ierr;
788:   EPS_LANCZOS            *lanczos = (EPS_LANCZOS*)eps->data;
789:   PetscBool              flg;
790:   EPSLanczosReorthogType reorthog;

793:   PetscOptionsHead("EPS Lanczos Options");
794:   PetscOptionsEnum("-eps_lanczos_reorthog","Lanczos reorthogonalization","EPSLanczosSetReorthog",EPSLanczosReorthogTypes,(PetscEnum)lanczos->reorthog,(PetscEnum*)&reorthog,&flg);
795:   if (flg) {
796:     EPSLanczosSetReorthog(eps,reorthog);
797:   }
798:   PetscOptionsTail();
799:   return(0);
800: }

804: static PetscErrorCode EPSLanczosSetReorthog_Lanczos(EPS eps,EPSLanczosReorthogType reorthog)
805: {
806:   EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;

809:   switch (reorthog) {
810:     case EPS_LANCZOS_REORTHOG_LOCAL:
811:     case EPS_LANCZOS_REORTHOG_FULL:
812:     case EPS_LANCZOS_REORTHOG_DELAYED:
813:     case EPS_LANCZOS_REORTHOG_SELECTIVE:
814:     case EPS_LANCZOS_REORTHOG_PERIODIC:
815:     case EPS_LANCZOS_REORTHOG_PARTIAL:
816:       lanczos->reorthog = reorthog;
817:       break;
818:     default:
819:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
820:   }
821:   return(0);
822: }

826: /*@
827:    EPSLanczosSetReorthog - Sets the type of reorthogonalization used during the Lanczos
828:    iteration.

830:    Logically Collective on EPS

832:    Input Parameters:
833: +  eps - the eigenproblem solver context
834: -  reorthog - the type of reorthogonalization

836:    Options Database Key:
837: .  -eps_lanczos_reorthog - Sets the reorthogonalization type (either 'local', 'selective',
838:                          'periodic', 'partial', 'full' or 'delayed')

840:    Level: advanced

842: .seealso: EPSLanczosGetReorthog(), EPSLanczosReorthogType
843: @*/
844: PetscErrorCode EPSLanczosSetReorthog(EPS eps,EPSLanczosReorthogType reorthog)
845: {

851:   PetscTryMethod(eps,"EPSLanczosSetReorthog_C",(EPS,EPSLanczosReorthogType),(eps,reorthog));
852:   return(0);
853: }

857: static PetscErrorCode EPSLanczosGetReorthog_Lanczos(EPS eps,EPSLanczosReorthogType *reorthog)
858: {
859:   EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;

862:   *reorthog = lanczos->reorthog;
863:   return(0);
864: }

868: /*@C
869:    EPSLanczosGetReorthog - Gets the type of reorthogonalization used during the Lanczos
870:    iteration.

872:    Not Collective

874:    Input Parameter:
875: .  eps - the eigenproblem solver context

877:    Output Parameter:
878: .  reorthog - the type of reorthogonalization

880:    Level: advanced

882: .seealso: EPSLanczosSetReorthog(), EPSLanczosReorthogType
883: @*/
884: PetscErrorCode EPSLanczosGetReorthog(EPS eps,EPSLanczosReorthogType *reorthog)
885: {

891:   PetscTryMethod(eps,"EPSLanczosGetReorthog_C",(EPS,EPSLanczosReorthogType*),(eps,reorthog));
892:   return(0);
893: }

897: PetscErrorCode EPSReset_Lanczos(EPS eps)
898: {
900:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS*)eps->data;

903:   VecDestroyVecs(eps->ncv,&lanczos->AV);
904:   EPSReset_Default(eps);
905:   return(0);
906: }

910: PetscErrorCode EPSDestroy_Lanczos(EPS eps)
911: {

915:   PetscFree(eps->data);
916:   PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",NULL);
917:   PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",NULL);
918:   return(0);
919: }

923: PetscErrorCode EPSView_Lanczos(EPS eps,PetscViewer viewer)
924: {
926:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS*)eps->data;
927:   PetscBool      isascii;

930:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
931:   if (isascii) {
932:     PetscViewerASCIIPrintf(viewer,"  Lanczos: %s reorthogonalization\n",EPSLanczosReorthogTypes[lanczos->reorthog]);
933:   }
934:   return(0);
935: }

939: PETSC_EXTERN PetscErrorCode EPSCreate_Lanczos(EPS eps)
940: {

944:   PetscNewLog(eps,EPS_LANCZOS,&eps->data);
945:   eps->ops->setup                = EPSSetUp_Lanczos;
946:   eps->ops->setfromoptions       = EPSSetFromOptions_Lanczos;
947:   eps->ops->destroy              = EPSDestroy_Lanczos;
948:   eps->ops->reset                = EPSReset_Lanczos;
949:   eps->ops->view                 = EPSView_Lanczos;
950:   eps->ops->backtransform        = EPSBackTransform_Default;
951:   eps->ops->computevectors       = EPSComputeVectors_Hermitian;
952:   PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",EPSLanczosSetReorthog_Lanczos);
953:   PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",EPSLanczosGetReorthog_Lanczos);
954:   return(0);
955: }