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-2009, Universidad 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 private/epsimpl.h
41: #include slepcblaslapack.h
43: typedef struct {
44: EPSLanczosReorthogType reorthog;
45: } EPS_LANCZOS;
49: PetscErrorCode EPSSetUp_LANCZOS(EPS eps)
50: {
51: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
53: PetscInt N;
56: VecGetSize(eps->vec_initial,&N);
57: if (eps->ncv) { /* ncv set */
58: if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
59: }
60: else if (eps->mpd) { /* mpd set */
61: eps->ncv = PetscMin(N,eps->nev+eps->mpd);
62: }
63: else { /* neither set: defaults depend on nev being small or large */
64: if (eps->nev<500) eps->ncv = PetscMin(N,PetscMax(2*eps->nev,eps->nev+15));
65: else { eps->mpd = 500; eps->ncv = PetscMin(N,eps->nev+eps->mpd); }
66: }
67: if (!eps->mpd) eps->mpd = eps->ncv;
68: if (eps->ncv>eps->nev+eps->mpd) SETERRQ(1,"The value of ncv must not be larger than nev+mpd");
69: if (!eps->max_it) eps->max_it = PetscMax(100,2*N/eps->ncv);
71: if (eps->solverclass==EPS_ONE_SIDE) {
72: if (eps->which == EPS_LARGEST_IMAGINARY || eps->which == EPS_SMALLEST_IMAGINARY)
73: SETERRQ(1,"Wrong value of eps->which");
74: if (!eps->ishermitian)
75: SETERRQ(PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
76: } else {
77: if (eps->which != EPS_LARGEST_MAGNITUDE)
78: SETERRQ(1,"Wrong value of eps->which");
79: }
80: if (!eps->extraction) {
81: EPSSetExtraction(eps,EPS_RITZ);
82: } else if (eps->extraction!=EPS_RITZ) {
83: SETERRQ(PETSC_ERR_SUP,"Unsupported extraction type\n");
84: }
86: EPSAllocateSolution(eps);
87: if (lanczos->reorthog == EPSLANCZOS_REORTHOG_SELECTIVE) {
88: VecDuplicateVecs(eps->vec_initial,eps->ncv,&eps->AV);
89: }
90: if (eps->solverclass==EPS_TWO_SIDE) {
91: PetscFree(eps->Tl);
92: PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->Tl);
93: }
94: EPSDefaultGetWork(eps,2);
95: return(0);
96: }
100: /*
101: EPSFullLanczos - Full reorthogonalization.
103: This is the simplest solution to the loss of orthogonality.
104: At each Lanczos step, the corresponding Lanczos vector is
105: orthogonalized with respect to all previous Lanczos vectors.
106: */
107: PetscErrorCode EPSFullLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscTruth *breakdown)
108: {
110: PetscInt j,m = *M;
111: PetscReal norm;
112: PetscScalar *swork,*hwork,lhwork[100];
115: if (eps->nds+m > 100) {
116: PetscMalloc((eps->nds+m)*sizeof(PetscScalar),&swork);
117: PetscMalloc((eps->nds+m)*sizeof(PetscScalar),&hwork);
118: } else {
119: swork = PETSC_NULL;
120: hwork = lhwork;
121: }
123: for (j=k;j<m-1;j++) {
124: STApply(eps->OP,V[j],V[j+1]);
125: IPOrthogonalize(eps->ip,eps->nds+j+1,PETSC_NULL,eps->DSV,V[j+1],hwork,&norm,breakdown,eps->work[0],swork);
126: alpha[j-k] = PetscRealPart(hwork[eps->nds+j]);
127: beta[j-k] = norm;
128: if (*breakdown) {
129: *M = j+1;
130: if (eps->nds+m > 100) {
131: PetscFree(swork);
132: PetscFree(hwork);
133: }
134: return(0);
135: } else {
136: VecScale(V[j+1],1.0/norm);
137: }
138: }
139: STApply(eps->OP,V[m-1],f);
140: IPOrthogonalize(eps->ip,eps->nds+m,PETSC_NULL,eps->DSV,f,hwork,&norm,PETSC_NULL,eps->work[0],swork);
141: alpha[m-1-k] = PetscRealPart(hwork[eps->nds+m-1]);
142: beta[m-1-k] = norm;
143:
144: if (eps->nds+m > 100) {
145: PetscFree(swork);
146: PetscFree(hwork);
147: }
148: return(0);
149: }
153: /*
154: EPSLocalLanczos - Local reorthogonalization.
156: This is the simplest variant. At each Lanczos step, the corresponding Lanczos vector
157: is orthogonalized with respect to the two previous Lanczos vectors, according to
158: the three term Lanczos recurrence. WARNING: This variant does not track the loss of
159: orthogonality that occurs in finite-precision arithmetic and, therefore, the
160: generated vectors are not guaranteed to be (semi-)orthogonal.
161: */
162: static PetscErrorCode EPSLocalLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscTruth *breakdown)
163: {
165: PetscInt i,j,m = *M;
166: PetscReal norm;
167: PetscTruth *which,lwhich[100];
168: PetscScalar *swork,*hwork,lhwork[100];
169:
171: if (eps->nds+m > 100) {
172: PetscMalloc(sizeof(PetscTruth)*(eps->nds+m),&which);
173: PetscMalloc((eps->nds+m)*sizeof(PetscScalar),&swork);
174: PetscMalloc((eps->nds+m)*sizeof(PetscScalar),&hwork);
175: } else {
176: which = lwhich;
177: swork = PETSC_NULL;
178: hwork = lhwork;
179: }
180: for (i=0;i<eps->nds+k;i++)
181: which[i] = PETSC_TRUE;
183: for (j=k;j<m-1;j++) {
184: STApply(eps->OP,V[j],V[j+1]);
185: which[eps->nds+j] = PETSC_TRUE;
186: if (j-2>=k) which[eps->nds+j-2] = PETSC_FALSE;
187: IPOrthogonalize(eps->ip,eps->nds+j+1,which,eps->DSV,V[j+1],hwork,&norm,breakdown,eps->work[0],swork);
188: alpha[j-k] = PetscRealPart(hwork[eps->nds+j]);
189: beta[j-k] = norm;
190: if (*breakdown) {
191: *M = j+1;
192: if (eps->nds+m > 100) {
193: PetscFree(which);
194: PetscFree(swork);
195: PetscFree(hwork);
196: }
197: return(0);
198: } else {
199: VecScale(V[j+1],1.0/norm);
200: }
201: }
202: STApply(eps->OP,V[m-1],f);
203: IPOrthogonalize(eps->ip,eps->nds+m,PETSC_NULL,eps->DSV,f,hwork,&norm,PETSC_NULL,eps->work[0],swork);
204: alpha[m-1-k] = PetscRealPart(hwork[eps->nds+m-1]);
205: beta[m-1-k] = norm;
207: if (eps->nds+m > 100) {
208: PetscFree(which);
209: PetscFree(swork);
210: PetscFree(hwork);
211: }
212: return(0);
213: }
217: /*
218: EPSSelectiveLanczos - Selective reorthogonalization.
219: */
220: static PetscErrorCode EPSSelectiveLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscTruth *breakdown,PetscReal anorm)
221: {
223: PetscInt i,j,m = *M,n,nritz=0,nritzo;
224: PetscReal *d,*e,*ritz,norm;
225: PetscScalar *Y,*swork,*hwork,lhwork[100];
226: PetscTruth *which,lwhich[100];
229: PetscMalloc(m*sizeof(PetscReal),&d);
230: PetscMalloc(m*sizeof(PetscReal),&e);
231: PetscMalloc(m*sizeof(PetscReal),&ritz);
232: PetscMalloc(m*m*sizeof(PetscScalar),&Y);
233: if (eps->nds+m > 100) {
234: PetscMalloc(sizeof(PetscTruth)*(eps->nds+m),&which);
235: PetscMalloc((eps->nds+m)*sizeof(PetscScalar),&swork);
236: PetscMalloc((eps->nds+m)*sizeof(PetscScalar),&hwork);
237: } else {
238: which = lwhich;
239: swork = PETSC_NULL;
240: hwork = lhwork;
241: }
242: for (i=0;i<eps->nds+k;i++)
243: which[i] = PETSC_TRUE;
245: for (j=k;j<m;j++) {
246: /* Lanczos step */
247: STApply(eps->OP,V[j],f);
248: which[eps->nds+j] = PETSC_TRUE;
249: if (j-2>=k) which[eps->nds+j-2] = PETSC_FALSE;
250: IPOrthogonalize(eps->ip,eps->nds+j+1,which,eps->DSV,f,hwork,&norm,breakdown,eps->work[0],swork);
251: alpha[j-k] = PetscRealPart(hwork[eps->nds+j]);
252: beta[j-k] = norm;
253: if (*breakdown) {
254: *M = j+1;
255: break;
256: }
258: /* Compute eigenvalues and eigenvectors Y of the tridiagonal block */
259: n = j-k+1;
260: for (i=0;i<n;i++) { d[i] = alpha[i]; e[i] = beta[i]; }
261: EPSDenseTridiagonal(n,d,e,ritz,Y);
262:
263: /* Estimate ||A|| */
264: for (i=0;i<n;i++)
265: if (PetscAbsReal(ritz[i]) > anorm) anorm = PetscAbsReal(ritz[i]);
267: /* Compute nearly converged Ritz vectors */
268: nritzo = 0;
269: for (i=0;i<n;i++)
270: if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm)
271: nritzo++;
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: VecSet(eps->AV[nritz],0.0);
278: VecMAXPY(eps->AV[nritz],n,Y+i*n,V+k);
279: nritz++;
280: }
281: }
282: }
284: if (nritz > 0) {
285: IPOrthogonalize(eps->ip,nritz,PETSC_NULL,eps->AV,f,hwork,&norm,breakdown,eps->work[0],swork);
286: if (*breakdown) {
287: *M = j+1;
288: break;
289: }
290: }
291:
292: if (j<m-1) {
293: VecScale(f,1.0 / norm);
294: VecCopy(f,V[j+1]);
295: }
296: }
297:
298: PetscFree(d);
299: PetscFree(e);
300: PetscFree(ritz);
301: PetscFree(Y);
302: if (eps->nds+m > 100) {
303: PetscFree(which);
304: PetscFree(swork);
305: PetscFree(hwork);
306: }
307: return(0);
308: }
312: static void update_omega(PetscReal *omega,PetscReal *omega_old,PetscInt j,PetscReal *alpha,PetscReal *beta,PetscReal eps1,PetscReal anorm)
313: {
314: PetscInt k;
315: PetscReal T,binv,temp;
318: /* Estimate of contribution to roundoff errors from A*v
319: fl(A*v) = A*v + f,
320: where ||f|| \approx eps1*||A||.
321: For a full matrix A, a rule-of-thumb estimate is eps1 = sqrt(n)*eps. */
322: T = eps1*anorm;
323: binv = 1.0/beta[j+1];
325: /* Update omega(1) using omega(0)==0. */
326: omega_old[0]= beta[1]*omega[1] + (alpha[0]-alpha[j])*omega[0] -
327: beta[j]*omega_old[0];
328: if (omega_old[0] > 0)
329: omega_old[0] = binv*(omega_old[0] + T);
330: else
331: omega_old[0] = binv*(omega_old[0] - T);
332:
333: /* Update remaining components. */
334: for (k=1;k<j-1;k++) {
335: omega_old[k] = beta[k+1]*omega[k+1] + (alpha[k]-alpha[j])*omega[k] +
336: beta[k]*omega[k-1] - beta[j]*omega_old[k];
337: if (omega_old[k] > 0)
338: omega_old[k] = binv*(omega_old[k] + T);
339: else
340: omega_old[k] = binv*(omega_old[k] - T);
341: }
342: omega_old[j-1] = binv*T;
343:
344: /* Swap omega and omega_old. */
345: for (k=0;k<j;k++) {
346: temp = omega[k];
347: omega[k] = omega_old[k];
348: omega_old[k] = omega[k];
349: }
350: omega[j] = eps1;
351: PetscFunctionReturnVoid();
352: }
356: static void compute_int(PetscTruth *which,PetscReal *mu,PetscInt j,PetscReal delta,PetscReal eta)
357: {
358: PetscInt i,k,maxpos;
359: PetscReal max;
360: PetscTruth found;
361:
363: /* initialize which */
364: found = PETSC_FALSE;
365: maxpos = 0;
366: max = 0.0;
367: for (i=0;i<j;i++) {
368: if (PetscAbsReal(mu[i]) >= delta) {
369: which[i] = PETSC_TRUE;
370: found = PETSC_TRUE;
371: } else which[i] = PETSC_FALSE;
372: if (PetscAbsReal(mu[i]) > max) {
373: maxpos = i;
374: max = PetscAbsReal(mu[i]);
375: }
376: }
377: if (!found) which[maxpos] = PETSC_TRUE;
378:
379: for (i=0;i<j;i++)
380: if (which[i]) {
381: /* find left interval */
382: for (k=i;k>=0;k--) {
383: if (PetscAbsReal(mu[k])<eta || which[k]) break;
384: else which[k] = PETSC_TRUE;
385: }
386: /* find right interval */
387: for (k=i+1;k<j;k++) {
388: if (PetscAbsReal(mu[k])<eta || which[k]) break;
389: else which[k] = PETSC_TRUE;
390: }
391: }
392: PetscFunctionReturnVoid();
393: }
397: /*
398: EPSPartialLanczos - Partial reorthogonalization.
399: */
400: static PetscErrorCode EPSPartialLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscTruth *breakdown,PetscReal anorm)
401: {
402: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
404: Mat A;
405: PetscInt i,j,m = *M;
406: PetscInt n;
407: PetscReal norm,*omega,lomega[100],*omega_old,lomega_old[100],eps1,delta,eta;
408: PetscTruth *which,lwhich[100],*which2,lwhich2[100],
409: reorth = PETSC_FALSE,force_reorth = PETSC_FALSE,fro = PETSC_FALSE,estimate_anorm = PETSC_FALSE;
410: PetscScalar *swork,*hwork,lhwork[100];
413: if (m>100) {
414: PetscMalloc(m*sizeof(PetscReal),&omega);
415: PetscMalloc(m*sizeof(PetscReal),&omega_old);
416: } else {
417: omega = lomega;
418: omega_old = lomega_old;
419: }
420: if (eps->nds+m > 100) {
421: PetscMalloc(sizeof(PetscTruth)*(eps->nds+m),&which);
422: PetscMalloc(sizeof(PetscTruth)*(eps->nds+m),&which2);
423: PetscMalloc((eps->nds+m)*sizeof(PetscScalar),&swork);
424: PetscMalloc((eps->nds+m)*sizeof(PetscScalar),&hwork);
425: } else {
426: which = lwhich;
427: which2 = lwhich2;
428: swork = PETSC_NULL;
429: hwork = lhwork;
430: }
432: STGetOperators(eps->OP,&A,PETSC_NULL);
433: MatGetSize(A,&n,PETSC_NULL);
434: eps1 = sqrt((PetscReal)n)*PETSC_MACHINE_EPSILON/2;
435: delta = PETSC_SQRT_MACHINE_EPSILON/sqrt((PetscReal)eps->ncv);
436: eta = pow(PETSC_MACHINE_EPSILON,3.0/4.0)/sqrt((PetscReal)eps->ncv);
437: if (anorm < 0.0) {
438: anorm = 1.0;
439: estimate_anorm = PETSC_TRUE;
440: }
441: for (i=0;i<m-k;i++)
442: omega[i] = omega_old[i] = 0.0;
443: for (i=0;i<eps->nds+k;i++)
444: which[i] = PETSC_TRUE;
445:
446: for (j=k;j<m;j++) {
447: STApply(eps->OP,V[j],f);
448: if (fro) {
449: /* Lanczos step with full reorthogonalization */
450: IPOrthogonalize(eps->ip,eps->nds+j+1,PETSC_NULL,eps->DSV,f,hwork,&norm,breakdown,eps->work[0],swork);
451: alpha[j-k] = PetscRealPart(hwork[eps->nds+j]);
452: } else {
453: /* Lanczos step */
454: which[eps->nds+j] = PETSC_TRUE;
455: if (j-2>=k) which[eps->nds+j-2] = PETSC_FALSE;
456: IPOrthogonalize(eps->ip,eps->nds+j+1,which,eps->DSV,f,hwork,&norm,breakdown,eps->work[0],swork);
457: alpha[j-k] = PetscRealPart(hwork[eps->nds+j]);
458: beta[j-k] = norm;
459:
460: /* Estimate ||A|| if needed */
461: if (estimate_anorm) {
462: if (j>k) anorm = PetscMax(anorm,PetscAbsReal(alpha[j-k])+norm+beta[j-k-1]);
463: else anorm = PetscMax(anorm,PetscAbsReal(alpha[j-k])+norm);
464: }
466: /* Check if reorthogonalization is needed */
467: reorth = PETSC_FALSE;
468: if (j>k) {
469: update_omega(omega,omega_old,j-k,alpha,beta-1,eps1,anorm);
470: for (i=0;i<j-k;i++)
471: if (PetscAbsScalar(omega[i]) > delta) reorth = PETSC_TRUE;
472: }
474: if (reorth || force_reorth) {
475: if (lanczos->reorthog == EPSLANCZOS_REORTHOG_PERIODIC) {
476: /* Periodic reorthogonalization */
477: if (force_reorth) force_reorth = PETSC_FALSE;
478: else force_reorth = PETSC_TRUE;
479: IPOrthogonalize(eps->ip,j-k,PETSC_NULL,V+k,f,hwork,&norm,breakdown,eps->work[0],swork);
480: for (i=0;i<j-k;i++)
481: omega[i] = eps1;
482: } else {
483: /* Partial reorthogonalization */
484: if (force_reorth) force_reorth = PETSC_FALSE;
485: else {
486: force_reorth = PETSC_TRUE;
487: compute_int(which2,omega,j-k,delta,eta);
488: for (i=0;i<j-k;i++)
489: if (which2[i]) omega[i] = eps1;
490: }
491: IPOrthogonalize(eps->ip,j-k,which2,V+k,f,hwork,&norm,breakdown,eps->work[0],swork);
492: }
493: }
494: }
495:
496: if (*breakdown || norm < n*anorm*PETSC_MACHINE_EPSILON) {
497: *M = j+1;
498: break;
499: }
500: if (!fro && norm*delta < anorm*eps1) {
501: fro = PETSC_TRUE;
502: PetscInfo1(eps,"Switching to full reorthogonalization at iteration %i\n",eps->its);
503: }
504:
505: beta[j-k] = norm;
506: if (j<m-1) {
507: VecScale(f,1.0/norm);
508: VecCopy(f,V[j+1]);
509: }
510: }
512: if (m>100) {
513: PetscFree(omega);
514: PetscFree(omega_old);
515: }
516: if (eps->nds+m > 100) {
517: PetscFree(which);
518: PetscFree(which2);
519: PetscFree(swork);
520: PetscFree(hwork);
521: }
522: return(0);
523: }
527: /*
528: EPSBasicLanczos - Computes an m-step Lanczos factorization. The first k
529: columns are assumed to be locked and therefore they are not modified. On
530: exit, the following relation is satisfied:
532: OP * V - V * T = f * e_m^T
534: where the columns of V are the Lanczos vectors, T is a tridiagonal matrix,
535: f is the residual vector and e_m is the m-th vector of the canonical basis.
536: The Lanczos vectors (together with vector f) are B-orthogonal (to working
537: accuracy) if full reorthogonalization is being used, otherwise they are
538: (B-)semi-orthogonal. On exit, beta contains the B-norm of f and the next
539: Lanczos vector can be computed as v_{m+1} = f / beta.
541: This function simply calls another function which depends on the selected
542: reorthogonalization strategy.
543: */
544: static PetscErrorCode EPSBasicLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *m,Vec f,PetscTruth *breakdown,PetscReal anorm)
545: {
546: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
548: IPOrthogonalizationRefinementType orthog_ref;
549: PetscScalar *T;
550: PetscInt i,n=*m;
551: PetscReal betam;
554: switch (lanczos->reorthog) {
555: case EPSLANCZOS_REORTHOG_LOCAL:
556: EPSLocalLanczos(eps,alpha,beta,V,k,m,f,breakdown);
557: break;
558: case EPSLANCZOS_REORTHOG_SELECTIVE:
559: EPSSelectiveLanczos(eps,alpha,beta,V,k,m,f,breakdown,anorm);
560: break;
561: case EPSLANCZOS_REORTHOG_FULL:
562: EPSFullLanczos(eps,alpha,beta,V,k,m,f,breakdown);
563: break;
564: case EPSLANCZOS_REORTHOG_PARTIAL:
565: case EPSLANCZOS_REORTHOG_PERIODIC:
566: EPSPartialLanczos(eps,alpha,beta,V,k,m,f,breakdown,anorm);
567: break;
568: case EPSLANCZOS_REORTHOG_DELAYED:
569: PetscMalloc(n*n*sizeof(PetscScalar),&T);
570: IPGetOrthogonalization(eps->ip,PETSC_NULL,&orthog_ref,PETSC_NULL);
571: if (orthog_ref == IP_ORTH_REFINE_NEVER) {
572: EPSDelayedArnoldi1(eps,T,n,V,k,m,f,&betam,breakdown);
573: } else {
574: EPSDelayedArnoldi(eps,T,n,V,k,m,f,&betam,breakdown);
575: }
576: for (i=k;i<n-1;i++) { alpha[i-k] = PetscRealPart(T[n*i+i]); beta[i-k] = PetscRealPart(T[n*i+i+1]); }
577: alpha[n-1] = PetscRealPart(T[n*(n-1)+n-1]);
578: beta[n-1] = betam;
579: PetscFree(T);
580: break;
581: default:
582: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
583: }
584: return(0);
585: }
589: PetscErrorCode EPSSolve_LANCZOS(EPS eps)
590: {
591: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
593: PetscInt nconv,i,j,k,l,x,n,m,*perm,restart,ncv=eps->ncv;
594: Vec w=eps->work[0],f=eps->work[1];
595: PetscScalar *Y,stmp;
596: PetscReal *d,*e,*ritz,*bnd,anorm,beta,norm,rtmp;
597: PetscTruth breakdown;
598: char *conv,ctmp;
601: PetscMalloc(ncv*sizeof(PetscReal),&d);
602: PetscMalloc(ncv*sizeof(PetscReal),&e);
603: PetscMalloc(ncv*sizeof(PetscReal),&ritz);
604: PetscMalloc(ncv*ncv*sizeof(PetscScalar),&Y);
605: PetscMalloc(ncv*sizeof(PetscReal),&bnd);
606: PetscMalloc(ncv*sizeof(PetscInt),&perm);
607: PetscMalloc(ncv*sizeof(char),&conv);
609: /* The first Lanczos vector is the normalized initial vector */
610: EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
611:
612: anorm = -1.0;
613: nconv = 0;
614:
615: /* Restart loop */
616: while (eps->reason == EPS_CONVERGED_ITERATING) {
617: eps->its++;
618: /* Compute an ncv-step Lanczos factorization */
619: m = PetscMin(nconv+eps->mpd,ncv);
620: EPSBasicLanczos(eps,d,e,eps->V,nconv,&m,f,&breakdown,anorm);
622: /* Compute eigenvalues and eigenvectors Y of the tridiagonal block */
623: n = m - nconv;
624: beta = e[n-1];
625: EPSDenseTridiagonal(n,d,e,ritz,Y);
626:
627: /* Estimate ||A|| */
628: for (i=0;i<n;i++)
629: if (PetscAbsReal(ritz[i]) > anorm) anorm = PetscAbsReal(ritz[i]);
630:
631: /* Compute residual norm estimates as beta*abs(Y(m,:)) + eps*||A|| */
632: for (i=0;i<n;i++) {
633: bnd[i] = beta*PetscAbsScalar(Y[i*n+n-1]) + PETSC_MACHINE_EPSILON*anorm;
634: bnd[i] = bnd[i] / PetscAbsScalar(ritz[i]);
635: if (bnd[i] < eps->tol) {
636: conv[i] = 'C';
637: } else {
638: conv[i] = 'N';
639: }
640: }
642: /* purge repeated ritz values */
643: if (lanczos->reorthog == EPSLANCZOS_REORTHOG_LOCAL)
644: for (i=1;i<n;i++)
645: if (conv[i] == 'C')
646: if (PetscAbsScalar((ritz[i]-ritz[i-1])/ritz[i]) < eps->tol)
647: conv[i] = 'R';
649: /* Compute restart vector */
650: if (breakdown) {
651: PetscInfo2(eps,"Breakdown in Lanczos method (it=%i norm=%g)\n",eps->its,beta);
652: } else {
653: restart = 0;
654: while (restart<n && conv[restart] != 'N') restart++;
655: if (restart >= n) {
656: breakdown = PETSC_TRUE;
657: } else {
658: switch (eps->which) {
659: case EPS_LARGEST_MAGNITUDE:
660: case EPS_SMALLEST_MAGNITUDE:
661: rtmp = PetscAbsReal(ritz[restart]);
662: break;
663: case EPS_LARGEST_REAL:
664: case EPS_SMALLEST_REAL:
665: rtmp = ritz[restart];
666: break;
667: default: SETERRQ(1,"Wrong value of which");
668: }
669: for (i=restart+1;i<n;i++)
670: if (conv[i] == 'N') {
671: switch (eps->which) {
672: case EPS_LARGEST_MAGNITUDE:
673: if (rtmp < PetscAbsReal(ritz[i])) { rtmp = PetscAbsReal(ritz[i]); restart = i; }
674: break;
675: case EPS_SMALLEST_MAGNITUDE:
676: if (rtmp > PetscAbsReal(ritz[i])) { rtmp = PetscAbsReal(ritz[i]); restart = i; }
677: break;
678: case EPS_LARGEST_REAL:
679: if (rtmp < ritz[i]) { rtmp = ritz[i]; restart = i; }
680: break;
681: case EPS_SMALLEST_REAL:
682: if (rtmp > ritz[i]) { rtmp = ritz[i]; restart = i; }
683: break;
684: default: SETERRQ(1,"Wrong value of which");
685: }
686: }
687: VecSet(f,0.0);
688: VecMAXPY(f,n,Y+restart*n,eps->V+nconv);
689: }
690: }
692: /* Count and put converged eigenvalues first */
693: for (i=0;i<n;i++) perm[i] = i;
694: for (k=0;k<n;k++)
695: if (conv[perm[k]] != 'C') {
696: j = k + 1;
697: while (j<n && conv[perm[j]] != 'C') j++;
698: if (j>=n) break;
699: l = perm[k]; perm[k] = perm[j]; perm[j] = l;
700: }
702: /* Sort eigenvectors according to permutation */
703: for (i=0;i<k;i++) {
704: x = perm[i];
705: if (x != i) {
706: j = i + 1;
707: while (perm[j] != i) j++;
708: /* swap eigenvalues i and j */
709: rtmp = ritz[x]; ritz[x] = ritz[i]; ritz[i] = rtmp;
710: rtmp = bnd[x]; bnd[x] = bnd[i]; bnd[i] = rtmp;
711: ctmp = conv[x]; conv[x] = conv[i]; conv[i] = ctmp;
712: perm[j] = x; perm[i] = i;
713: /* swap eigenvectors i and j */
714: for (l=0;l<n;l++) {
715: stmp = Y[x*n+l]; Y[x*n+l] = Y[i*n+l]; Y[i*n+l] = stmp;
716: }
717: }
718: }
719:
720: /* compute converged eigenvectors */
721: SlepcUpdateVectors(n,eps->V+nconv,0,k,Y,n,PETSC_FALSE);
722:
723: /* purge spurious ritz values */
724: if (lanczos->reorthog == EPSLANCZOS_REORTHOG_LOCAL) {
725: for (i=0;i<k;i++) {
726: VecNorm(eps->V[nconv+i],NORM_2,&norm);
727: VecScale(eps->V[nconv+i],1.0/norm);
728: STApply(eps->OP,eps->V[nconv+i],w);
729: VecAXPY(w,-ritz[i],eps->V[nconv+i]);
730: VecNorm(w,NORM_2,&norm);
731: bnd[i] = norm / PetscAbsScalar(ritz[i]);
732: if (bnd[i] >= eps->tol) conv[i] = 'S';
733: }
734: for (i=0;i<k;i++)
735: if (conv[i] != 'C') {
736: j = i + 1;
737: while (j<k && conv[j] != 'C') j++;
738: if (j>=k) break;
739: /* swap eigenvalues i and j */
740: rtmp = ritz[j]; ritz[j] = ritz[i]; ritz[i] = rtmp;
741: rtmp = bnd[j]; bnd[j] = bnd[i]; bnd[i] = rtmp;
742: ctmp = conv[j]; conv[j] = conv[i]; conv[i] = ctmp;
743: /* swap eigenvectors i and j */
744: VecSwap(eps->V[nconv+i],eps->V[nconv+j]);
745: }
746: k = i;
747: }
748:
749: /* store ritz values and estimated errors */
750: for (i=0;i<n;i++) {
751: eps->eigr[nconv+i] = ritz[i];
752: eps->errest[nconv+i] = bnd[i];
753: }
754: EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nconv+n);
755: nconv = nconv+k;
756: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
757: if (nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
758:
759: if (eps->reason == EPS_CONVERGED_ITERATING) { /* copy restart vector */
760: if (lanczos->reorthog == EPSLANCZOS_REORTHOG_LOCAL && !breakdown) {
761: /* Reorthonormalize restart vector */
762: IPOrthogonalize(eps->ip,eps->nds+nconv,PETSC_NULL,eps->DSV,f,PETSC_NULL,&norm,&breakdown,w,PETSC_NULL);
763: VecScale(f,1.0/norm);
764: }
765: if (breakdown) {
766: /* Use random vector for restarting */
767: PetscInfo(eps,"Using random vector for restart\n");
768: EPSGetStartVector(eps,nconv,f,&breakdown);
769: }
770: if (breakdown) { /* give up */
771: eps->reason = EPS_DIVERGED_BREAKDOWN;
772: PetscInfo(eps,"Unable to generate more start vectors\n");
773: } else {
774: VecCopy(f,eps->V[nconv]);
775: }
776: }
777: }
778:
779: eps->nconv = nconv;
781: PetscFree(d);
782: PetscFree(e);
783: PetscFree(ritz);
784: PetscFree(Y);
785: PetscFree(bnd);
786: PetscFree(perm);
787: PetscFree(conv);
788: return(0);
789: }
791: static const char *lanczoslist[6] = { "local", "full", "selective", "periodic", "partial" , "delayed" };
795: PetscErrorCode EPSSetFromOptions_LANCZOS(EPS eps)
796: {
798: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
799: PetscTruth flg;
800: PetscInt i;
803: PetscOptionsHead("LANCZOS options");
804: PetscOptionsEList("-eps_lanczos_reorthog","Lanczos reorthogonalization","EPSLanczosSetReorthog",lanczoslist,6,lanczoslist[lanczos->reorthog],&i,&flg);
805: if (flg) lanczos->reorthog = (EPSLanczosReorthogType)i;
806: PetscOptionsTail();
807: return(0);
808: }
813: PetscErrorCode EPSLanczosSetReorthog_LANCZOS(EPS eps,EPSLanczosReorthogType reorthog)
814: {
815: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
818: switch (reorthog) {
819: case EPSLANCZOS_REORTHOG_LOCAL:
820: case EPSLANCZOS_REORTHOG_FULL:
821: case EPSLANCZOS_REORTHOG_DELAYED:
822: case EPSLANCZOS_REORTHOG_SELECTIVE:
823: case EPSLANCZOS_REORTHOG_PERIODIC:
824: case EPSLANCZOS_REORTHOG_PARTIAL:
825: lanczos->reorthog = reorthog;
826: break;
827: default:
828: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
829: }
830: return(0);
831: }
836: /*@
837: EPSLanczosSetReorthog - Sets the type of reorthogonalization used during the Lanczos
838: iteration.
840: Collective on EPS
842: Input Parameters:
843: + eps - the eigenproblem solver context
844: - reorthog - the type of reorthogonalization
846: Options Database Key:
847: . -eps_lanczos_reorthog - Sets the reorthogonalization type (either 'local', 'selective',
848: 'periodic', 'partial', 'full' or 'delayed')
849:
850: Level: advanced
852: .seealso: EPSLanczosGetReorthog(), EPSLanczosReorthogType
853: @*/
854: PetscErrorCode EPSLanczosSetReorthog(EPS eps,EPSLanczosReorthogType reorthog)
855: {
856: PetscErrorCode ierr, (*f)(EPS,EPSLanczosReorthogType);
860: PetscObjectQueryFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",(void (**)())&f);
861: if (f) {
862: (*f)(eps,reorthog);
863: }
864: return(0);
865: }
870: PetscErrorCode EPSLanczosGetReorthog_LANCZOS(EPS eps,EPSLanczosReorthogType *reorthog)
871: {
872: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
874: *reorthog = lanczos->reorthog;
875: return(0);
876: }
881: /*@C
882: EPSLanczosGetReorthog - Gets the type of reorthogonalization used during the Lanczos
883: iteration.
885: Collective on EPS
887: Input Parameter:
888: . eps - the eigenproblem solver context
890: Input Parameter:
891: . reorthog - the type of reorthogonalization
893: Level: advanced
895: .seealso: EPSLanczosSetReorthog(), EPSLanczosReorthogType
896: @*/
897: PetscErrorCode EPSLanczosGetReorthog(EPS eps,EPSLanczosReorthogType *reorthog)
898: {
899: PetscErrorCode ierr, (*f)(EPS,EPSLanczosReorthogType*);
903: PetscObjectQueryFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",(void (**)())&f);
904: if (f) {
905: (*f)(eps,reorthog);
906: }
907: return(0);
908: }
912: PetscErrorCode EPSView_LANCZOS(EPS eps,PetscViewer viewer)
913: {
915: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
916: PetscTruth isascii;
919: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
920: if (!isascii) {
921: SETERRQ1(1,"Viewer type %s not supported for EPSLANCZOS",((PetscObject)viewer)->type_name);
922: }
923: PetscViewerASCIIPrintf(viewer,"reorthogonalization: %s\n",lanczoslist[lanczos->reorthog]);
924: return(0);
925: }
927: /*
928: EXTERN PetscErrorCode EPSSolve_TS_LANCZOS(EPS);
929: */
934: PetscErrorCode EPSCreate_LANCZOS(EPS eps)
935: {
937: EPS_LANCZOS *lanczos;
940: PetscNew(EPS_LANCZOS,&lanczos);
941: PetscLogObjectMemory(eps,sizeof(EPS_LANCZOS));
942: eps->data = (void *) lanczos;
943: eps->ops->solve = EPSSolve_LANCZOS;
944: /* eps->ops->solvets = EPSSolve_TS_LANCZOS;*/
945: eps->ops->setup = EPSSetUp_LANCZOS;
946: eps->ops->setfromoptions = EPSSetFromOptions_LANCZOS;
947: eps->ops->destroy = EPSDestroy_Default;
948: eps->ops->view = EPSView_LANCZOS;
949: eps->ops->backtransform = EPSBackTransform_Default;
950: /*if (eps->solverclass==EPS_TWO_SIDE)
951: eps->ops->computevectors = EPSComputeVectors_Schur;
952: else*/ eps->ops->computevectors = EPSComputeVectors_Hermitian;
953: lanczos->reorthog = EPSLANCZOS_REORTHOG_LOCAL;
954: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosSetReorthog_C","EPSLanczosSetReorthog_LANCZOS",EPSLanczosSetReorthog_LANCZOS);
955: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosGetReorthog_C","EPSLanczosGetReorthog_LANCZOS",EPSLanczosGetReorthog_LANCZOS);
956: return(0);
957: }