Actual source code: lanczos.c
slepc-3.5.2 2014-10-10
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: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
19: SLEPc - Scalable Library for Eigenvalue Problem Computations
20: Copyright (c) 2002-2014, 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: PetscReal eta;
55: PetscErrorCode ierr;
58: EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
59: if (eps->ncv>eps->nev+eps->mpd) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must not be larger than nev+mpd");
60: if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
61: if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
62: switch (eps->which) {
63: case EPS_LARGEST_IMAGINARY:
64: case EPS_SMALLEST_IMAGINARY:
65: case EPS_TARGET_IMAGINARY:
66: SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
67: default: ; /* default case to remove warning */
68: }
69: if (!eps->extraction) {
70: EPSSetExtraction(eps,EPS_RITZ);
71: } else if (eps->extraction!=EPS_RITZ) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
72: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
74: EPSAllocateSolution(eps,1);
75: EPS_SetInnerProduct(eps);
76: if (lanczos->reorthog != EPS_LANCZOS_REORTHOG_FULL) {
77: BVGetOrthogonalization(eps->V,NULL,&refine,&eta);
78: BVSetOrthogonalization(eps->V,BV_ORTHOG_MGS,refine,eta);
79: PetscInfo(eps,"Switching to MGS orthogonalization\n");
80: }
81: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_SELECTIVE) {
82: BVDuplicate(eps->V,&lanczos->AV);
83: }
85: DSSetType(eps->ds,DSHEP);
86: DSSetCompact(eps->ds,PETSC_TRUE);
87: DSAllocate(eps->ds,eps->ncv+1);
88: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
89: EPSSetWorkVecs(eps,1);
90: }
92: /* dispatch solve method */
93: if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
94: if (eps->isgeneralized && eps->ishermitian && !eps->ispositive) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method does not work for indefinite problems");
95: eps->ops->solve = EPSSolve_Lanczos;
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,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *M,PetscBool *breakdown)
111: {
113: PetscInt i,j,m = *M;
114: Vec vj,vj1;
115: PetscBool *which,lwhich[100];
116: PetscScalar *hwork,lhwork[100];
119: if (m > 100) {
120: PetscMalloc2(m,&which,m,&hwork);
121: } else {
122: which = lwhich;
123: hwork = lhwork;
124: }
125: for (i=0;i<k;i++) which[i] = PETSC_TRUE;
127: BVSetActiveColumns(eps->V,0,m);
128: for (j=k;j<m;j++) {
129: BVGetColumn(eps->V,j,&vj);
130: BVGetColumn(eps->V,j+1,&vj1);
131: STApply(eps->st,vj,vj1);
132: BVRestoreColumn(eps->V,j,&vj);
133: BVRestoreColumn(eps->V,j+1,&vj1);
134: which[j] = PETSC_TRUE;
135: if (j-2>=k) which[j-2] = PETSC_FALSE;
136: BVOrthogonalizeSomeColumn(eps->V,j+1,which,hwork,beta+j,breakdown);
137: alpha[j] = PetscRealPart(hwork[j]);
138: if (*breakdown) {
139: *M = j+1;
140: break;
141: } else {
142: BVScaleColumn(eps->V,j+1,1/beta[j]);
143: }
144: }
145: if (m > 100) {
146: PetscFree2(which,hwork);
147: }
148: return(0);
149: }
153: /*
154: DenseTridiagonal - Solves a real tridiagonal Hermitian Eigenvalue Problem.
156: Input Parameters:
157: + n - dimension of the eigenproblem
158: . D - pointer to the array containing the diagonal elements
159: - E - pointer to the array containing the off-diagonal elements
161: Output Parameters:
162: + w - pointer to the array to store the computed eigenvalues
163: - V - pointer to the array to store the eigenvectors
165: Notes:
166: If V is NULL then the eigenvectors are not computed.
168: This routine use LAPACK routines xSTEVR.
169: */
170: static PetscErrorCode DenseTridiagonal(PetscInt n_,PetscReal *D,PetscReal *E,PetscReal *w,PetscScalar *V)
171: {
172: #if defined(SLEPC_MISSING_LAPACK_STEVR)
174: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"STEVR - Lapack routine is unavailable");
175: #else
177: PetscReal abstol = 0.0,vl,vu,*work;
178: PetscBLASInt il,iu,m,*isuppz,n,lwork,*iwork,liwork,info;
179: const char *jobz;
180: #if defined(PETSC_USE_COMPLEX)
181: PetscInt i,j;
182: PetscReal *VV;
183: #endif
186: PetscBLASIntCast(n_,&n);
187: PetscBLASIntCast(20*n_,&lwork);
188: PetscBLASIntCast(10*n_,&liwork);
189: if (V) {
190: jobz = "V";
191: #if defined(PETSC_USE_COMPLEX)
192: PetscMalloc1(n*n,&VV);
193: #endif
194: } else jobz = "N";
195: PetscMalloc3(2*n,&isuppz,lwork,&work,liwork,&iwork);
196: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
197: #if defined(PETSC_USE_COMPLEX)
198: PetscStackCallBLAS("LAPACKstevr",LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,VV,&n,isuppz,work,&lwork,iwork,&liwork,&info));
199: #else
200: PetscStackCallBLAS("LAPACKstevr",LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info));
201: #endif
202: PetscFPTrapPop();
203: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack DSTEVR %d",info);
204: #if defined(PETSC_USE_COMPLEX)
205: if (V) {
206: for (i=0;i<n;i++)
207: for (j=0;j<n;j++)
208: V[i*n+j] = VV[i*n+j];
209: PetscFree(VV);
210: }
211: #endif
212: PetscFree3(isuppz,work,iwork);
213: return(0);
214: #endif
215: }
219: /*
220: EPSSelectiveLanczos - Selective reorthogonalization.
221: */
222: static PetscErrorCode EPSSelectiveLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscReal anorm)
223: {
225: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
226: PetscInt i,j,m = *M,n,nritz=0,nritzo;
227: Vec vj,vj1,av;
228: PetscReal *d,*e,*ritz,norm;
229: PetscScalar *Y,*hwork;
230: PetscBool *which;
233: PetscCalloc6(m+1,&d,m,&e,m,&ritz,m*m,&Y,m,&which,m,&hwork);
234: for (i=0;i<k;i++) which[i] = PETSC_TRUE;
236: for (j=k;j<m;j++) {
237: BVSetActiveColumns(eps->V,0,m);
239: /* Lanczos step */
240: BVGetColumn(eps->V,j,&vj);
241: BVGetColumn(eps->V,j+1,&vj1);
242: STApply(eps->st,vj,vj1);
243: BVRestoreColumn(eps->V,j,&vj);
244: BVRestoreColumn(eps->V,j+1,&vj1);
245: which[j] = PETSC_TRUE;
246: if (j-2>=k) which[j-2] = PETSC_FALSE;
247: BVOrthogonalizeSomeColumn(eps->V,j+1,which,hwork,&norm,breakdown);
248: alpha[j] = PetscRealPart(hwork[j]);
249: beta[j] = norm;
250: if (*breakdown) {
251: *M = j+1;
252: break;
253: }
255: /* Compute eigenvalues and eigenvectors Y of the tridiagonal block */
256: n = j-k+1;
257: for (i=0;i<n;i++) {
258: d[i] = alpha[i+k];
259: e[i] = beta[i+k];
260: }
261: DenseTridiagonal(n,d,e,ritz,Y);
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) nritzo++;
271: }
272: if (nritzo>nritz) {
273: nritz = 0;
274: for (i=0;i<n;i++) {
275: if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm) {
276: BVSetActiveColumns(eps->V,k,k+n);
277: BVGetColumn(lanczos->AV,nritz,&av);
278: BVMultVec(eps->V,1.0,0.0,av,Y+i*n);
279: BVRestoreColumn(lanczos->AV,nritz,&av);
280: nritz++;
281: }
282: }
283: }
284: if (nritz > 0) {
285: BVGetColumn(eps->V,j+1,&vj1);
286: BVSetActiveColumns(lanczos->AV,0,nritz);
287: BVOrthogonalizeVec(lanczos->AV,vj1,hwork,&norm,breakdown);
288: BVRestoreColumn(eps->V,j+1,&vj1);
289: if (*breakdown) {
290: *M = j+1;
291: break;
292: }
293: }
294: BVScaleColumn(eps->V,j+1,1.0/norm);
295: }
297: PetscFree6(d,e,ritz,Y,which,hwork);
298: return(0);
299: }
303: static void update_omega(PetscReal *omega,PetscReal *omega_old,PetscInt j,PetscReal *alpha,PetscReal *beta,PetscReal eps1,PetscReal anorm)
304: {
305: PetscInt k;
306: PetscReal T,binv;
309: /* Estimate of contribution to roundoff errors from A*v
310: fl(A*v) = A*v + f,
311: where ||f|| \approx eps1*||A||.
312: For a full matrix A, a rule-of-thumb estimate is eps1 = sqrt(n)*eps */
313: T = eps1*anorm;
314: binv = 1.0/beta[j+1];
316: /* Update omega(1) using omega(0)==0 */
317: omega_old[0]= beta[1]*omega[1] + (alpha[0]-alpha[j])*omega[0] - beta[j]*omega_old[0];
318: if (omega_old[0] > 0) omega_old[0] = binv*(omega_old[0] + T);
319: else omega_old[0] = binv*(omega_old[0] - T);
321: /* Update remaining components */
322: for (k=1;k<j-1;k++) {
323: 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];
324: if (omega_old[k] > 0) omega_old[k] = binv*(omega_old[k] + T);
325: else omega_old[k] = binv*(omega_old[k] - T);
326: }
327: omega_old[j-1] = binv*T;
329: /* Swap omega and omega_old */
330: for (k=0;k<j;k++) {
331: omega[k] = omega_old[k];
332: omega_old[k] = omega[k];
333: }
334: omega[j] = eps1;
335: PetscFunctionReturnVoid();
336: }
340: static void compute_int(PetscBool *which,PetscReal *mu,PetscInt j,PetscReal delta,PetscReal eta)
341: {
342: PetscInt i,k,maxpos;
343: PetscReal max;
344: PetscBool found;
347: /* initialize which */
348: found = PETSC_FALSE;
349: maxpos = 0;
350: max = 0.0;
351: for (i=0;i<j;i++) {
352: if (PetscAbsReal(mu[i]) >= delta) {
353: which[i] = PETSC_TRUE;
354: found = PETSC_TRUE;
355: } else which[i] = PETSC_FALSE;
356: if (PetscAbsReal(mu[i]) > max) {
357: maxpos = i;
358: max = PetscAbsReal(mu[i]);
359: }
360: }
361: if (!found) which[maxpos] = PETSC_TRUE;
363: for (i=0;i<j;i++) {
364: if (which[i]) {
365: /* find left interval */
366: for (k=i;k>=0;k--) {
367: if (PetscAbsReal(mu[k])<eta || which[k]) break;
368: else which[k] = PETSC_TRUE;
369: }
370: /* find right interval */
371: for (k=i+1;k<j;k++) {
372: if (PetscAbsReal(mu[k])<eta || which[k]) break;
373: else which[k] = PETSC_TRUE;
374: }
375: }
376: }
377: PetscFunctionReturnVoid();
378: }
382: /*
383: EPSPartialLanczos - Partial reorthogonalization.
384: */
385: static PetscErrorCode EPSPartialLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscReal anorm)
386: {
388: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
389: PetscInt i,j,m = *M;
390: Vec vj,vj1;
391: PetscReal norm,*omega,lomega[100],*omega_old,lomega_old[100],eps1,delta,eta;
392: PetscBool *which,lwhich[100],*which2,lwhich2[100];
393: PetscBool reorth = PETSC_FALSE,force_reorth = PETSC_FALSE;
394: PetscBool fro = PETSC_FALSE,estimate_anorm = PETSC_FALSE;
395: PetscScalar *hwork,lhwork[100];
398: if (m>100) {
399: PetscMalloc5(m,&omega,m,&omega_old,m,&which,m,&which2,m,&hwork);
400: } else {
401: omega = lomega;
402: omega_old = lomega_old;
403: which = lwhich;
404: which2 = lwhich2;
405: hwork = lhwork;
406: }
408: eps1 = PetscSqrtReal((PetscReal)eps->n)*PETSC_MACHINE_EPSILON/2;
409: delta = PETSC_SQRT_MACHINE_EPSILON/PetscSqrtReal((PetscReal)eps->ncv);
410: eta = PetscPowReal(PETSC_MACHINE_EPSILON,3.0/4.0)/PetscSqrtReal((PetscReal)eps->ncv);
411: if (anorm < 0.0) {
412: anorm = 1.0;
413: estimate_anorm = PETSC_TRUE;
414: }
415: for (i=0;i<m-k;i++) omega[i] = omega_old[i] = 0.0;
416: for (i=0;i<k;i++) which[i] = PETSC_TRUE;
418: BVSetActiveColumns(eps->V,0,m);
419: for (j=k;j<m;j++) {
420: BVGetColumn(eps->V,j,&vj);
421: BVGetColumn(eps->V,j+1,&vj1);
422: STApply(eps->st,vj,vj1);
423: BVRestoreColumn(eps->V,j,&vj);
424: BVRestoreColumn(eps->V,j+1,&vj1);
425: if (fro) {
426: /* Lanczos step with full reorthogonalization */
427: BVOrthogonalizeColumn(eps->V,j+1,hwork,&norm,breakdown);
428: alpha[j] = PetscRealPart(hwork[j]);
429: } else {
430: /* Lanczos step */
431: which[j] = PETSC_TRUE;
432: if (j-2>=k) which[j-2] = PETSC_FALSE;
433: BVOrthogonalizeSomeColumn(eps->V,j+1,which,hwork,&norm,breakdown);
434: alpha[j] = PetscRealPart(hwork[j]);
435: beta[j] = norm;
437: /* Estimate ||A|| if needed */
438: if (estimate_anorm) {
439: if (j>k) anorm = PetscMax(anorm,PetscAbsReal(alpha[j])+norm+beta[j-1]);
440: else anorm = PetscMax(anorm,PetscAbsReal(alpha[j])+norm);
441: }
443: /* Check if reorthogonalization is needed */
444: reorth = PETSC_FALSE;
445: if (j>k) {
446: update_omega(omega,omega_old,j,alpha,beta-1,eps1,anorm);
447: for (i=0;i<j-k;i++) {
448: if (PetscAbsScalar(omega[i]) > delta) reorth = PETSC_TRUE;
449: }
450: }
451: if (reorth || force_reorth) {
452: for (i=0;i<k;i++) which2[i] = PETSC_FALSE;
453: for (i=k;i<=j;i++) which2[i] = PETSC_TRUE;
454: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_PERIODIC) {
455: /* Periodic reorthogonalization */
456: if (force_reorth) force_reorth = PETSC_FALSE;
457: else force_reorth = PETSC_TRUE;
458: for (i=0;i<j-k;i++) omega[i] = eps1;
459: } else {
460: /* Partial reorthogonalization */
461: if (force_reorth) force_reorth = PETSC_FALSE;
462: else {
463: force_reorth = PETSC_TRUE;
464: compute_int(which2+k,omega,j-k,delta,eta);
465: for (i=0;i<j-k;i++) {
466: if (which2[i+k]) omega[i] = eps1;
467: }
468: }
469: }
470: BVOrthogonalizeSomeColumn(eps->V,j+1,which2,hwork,&norm,breakdown);
471: }
472: }
474: if (*breakdown || norm < eps->n*anorm*PETSC_MACHINE_EPSILON) {
475: *M = j+1;
476: break;
477: }
478: if (!fro && norm*delta < anorm*eps1) {
479: fro = PETSC_TRUE;
480: PetscInfo1(eps,"Switching to full reorthogonalization at iteration %D\n",eps->its);
481: }
482: beta[j] = norm;
483: BVScaleColumn(eps->V,j+1,1.0/norm);
484: }
486: if (m>100) {
487: PetscFree5(omega,omega_old,which,which2,hwork);
488: }
489: return(0);
490: }
494: /*
495: EPSBasicLanczos - Computes an m-step Lanczos factorization. The first k
496: columns are assumed to be locked and therefore they are not modified. On
497: exit, the following relation is satisfied:
499: OP * V - V * T = f * e_m^T
501: where the columns of V are the Lanczos vectors, T is a tridiagonal matrix,
502: f is the residual vector and e_m is the m-th vector of the canonical basis.
503: The Lanczos vectors (together with vector f) are B-orthogonal (to working
504: accuracy) if full reorthogonalization is being used, otherwise they are
505: (B-)semi-orthogonal. On exit, beta contains the B-norm of f and the next
506: Lanczos vector can be computed as v_{m+1} = f / beta.
508: This function simply calls another function which depends on the selected
509: reorthogonalization strategy.
510: */
511: static PetscErrorCode EPSBasicLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *m,PetscBool *breakdown,PetscReal anorm)
512: {
513: PetscErrorCode ierr;
514: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
515: /*
516: PetscScalar *T;
517: PetscInt i,n=*m;
518: PetscReal betam;
519: BVOrthogRefineType orthog_ref;
520: */
523: switch (lanczos->reorthog) {
524: case EPS_LANCZOS_REORTHOG_LOCAL:
525: EPSLocalLanczos(eps,alpha,beta,k,m,breakdown);
526: break;
527: case EPS_LANCZOS_REORTHOG_SELECTIVE:
528: EPSSelectiveLanczos(eps,alpha,beta,k,m,breakdown,anorm);
529: break;
530: case EPS_LANCZOS_REORTHOG_FULL:
531: EPSFullLanczos(eps,alpha,beta,k,m,breakdown);
532: break;
533: case EPS_LANCZOS_REORTHOG_PARTIAL:
534: case EPS_LANCZOS_REORTHOG_PERIODIC:
535: EPSPartialLanczos(eps,alpha,beta,k,m,breakdown,anorm);
536: break;
537: case EPS_LANCZOS_REORTHOG_DELAYED:
538: SETERRQ(PetscObjectComm((PetscObject)eps),1,"Not implemented");
539: /*
540: PetscMalloc1(n*n,&T);
541: BVGetOrthogonalization(eps->ip,NULL,&orthog_ref,NULL);
542: if (orthog_ref == BV_ORTHOG_REFINE_NEVER) {
543: EPSDelayedArnoldi1(eps,T,n,V,k,m,f,&betam,breakdown);
544: } else {
545: EPSDelayedArnoldi(eps,T,n,V,k,m,f,&betam,breakdown);
546: }
547: for (i=k;i<n-1;i++) {
548: alpha[i] = PetscRealPart(T[n*i+i]);
549: beta[i] = PetscRealPart(T[n*i+i+1]);
550: }
551: alpha[n-1] = PetscRealPart(T[n*(n-1)+n-1]);
552: beta[n-1] = betam;
553: PetscFree(T);
554: break;
555: */
556: default:
557: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
558: }
559: return(0);
560: }
564: PetscErrorCode EPSSolve_Lanczos(EPS eps)
565: {
566: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
568: PetscInt nconv,i,j,k,l,x,n,*perm,restart,ncv=eps->ncv,r,ld;
569: Vec vi,vj,w;
570: Mat U;
571: PetscScalar *Y,*ritz,stmp;
572: PetscReal *d,*e,*bnd,anorm,beta,norm,rtmp,resnorm;
573: PetscBool breakdown;
574: char *conv,ctmp;
577: DSGetLeadingDimension(eps->ds,&ld);
578: PetscMalloc4(ncv,&ritz,ncv,&bnd,ncv,&perm,ncv,&conv);
580: /* The first Lanczos vector is the normalized initial vector */
581: EPSGetStartVector(eps,0,NULL);
583: anorm = -1.0;
584: nconv = 0;
586: /* Restart loop */
587: while (eps->reason == EPS_CONVERGED_ITERATING) {
588: eps->its++;
590: /* Compute an ncv-step Lanczos factorization */
591: n = PetscMin(nconv+eps->mpd,ncv);
592: DSGetArrayReal(eps->ds,DS_MAT_T,&d);
593: e = d + ld;
594: EPSBasicLanczos(eps,d,e,nconv,&n,&breakdown,anorm);
595: beta = e[n-1];
596: DSRestoreArrayReal(eps->ds,DS_MAT_T,&d);
597: DSSetDimensions(eps->ds,n,0,nconv,0);
598: DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
599: BVSetActiveColumns(eps->V,nconv,n);
601: /* Solve projected problem */
602: DSSolve(eps->ds,ritz,NULL);
603: DSSort(eps->ds,ritz,NULL,NULL,NULL,NULL);
605: /* Estimate ||A|| */
606: for (i=nconv;i<n;i++)
607: anorm = PetscMax(anorm,PetscAbsReal(PetscRealPart(ritz[i])));
609: /* Compute residual norm estimates as beta*abs(Y(m,:)) + eps*||A|| */
610: DSGetArray(eps->ds,DS_MAT_Q,&Y);
611: for (i=nconv;i<n;i++) {
612: resnorm = beta*PetscAbsScalar(Y[n-1+i*ld]) + PETSC_MACHINE_EPSILON*anorm;
613: (*eps->converged)(eps,ritz[i],eps->eigi[i],resnorm,&bnd[i],eps->convergedctx);
614: if (bnd[i]<eps->tol) conv[i] = 'C';
615: else conv[i] = 'N';
616: }
617: DSRestoreArray(eps->ds,DS_MAT_Q,&Y);
619: /* purge repeated ritz values */
620: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
621: for (i=nconv+1;i<n;i++) {
622: if (conv[i] == 'C' && PetscAbsScalar((ritz[i]-ritz[i-1])/ritz[i]) < eps->tol) conv[i] = 'R';
623: }
624: }
626: /* Compute restart vector */
627: if (breakdown) {
628: PetscInfo2(eps,"Breakdown in Lanczos method (it=%D norm=%g)\n",eps->its,(double)beta);
629: } else {
630: restart = nconv;
631: while (restart<n && conv[restart] != 'N') restart++;
632: if (restart >= n) {
633: breakdown = PETSC_TRUE;
634: } else {
635: for (i=restart+1;i<n;i++) {
636: if (conv[i] == 'N') {
637: SlepcSCCompare(eps->sc,ritz[restart],0.0,ritz[i],0.0,&r);
638: if (r>0) restart = i;
639: }
640: }
641: DSGetArray(eps->ds,DS_MAT_Q,&Y);
642: BVMultColumn(eps->V,1.0,0.0,n,Y+restart*ld+nconv);
643: DSRestoreArray(eps->ds,DS_MAT_Q,&Y);
644: }
645: }
647: /* Count and put converged eigenvalues first */
648: for (i=nconv;i<n;i++) perm[i] = i;
649: for (k=nconv;k<n;k++) {
650: if (conv[perm[k]] != 'C') {
651: j = k + 1;
652: while (j<n && conv[perm[j]] != 'C') j++;
653: if (j>=n) break;
654: l = perm[k]; perm[k] = perm[j]; perm[j] = l;
655: }
656: }
658: /* Sort eigenvectors according to permutation */
659: DSGetArray(eps->ds,DS_MAT_Q,&Y);
660: for (i=nconv;i<k;i++) {
661: x = perm[i];
662: if (x != i) {
663: j = i + 1;
664: while (perm[j] != i) j++;
665: /* swap eigenvalues i and j */
666: stmp = ritz[x]; ritz[x] = ritz[i]; ritz[i] = stmp;
667: rtmp = bnd[x]; bnd[x] = bnd[i]; bnd[i] = rtmp;
668: ctmp = conv[x]; conv[x] = conv[i]; conv[i] = ctmp;
669: perm[j] = x; perm[i] = i;
670: /* swap eigenvectors i and j */
671: for (l=0;l<n;l++) {
672: stmp = Y[l+x*ld]; Y[l+x*ld] = Y[l+i*ld]; Y[l+i*ld] = stmp;
673: }
674: }
675: }
676: DSRestoreArray(eps->ds,DS_MAT_Q,&Y);
678: /* compute converged eigenvectors */
679: DSGetMat(eps->ds,DS_MAT_Q,&U);
680: BVMultInPlace(eps->V,U,nconv,k);
681: MatDestroy(&U);
683: /* purge spurious ritz values */
684: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
685: for (i=nconv;i<k;i++) {
686: BVGetColumn(eps->V,i,&vi);
687: VecNorm(vi,NORM_2,&norm);
688: VecScale(vi,1.0/norm);
689: w = eps->work[0];
690: STApply(eps->st,vi,w);
691: VecAXPY(w,-ritz[i],vi);
692: BVRestoreColumn(eps->V,i,&vi);
693: VecNorm(w,NORM_2,&norm);
694: (*eps->converged)(eps,ritz[i],eps->eigi[i],norm,&bnd[i],eps->convergedctx);
695: if (bnd[i]>=eps->tol) conv[i] = 'S';
696: }
697: for (i=nconv;i<k;i++) {
698: if (conv[i] != 'C') {
699: j = i + 1;
700: while (j<k && conv[j] != 'C') j++;
701: if (j>=k) break;
702: /* swap eigenvalues i and j */
703: stmp = ritz[j]; ritz[j] = ritz[i]; ritz[i] = stmp;
704: rtmp = bnd[j]; bnd[j] = bnd[i]; bnd[i] = rtmp;
705: ctmp = conv[j]; conv[j] = conv[i]; conv[i] = ctmp;
706: /* swap eigenvectors i and j */
707: BVGetColumn(eps->V,i,&vi);
708: BVGetColumn(eps->V,j,&vj);
709: VecSwap(vi,vj);
710: BVRestoreColumn(eps->V,i,&vi);
711: BVRestoreColumn(eps->V,j,&vj);
712: }
713: }
714: k = i;
715: }
717: /* store ritz values and estimated errors */
718: for (i=nconv;i<n;i++) {
719: eps->eigr[i] = ritz[i];
720: eps->errest[i] = bnd[i];
721: }
722: EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,n);
723: nconv = k;
724: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
725: if (nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
727: if (eps->reason == EPS_CONVERGED_ITERATING) { /* copy restart vector */
728: BVCopyColumn(eps->V,n,nconv);
729: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL && !breakdown) {
730: /* Reorthonormalize restart vector */
731: BVOrthogonalizeColumn(eps->V,nconv,NULL,&norm,&breakdown);
732: BVScaleColumn(eps->V,nconv,1.0/norm);
733: }
734: if (breakdown) {
735: /* Use random vector for restarting */
736: PetscInfo(eps,"Using random vector for restart\n");
737: EPSGetStartVector(eps,nconv,&breakdown);
738: }
739: if (breakdown) { /* give up */
740: eps->reason = EPS_DIVERGED_BREAKDOWN;
741: PetscInfo(eps,"Unable to generate more start vectors\n");
742: }
743: }
744: }
745: eps->nconv = nconv;
747: PetscFree4(ritz,bnd,perm,conv);
748: return(0);
749: }
753: PetscErrorCode EPSSetFromOptions_Lanczos(EPS eps)
754: {
755: PetscErrorCode ierr;
756: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
757: PetscBool flg;
758: EPSLanczosReorthogType reorthog;
761: PetscOptionsHead("EPS Lanczos Options");
762: PetscOptionsEnum("-eps_lanczos_reorthog","Lanczos reorthogonalization","EPSLanczosSetReorthog",EPSLanczosReorthogTypes,(PetscEnum)lanczos->reorthog,(PetscEnum*)&reorthog,&flg);
763: if (flg) {
764: EPSLanczosSetReorthog(eps,reorthog);
765: }
766: PetscOptionsTail();
767: return(0);
768: }
772: static PetscErrorCode EPSLanczosSetReorthog_Lanczos(EPS eps,EPSLanczosReorthogType reorthog)
773: {
774: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
777: switch (reorthog) {
778: case EPS_LANCZOS_REORTHOG_LOCAL:
779: case EPS_LANCZOS_REORTHOG_FULL:
780: case EPS_LANCZOS_REORTHOG_DELAYED:
781: case EPS_LANCZOS_REORTHOG_SELECTIVE:
782: case EPS_LANCZOS_REORTHOG_PERIODIC:
783: case EPS_LANCZOS_REORTHOG_PARTIAL:
784: lanczos->reorthog = reorthog;
785: break;
786: default:
787: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
788: }
789: return(0);
790: }
794: /*@
795: EPSLanczosSetReorthog - Sets the type of reorthogonalization used during the Lanczos
796: iteration.
798: Logically Collective on EPS
800: Input Parameters:
801: + eps - the eigenproblem solver context
802: - reorthog - the type of reorthogonalization
804: Options Database Key:
805: . -eps_lanczos_reorthog - Sets the reorthogonalization type (either 'local', 'selective',
806: 'periodic', 'partial', 'full' or 'delayed')
808: Level: advanced
810: .seealso: EPSLanczosGetReorthog(), EPSLanczosReorthogType
811: @*/
812: PetscErrorCode EPSLanczosSetReorthog(EPS eps,EPSLanczosReorthogType reorthog)
813: {
819: PetscTryMethod(eps,"EPSLanczosSetReorthog_C",(EPS,EPSLanczosReorthogType),(eps,reorthog));
820: return(0);
821: }
825: static PetscErrorCode EPSLanczosGetReorthog_Lanczos(EPS eps,EPSLanczosReorthogType *reorthog)
826: {
827: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
830: *reorthog = lanczos->reorthog;
831: return(0);
832: }
836: /*@C
837: EPSLanczosGetReorthog - Gets the type of reorthogonalization used during the Lanczos
838: iteration.
840: Not Collective
842: Input Parameter:
843: . eps - the eigenproblem solver context
845: Output Parameter:
846: . reorthog - the type of reorthogonalization
848: Level: advanced
850: .seealso: EPSLanczosSetReorthog(), EPSLanczosReorthogType
851: @*/
852: PetscErrorCode EPSLanczosGetReorthog(EPS eps,EPSLanczosReorthogType *reorthog)
853: {
859: PetscTryMethod(eps,"EPSLanczosGetReorthog_C",(EPS,EPSLanczosReorthogType*),(eps,reorthog));
860: return(0);
861: }
865: PetscErrorCode EPSReset_Lanczos(EPS eps)
866: {
868: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
871: BVDestroy(&lanczos->AV);
872: return(0);
873: }
877: PetscErrorCode EPSDestroy_Lanczos(EPS eps)
878: {
882: PetscFree(eps->data);
883: PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",NULL);
884: PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",NULL);
885: return(0);
886: }
890: PetscErrorCode EPSView_Lanczos(EPS eps,PetscViewer viewer)
891: {
893: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
894: PetscBool isascii;
897: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
898: if (isascii) {
899: PetscViewerASCIIPrintf(viewer," Lanczos: %s reorthogonalization\n",EPSLanczosReorthogTypes[lanczos->reorthog]);
900: }
901: return(0);
902: }
906: PETSC_EXTERN PetscErrorCode EPSCreate_Lanczos(EPS eps)
907: {
908: EPS_LANCZOS *ctx;
912: PetscNewLog(eps,&ctx);
913: eps->data = (void*)ctx;
915: eps->ops->setup = EPSSetUp_Lanczos;
916: eps->ops->setfromoptions = EPSSetFromOptions_Lanczos;
917: eps->ops->destroy = EPSDestroy_Lanczos;
918: eps->ops->reset = EPSReset_Lanczos;
919: eps->ops->view = EPSView_Lanczos;
920: eps->ops->backtransform = EPSBackTransform_Default;
921: eps->ops->computevectors = EPSComputeVectors_Hermitian;
922: PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",EPSLanczosSetReorthog_Lanczos);
923: PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",EPSLanczosGetReorthog_Lanczos);
924: return(0);
925: }