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-2012, Universitat Politecnica de Valencia, Spain
24: This file is part of SLEPc.
25:
26: SLEPc is free software: you can redistribute it and/or modify it under the
27: terms of version 3 of the GNU Lesser General Public License as published by
28: the Free Software Foundation.
30: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
31: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
32: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
33: more details.
35: You should have received a copy of the GNU Lesser General Public License
36: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
37: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
38: */
40: #include <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(((PetscObject)eps)->comm,1,"The value of ncv must be at least nev");
60: }
61: else if (eps->mpd) { /* mpd set */
62: eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd);
63: }
64: else { /* neither set: defaults depend on nev being small or large */
65: if (eps->nev<500) eps->ncv = PetscMin(eps->n,PetscMax(2*eps->nev,eps->nev+15));
66: else { eps->mpd = 500; eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd); }
67: }
68: if (!eps->mpd) eps->mpd = eps->ncv;
69: if (eps->ncv>eps->nev+eps->mpd) SETERRQ(((PetscObject)eps)->comm,1,"The value of ncv must not be larger than nev+mpd");
70: if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
72: if (!eps->which) { EPSDefaultSetWhich(eps); }
73: switch (eps->which) {
74: case EPS_LARGEST_IMAGINARY:
75: case EPS_SMALLEST_IMAGINARY:
76: case EPS_TARGET_IMAGINARY:
77: SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
78: default: ; /* default case to remove warning */
79: }
80: if (!eps->ishermitian) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
81: if (!eps->extraction) {
82: EPSSetExtraction(eps,EPS_RITZ);
83: } else if (eps->extraction!=EPS_RITZ) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
84: if (eps->arbit_func) SETERRQ(((PetscObject)eps)->comm,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: }
90: DSSetType(eps->ds,DSHEP);
91: DSSetCompact(eps->ds,PETSC_TRUE);
92: DSAllocate(eps->ds,eps->ncv);
93: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
94: EPSDefaultGetWork(eps,2);
95: } else {
96: EPSDefaultGetWork(eps,1);
97: }
99: /* dispatch solve method */
100: if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
101: eps->ops->solve = EPSSolve_Lanczos;
102: return(0);
103: }
107: /*
108: EPSLocalLanczos - Local reorthogonalization.
110: This is the simplest variant. At each Lanczos step, the corresponding Lanczos vector
111: is orthogonalized with respect to the two previous Lanczos vectors, according to
112: the three term Lanczos recurrence. WARNING: This variant does not track the loss of
113: orthogonality that occurs in finite-precision arithmetic and, therefore, the
114: generated vectors are not guaranteed to be (semi-)orthogonal.
115: */
116: static PetscErrorCode EPSLocalLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscBool *breakdown)
117: {
119: PetscInt i,j,m = *M;
120: PetscReal norm;
121: PetscBool *which,lwhich[100];
122: PetscScalar *hwork,lhwork[100];
123:
125: if (m > 100) {
126: PetscMalloc(sizeof(PetscBool)*m,&which);
127: PetscMalloc(m*sizeof(PetscScalar),&hwork);
128: } else {
129: which = lwhich;
130: hwork = lhwork;
131: }
132: for (i=0;i<k;i++)
133: which[i] = PETSC_TRUE;
135: for (j=k;j<m-1;j++) {
136: STApply(eps->OP,V[j],V[j+1]);
137: which[j] = PETSC_TRUE;
138: if (j-2>=k) which[j-2] = PETSC_FALSE;
139: IPOrthogonalize(eps->ip,eps->nds,eps->defl,j+1,which,V,V[j+1],hwork,&norm,breakdown);
140: alpha[j] = PetscRealPart(hwork[j]);
141: beta[j] = norm;
142: if (*breakdown) {
143: *M = j+1;
144: if (m > 100) {
145: PetscFree(which);
146: PetscFree(hwork);
147: }
148: return(0);
149: } else {
150: VecScale(V[j+1],1.0/norm);
151: }
152: }
153: STApply(eps->OP,V[m-1],f);
154: IPOrthogonalize(eps->ip,eps->nds,eps->defl,m,PETSC_NULL,V,f,hwork,&norm,PETSC_NULL);
155: alpha[m-1] = PetscRealPart(hwork[m-1]);
156: beta[m-1] = norm;
158: if (m > 100) {
159: PetscFree(which);
160: PetscFree(hwork);
161: }
162: return(0);
163: }
167: /*
168: DenseTridiagonal - Solves a real tridiagonal Hermitian Eigenvalue Problem.
170: Input Parameters:
171: + n - dimension of the eigenproblem
172: . D - pointer to the array containing the diagonal elements
173: - E - pointer to the array containing the off-diagonal elements
175: Output Parameters:
176: + w - pointer to the array to store the computed eigenvalues
177: - V - pointer to the array to store the eigenvectors
179: Notes:
180: If V is PETSC_NULL then the eigenvectors are not computed.
182: This routine use LAPACK routines xSTEVR.
184: */
185: static PetscErrorCode DenseTridiagonal(PetscInt n_,PetscReal *D,PetscReal *E,PetscReal *w,PetscScalar *V)
186: {
187: #if defined(SLEPC_MISSING_LAPACK_STEVR)
189: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"STEVR - Lapack routine is unavailable");
190: #else
192: PetscReal abstol = 0.0,vl,vu,*work;
193: PetscBLASInt il,iu,m,*isuppz,n,lwork,*iwork,liwork,info;
194: const char *jobz;
195: #if defined(PETSC_USE_COMPLEX)
196: PetscInt i,j;
197: PetscReal *VV;
198: #endif
199:
201: n = PetscBLASIntCast(n_);
202: lwork = PetscBLASIntCast(20*n_);
203: liwork = PetscBLASIntCast(10*n_);
204: if (V) {
205: jobz = "V";
206: #if defined(PETSC_USE_COMPLEX)
207: PetscMalloc(n*n*sizeof(PetscReal),&VV);
208: #endif
209: } else jobz = "N";
210: PetscMalloc(2*n*sizeof(PetscBLASInt),&isuppz);
211: PetscMalloc(lwork*sizeof(PetscReal),&work);
212: PetscMalloc(liwork*sizeof(PetscBLASInt),&iwork);
213: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
214: #if defined(PETSC_USE_COMPLEX)
215: LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,VV,&n,isuppz,work,&lwork,iwork,&liwork,&info);
216: #else
217: LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info);
218: #endif
219: PetscFPTrapPop();
220: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack DSTEVR %d",info);
221: #if defined(PETSC_USE_COMPLEX)
222: if (V) {
223: for (i=0;i<n;i++)
224: for (j=0;j<n;j++)
225: V[i*n+j] = VV[i*n+j];
226: PetscFree(VV);
227: }
228: #endif
229: PetscFree(isuppz);
230: PetscFree(work);
231: PetscFree(iwork);
232: return(0);
233: #endif
234: }
238: /*
239: EPSSelectiveLanczos - Selective reorthogonalization.
240: */
241: static PetscErrorCode EPSSelectiveLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscBool *breakdown,PetscReal anorm)
242: {
244: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
245: PetscInt i,j,m = *M,n,nritz=0,nritzo;
246: PetscReal *d,*e,*ritz,norm;
247: PetscScalar *Y,*hwork,lhwork[100];
248: PetscBool *which,lwhich[100];
251: PetscMalloc(m*sizeof(PetscReal),&d);
252: PetscMalloc(m*sizeof(PetscReal),&e);
253: PetscMalloc(m*sizeof(PetscReal),&ritz);
254: PetscMalloc(m*m*sizeof(PetscScalar),&Y);
255: if (m > 100) {
256: PetscMalloc(sizeof(PetscBool)*m,&which);
257: PetscMalloc(m*sizeof(PetscScalar),&hwork);
258: } else {
259: which = lwhich;
260: hwork = lhwork;
261: }
262: for (i=0;i<k;i++)
263: which[i] = PETSC_TRUE;
265: for (j=k;j<m;j++) {
266: /* Lanczos step */
267: STApply(eps->OP,V[j],f);
268: which[j] = PETSC_TRUE;
269: if (j-2>=k) which[j-2] = PETSC_FALSE;
270: IPOrthogonalize(eps->ip,eps->nds,eps->defl,j+1,which,V,f,hwork,&norm,breakdown);
271: alpha[j] = PetscRealPart(hwork[j]);
272: beta[j] = norm;
273: if (*breakdown) {
274: *M = j+1;
275: break;
276: }
278: /* Compute eigenvalues and eigenvectors Y of the tridiagonal block */
279: n = j-k+1;
280: for (i=0;i<n;i++) { d[i] = alpha[i+k]; e[i] = beta[i+k]; }
281: DenseTridiagonal(n,d,e,ritz,Y);
282:
283: /* Estimate ||A|| */
284: for (i=0;i<n;i++)
285: if (PetscAbsReal(ritz[i]) > anorm) anorm = PetscAbsReal(ritz[i]);
287: /* Compute nearly converged Ritz vectors */
288: nritzo = 0;
289: for (i=0;i<n;i++)
290: if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm)
291: nritzo++;
293: if (nritzo>nritz) {
294: nritz = 0;
295: for (i=0;i<n;i++) {
296: if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm) {
297: SlepcVecMAXPBY(lanczos->AV[nritz],0.0,1.0,n,Y+i*n,V+k);
298: nritz++;
299: }
300: }
301: }
303: if (nritz > 0) {
304: IPOrthogonalize(eps->ip,0,PETSC_NULL,nritz,PETSC_NULL,lanczos->AV,f,hwork,&norm,breakdown);
305: if (*breakdown) {
306: *M = j+1;
307: break;
308: }
309: }
310:
311: if (j<m-1) {
312: VecScale(f,1.0 / norm);
313: VecCopy(f,V[j+1]);
314: }
315: }
316:
317: PetscFree(d);
318: PetscFree(e);
319: PetscFree(ritz);
320: PetscFree(Y);
321: if (m > 100) {
322: PetscFree(which);
323: PetscFree(hwork);
324: }
325: return(0);
326: }
330: static void update_omega(PetscReal *omega,PetscReal *omega_old,PetscInt j,PetscReal *alpha,PetscReal *beta,PetscReal eps1,PetscReal anorm)
331: {
332: PetscInt k;
333: PetscReal T,binv;
336: /* Estimate of contribution to roundoff errors from A*v
337: fl(A*v) = A*v + f,
338: where ||f|| \approx eps1*||A||.
339: For a full matrix A, a rule-of-thumb estimate is eps1 = sqrt(n)*eps. */
340: T = eps1*anorm;
341: binv = 1.0/beta[j+1];
343: /* Update omega(1) using omega(0)==0. */
344: omega_old[0]= beta[1]*omega[1] + (alpha[0]-alpha[j])*omega[0] -
345: beta[j]*omega_old[0];
346: if (omega_old[0] > 0)
347: omega_old[0] = binv*(omega_old[0] + T);
348: else
349: omega_old[0] = binv*(omega_old[0] - T);
350:
351: /* Update remaining components. */
352: for (k=1;k<j-1;k++) {
353: omega_old[k] = beta[k+1]*omega[k+1] + (alpha[k]-alpha[j])*omega[k] +
354: beta[k]*omega[k-1] - beta[j]*omega_old[k];
355: if (omega_old[k] > 0)
356: omega_old[k] = binv*(omega_old[k] + T);
357: else
358: omega_old[k] = binv*(omega_old[k] - T);
359: }
360: omega_old[j-1] = binv*T;
361:
362: /* Swap omega and omega_old. */
363: for (k=0;k<j;k++) {
364: omega[k] = omega_old[k];
365: omega_old[k] = omega[k];
366: }
367: omega[j] = eps1;
368: PetscFunctionReturnVoid();
369: }
373: static void compute_int(PetscBool *which,PetscReal *mu,PetscInt j,PetscReal delta,PetscReal eta)
374: {
375: PetscInt i,k,maxpos;
376: PetscReal max;
377: PetscBool found;
378:
380: /* initialize which */
381: found = PETSC_FALSE;
382: maxpos = 0;
383: max = 0.0;
384: for (i=0;i<j;i++) {
385: if (PetscAbsReal(mu[i]) >= delta) {
386: which[i] = PETSC_TRUE;
387: found = PETSC_TRUE;
388: } else which[i] = PETSC_FALSE;
389: if (PetscAbsReal(mu[i]) > max) {
390: maxpos = i;
391: max = PetscAbsReal(mu[i]);
392: }
393: }
394: if (!found) which[maxpos] = PETSC_TRUE;
395:
396: for (i=0;i<j;i++)
397: if (which[i]) {
398: /* find left interval */
399: for (k=i;k>=0;k--) {
400: if (PetscAbsReal(mu[k])<eta || which[k]) break;
401: else which[k] = PETSC_TRUE;
402: }
403: /* find right interval */
404: for (k=i+1;k<j;k++) {
405: if (PetscAbsReal(mu[k])<eta || which[k]) break;
406: else which[k] = PETSC_TRUE;
407: }
408: }
409: PetscFunctionReturnVoid();
410: }
414: /*
415: EPSPartialLanczos - Partial reorthogonalization.
416: */
417: static PetscErrorCode EPSPartialLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscBool *breakdown,PetscReal anorm)
418: {
419: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
421: PetscInt i,j,m = *M;
422: PetscReal norm,*omega,lomega[100],*omega_old,lomega_old[100],eps1,delta,eta;
423: PetscBool *which,lwhich[100],*which2,lwhich2[100],
424: reorth = PETSC_FALSE,force_reorth = PETSC_FALSE,
425: fro = PETSC_FALSE,estimate_anorm = PETSC_FALSE;
426: PetscScalar *hwork,lhwork[100];
429: if (m>100) {
430: PetscMalloc(m*sizeof(PetscReal),&omega);
431: PetscMalloc(m*sizeof(PetscReal),&omega_old);
432: } else {
433: omega = lomega;
434: omega_old = lomega_old;
435: }
436: if (m > 100) {
437: PetscMalloc(sizeof(PetscBool)*m,&which);
438: PetscMalloc(sizeof(PetscBool)*m,&which2);
439: PetscMalloc(m*sizeof(PetscScalar),&hwork);
440: } else {
441: which = lwhich;
442: which2 = lwhich2;
443: hwork = lhwork;
444: }
446: eps1 = PetscSqrtReal((PetscReal)eps->n)*PETSC_MACHINE_EPSILON/2;
447: delta = PETSC_SQRT_MACHINE_EPSILON/PetscSqrtReal((PetscReal)eps->ncv);
448: eta = pow(PETSC_MACHINE_EPSILON,3.0/4.0)/PetscSqrtReal((PetscReal)eps->ncv);
449: if (anorm < 0.0) {
450: anorm = 1.0;
451: estimate_anorm = PETSC_TRUE;
452: }
453: for (i=0;i<m-k;i++)
454: omega[i] = omega_old[i] = 0.0;
455: for (i=0;i<k;i++)
456: which[i] = PETSC_TRUE;
457:
458: for (j=k;j<m;j++) {
459: STApply(eps->OP,V[j],f);
460: if (fro) {
461: /* Lanczos step with full reorthogonalization */
462: IPOrthogonalize(eps->ip,eps->nds,eps->defl,j+1,PETSC_NULL,V,f,hwork,&norm,breakdown);
463: alpha[j] = PetscRealPart(hwork[j]);
464: } else {
465: /* Lanczos step */
466: which[j] = PETSC_TRUE;
467: if (j-2>=k) which[j-2] = PETSC_FALSE;
468: IPOrthogonalize(eps->ip,eps->nds,eps->defl,j+1,which,V,f,hwork,&norm,breakdown);
469: alpha[j] = PetscRealPart(hwork[j]);
470: beta[j] = norm;
471:
472: /* Estimate ||A|| if needed */
473: if (estimate_anorm) {
474: if (j>k) anorm = PetscMax(anorm,PetscAbsReal(alpha[j])+norm+beta[j-1]);
475: else anorm = PetscMax(anorm,PetscAbsReal(alpha[j])+norm);
476: }
478: /* Check if reorthogonalization is needed */
479: reorth = PETSC_FALSE;
480: if (j>k) {
481: update_omega(omega,omega_old,j,alpha,beta-1,eps1,anorm);
482: for (i=0;i<j-k;i++)
483: if (PetscAbsScalar(omega[i]) > delta) reorth = PETSC_TRUE;
484: }
486: if (reorth || force_reorth) {
487: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_PERIODIC) {
488: /* Periodic reorthogonalization */
489: if (force_reorth) force_reorth = PETSC_FALSE;
490: else force_reorth = PETSC_TRUE;
491: IPOrthogonalize(eps->ip,0,PETSC_NULL,j-k,PETSC_NULL,V+k,f,hwork,&norm,breakdown);
492: for (i=0;i<j-k;i++)
493: omega[i] = eps1;
494: } else {
495: /* Partial reorthogonalization */
496: if (force_reorth) force_reorth = PETSC_FALSE;
497: else {
498: force_reorth = PETSC_TRUE;
499: compute_int(which2,omega,j-k,delta,eta);
500: for (i=0;i<j-k;i++)
501: if (which2[i]) omega[i] = eps1;
502: }
503: IPOrthogonalize(eps->ip,0,PETSC_NULL,j-k,which2,V+k,f,hwork,&norm,breakdown);
504: }
505: }
506: }
507:
508: if (*breakdown || norm < eps->n*anorm*PETSC_MACHINE_EPSILON) {
509: *M = j+1;
510: break;
511: }
512: if (!fro && norm*delta < anorm*eps1) {
513: fro = PETSC_TRUE;
514: PetscInfo1(eps,"Switching to full reorthogonalization at iteration %D\n",eps->its);
515: }
516:
517: beta[j] = norm;
518: if (j<m-1) {
519: VecScale(f,1.0/norm);
520: VecCopy(f,V[j+1]);
521: }
522: }
524: if (m>100) {
525: PetscFree(omega);
526: PetscFree(omega_old);
527: PetscFree(which);
528: PetscFree(which2);
529: PetscFree(hwork);
530: }
531: return(0);
532: }
536: /*
537: EPSBasicLanczos - Computes an m-step Lanczos factorization. The first k
538: columns are assumed to be locked and therefore they are not modified. On
539: exit, the following relation is satisfied:
541: OP * V - V * T = f * e_m^T
543: where the columns of V are the Lanczos vectors, T is a tridiagonal matrix,
544: f is the residual vector and e_m is the m-th vector of the canonical basis.
545: The Lanczos vectors (together with vector f) are B-orthogonal (to working
546: accuracy) if full reorthogonalization is being used, otherwise they are
547: (B-)semi-orthogonal. On exit, beta contains the B-norm of f and the next
548: Lanczos vector can be computed as v_{m+1} = f / beta.
550: This function simply calls another function which depends on the selected
551: reorthogonalization strategy.
552: */
553: static PetscErrorCode EPSBasicLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *m,Vec f,PetscBool *breakdown,PetscReal anorm)
554: {
555: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
556: PetscScalar *T;
557: PetscInt i,n=*m;
558: PetscReal betam;
559: PetscErrorCode ierr;
560: IPOrthogRefineType orthog_ref;
563: switch (lanczos->reorthog) {
564: case EPS_LANCZOS_REORTHOG_LOCAL:
565: EPSLocalLanczos(eps,alpha,beta,V,k,m,f,breakdown);
566: break;
567: case EPS_LANCZOS_REORTHOG_SELECTIVE:
568: EPSSelectiveLanczos(eps,alpha,beta,V,k,m,f,breakdown,anorm);
569: break;
570: case EPS_LANCZOS_REORTHOG_FULL:
571: EPSFullLanczos(eps,alpha,beta,V,k,m,f,breakdown);
572: break;
573: case EPS_LANCZOS_REORTHOG_PARTIAL:
574: case EPS_LANCZOS_REORTHOG_PERIODIC:
575: EPSPartialLanczos(eps,alpha,beta,V,k,m,f,breakdown,anorm);
576: break;
577: case EPS_LANCZOS_REORTHOG_DELAYED:
578: PetscMalloc(n*n*sizeof(PetscScalar),&T);
579: IPGetOrthogonalization(eps->ip,PETSC_NULL,&orthog_ref,PETSC_NULL);
580: if (orthog_ref == IP_ORTHOG_REFINE_NEVER) {
581: EPSDelayedArnoldi1(eps,T,n,V,k,m,f,&betam,breakdown);
582: } else {
583: EPSDelayedArnoldi(eps,T,n,V,k,m,f,&betam,breakdown);
584: }
585: for (i=k;i<n-1;i++) { alpha[i] = PetscRealPart(T[n*i+i]); beta[i] = PetscRealPart(T[n*i+i+1]); }
586: alpha[n-1] = PetscRealPart(T[n*(n-1)+n-1]);
587: beta[n-1] = betam;
588: PetscFree(T);
589: break;
590: default:
591: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
592: }
593: return(0);
594: }
598: PetscErrorCode EPSSolve_Lanczos(EPS eps)
599: {
600: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
602: PetscInt nconv,i,j,k,l,x,n,*perm,restart,ncv=eps->ncv,r,ld;
603: Vec w=eps->work[1],f=eps->work[0];
604: PetscScalar *Y,*ritz,stmp;
605: PetscReal *d,*e,*bnd,anorm,beta,norm,rtmp,resnorm;
606: PetscBool breakdown;
607: char *conv,ctmp;
610: DSGetLeadingDimension(eps->ds,&ld);
611: PetscMalloc(ncv*sizeof(PetscScalar),&ritz);
612: PetscMalloc(ncv*sizeof(PetscReal),&bnd);
613: PetscMalloc(ncv*sizeof(PetscInt),&perm);
614: PetscMalloc(ncv*sizeof(char),&conv);
616: /* The first Lanczos vector is the normalized initial vector */
617: EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
618:
619: anorm = -1.0;
620: nconv = 0;
621:
622: /* Restart loop */
623: while (eps->reason == EPS_CONVERGED_ITERATING) {
624: eps->its++;
626: /* Compute an ncv-step Lanczos factorization */
627: n = PetscMin(nconv+eps->mpd,ncv);
628: DSGetArrayReal(eps->ds,DS_MAT_T,&d);
629: e = d + ld;
630: EPSBasicLanczos(eps,d,e,eps->V,nconv,&n,f,&breakdown,anorm);
631: beta = e[n-1];
632: DSRestoreArrayReal(eps->ds,DS_MAT_T,&d);
633: DSSetDimensions(eps->ds,n,PETSC_IGNORE,nconv,0);
634: DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
636: /* Solve projected problem */
637: DSSolve(eps->ds,ritz,PETSC_NULL);
638: DSSort(eps->ds,ritz,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);
639:
640: /* Estimate ||A|| */
641: for (i=nconv;i<n;i++)
642: anorm = PetscMax(anorm,PetscAbsReal(PetscRealPart(ritz[i])));
643:
644: /* Compute residual norm estimates as beta*abs(Y(m,:)) + eps*||A|| */
645: DSGetArray(eps->ds,DS_MAT_Q,&Y);
646: for (i=nconv;i<n;i++) {
647: resnorm = beta*PetscAbsScalar(Y[n-1+i*ld]) + PETSC_MACHINE_EPSILON*anorm;
648: (*eps->conv_func)(eps,ritz[i],eps->eigi[i],resnorm,&bnd[i],eps->conv_ctx);
649: if (bnd[i]<eps->tol) {
650: conv[i] = 'C';
651: } else {
652: conv[i] = 'N';
653: }
654: }
655: DSRestoreArray(eps->ds,DS_MAT_Q,&Y);
657: /* purge repeated ritz values */
658: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL)
659: for (i=nconv+1;i<n;i++)
660: if (conv[i] == 'C')
661: if (PetscAbsScalar((ritz[i]-ritz[i-1])/ritz[i]) < eps->tol)
662: conv[i] = 'R';
664: /* Compute restart vector */
665: if (breakdown) {
666: PetscInfo2(eps,"Breakdown in Lanczos method (it=%D norm=%G)\n",eps->its,beta);
667: } else {
668: restart = nconv;
669: while (restart<n && conv[restart] != 'N') restart++;
670: if (restart >= n) {
671: breakdown = PETSC_TRUE;
672: } else {
673: for (i=restart+1;i<n;i++)
674: if (conv[i] == 'N') {
675: (*eps->which_func)(ritz[restart],0.0,ritz[i],0.0,&r,eps->which_ctx);
676: if (r>0) restart = i;
677: }
678: DSGetArray(eps->ds,DS_MAT_Q,&Y);
679: SlepcVecMAXPBY(f,0.0,1.0,n-nconv,Y+restart*ld+nconv,eps->V+nconv);
680: DSRestoreArray(eps->ds,DS_MAT_Q,&Y);
681: }
682: }
684: /* Count and put converged eigenvalues first */
685: for (i=nconv;i<n;i++) perm[i] = i;
686: for (k=nconv;k<n;k++)
687: if (conv[perm[k]] != 'C') {
688: j = k + 1;
689: while (j<n && conv[perm[j]] != 'C') j++;
690: if (j>=n) break;
691: l = perm[k]; perm[k] = perm[j]; perm[j] = l;
692: }
694: /* Sort eigenvectors according to permutation */
695: DSGetArray(eps->ds,DS_MAT_Q,&Y);
696: for (i=nconv;i<k;i++) {
697: x = perm[i];
698: if (x != i) {
699: j = i + 1;
700: while (perm[j] != i) j++;
701: /* swap eigenvalues i and j */
702: stmp = ritz[x]; ritz[x] = ritz[i]; ritz[i] = stmp;
703: rtmp = bnd[x]; bnd[x] = bnd[i]; bnd[i] = rtmp;
704: ctmp = conv[x]; conv[x] = conv[i]; conv[i] = ctmp;
705: perm[j] = x; perm[i] = i;
706: /* swap eigenvectors i and j */
707: for (l=0;l<n;l++) {
708: stmp = Y[l+x*ld]; Y[l+x*ld] = Y[l+i*ld]; Y[l+i*ld] = stmp;
709: }
710: }
711: }
712:
713: /* compute converged eigenvectors */
714: SlepcUpdateVectors(n,eps->V,nconv,nconv+k,Y,ld,PETSC_FALSE);
715: DSRestoreArray(eps->ds,DS_MAT_Q,&Y);
716:
717: /* purge spurious ritz values */
718: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
719: for (i=nconv;i<k;i++) {
720: VecNorm(eps->V[i],NORM_2,&norm);
721: VecScale(eps->V[i],1.0/norm);
722: STApply(eps->OP,eps->V[i],w);
723: VecAXPY(w,-ritz[i],eps->V[i]);
724: VecNorm(w,NORM_2,&norm);
725: (*eps->conv_func)(eps,ritz[i],eps->eigi[i],norm,&bnd[i],eps->conv_ctx);
726: if (bnd[i]>=eps->tol) conv[i] = 'S';
727: }
728: for (i=nconv;i<k;i++)
729: if (conv[i] != 'C') {
730: j = i + 1;
731: while (j<k && conv[j] != 'C') j++;
732: if (j>=k) break;
733: /* swap eigenvalues i and j */
734: stmp = ritz[j]; ritz[j] = ritz[i]; ritz[i] = stmp;
735: rtmp = bnd[j]; bnd[j] = bnd[i]; bnd[i] = rtmp;
736: ctmp = conv[j]; conv[j] = conv[i]; conv[i] = ctmp;
737: /* swap eigenvectors i and j */
738: VecSwap(eps->V[i],eps->V[j]);
739: }
740: k = i;
741: }
742:
743: /* store ritz values and estimated errors */
744: for (i=nconv;i<n;i++) {
745: eps->eigr[i] = ritz[i];
746: eps->errest[i] = bnd[i];
747: }
748: EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,n);
749: nconv = k;
750: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
751: if (nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
752:
753: if (eps->reason == EPS_CONVERGED_ITERATING) { /* copy restart vector */
754: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL && !breakdown) {
755: /* Reorthonormalize restart vector */
756: IPOrthogonalize(eps->ip,eps->nds,eps->defl,nconv,PETSC_NULL,eps->V,f,PETSC_NULL,&norm,&breakdown);
757: VecScale(f,1.0/norm);
758: }
759: if (breakdown) {
760: /* Use random vector for restarting */
761: PetscInfo(eps,"Using random vector for restart\n");
762: EPSGetStartVector(eps,nconv,f,&breakdown);
763: }
764: if (breakdown) { /* give up */
765: eps->reason = EPS_DIVERGED_BREAKDOWN;
766: PetscInfo(eps,"Unable to generate more start vectors\n");
767: } else {
768: VecCopy(f,eps->V[nconv]);
769: }
770: }
771: }
772:
773: eps->nconv = nconv;
775: PetscFree(ritz);
776: PetscFree(bnd);
777: PetscFree(perm);
778: PetscFree(conv);
779: return(0);
780: }
784: PetscErrorCode EPSSetFromOptions_Lanczos(EPS eps)
785: {
786: PetscErrorCode ierr;
787: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
788: PetscBool flg;
789: EPSLanczosReorthogType reorthog;
792: PetscOptionsHead("EPS Lanczos Options");
793: PetscOptionsEnum("-eps_lanczos_reorthog","Lanczos reorthogonalization","EPSLanczosSetReorthog",EPSLanczosReorthogTypes,(PetscEnum)lanczos->reorthog,(PetscEnum*)&reorthog,&flg);
794: if (flg) { EPSLanczosSetReorthog(eps,reorthog); }
795: PetscOptionsTail();
796: return(0);
797: }
799: EXTERN_C_BEGIN
802: PetscErrorCode EPSLanczosSetReorthog_Lanczos(EPS eps,EPSLanczosReorthogType reorthog)
803: {
804: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
807: switch (reorthog) {
808: case EPS_LANCZOS_REORTHOG_LOCAL:
809: case EPS_LANCZOS_REORTHOG_FULL:
810: case EPS_LANCZOS_REORTHOG_DELAYED:
811: case EPS_LANCZOS_REORTHOG_SELECTIVE:
812: case EPS_LANCZOS_REORTHOG_PERIODIC:
813: case EPS_LANCZOS_REORTHOG_PARTIAL:
814: lanczos->reorthog = reorthog;
815: break;
816: default:
817: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
818: }
819: return(0);
820: }
821: EXTERN_C_END
825: /*@
826: EPSLanczosSetReorthog - Sets the type of reorthogonalization used during the Lanczos
827: iteration.
829: Logically Collective on EPS
831: Input Parameters:
832: + eps - the eigenproblem solver context
833: - reorthog - the type of reorthogonalization
835: Options Database Key:
836: . -eps_lanczos_reorthog - Sets the reorthogonalization type (either 'local', 'selective',
837: 'periodic', 'partial', 'full' or 'delayed')
838:
839: Level: advanced
841: .seealso: EPSLanczosGetReorthog(), EPSLanczosReorthogType
842: @*/
843: PetscErrorCode EPSLanczosSetReorthog(EPS eps,EPSLanczosReorthogType reorthog)
844: {
850: PetscTryMethod(eps,"EPSLanczosSetReorthog_C",(EPS,EPSLanczosReorthogType),(eps,reorthog));
851: return(0);
852: }
854: EXTERN_C_BEGIN
857: PetscErrorCode EPSLanczosGetReorthog_Lanczos(EPS eps,EPSLanczosReorthogType *reorthog)
858: {
859: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
862: *reorthog = lanczos->reorthog;
863: return(0);
864: }
865: EXTERN_C_END
869: /*@C
870: EPSLanczosGetReorthog - Gets the type of reorthogonalization used during the Lanczos
871: iteration.
873: Not Collective
875: Input Parameter:
876: . eps - the eigenproblem solver context
878: Output Parameter:
879: . reorthog - the type of reorthogonalization
881: Level: advanced
883: .seealso: EPSLanczosSetReorthog(), EPSLanczosReorthogType
884: @*/
885: PetscErrorCode EPSLanczosGetReorthog(EPS eps,EPSLanczosReorthogType *reorthog)
886: {
892: PetscTryMethod(eps,"EPSLanczosGetReorthog_C",(EPS,EPSLanczosReorthogType*),(eps,reorthog));
893: return(0);
894: }
898: PetscErrorCode EPSReset_Lanczos(EPS eps)
899: {
901: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
904: VecDestroyVecs(eps->ncv,&lanczos->AV);
905: EPSReset_Default(eps);
906: return(0);
907: }
911: PetscErrorCode EPSDestroy_Lanczos(EPS eps)
912: {
916: PetscFree(eps->data);
917: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosSetReorthog_C","",PETSC_NULL);
918: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosGetReorthog_C","",PETSC_NULL);
919: return(0);
920: }
924: PetscErrorCode EPSView_Lanczos(EPS eps,PetscViewer viewer)
925: {
927: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
928: PetscBool isascii;
931: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
932: if (!isascii) SETERRQ1(((PetscObject)eps)->comm,1,"Viewer type %s not supported for EPS Lanczos",((PetscObject)viewer)->type_name);
933: PetscViewerASCIIPrintf(viewer," Lanczos: %s reorthogonalization\n",EPSLanczosReorthogTypes[lanczos->reorthog]);
934: return(0);
935: }
937: EXTERN_C_BEGIN
940: PetscErrorCode EPSCreate_Lanczos(EPS eps)
941: {
945: PetscNewLog(eps,EPS_LANCZOS,&eps->data);
946: eps->ops->setup = EPSSetUp_Lanczos;
947: eps->ops->setfromoptions = EPSSetFromOptions_Lanczos;
948: eps->ops->destroy = EPSDestroy_Lanczos;
949: eps->ops->reset = EPSReset_Lanczos;
950: eps->ops->view = EPSView_Lanczos;
951: eps->ops->backtransform = EPSBackTransform_Default;
952: eps->ops->computevectors = EPSComputeVectors_Hermitian;
953: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosSetReorthog_C","EPSLanczosSetReorthog_Lanczos",EPSLanczosSetReorthog_Lanczos);
954: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosGetReorthog_C","EPSLanczosGetReorthog_Lanczos",EPSLanczosGetReorthog_Lanczos);
955: return(0);
956: }
957: EXTERN_C_END