Actual source code: lanczos.c

slepc-3.7.0 2016-05-16
Report Typos and Errors
  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://slepc.upv.es.

 18:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 19:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 20:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

 22:    This file is part of SLEPc.

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

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

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

 38: #include <slepc/private/epsimpl.h>                /*I "slepceps.h" I*/
 39: #include <slepcblaslapack.h>

 41: PetscErrorCode EPSSolve_Lanczos(EPS);

 43: typedef struct {
 44:   EPSLanczosReorthogType reorthog;
 45:   BV                     AV;
 46: } EPS_LANCZOS;

 50: PetscErrorCode EPSSetUp_Lanczos(EPS eps)
 51: {
 52:   EPS_LANCZOS        *lanczos = (EPS_LANCZOS*)eps->data;
 53:   BVOrthogRefineType refine;
 54:   BVOrthogBlockType  btype;
 55:   PetscReal          eta;
 56:   PetscErrorCode     ierr;

 59:   EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
 60:   if (eps->ncv>eps->nev+eps->mpd) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must not be larger than nev+mpd");
 61:   if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
 62:   if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
 63:   switch (eps->which) {
 64:     case EPS_LARGEST_IMAGINARY:
 65:     case EPS_SMALLEST_IMAGINARY:
 66:     case EPS_TARGET_IMAGINARY:
 67:       SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
 68:     default: ; /* default case to remove warning */
 69:   }
 70:   if (!eps->extraction) {
 71:     EPSSetExtraction(eps,EPS_RITZ);
 72:   } else if (eps->extraction!=EPS_RITZ) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
 73:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");

 75:   EPSAllocateSolution(eps,1);
 76:   EPS_SetInnerProduct(eps);
 77:   if (lanczos->reorthog != EPS_LANCZOS_REORTHOG_FULL) {
 78:     BVGetOrthogonalization(eps->V,NULL,&refine,&eta,&btype);
 79:     BVSetOrthogonalization(eps->V,BV_ORTHOG_MGS,refine,eta,btype);
 80:     PetscInfo(eps,"Switching to MGS orthogonalization\n");
 81:   }
 82:   if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_SELECTIVE) {
 83:     BVDuplicate(eps->V,&lanczos->AV);
 84:   }

 86:   DSSetType(eps->ds,DSHEP);
 87:   DSSetCompact(eps->ds,PETSC_TRUE);
 88:   DSAllocate(eps->ds,eps->ncv+1);
 89:   if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
 90:     EPSSetWorkVecs(eps,1);
 91:   }

 93:   /* dispatch solve method */
 94:   if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
 95:   if (eps->isgeneralized && eps->ishermitian && !eps->ispositive) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method does not work for indefinite problems");
 96:   eps->ops->solve = EPSSolve_Lanczos;
 97:   return(0);
 98: }

102: /*
103:    EPSLocalLanczos - Local reorthogonalization.

105:    This is the simplest variant. At each Lanczos step, the corresponding Lanczos vector
106:    is orthogonalized with respect to the two previous Lanczos vectors, according to
107:    the three term Lanczos recurrence. WARNING: This variant does not track the loss of
108:    orthogonality that occurs in finite-precision arithmetic and, therefore, the
109:    generated vectors are not guaranteed to be (semi-)orthogonal.
110: */
111: static PetscErrorCode EPSLocalLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *M,PetscBool *breakdown)
112: {
114:   PetscInt       i,j,m = *M;
115:   Vec            vj,vj1;
116:   PetscBool      *which,lwhich[100];
117:   PetscScalar    *hwork,lhwork[100];

120:   if (m > 100) {
121:     PetscMalloc2(m,&which,m,&hwork);
122:   } else {
123:     which = lwhich;
124:     hwork = lhwork;
125:   }
126:   for (i=0;i<k;i++) which[i] = PETSC_TRUE;

128:   BVSetActiveColumns(eps->V,0,m);
129:   for (j=k;j<m;j++) {
130:     BVGetColumn(eps->V,j,&vj);
131:     BVGetColumn(eps->V,j+1,&vj1);
132:     STApply(eps->st,vj,vj1);
133:     BVRestoreColumn(eps->V,j,&vj);
134:     BVRestoreColumn(eps->V,j+1,&vj1);
135:     which[j] = PETSC_TRUE;
136:     if (j-2>=k) which[j-2] = PETSC_FALSE;
137:     BVOrthogonalizeSomeColumn(eps->V,j+1,which,hwork,beta+j,breakdown);
138:     alpha[j] = PetscRealPart(hwork[j]);
139:     if (*breakdown) {
140:       *M = j+1;
141:       break;
142:     } else {
143:       BVScaleColumn(eps->V,j+1,1/beta[j]);
144:     }
145:   }
146:   if (m > 100) {
147:     PetscFree2(which,hwork);
148:   }
149:   return(0);
150: }

154: /*
155:    DenseTridiagonal - Solves a real tridiagonal Hermitian Eigenvalue Problem.

157:    Input Parameters:
158: +  n   - dimension of the eigenproblem
159: .  D   - pointer to the array containing the diagonal elements
160: -  E   - pointer to the array containing the off-diagonal elements

162:    Output Parameters:
163: +  w  - pointer to the array to store the computed eigenvalues
164: -  V  - pointer to the array to store the eigenvectors

166:    Notes:
167:    If V is NULL then the eigenvectors are not computed.

169:    This routine use LAPACK routines xSTEVR.
170: */
171: static PetscErrorCode DenseTridiagonal(PetscInt n_,PetscReal *D,PetscReal *E,PetscReal *w,PetscScalar *V)
172: {
173: #if defined(SLEPC_MISSING_LAPACK_STEVR)
175:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"STEVR - Lapack routine is unavailable");
176: #else
178:   PetscReal      abstol = 0.0,vl,vu,*work;
179:   PetscBLASInt   il,iu,m,*isuppz,n,lwork,*iwork,liwork,info;
180:   const char     *jobz;
181: #if defined(PETSC_USE_COMPLEX)
182:   PetscInt       i,j;
183:   PetscReal      *VV;
184: #endif

187:   PetscBLASIntCast(n_,&n);
188:   PetscBLASIntCast(20*n_,&lwork);
189:   PetscBLASIntCast(10*n_,&liwork);
190:   if (V) {
191:     jobz = "V";
192: #if defined(PETSC_USE_COMPLEX)
193:     PetscMalloc1(n*n,&VV);
194: #endif
195:   } else jobz = "N";
196:   PetscMalloc3(2*n,&isuppz,lwork,&work,liwork,&iwork);
197:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
198: #if defined(PETSC_USE_COMPLEX)
199:   PetscStackCallBLAS("LAPACKstevr",LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,VV,&n,isuppz,work,&lwork,iwork,&liwork,&info));
200: #else
201:   PetscStackCallBLAS("LAPACKstevr",LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info));
202: #endif
203:   PetscFPTrapPop();
204:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack DSTEVR %d",info);
205: #if defined(PETSC_USE_COMPLEX)
206:   if (V) {
207:     for (i=0;i<n;i++)
208:       for (j=0;j<n;j++)
209:         V[i*n+j] = VV[i*n+j];
210:     PetscFree(VV);
211:   }
212: #endif
213:   PetscFree3(isuppz,work,iwork);
214:   return(0);
215: #endif
216: }

220: /*
221:    EPSSelectiveLanczos - Selective reorthogonalization.
222: */
223: static PetscErrorCode EPSSelectiveLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscReal anorm)
224: {
226:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS*)eps->data;
227:   PetscInt       i,j,m = *M,n,nritz=0,nritzo;
228:   Vec            vj,vj1,av;
229:   PetscReal      *d,*e,*ritz,norm;
230:   PetscScalar    *Y,*hwork;
231:   PetscBool      *which;

234:   PetscCalloc6(m+1,&d,m,&e,m,&ritz,m*m,&Y,m,&which,m,&hwork);
235:   for (i=0;i<k;i++) which[i] = PETSC_TRUE;

237:   for (j=k;j<m;j++) {
238:     BVSetActiveColumns(eps->V,0,m);

240:     /* Lanczos step */
241:     BVGetColumn(eps->V,j,&vj);
242:     BVGetColumn(eps->V,j+1,&vj1);
243:     STApply(eps->st,vj,vj1);
244:     BVRestoreColumn(eps->V,j,&vj);
245:     BVRestoreColumn(eps->V,j+1,&vj1);
246:     which[j] = PETSC_TRUE;
247:     if (j-2>=k) which[j-2] = PETSC_FALSE;
248:     BVOrthogonalizeSomeColumn(eps->V,j+1,which,hwork,&norm,breakdown);
249:     alpha[j] = PetscRealPart(hwork[j]);
250:     beta[j] = norm;
251:     if (*breakdown) {
252:       *M = j+1;
253:       break;
254:     }

256:     /* Compute eigenvalues and eigenvectors Y of the tridiagonal block */
257:     n = j-k+1;
258:     for (i=0;i<n;i++) {
259:       d[i] = alpha[i+k];
260:       e[i] = beta[i+k];
261:     }
262:     DenseTridiagonal(n,d,e,ritz,Y);

264:     /* Estimate ||A|| */
265:     for (i=0;i<n;i++)
266:       if (PetscAbsReal(ritz[i]) > anorm) anorm = PetscAbsReal(ritz[i]);

268:     /* Compute nearly converged Ritz vectors */
269:     nritzo = 0;
270:     for (i=0;i<n;i++) {
271:       if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm) nritzo++;
272:     }
273:     if (nritzo>nritz) {
274:       nritz = 0;
275:       for (i=0;i<n;i++) {
276:         if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm) {
277:           BVSetActiveColumns(eps->V,k,k+n);
278:           BVGetColumn(lanczos->AV,nritz,&av);
279:           BVMultVec(eps->V,1.0,0.0,av,Y+i*n);
280:           BVRestoreColumn(lanczos->AV,nritz,&av);
281:           nritz++;
282:         }
283:       }
284:     }
285:     if (nritz > 0) {
286:       BVGetColumn(eps->V,j+1,&vj1);
287:       BVSetActiveColumns(lanczos->AV,0,nritz);
288:       BVOrthogonalizeVec(lanczos->AV,vj1,hwork,&norm,breakdown);
289:       BVRestoreColumn(eps->V,j+1,&vj1);
290:       if (*breakdown) {
291:         *M = j+1;
292:         break;
293:       }
294:     }
295:     BVScaleColumn(eps->V,j+1,1.0/norm);
296:   }

298:   PetscFree6(d,e,ritz,Y,which,hwork);
299:   return(0);
300: }

304: static void update_omega(PetscReal *omega,PetscReal *omega_old,PetscInt j,PetscReal *alpha,PetscReal *beta,PetscReal eps1,PetscReal anorm)
305: {
306:   PetscInt  k;
307:   PetscReal T,binv;

310:   /* Estimate of contribution to roundoff errors from A*v
311:        fl(A*v) = A*v + f,
312:      where ||f|| \approx eps1*||A||.
313:      For a full matrix A, a rule-of-thumb estimate is eps1 = sqrt(n)*eps */
314:   T = eps1*anorm;
315:   binv = 1.0/beta[j+1];

317:   /* Update omega(1) using omega(0)==0 */
318:   omega_old[0]= beta[1]*omega[1] + (alpha[0]-alpha[j])*omega[0] - beta[j]*omega_old[0];
319:   if (omega_old[0] > 0) omega_old[0] = binv*(omega_old[0] + T);
320:   else omega_old[0] = binv*(omega_old[0] - T);

322:   /* Update remaining components */
323:   for (k=1;k<j-1;k++) {
324:     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];
325:     if (omega_old[k] > 0) omega_old[k] = binv*(omega_old[k] + T);
326:     else omega_old[k] = binv*(omega_old[k] - T);
327:   }
328:   omega_old[j-1] = binv*T;

330:   /* Swap omega and omega_old */
331:   for (k=0;k<j;k++) {
332:     omega[k] = omega_old[k];
333:     omega_old[k] = omega[k];
334:   }
335:   omega[j] = eps1;
336:   PetscFunctionReturnVoid();
337: }

341: static void compute_int(PetscBool *which,PetscReal *mu,PetscInt j,PetscReal delta,PetscReal eta)
342: {
343:   PetscInt  i,k,maxpos;
344:   PetscReal max;
345:   PetscBool found;

348:   /* initialize which */
349:   found = PETSC_FALSE;
350:   maxpos = 0;
351:   max = 0.0;
352:   for (i=0;i<j;i++) {
353:     if (PetscAbsReal(mu[i]) >= delta) {
354:       which[i] = PETSC_TRUE;
355:       found = PETSC_TRUE;
356:     } else which[i] = PETSC_FALSE;
357:     if (PetscAbsReal(mu[i]) > max) {
358:       maxpos = i;
359:       max = PetscAbsReal(mu[i]);
360:     }
361:   }
362:   if (!found) which[maxpos] = PETSC_TRUE;

364:   for (i=0;i<j;i++) {
365:     if (which[i]) {
366:       /* find left interval */
367:       for (k=i;k>=0;k--) {
368:         if (PetscAbsReal(mu[k])<eta || which[k]) break;
369:         else which[k] = PETSC_TRUE;
370:       }
371:       /* find right interval */
372:       for (k=i+1;k<j;k++) {
373:         if (PetscAbsReal(mu[k])<eta || which[k]) break;
374:         else which[k] = PETSC_TRUE;
375:       }
376:     }
377:   }
378:   PetscFunctionReturnVoid();
379: }

383: /*
384:    EPSPartialLanczos - Partial reorthogonalization.
385: */
386: static PetscErrorCode EPSPartialLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscReal anorm)
387: {
389:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS*)eps->data;
390:   PetscInt       i,j,m = *M;
391:   Vec            vj,vj1;
392:   PetscReal      norm,*omega,lomega[100],*omega_old,lomega_old[100],eps1,delta,eta;
393:   PetscBool      *which,lwhich[100],*which2,lwhich2[100];
394:   PetscBool      reorth = PETSC_FALSE,force_reorth = PETSC_FALSE;
395:   PetscBool      fro = PETSC_FALSE,estimate_anorm = PETSC_FALSE;
396:   PetscScalar    *hwork,lhwork[100];

399:   if (m>100) {
400:     PetscMalloc5(m,&omega,m,&omega_old,m,&which,m,&which2,m,&hwork);
401:   } else {
402:     omega     = lomega;
403:     omega_old = lomega_old;
404:     which     = lwhich;
405:     which2    = lwhich2;
406:     hwork     = lhwork;
407:   }

409:   eps1 = PetscSqrtReal((PetscReal)eps->n)*PETSC_MACHINE_EPSILON/2;
410:   delta = PETSC_SQRT_MACHINE_EPSILON/PetscSqrtReal((PetscReal)eps->ncv);
411:   eta = PetscPowReal(PETSC_MACHINE_EPSILON,3.0/4.0)/PetscSqrtReal((PetscReal)eps->ncv);
412:   if (anorm < 0.0) {
413:     anorm = 1.0;
414:     estimate_anorm = PETSC_TRUE;
415:   }
416:   for (i=0;i<m-k;i++) omega[i] = omega_old[i] = 0.0;
417:   for (i=0;i<k;i++) which[i] = PETSC_TRUE;

419:   BVSetActiveColumns(eps->V,0,m);
420:   for (j=k;j<m;j++) {
421:     BVGetColumn(eps->V,j,&vj);
422:     BVGetColumn(eps->V,j+1,&vj1);
423:     STApply(eps->st,vj,vj1);
424:     BVRestoreColumn(eps->V,j,&vj);
425:     BVRestoreColumn(eps->V,j+1,&vj1);
426:     if (fro) {
427:       /* Lanczos step with full reorthogonalization */
428:       BVOrthogonalizeColumn(eps->V,j+1,hwork,&norm,breakdown);
429:       alpha[j] = PetscRealPart(hwork[j]);
430:     } else {
431:       /* Lanczos step */
432:       which[j] = PETSC_TRUE;
433:       if (j-2>=k) which[j-2] = PETSC_FALSE;
434:       BVOrthogonalizeSomeColumn(eps->V,j+1,which,hwork,&norm,breakdown);
435:       alpha[j] = PetscRealPart(hwork[j]);
436:       beta[j] = norm;

438:       /* Estimate ||A|| if needed */
439:       if (estimate_anorm) {
440:         if (j>k) anorm = PetscMax(anorm,PetscAbsReal(alpha[j])+norm+beta[j-1]);
441:         else anorm = PetscMax(anorm,PetscAbsReal(alpha[j])+norm);
442:       }

444:       /* Check if reorthogonalization is needed */
445:       reorth = PETSC_FALSE;
446:       if (j>k) {
447:         update_omega(omega,omega_old,j,alpha,beta-1,eps1,anorm);
448:         for (i=0;i<j-k;i++) {
449:           if (PetscAbsReal(omega[i]) > delta) reorth = PETSC_TRUE;
450:         }
451:       }
452:       if (reorth || force_reorth) {
453:         for (i=0;i<k;i++) which2[i] = PETSC_FALSE;
454:         for (i=k;i<=j;i++) which2[i] = PETSC_TRUE;
455:         if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_PERIODIC) {
456:           /* Periodic reorthogonalization */
457:           if (force_reorth) force_reorth = PETSC_FALSE;
458:           else force_reorth = PETSC_TRUE;
459:           for (i=0;i<j-k;i++) omega[i] = eps1;
460:         } else {
461:           /* Partial reorthogonalization */
462:           if (force_reorth) force_reorth = PETSC_FALSE;
463:           else {
464:             force_reorth = PETSC_TRUE;
465:             compute_int(which2+k,omega,j-k,delta,eta);
466:             for (i=0;i<j-k;i++) {
467:               if (which2[i+k]) omega[i] = eps1;
468:             }
469:           }
470:         }
471:         BVOrthogonalizeSomeColumn(eps->V,j+1,which2,hwork,&norm,breakdown);
472:       }
473:     }

475:     if (*breakdown || norm < eps->n*anorm*PETSC_MACHINE_EPSILON) {
476:       *M = j+1;
477:       break;
478:     }
479:     if (!fro && norm*delta < anorm*eps1) {
480:       fro = PETSC_TRUE;
481:       PetscInfo1(eps,"Switching to full reorthogonalization at iteration %D\n",eps->its);
482:     }
483:     beta[j] = norm;
484:     BVScaleColumn(eps->V,j+1,1.0/norm);
485:   }

487:   if (m>100) {
488:     PetscFree5(omega,omega_old,which,which2,hwork);
489:   }
490:   return(0);
491: }

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

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

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

509:    This function simply calls another function which depends on the selected
510:    reorthogonalization strategy.
511: */
512: static PetscErrorCode EPSBasicLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *m,PetscBool *breakdown,PetscReal anorm)
513: {
514:   PetscErrorCode     ierr;
515:   EPS_LANCZOS        *lanczos = (EPS_LANCZOS*)eps->data;
516:   PetscScalar        *T;
517:   PetscInt           i,n=*m;
518:   PetscReal          betam;
519:   BVOrthogRefineType orthog_ref;

522:   switch (lanczos->reorthog) {
523:     case EPS_LANCZOS_REORTHOG_LOCAL:
524:       EPSLocalLanczos(eps,alpha,beta,k,m,breakdown);
525:       break;
526:     case EPS_LANCZOS_REORTHOG_SELECTIVE:
527:       EPSSelectiveLanczos(eps,alpha,beta,k,m,breakdown,anorm);
528:       break;
529:     case EPS_LANCZOS_REORTHOG_FULL:
530:       EPSFullLanczos(eps,alpha,beta,k,m,breakdown);
531:       break;
532:     case EPS_LANCZOS_REORTHOG_PARTIAL:
533:     case EPS_LANCZOS_REORTHOG_PERIODIC:
534:       EPSPartialLanczos(eps,alpha,beta,k,m,breakdown,anorm);
535:       break;
536:     case EPS_LANCZOS_REORTHOG_DELAYED:
537:       PetscMalloc1(n*n,&T);
538:       BVGetOrthogonalization(eps->V,NULL,&orthog_ref,NULL,NULL);
539:       if (orthog_ref == BV_ORTHOG_REFINE_NEVER) {
540:         EPSDelayedArnoldi1(eps,T,n,k,m,&betam,breakdown);
541:       } else {
542:         EPSDelayedArnoldi(eps,T,n,k,m,&betam,breakdown);
543:       }
544:       for (i=k;i<n-1;i++) {
545:         alpha[i] = PetscRealPart(T[n*i+i]);
546:         beta[i] = PetscRealPart(T[n*i+i+1]);
547:       }
548:       alpha[n-1] = PetscRealPart(T[n*(n-1)+n-1]);
549:       beta[n-1] = betam;
550:       PetscFree(T);
551:       break;
552:     default:
553:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
554:   }
555:   return(0);
556: }

560: PetscErrorCode EPSSolve_Lanczos(EPS eps)
561: {
562:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS*)eps->data;
564:   PetscInt       nconv,i,j,k,l,x,n,*perm,restart,ncv=eps->ncv,r,ld;
565:   Vec            vi,vj,w;
566:   Mat            U;
567:   PetscScalar    *Y,*ritz,stmp;
568:   PetscReal      *d,*e,*bnd,anorm,beta,norm,rtmp,resnorm;
569:   PetscBool      breakdown;
570:   char           *conv,ctmp;

573:   DSGetLeadingDimension(eps->ds,&ld);
574:   PetscMalloc4(ncv,&ritz,ncv,&bnd,ncv,&perm,ncv,&conv);

576:   /* The first Lanczos vector is the normalized initial vector */
577:   EPSGetStartVector(eps,0,NULL);

579:   anorm = -1.0;
580:   nconv = 0;

582:   /* Restart loop */
583:   while (eps->reason == EPS_CONVERGED_ITERATING) {
584:     eps->its++;

586:     /* Compute an ncv-step Lanczos factorization */
587:     n = PetscMin(nconv+eps->mpd,ncv);
588:     DSGetArrayReal(eps->ds,DS_MAT_T,&d);
589:     e = d + ld;
590:     EPSBasicLanczos(eps,d,e,nconv,&n,&breakdown,anorm);
591:     beta = e[n-1];
592:     DSRestoreArrayReal(eps->ds,DS_MAT_T,&d);
593:     DSSetDimensions(eps->ds,n,0,nconv,0);
594:     DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
595:     BVSetActiveColumns(eps->V,nconv,n);

597:     /* Solve projected problem */
598:     DSSolve(eps->ds,ritz,NULL);
599:     DSSort(eps->ds,ritz,NULL,NULL,NULL,NULL);

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

605:     /* Compute residual norm estimates as beta*abs(Y(m,:)) + eps*||A|| */
606:     DSGetArray(eps->ds,DS_MAT_Q,&Y);
607:     for (i=nconv;i<n;i++) {
608:       resnorm = beta*PetscAbsScalar(Y[n-1+i*ld]) + PETSC_MACHINE_EPSILON*anorm;
609:       (*eps->converged)(eps,ritz[i],eps->eigi[i],resnorm,&bnd[i],eps->convergedctx);
610:       if (bnd[i]<eps->tol) conv[i] = 'C';
611:       else conv[i] = 'N';
612:     }
613:     DSRestoreArray(eps->ds,DS_MAT_Q,&Y);

615:     /* purge repeated ritz values */
616:     if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
617:       for (i=nconv+1;i<n;i++) {
618:         if (conv[i] == 'C' && PetscAbsScalar((ritz[i]-ritz[i-1])/ritz[i]) < eps->tol) conv[i] = 'R';
619:       }
620:     }

622:     /* Compute restart vector */
623:     if (breakdown) {
624:       PetscInfo2(eps,"Breakdown in Lanczos method (it=%D norm=%g)\n",eps->its,(double)beta);
625:     } else {
626:       restart = nconv;
627:       while (restart<n && conv[restart] != 'N') restart++;
628:       if (restart >= n) {
629:         breakdown = PETSC_TRUE;
630:       } else {
631:         for (i=restart+1;i<n;i++) {
632:           if (conv[i] == 'N') {
633:             SlepcSCCompare(eps->sc,ritz[restart],0.0,ritz[i],0.0,&r);
634:             if (r>0) restart = i;
635:           }
636:         }
637:         DSGetArray(eps->ds,DS_MAT_Q,&Y);
638:         BVMultColumn(eps->V,1.0,0.0,n,Y+restart*ld+nconv);
639:         DSRestoreArray(eps->ds,DS_MAT_Q,&Y);
640:       }
641:     }

643:     /* Count and put converged eigenvalues first */
644:     for (i=nconv;i<n;i++) perm[i] = i;
645:     for (k=nconv;k<n;k++) {
646:       if (conv[perm[k]] != 'C') {
647:         j = k + 1;
648:         while (j<n && conv[perm[j]] != 'C') j++;
649:         if (j>=n) break;
650:         l = perm[k]; perm[k] = perm[j]; perm[j] = l;
651:       }
652:     }

654:     /* Sort eigenvectors according to permutation */
655:     DSGetArray(eps->ds,DS_MAT_Q,&Y);
656:     for (i=nconv;i<k;i++) {
657:       x = perm[i];
658:       if (x != i) {
659:         j = i + 1;
660:         while (perm[j] != i) j++;
661:         /* swap eigenvalues i and j */
662:         stmp = ritz[x]; ritz[x] = ritz[i]; ritz[i] = stmp;
663:         rtmp = bnd[x]; bnd[x] = bnd[i]; bnd[i] = rtmp;
664:         ctmp = conv[x]; conv[x] = conv[i]; conv[i] = ctmp;
665:         perm[j] = x; perm[i] = i;
666:         /* swap eigenvectors i and j */
667:         for (l=0;l<n;l++) {
668:           stmp = Y[l+x*ld]; Y[l+x*ld] = Y[l+i*ld]; Y[l+i*ld] = stmp;
669:         }
670:       }
671:     }
672:     DSRestoreArray(eps->ds,DS_MAT_Q,&Y);

674:     /* compute converged eigenvectors */
675:     DSGetMat(eps->ds,DS_MAT_Q,&U);
676:     BVMultInPlace(eps->V,U,nconv,k);
677:     MatDestroy(&U);

679:     /* purge spurious ritz values */
680:     if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
681:       for (i=nconv;i<k;i++) {
682:         BVGetColumn(eps->V,i,&vi);
683:         VecNorm(vi,NORM_2,&norm);
684:         VecScale(vi,1.0/norm);
685:         w = eps->work[0];
686:         STApply(eps->st,vi,w);
687:         VecAXPY(w,-ritz[i],vi);
688:         BVRestoreColumn(eps->V,i,&vi);
689:         VecNorm(w,NORM_2,&norm);
690:         (*eps->converged)(eps,ritz[i],eps->eigi[i],norm,&bnd[i],eps->convergedctx);
691:         if (bnd[i]>=eps->tol) conv[i] = 'S';
692:       }
693:       for (i=nconv;i<k;i++) {
694:         if (conv[i] != 'C') {
695:           j = i + 1;
696:           while (j<k && conv[j] != 'C') j++;
697:           if (j>=k) break;
698:           /* swap eigenvalues i and j */
699:           stmp = ritz[j]; ritz[j] = ritz[i]; ritz[i] = stmp;
700:           rtmp = bnd[j]; bnd[j] = bnd[i]; bnd[i] = rtmp;
701:           ctmp = conv[j]; conv[j] = conv[i]; conv[i] = ctmp;
702:           /* swap eigenvectors i and j */
703:           BVGetColumn(eps->V,i,&vi);
704:           BVGetColumn(eps->V,j,&vj);
705:           VecSwap(vi,vj);
706:           BVRestoreColumn(eps->V,i,&vi);
707:           BVRestoreColumn(eps->V,j,&vj);
708:         }
709:       }
710:       k = i;
711:     }

713:     /* store ritz values and estimated errors */
714:     for (i=nconv;i<n;i++) {
715:       eps->eigr[i] = ritz[i];
716:       eps->errest[i] = bnd[i];
717:     }
718:     nconv = k;
719:     EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,n);
720:     (*eps->stopping)(eps,eps->its,eps->max_it,nconv,eps->nev,&eps->reason,eps->stoppingctx);

722:     if (eps->reason == EPS_CONVERGED_ITERATING) { /* copy restart vector */
723:       BVCopyColumn(eps->V,n,nconv);
724:       if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL && !breakdown) {
725:         /* Reorthonormalize restart vector */
726:         BVOrthogonalizeColumn(eps->V,nconv,NULL,&norm,&breakdown);
727:         BVScaleColumn(eps->V,nconv,1.0/norm);
728:       }
729:       if (breakdown) {
730:         /* Use random vector for restarting */
731:         PetscInfo(eps,"Using random vector for restart\n");
732:         EPSGetStartVector(eps,nconv,&breakdown);
733:       }
734:       if (breakdown) { /* give up */
735:         eps->reason = EPS_DIVERGED_BREAKDOWN;
736:         PetscInfo(eps,"Unable to generate more start vectors\n");
737:       }
738:     }
739:   }
740:   eps->nconv = nconv;

742:   PetscFree4(ritz,bnd,perm,conv);
743:   return(0);
744: }

748: PetscErrorCode EPSSetFromOptions_Lanczos(PetscOptionItems *PetscOptionsObject,EPS eps)
749: {
750:   PetscErrorCode         ierr;
751:   EPS_LANCZOS            *lanczos = (EPS_LANCZOS*)eps->data;
752:   PetscBool              flg;
753:   EPSLanczosReorthogType reorthog;

756:   PetscOptionsHead(PetscOptionsObject,"EPS Lanczos Options");
757:   PetscOptionsEnum("-eps_lanczos_reorthog","Lanczos reorthogonalization","EPSLanczosSetReorthog",EPSLanczosReorthogTypes,(PetscEnum)lanczos->reorthog,(PetscEnum*)&reorthog,&flg);
758:   if (flg) {
759:     EPSLanczosSetReorthog(eps,reorthog);
760:   }
761:   PetscOptionsTail();
762:   return(0);
763: }

767: static PetscErrorCode EPSLanczosSetReorthog_Lanczos(EPS eps,EPSLanczosReorthogType reorthog)
768: {
769:   EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;

772:   switch (reorthog) {
773:     case EPS_LANCZOS_REORTHOG_LOCAL:
774:     case EPS_LANCZOS_REORTHOG_FULL:
775:     case EPS_LANCZOS_REORTHOG_DELAYED:
776:     case EPS_LANCZOS_REORTHOG_SELECTIVE:
777:     case EPS_LANCZOS_REORTHOG_PERIODIC:
778:     case EPS_LANCZOS_REORTHOG_PARTIAL:
779:       lanczos->reorthog = reorthog;
780:       break;
781:     default:
782:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
783:   }
784:   return(0);
785: }

789: /*@
790:    EPSLanczosSetReorthog - Sets the type of reorthogonalization used during the Lanczos
791:    iteration.

793:    Logically Collective on EPS

795:    Input Parameters:
796: +  eps - the eigenproblem solver context
797: -  reorthog - the type of reorthogonalization

799:    Options Database Key:
800: .  -eps_lanczos_reorthog - Sets the reorthogonalization type (either 'local', 'selective',
801:                          'periodic', 'partial', 'full' or 'delayed')

803:    Level: advanced

805: .seealso: EPSLanczosGetReorthog(), EPSLanczosReorthogType
806: @*/
807: PetscErrorCode EPSLanczosSetReorthog(EPS eps,EPSLanczosReorthogType reorthog)
808: {

814:   PetscTryMethod(eps,"EPSLanczosSetReorthog_C",(EPS,EPSLanczosReorthogType),(eps,reorthog));
815:   return(0);
816: }

820: static PetscErrorCode EPSLanczosGetReorthog_Lanczos(EPS eps,EPSLanczosReorthogType *reorthog)
821: {
822:   EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;

825:   *reorthog = lanczos->reorthog;
826:   return(0);
827: }

831: /*@
832:    EPSLanczosGetReorthog - Gets the type of reorthogonalization used during
833:    the Lanczos iteration.

835:    Not Collective

837:    Input Parameter:
838: .  eps - the eigenproblem solver context

840:    Output Parameter:
841: .  reorthog - the type of reorthogonalization

843:    Level: advanced

845: .seealso: EPSLanczosSetReorthog(), EPSLanczosReorthogType
846: @*/
847: PetscErrorCode EPSLanczosGetReorthog(EPS eps,EPSLanczosReorthogType *reorthog)
848: {

854:   PetscUseMethod(eps,"EPSLanczosGetReorthog_C",(EPS,EPSLanczosReorthogType*),(eps,reorthog));
855:   return(0);
856: }

860: PetscErrorCode EPSReset_Lanczos(EPS eps)
861: {
863:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS*)eps->data;

866:   BVDestroy(&lanczos->AV);
867:   return(0);
868: }

872: PetscErrorCode EPSDestroy_Lanczos(EPS eps)
873: {

877:   PetscFree(eps->data);
878:   PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",NULL);
879:   PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",NULL);
880:   return(0);
881: }

885: PetscErrorCode EPSView_Lanczos(EPS eps,PetscViewer viewer)
886: {
888:   EPS_LANCZOS    *lanczos = (EPS_LANCZOS*)eps->data;
889:   PetscBool      isascii;

892:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
893:   if (isascii) {
894:     PetscViewerASCIIPrintf(viewer,"  Lanczos: %s reorthogonalization\n",EPSLanczosReorthogTypes[lanczos->reorthog]);
895:   }
896:   return(0);
897: }

901: PETSC_EXTERN PetscErrorCode EPSCreate_Lanczos(EPS eps)
902: {
903:   EPS_LANCZOS    *ctx;

907:   PetscNewLog(eps,&ctx);
908:   eps->data = (void*)ctx;

910:   eps->ops->setup                = EPSSetUp_Lanczos;
911:   eps->ops->setfromoptions       = EPSSetFromOptions_Lanczos;
912:   eps->ops->destroy              = EPSDestroy_Lanczos;
913:   eps->ops->reset                = EPSReset_Lanczos;
914:   eps->ops->view                 = EPSView_Lanczos;
915:   eps->ops->backtransform        = EPSBackTransform_Default;
916:   eps->ops->computevectors       = EPSComputeVectors_Hermitian;
917:   PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",EPSLanczosSetReorthog_Lanczos);
918:   PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",EPSLanczosGetReorthog_Lanczos);
919:   return(0);
920: }