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: June 2005
20: */
21: #include src/eps/epsimpl.h
22: #include slepcblaslapack.h
24: typedef struct {
25: EPSLanczosReorthogType reorthog;
26: } EPS_LANCZOS;
30: PetscErrorCode EPSSetUp_LANCZOS(EPS eps)
31: {
33: PetscInt N;
36: VecGetSize(eps->vec_initial,&N);
37: if (eps->nev > N) eps->nev = N;
38: if (eps->ncv) {
39: if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
40: }
41: else eps->ncv = PetscMin(N,PetscMax(2*eps->nev,eps->nev+15));
42: if (!eps->max_it) eps->max_it = PetscMax(100,N);
43: if (!eps->tol) eps->tol = 1.e-7;
44: if (eps->solverclass==EPS_ONE_SIDE) {
45: if (eps->which == EPS_LARGEST_IMAGINARY || eps->which == EPS_SMALLEST_IMAGINARY)
46: SETERRQ(1,"Wrong value of eps->which");
47: if (!eps->ishermitian)
48: SETERRQ(PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
49: } else {
50: if (eps->which != EPS_LARGEST_MAGNITUDE)
51: SETERRQ(1,"Wrong value of eps->which");
52: }
53: EPSAllocateSolution(eps);
54: PetscFree(eps->T);
55: PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);
56: if (eps->solverclass==EPS_TWO_SIDE) {
57: PetscFree(eps->Tl);
58: PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->Tl);
59: EPSDefaultGetWork(eps,2);
60: }
61: else { EPSDefaultGetWork(eps,1); }
62: return(0);
63: }
67: /*
68: EPSFullLanczos - Full reorthogonalization.
70: In this variant, at each Lanczos step, the corresponding Lanczos vector
71: is orthogonalized with respect to all the previous Lanczos vectors.
72: */
73: static PetscErrorCode EPSFullLanczos(EPS eps,PetscScalar *T,Vec *V,int k,int *M,Vec f,PetscReal *beta,PetscTruth *breakdown)
74: {
76: int j,m = *M;
77: PetscReal norm;
80: for (j=k;j<m;j++) {
81: STApply(eps->OP,V[j],f);
82: eps->its++;
83: EPSOrthogonalize(eps,eps->nds,eps->DS,f,PETSC_NULL,PETSC_NULL,PETSC_NULL);
84: EPSOrthogonalize(eps,j+1,V,f,T+m*j,&norm,breakdown);
85: if (*breakdown) {
86: *M = j+1;
87: break;
88: }
89: if (j<m-1) {
90: T[m*j+j+1] = norm;
91: VecScale(f,1.0/norm);
92: VecCopy(f,V[j+1]);
93: }
94: }
95: *beta = norm;
96: return(0);
97: }
101: /*
102: EPSLocalLanczos - Local reorthogonalization.
104: This is the simplest variant. At each Lanczos step, the corresponding Lanczos vector
105: is orthogonalized with respect to the two previous Lanczos vectors, according to
106: the three term Lanczos recurrence. WARNING: This variant does not track the loss of
107: orthogonality that occurs in finite-precision arithmetic and, therefore, the
108: generated vectors are not guaranteed to be (semi-)orthogonal.
109: */
110: static PetscErrorCode EPSLocalLanczos(EPS eps,PetscScalar *T,Vec *V,int k,int *M,Vec f,PetscReal *beta,PetscTruth *breakdown)
111: {
113: int j,m = *M;
114: PetscReal norm;
117: for (j=k;j<m;j++) {
118: STApply(eps->OP,V[j],f);
119: eps->its++;
120: EPSOrthogonalize(eps,eps->nds+k,eps->DSV,f,PETSC_NULL,PETSC_NULL,PETSC_NULL);
121: if (j == k) {
122: EPSOrthogonalize(eps,1,V+j,f,T+m*j+j,&norm,breakdown);
123: } else {
124: EPSOrthogonalize(eps,2,V+j-1,f,T+m*j+j-1,&norm,breakdown);
125: }
126: if (*breakdown) {
127: *M = j+1;
128: break;
129: }
130: if (j<m-1) {
131: T[m*j+j+1] = norm;
132: VecScale(f,1.0/norm);
133: VecCopy(f,V[j+1]);
134: }
135: }
136: *beta = norm;
137: return(0);
138: }
142: /*
143: EPSSelectiveLanczos - Selective reorthogonalization.
144: */
145: static PetscErrorCode EPSSelectiveLanczos(EPS eps,PetscScalar *T,Vec *V,int k,int *M,Vec f,PetscReal *beta,PetscTruth *breakdown,PetscReal anorm)
146: {
147: #if defined(SLEPC_MISSING_LAPACK_DSTEVR) || defined(SLEPC_MISSING_LAPACK_DLAMCH)
149: SETERRQ(PETSC_ERR_SUP,"DSTEVR - Lapack routine is unavailable.");
150: #else
152: int i,j,m = *M,n,il,iu,mout,*isuppz,*iwork,lwork,liwork,info;
153: PetscReal *D,*E,*ritz,*Y,*work,abstol,vl,vu,norm;
154: PetscTruth conv;
157: PetscMalloc(m*sizeof(PetscReal),&D);
158: PetscMalloc(m*sizeof(PetscReal),&E);
159: PetscMalloc(m*sizeof(PetscReal),&ritz);
160: PetscMalloc(m*m*sizeof(PetscReal),&Y);
161: PetscMalloc(2*m*sizeof(int),&isuppz);
162: lwork = 20*m;
163: PetscMalloc(lwork*sizeof(PetscReal),&work);
164: liwork = 10*m;
165: PetscMalloc(liwork*sizeof(int),&iwork);
166: abstol = 2.0*LAPACKlamch_("S",1);
168: for (j=k;j<m;j++) {
169: STApply(eps->OP,V[j],f);
170: eps->its++;
171: EPSOrthogonalize(eps,eps->nds+k,eps->DSV,f,PETSC_NULL,PETSC_NULL,PETSC_NULL);
172: if (j == k) {
173: EPSOrthogonalize(eps,1,V+j,f,T+m*j+j,&norm,breakdown);
174: } else {
175: EPSOrthogonalize(eps,2,V+j-1,f,T+m*j+j-1,&norm,breakdown);
176: }
177:
178: if (*breakdown) {
179: *M = j+1;
180: break;
181: }
183: n = j-k+1;
184: for (i=0;i<n-1;i++) {
185: ritz[i] = PetscRealPart(T[(i+k)*(m+1)]);
186: E[i] = PetscRealPart(T[(i+k)*(m+1)+1]);
187: }
188: ritz[n-1] = PetscRealPart(T[(n-1+k)*(m+1)]);
190: /* Compute eigenvalues and eigenvectors Y of the tridiagonal block */
191: LAPACKstevr_("V","A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&mout,ritz,Y,&n,isuppz,work,&lwork,iwork,&liwork,&info,1,1);
192: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xSTEQR %i",info);
193:
194: /* Estimate ||A|| */
195: for (i=0;i<n;i++)
196: if (PetscAbsReal(ritz[i]) > anorm) anorm = PetscAbsReal(ritz[i]);
197:
198: /* Exit if residual norm is small [Parlett, page 300] */
199: conv = PETSC_FALSE;
200: for (i=0;i<n && !conv;i++)
201: if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm)
202: conv = PETSC_TRUE;
203:
204: if (conv) {
205: *M = j+1;
206: break;
207: }
209: if (j<m-1) {
210: T[m*j+j+1] = norm;
211: VecScale(f,1.0 / norm);
212: VecCopy(f,V[j+1]);
213: }
214: }
215: *beta = norm;
216: PetscFree(D);
217: PetscFree(E);
218: PetscFree(ritz);
219: PetscFree(Y);
220: PetscFree(isuppz);
221: PetscFree(work);
222: PetscFree(iwork);
223: return(0);
224: #endif
225: }
229: /*
230: EPSPeriodicLanczos - Periodic reorthogonalization.
231: */
232: static PetscErrorCode EPSPeriodicLanczos(EPS eps,PetscScalar *T,Vec *V,int k,int *M,Vec f,PetscReal *beta, PetscTruth *breakdown)
233: {
235: int i,j,m = *M;
236: PetscReal norm;
237: PetscScalar *omega;
238: PetscTruth reorthog;
241: reorthog = PETSC_FALSE;
242: PetscMalloc(m*sizeof(PetscScalar),&omega);
243: for (j=k;j<m;j++) {
244: STApply(eps->OP,V[j],f);
245: eps->its++;
246: EPSOrthogonalize(eps,eps->nds,eps->DS,f,PETSC_NULL,PETSC_NULL,PETSC_NULL);
248: if (reorthog) {
249: EPSOrthogonalize(eps,j+1,V,f,T+m*j,&norm,breakdown);
250: } else {
251: for (i=0;i<j-1;i++) T[m*j+i] = 0.0;
252: if (j == k) {
253: EPSOrthogonalize(eps,1,V+j,f,T+m*j+j,&norm,breakdown);
254: } else {
255: EPSOrthogonalize(eps,2,V+j-1,f,T+m*j+j-1,&norm,breakdown);
256: }
257: }
259: if (*breakdown) {
260: *M = j+1;
261: break;
262: }
264: if (reorthog) {
265: reorthog = PETSC_FALSE;
266: } else if (j>1) {
267: VecMDot(f,j-1,eps->V,omega);
268: for (i=0;i<j-1 && !reorthog;i++)
269: if (PetscAbsScalar(omega[i]) > PETSC_SQRT_MACHINE_EPSILON/sqrt((PetscReal)j)) {
270: reorthog = PETSC_TRUE;
271: }
272: if (reorthog) {
273: EPSOrthogonalize(eps,j-1,V,f,T+m*j,&norm,PETSC_NULL);
274: }
275: }
277: if (j<m-1) {
278: T[m*j+j+1] = norm;
279: VecScale(f,1.0/norm);
280: VecCopy(f,V[j+1]);
281: }
282: }
283: *beta = norm;
284: PetscFree(omega);
285: return(0);
286: }
290: /*
291: EPSPartialLanczos - Partial reorthogonalization.
292: */
293: static PetscErrorCode EPSPartialLanczos(EPS eps,PetscScalar *T,Vec *V,int k,int *M,Vec f,PetscReal *beta, PetscTruth *breakdown)
294: {
296: int i,j,l,m = *M;
297: PetscReal norm,nu;
298: PetscScalar *omega;
299: PetscTruth reorthog,*which;
302: nu = sqrt(PETSC_MACHINE_EPSILON*PETSC_SQRT_MACHINE_EPSILON);
303: reorthog = PETSC_FALSE;
304: PetscMalloc(m*sizeof(PetscScalar),&omega);
305: PetscMalloc(m*sizeof(PetscTruth),&which);
306: for (j=k;j<m;j++) {
307: STApply(eps->OP,V[j],f);
308: eps->its++;
309: EPSOrthogonalize(eps,eps->nds,eps->DS,f,PETSC_NULL,PETSC_NULL,PETSC_NULL);
311: if (j == k) {
312: EPSOrthogonalize(eps,1,V+j,f,T+m*j+j,&norm,breakdown);
313: } else {
314: EPSOrthogonalize(eps,2,V+j-1,f,T+m*j+j-1,&norm,breakdown);
315: }
317: if (*breakdown) {
318: *M = j+1;
319: break;
320: }
322: if (reorthog) {
323: for (i=0;i<j-1;i++)
324: if (which[i]) {
325: EPSOrthogonalize(eps,1,V+i,f,PETSC_NULL,&norm,PETSC_NULL);
326: }
327: reorthog = PETSC_FALSE;
328: } else if (j>1) {
329: VecMDot(f,j-1,eps->V,omega);
330: for (i=0;i<j-1;i++) {
331: omega[i] /= norm;
332: which[i] = PETSC_FALSE;
333: }
334: for (i=0;i<j-1;i++)
335: if (PetscAbsScalar(omega[i]) > PETSC_SQRT_MACHINE_EPSILON) {
336: reorthog = PETSC_TRUE;
337: which[i] = PETSC_TRUE;
338: for (l=i-1;l>0 && PetscAbsScalar(omega[l]) > nu;l--) which[l] = PETSC_TRUE;
339: for (l=i+1;l<j-1 && PetscAbsScalar(omega[l]) > nu;l++) which[l] = PETSC_TRUE;
340: }
341: if (reorthog)
342: for (i=0;i<j-1;i++)
343: if (which[i]) {
344: EPSOrthogonalize(eps,1,V+i,f,PETSC_NULL,&norm,PETSC_NULL);
345: }
346: }
348: if (j<m-1) {
349: T[m*j+j+1] = norm;
350: VecScale(f,1.0/norm);
351: VecCopy(f,V[j+1]);
352: }
353: }
354: *beta = norm;
355: PetscFree(omega);
356: PetscFree(which);
357: return(0);
358: }
362: /*
363: EPSBasicLanczos - Computes an m-step Lanczos factorization. The first k
364: columns are assumed to be locked and therefore they are not modified. On
365: exit, the following relation is satisfied:
367: OP * V - V * T = f * e_m^T
369: where the columns of V are the Lanczos vectors, T is a tridiagonal matrix,
370: f is the residual vector and e_m is the m-th vector of the canonical basis.
371: The Lanczos vectors (together with vector f) are B-orthogonal (to working
372: accuracy) if full reorthogonalization is being used, otherwise they are
373: (B-)semi-orthogonal. On exit, beta contains the B-norm of f and the next
374: Lanczos vector can be computed as v_{m+1} = f / beta.
376: This function simply calls another function which depends on the selected
377: reorthogonalization strategy.
378: */
379: static PetscErrorCode EPSBasicLanczos(EPS eps,PetscScalar *T,Vec *V,int k,int *m,Vec f,PetscReal *beta,PetscTruth *breakdown,PetscReal anorm)
380: {
381: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
385: switch (lanczos->reorthog) {
386: case EPSLANCZOS_REORTHOG_LOCAL:
387: EPSLocalLanczos(eps,T,V,k,m,f,beta,breakdown);
388: break;
389: case EPSLANCZOS_REORTHOG_SELECTIVE:
390: EPSSelectiveLanczos(eps,T,V,k,m,f,beta,breakdown,anorm);
391: break;
392: case EPSLANCZOS_REORTHOG_PARTIAL:
393: EPSPartialLanczos(eps,T,V,k,m,f,beta,breakdown);
394: break;
395: case EPSLANCZOS_REORTHOG_PERIODIC:
396: EPSPeriodicLanczos(eps,T,V,k,m,f,beta,breakdown);
397: break;
398: case EPSLANCZOS_REORTHOG_FULL:
399: EPSFullLanczos(eps,T,V,k,m,f,beta,breakdown);
400: break;
401: }
402: return(0);
403: }
407: PetscErrorCode EPSSolve_LANCZOS(EPS eps)
408: {
409: #if defined(SLEPC_MISSING_LAPACK_DSTEVR) || defined(SLEPC_MISSING_LAPACK_DLAMCH)
411: SETERRQ(PETSC_ERR_SUP,"DSTEVR - Lapack routine is unavailable.");
412: #else
413: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
415: int nconv,i,j,k,n,m,*perm,restart,
416: *isuppz,*iwork,mout,lwork,liwork,il,iu,info,
417: ncv=eps->ncv;
418: Vec f=eps->work[0];
419: PetscScalar *T=eps->T;
420: PetscReal *ritz,*Y,*bnd,*D,*E,*work,anorm,beta,norm,abstol,vl,vu,restart_ritz;
421: PetscTruth breakdown;
422: char *conv;
425: PetscMalloc(ncv*sizeof(PetscReal),&D);
426: PetscMalloc(ncv*sizeof(PetscReal),&E);
427: PetscMalloc(ncv*sizeof(PetscReal),&ritz);
428: PetscMalloc(ncv*ncv*sizeof(PetscReal),&Y);
429: PetscMalloc(2*ncv*sizeof(int),&isuppz);
430: lwork = 20*ncv;
431: PetscMalloc(lwork*sizeof(PetscReal),&work);
432: liwork = 10*ncv;
433: PetscMalloc(liwork*sizeof(int),&iwork);
435: PetscMalloc(ncv*sizeof(PetscReal),&bnd);
436: PetscMalloc(ncv*sizeof(int),&perm);
437: PetscMalloc(ncv*sizeof(char),&conv);
439: /* The first Lanczos vector is the normalized initial vector */
440: EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
441:
442: abstol = 2.0*LAPACKlamch_("S",1);
443: anorm = 1.0;
444: nconv = 0;
445: eps->its = 0;
446: restart_ritz = 0.0;
447: for (i=0;i<eps->ncv;i++) eps->eigi[i]=0.0;
448: EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,ncv);
449:
450: /* Restart loop */
451: while (eps->reason == EPS_CONVERGED_ITERATING) {
452: /* Compute an ncv-step Lanczos factorization */
453: m = ncv;
454: EPSBasicLanczos(eps,T,eps->V,nconv,&m,f,&beta,&breakdown,anorm);
456: /* Extract the tridiagonal block from T. Store the diagonal elements in D
457: and the off-diagonal elements in E */
458: n = m - nconv;
459: for (i=0;i<n-1;i++) {
460: D[i] = PetscRealPart(T[(i+nconv)*(ncv+1)]);
461: E[i] = PetscRealPart(T[(i+nconv)*(ncv+1)+1]);
462: }
463: D[n-1] = PetscRealPart(T[(n-1+nconv)*(ncv+1)]);
465: /* Compute eigenvalues and eigenvectors Y of the tridiagonal block */
466: LAPACKstevr_("V","A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&mout,ritz,Y,&n,isuppz,work,&lwork,iwork,&liwork,&info,1,1);
467: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xSTEVR %i",info);
469: /* Estimate ||A|| */
470: for (i=0;i<n;i++)
471: if (PetscAbsReal(ritz[i]) > anorm) anorm = PetscAbsReal(ritz[i]);
472:
473: /* Compute residual norm estimates as beta*abs(Y(m,:)) + eps*||A|| */
474: for (i=0;i<n;i++)
475: bnd[i] = beta*PetscAbsReal(Y[i*n+n-1]) + PETSC_MACHINE_EPSILON*anorm;
477: #ifdef PETSC_USE_COMPLEX
478: for (i=0;i<n;i++) {
479: eps->eigr[i+nconv] = ritz[i];
480: }
481: EPSSortEigenvalues(n,eps->eigr+nconv,eps->eigi,eps->which,n,perm);
482: #else
483: EPSSortEigenvalues(n,ritz,eps->eigi,eps->which,n,perm);
484: #endif
486: /* Look for converged eigenpairs */
487: k = nconv;
488: for (i=0;i<n;i++) {
489: eps->eigr[k] = ritz[perm[i]];
490: eps->errest[k] = bnd[perm[i]] / PetscAbsScalar(eps->eigr[k]);
491: if (eps->errest[k] < eps->tol) {
492:
493: if (lanczos->reorthog == EPSLANCZOS_REORTHOG_LOCAL) {
494: if (i>0 && PetscAbsScalar((eps->eigr[k]-ritz[perm[i-1]])/eps->eigr[k]) < eps->tol) {
495: /* Discard repeated eigenvalues */
496: conv[i] = 'R';
497: continue;
498: }
499: }
500:
501: VecSet(eps->AV[k],0.0);
502: #ifndef PETSC_USE_COMPLEX
503: VecMAXPY(eps->AV[k],n,Y+perm[i]*n,eps->V+nconv);
504: #else
505: for (j=0;j<n;j++) {
506: VecAXPY(eps->AV[k],Y[perm[i]*n+j],eps->V[nconv+j]);
507: }
508: #endif
509:
510: if (lanczos->reorthog == EPSLANCZOS_REORTHOG_LOCAL) {
511: VecNorm(eps->AV[k],NORM_2,&norm);
512: VecScale(eps->AV[k],1.0/norm);
513: eps->errest[k] = eps->errest[k] / norm;
514: if (i<n-1 && PetscAbsScalar((eps->eigr[k]-ritz[perm[i+1]])/eps->eigr[k]) >= eps->tol) {
515: PetscReal res;
516: STApply(eps->OP,eps->AV[k],f);
517: VecAXPY(f,-eps->eigr[k],eps->AV[k]);
518: VecNorm(f,NORM_2,&res);
519: eps->errest[k] = res / PetscAbsScalar(eps->eigr[k]);
520: }
521: if (eps->errest[k] >= eps->tol) {
522: conv[i] = 'S';
523: continue;
524: }
525: }
526:
527: conv[i] = 'C';
528: k++;
529: } else conv[i] = 'N';
530: }
532: /* Look for non-converged eigenpairs */
533: j = k;
534: restart = -1;
535: for (i=0;i<n;i++) {
536: if (conv[i] != 'C') {
537: if (restart == -1 && conv[i] == 'N') restart = i;
538: eps->eigr[j] = ritz[perm[i]];
539: eps->errest[j] = bnd[perm[i]] / ritz[perm[i]];
540: j++;
541: }
542: }
544: if (breakdown) {
545: restart = -1;
546: PetscInfo1(eps,"Breakdown in Lanczos method (norm=%g)\n",beta);
547: eps->count_breakdown++;
548: }
549:
550: if (k<eps->nev) {
551: if (lanczos->reorthog == EPSLANCZOS_REORTHOG_SELECTIVE && restart != -1) {
552: /* Avoid stagnation in selective reorthogonalization */
553: if (PetscAbsScalar(restart_ritz - ritz[perm[restart]]) < PETSC_MACHINE_EPSILON) {
554: restart = -1;
555: restart_ritz = 0;
556: } else restart_ritz = ritz[perm[restart]];
557: }
558: if (restart != -1) {
559: /* Use first non-converged vector for restarting */
560: VecSet(eps->AV[k],0.0);
561: #ifndef PETSC_USE_COMPLEX
562: VecMAXPY(eps->AV[k],n,Y+perm[restart]*n,eps->V+nconv);
563: #else
564: for (j=0;j<n;j++) {
565: VecAXPY(eps->AV[k],Y[perm[restart]*n+j],eps->V[nconv+j]);
566: }
567: #endif
568: VecCopy(eps->AV[k],eps->V[k]);
569: }
570: }
571:
572: /* Copy converged vectors to V */
573: for (i=nconv;i<k;i++) {
574: VecCopy(eps->AV[i],eps->V[i]);
575: }
577: if (k<eps->nev) {
578: if (restart == -1) {
579: /* Use random vector for restarting */
580: PetscInfo(eps,"Using random vector for restart\n");
581: EPSGetStartVector(eps,k,eps->V[k],&breakdown);
582: } else switch (lanczos->reorthog) {
583: case EPSLANCZOS_REORTHOG_LOCAL:
584: case EPSLANCZOS_REORTHOG_PERIODIC:
585: case EPSLANCZOS_REORTHOG_PARTIAL:
586: /* Reorthonormalize restart vector */
587: EPSOrthogonalize(eps,eps->nds+k,eps->DSV,eps->V[k],PETSC_NULL,&norm,&breakdown);
588: VecScale(eps->V[k],1.0/norm);
589: break;
590: default:
591: breakdown = PETSC_FALSE;
592: }
593: if (breakdown) {
594: eps->reason = EPS_DIVERGED_BREAKDOWN;
595: PetscInfo(eps,"Unable to generate more start vectors\n");
596: }
597: }
599: nconv = k;
600: EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nconv+n);
601: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
602: if (nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
603: }
604:
605: eps->nconv = nconv;
607: PetscFree(D);
608: PetscFree(E);
609: PetscFree(ritz);
610: PetscFree(Y);
611: PetscFree(isuppz);
612: PetscFree(work);
613: PetscFree(iwork);
615: PetscFree(bnd);
616: PetscFree(perm);
617: PetscFree(conv);
618: return(0);
619: #endif
620: }
622: static const char *lanczoslist[5] = { "local", "full", "selective", "periodic", "partial" };
626: PetscErrorCode EPSSetFromOptions_LANCZOS(EPS eps)
627: {
629: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
630: PetscTruth flg;
631: PetscInt i;
634: PetscOptionsHead("LANCZOS options");
635: PetscOptionsEList("-eps_lanczos_reorthog","Lanczos reorthogonalization","EPSLanczosSetReorthog",lanczoslist,5,lanczoslist[lanczos->reorthog],&i,&flg);
636: if (flg) lanczos->reorthog = (EPSLanczosReorthogType)i;
637: PetscOptionsTail();
638: return(0);
639: }
644: PetscErrorCode EPSLanczosSetReorthog_LANCZOS(EPS eps,EPSLanczosReorthogType reorthog)
645: {
646: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
649: switch (reorthog) {
650: case EPSLANCZOS_REORTHOG_LOCAL:
651: case EPSLANCZOS_REORTHOG_FULL:
652: case EPSLANCZOS_REORTHOG_SELECTIVE:
653: case EPSLANCZOS_REORTHOG_PERIODIC:
654: case EPSLANCZOS_REORTHOG_PARTIAL:
655: lanczos->reorthog = reorthog;
656: break;
657: default:
658: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
659: }
660: return(0);
661: }
666: /*@
667: EPSLanczosSetReorthog - Sets the type of reorthogonalization used during the Lanczos
668: iteration.
670: Collective on EPS
672: Input Parameters:
673: + eps - the eigenproblem solver context
674: - reorthog - the type of reorthogonalization
676: Options Database Key:
677: . -eps_lanczos_orthog - Sets the reorthogonalization type (either 'local', 'selective',
678: 'periodic', 'partial' or 'full')
679:
680: Level: advanced
682: .seealso: EPSLanczosGetReorthog()
683: @*/
684: PetscErrorCode EPSLanczosSetReorthog(EPS eps,EPSLanczosReorthogType reorthog)
685: {
686: PetscErrorCode ierr, (*f)(EPS,EPSLanczosReorthogType);
690: PetscObjectQueryFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",(void (**)())&f);
691: if (f) {
692: (*f)(eps,reorthog);
693: }
694: return(0);
695: }
700: PetscErrorCode EPSLanczosGetReorthog_LANCZOS(EPS eps,EPSLanczosReorthogType *reorthog)
701: {
702: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
704: *reorthog = lanczos->reorthog;
705: return(0);
706: }
711: /*@C
712: EPSLanczosGetReorthog - Gets the type of reorthogonalization used during the Lanczos
713: iteration.
715: Collective on EPS
717: Input Parameter:
718: . eps - the eigenproblem solver context
720: Input Parameter:
721: . reorthog - the type of reorthogonalization
723: Level: advanced
725: .seealso: EPSLanczosSetReorthog()
726: @*/
727: PetscErrorCode EPSLanczosGetReorthog(EPS eps,EPSLanczosReorthogType *reorthog)
728: {
729: PetscErrorCode ierr, (*f)(EPS,EPSLanczosReorthogType*);
733: PetscObjectQueryFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",(void (**)())&f);
734: if (f) {
735: (*f)(eps,reorthog);
736: }
737: return(0);
738: }
742: PetscErrorCode EPSView_LANCZOS(EPS eps,PetscViewer viewer)
743: {
745: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
746: PetscTruth isascii;
749: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
750: if (!isascii) {
751: SETERRQ1(1,"Viewer type %s not supported for EPSLANCZOS",((PetscObject)viewer)->type_name);
752: }
753: PetscViewerASCIIPrintf(viewer,"reorthogonalization: %s\n",lanczoslist[lanczos->reorthog]);
754: return(0);
755: }
757: /*
758: EXTERN PetscErrorCode EPSSolve_TS_LANCZOS(EPS);
759: */
764: PetscErrorCode EPSCreate_LANCZOS(EPS eps)
765: {
767: EPS_LANCZOS *lanczos;
770: PetscNew(EPS_LANCZOS,&lanczos);
771: PetscLogObjectMemory(eps,sizeof(EPS_LANCZOS));
772: eps->data = (void *) lanczos;
773: eps->ops->solve = EPSSolve_LANCZOS;
774: /* eps->ops->solvets = EPSSolve_TS_LANCZOS;*/
775: eps->ops->setup = EPSSetUp_LANCZOS;
776: eps->ops->setfromoptions = EPSSetFromOptions_LANCZOS;
777: eps->ops->destroy = EPSDestroy_Default;
778: eps->ops->view = EPSView_LANCZOS;
779: eps->ops->backtransform = EPSBackTransform_Default;
780: /*if (eps->solverclass==EPS_TWO_SIDE)
781: eps->ops->computevectors = EPSComputeVectors_Schur;
782: else*/ eps->ops->computevectors = EPSComputeVectors_Default;
783: lanczos->reorthog = EPSLANCZOS_REORTHOG_LOCAL;
784: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosSetReorthog_C","EPSLanczosSetReorthog_LANCZOS",EPSLanczosSetReorthog_LANCZOS);
785: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosGetReorthog_C","EPSLanczosGetReorthog_LANCZOS",EPSLanczosGetReorthog_LANCZOS);
786: return(0);
787: }