Actual source code: lanczos.c
slepc-3.8.2 2017-12-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2017, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: SLEPc eigensolver: "lanczos"
13: Method: Explicitly Restarted Symmetric/Hermitian Lanczos
15: Algorithm:
17: Lanczos method for symmetric (Hermitian) problems, with explicit
18: restart and deflation. Several reorthogonalization strategies can
19: be selected.
21: References:
23: [1] "Lanczos Methods in SLEPc", SLEPc Technical Report STR-5,
24: available at http://slepc.upv.es.
25: */
27: #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
28: #include <slepcblaslapack.h>
30: typedef struct {
31: EPSLanczosReorthogType reorthog; /* user-provided reorthogonalization parameter */
32: PetscInt allocsize; /* number of columns of work BV's allocated at setup */
33: BV AV; /* work BV used in selective reorthogonalization */
34: } EPS_LANCZOS;
36: PetscErrorCode EPSSetUp_Lanczos(EPS eps)
37: {
38: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
39: BVOrthogRefineType refine;
40: BVOrthogBlockType btype;
41: PetscReal eta;
42: PetscErrorCode ierr;
45: EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
46: if (eps->ncv>eps->nev+eps->mpd) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must not be larger than nev+mpd");
47: if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
48: if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
49: switch (eps->which) {
50: case EPS_LARGEST_IMAGINARY:
51: case EPS_SMALLEST_IMAGINARY:
52: case EPS_TARGET_IMAGINARY:
53: SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
54: default: ; /* default case to remove warning */
55: }
56: if (!eps->extraction) {
57: EPSSetExtraction(eps,EPS_RITZ);
58: } else if (eps->extraction!=EPS_RITZ) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
59: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
61: EPSAllocateSolution(eps,1);
62: EPS_SetInnerProduct(eps);
63: if (lanczos->reorthog != EPS_LANCZOS_REORTHOG_FULL) {
64: BVGetOrthogonalization(eps->V,NULL,&refine,&eta,&btype);
65: BVSetOrthogonalization(eps->V,BV_ORTHOG_MGS,refine,eta,btype);
66: PetscInfo(eps,"Switching to MGS orthogonalization\n");
67: }
68: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_SELECTIVE) {
69: if (!lanczos->allocsize) {
70: BVDuplicate(eps->V,&lanczos->AV);
71: BVGetSizes(lanczos->AV,NULL,NULL,&lanczos->allocsize);
72: } else { /* make sure V and AV have the same size */
73: BVGetSizes(eps->V,NULL,NULL,&lanczos->allocsize);
74: BVResize(lanczos->AV,lanczos->allocsize,PETSC_FALSE);
75: }
76: }
78: DSSetType(eps->ds,DSHEP);
79: DSSetCompact(eps->ds,PETSC_TRUE);
80: DSAllocate(eps->ds,eps->ncv+1);
81: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
82: EPSSetWorkVecs(eps,1);
83: }
85: if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
86: if (eps->isgeneralized && eps->ishermitian && !eps->ispositive) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method does not work for indefinite problems");
87: return(0);
88: }
90: /*
91: EPSLocalLanczos - Local reorthogonalization.
93: This is the simplest variant. At each Lanczos step, the corresponding Lanczos vector
94: is orthogonalized with respect to the two previous Lanczos vectors, according to
95: the three term Lanczos recurrence. WARNING: This variant does not track the loss of
96: orthogonality that occurs in finite-precision arithmetic and, therefore, the
97: generated vectors are not guaranteed to be (semi-)orthogonal.
98: */
99: static PetscErrorCode EPSLocalLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *M,PetscBool *breakdown)
100: {
102: PetscInt i,j,m = *M;
103: Vec vj,vj1;
104: PetscBool *which,lwhich[100];
105: PetscScalar *hwork,lhwork[100];
108: if (m > 100) {
109: PetscMalloc2(m,&which,m,&hwork);
110: } else {
111: which = lwhich;
112: hwork = lhwork;
113: }
114: for (i=0;i<k;i++) which[i] = PETSC_TRUE;
116: BVSetActiveColumns(eps->V,0,m);
117: for (j=k;j<m;j++) {
118: BVGetColumn(eps->V,j,&vj);
119: BVGetColumn(eps->V,j+1,&vj1);
120: STApply(eps->st,vj,vj1);
121: BVRestoreColumn(eps->V,j,&vj);
122: BVRestoreColumn(eps->V,j+1,&vj1);
123: which[j] = PETSC_TRUE;
124: if (j-2>=k) which[j-2] = PETSC_FALSE;
125: BVOrthogonalizeSomeColumn(eps->V,j+1,which,hwork,beta+j,breakdown);
126: alpha[j] = PetscRealPart(hwork[j]);
127: if (*breakdown) {
128: *M = j+1;
129: break;
130: } else {
131: BVScaleColumn(eps->V,j+1,1/beta[j]);
132: }
133: }
134: if (m > 100) {
135: PetscFree2(which,hwork);
136: }
137: return(0);
138: }
140: /*
141: DenseTridiagonal - Solves a real tridiagonal Hermitian Eigenvalue Problem.
143: Input Parameters:
144: + n - dimension of the eigenproblem
145: . D - pointer to the array containing the diagonal elements
146: - E - pointer to the array containing the off-diagonal elements
148: Output Parameters:
149: + w - pointer to the array to store the computed eigenvalues
150: - V - pointer to the array to store the eigenvectors
152: Notes:
153: If V is NULL then the eigenvectors are not computed.
155: This routine use LAPACK routines xSTEVR.
156: */
157: static PetscErrorCode DenseTridiagonal(PetscInt n_,PetscReal *D,PetscReal *E,PetscReal *w,PetscScalar *V)
158: {
159: #if defined(SLEPC_MISSING_LAPACK_STEVR)
161: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"STEVR - Lapack routine is unavailable");
162: #else
164: PetscReal abstol = 0.0,vl,vu,*work;
165: PetscBLASInt il,iu,m,*isuppz,n,lwork,*iwork,liwork,info;
166: const char *jobz;
167: #if defined(PETSC_USE_COMPLEX)
168: PetscInt i,j;
169: PetscReal *VV;
170: #endif
173: PetscBLASIntCast(n_,&n);
174: PetscBLASIntCast(20*n_,&lwork);
175: PetscBLASIntCast(10*n_,&liwork);
176: if (V) {
177: jobz = "V";
178: #if defined(PETSC_USE_COMPLEX)
179: PetscMalloc1(n*n,&VV);
180: #endif
181: } else jobz = "N";
182: PetscMalloc3(2*n,&isuppz,lwork,&work,liwork,&iwork);
183: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
184: #if defined(PETSC_USE_COMPLEX)
185: PetscStackCallBLAS("LAPACKstevr",LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,VV,&n,isuppz,work,&lwork,iwork,&liwork,&info));
186: #else
187: PetscStackCallBLAS("LAPACKstevr",LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info));
188: #endif
189: PetscFPTrapPop();
190: SlepcCheckLapackInfo("stevr",info);
191: #if defined(PETSC_USE_COMPLEX)
192: if (V) {
193: for (i=0;i<n;i++)
194: for (j=0;j<n;j++)
195: V[i*n+j] = VV[i*n+j];
196: PetscFree(VV);
197: }
198: #endif
199: PetscFree3(isuppz,work,iwork);
200: return(0);
201: #endif
202: }
204: /*
205: EPSSelectiveLanczos - Selective reorthogonalization.
206: */
207: static PetscErrorCode EPSSelectiveLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscReal anorm)
208: {
210: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
211: PetscInt i,j,m = *M,n,nritz=0,nritzo;
212: Vec vj,vj1,av;
213: PetscReal *d,*e,*ritz,norm;
214: PetscScalar *Y,*hwork;
215: PetscBool *which;
218: PetscCalloc6(m+1,&d,m,&e,m,&ritz,m*m,&Y,m,&which,m,&hwork);
219: for (i=0;i<k;i++) which[i] = PETSC_TRUE;
221: for (j=k;j<m;j++) {
222: BVSetActiveColumns(eps->V,0,m);
224: /* Lanczos step */
225: BVGetColumn(eps->V,j,&vj);
226: BVGetColumn(eps->V,j+1,&vj1);
227: STApply(eps->st,vj,vj1);
228: BVRestoreColumn(eps->V,j,&vj);
229: BVRestoreColumn(eps->V,j+1,&vj1);
230: which[j] = PETSC_TRUE;
231: if (j-2>=k) which[j-2] = PETSC_FALSE;
232: BVOrthogonalizeSomeColumn(eps->V,j+1,which,hwork,&norm,breakdown);
233: alpha[j] = PetscRealPart(hwork[j]);
234: beta[j] = norm;
235: if (*breakdown) {
236: *M = j+1;
237: break;
238: }
240: /* Compute eigenvalues and eigenvectors Y of the tridiagonal block */
241: n = j-k+1;
242: for (i=0;i<n;i++) {
243: d[i] = alpha[i+k];
244: e[i] = beta[i+k];
245: }
246: DenseTridiagonal(n,d,e,ritz,Y);
248: /* Estimate ||A|| */
249: for (i=0;i<n;i++)
250: if (PetscAbsReal(ritz[i]) > anorm) anorm = PetscAbsReal(ritz[i]);
252: /* Compute nearly converged Ritz vectors */
253: nritzo = 0;
254: for (i=0;i<n;i++) {
255: if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm) nritzo++;
256: }
257: if (nritzo>nritz) {
258: nritz = 0;
259: for (i=0;i<n;i++) {
260: if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm) {
261: BVSetActiveColumns(eps->V,k,k+n);
262: BVGetColumn(lanczos->AV,nritz,&av);
263: BVMultVec(eps->V,1.0,0.0,av,Y+i*n);
264: BVRestoreColumn(lanczos->AV,nritz,&av);
265: nritz++;
266: }
267: }
268: }
269: if (nritz > 0) {
270: BVGetColumn(eps->V,j+1,&vj1);
271: BVSetActiveColumns(lanczos->AV,0,nritz);
272: BVOrthogonalizeVec(lanczos->AV,vj1,hwork,&norm,breakdown);
273: BVRestoreColumn(eps->V,j+1,&vj1);
274: if (*breakdown) {
275: *M = j+1;
276: break;
277: }
278: }
279: BVScaleColumn(eps->V,j+1,1.0/norm);
280: }
282: PetscFree6(d,e,ritz,Y,which,hwork);
283: return(0);
284: }
286: static void update_omega(PetscReal *omega,PetscReal *omega_old,PetscInt j,PetscReal *alpha,PetscReal *beta,PetscReal eps1,PetscReal anorm)
287: {
288: PetscInt k;
289: PetscReal T,binv;
292: /* Estimate of contribution to roundoff errors from A*v
293: fl(A*v) = A*v + f,
294: where ||f|| \approx eps1*||A||.
295: For a full matrix A, a rule-of-thumb estimate is eps1 = sqrt(n)*eps */
296: T = eps1*anorm;
297: binv = 1.0/beta[j+1];
299: /* Update omega(1) using omega(0)==0 */
300: omega_old[0]= beta[1]*omega[1] + (alpha[0]-alpha[j])*omega[0] - beta[j]*omega_old[0];
301: if (omega_old[0] > 0) omega_old[0] = binv*(omega_old[0] + T);
302: else omega_old[0] = binv*(omega_old[0] - T);
304: /* Update remaining components */
305: for (k=1;k<j-1;k++) {
306: 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];
307: if (omega_old[k] > 0) omega_old[k] = binv*(omega_old[k] + T);
308: else omega_old[k] = binv*(omega_old[k] - T);
309: }
310: omega_old[j-1] = binv*T;
312: /* Swap omega and omega_old */
313: for (k=0;k<j;k++) {
314: omega[k] = omega_old[k];
315: omega_old[k] = omega[k];
316: }
317: omega[j] = eps1;
318: PetscFunctionReturnVoid();
319: }
321: static void compute_int(PetscBool *which,PetscReal *mu,PetscInt j,PetscReal delta,PetscReal eta)
322: {
323: PetscInt i,k,maxpos;
324: PetscReal max;
325: PetscBool found;
328: /* initialize which */
329: found = PETSC_FALSE;
330: maxpos = 0;
331: max = 0.0;
332: for (i=0;i<j;i++) {
333: if (PetscAbsReal(mu[i]) >= delta) {
334: which[i] = PETSC_TRUE;
335: found = PETSC_TRUE;
336: } else which[i] = PETSC_FALSE;
337: if (PetscAbsReal(mu[i]) > max) {
338: maxpos = i;
339: max = PetscAbsReal(mu[i]);
340: }
341: }
342: if (!found) which[maxpos] = PETSC_TRUE;
344: for (i=0;i<j;i++) {
345: if (which[i]) {
346: /* find left interval */
347: for (k=i;k>=0;k--) {
348: if (PetscAbsReal(mu[k])<eta || which[k]) break;
349: else which[k] = PETSC_TRUE;
350: }
351: /* find right interval */
352: for (k=i+1;k<j;k++) {
353: if (PetscAbsReal(mu[k])<eta || which[k]) break;
354: else which[k] = PETSC_TRUE;
355: }
356: }
357: }
358: PetscFunctionReturnVoid();
359: }
361: /*
362: EPSPartialLanczos - Partial reorthogonalization.
363: */
364: static PetscErrorCode EPSPartialLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscReal anorm)
365: {
367: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
368: PetscInt i,j,m = *M;
369: Vec vj,vj1;
370: PetscReal norm,*omega,lomega[100],*omega_old,lomega_old[100],eps1,delta,eta;
371: PetscBool *which,lwhich[100],*which2,lwhich2[100];
372: PetscBool reorth = PETSC_FALSE,force_reorth = PETSC_FALSE;
373: PetscBool fro = PETSC_FALSE,estimate_anorm = PETSC_FALSE;
374: PetscScalar *hwork,lhwork[100];
377: if (m>100) {
378: PetscMalloc5(m,&omega,m,&omega_old,m,&which,m,&which2,m,&hwork);
379: } else {
380: omega = lomega;
381: omega_old = lomega_old;
382: which = lwhich;
383: which2 = lwhich2;
384: hwork = lhwork;
385: }
387: eps1 = PetscSqrtReal((PetscReal)eps->n)*PETSC_MACHINE_EPSILON/2;
388: delta = PETSC_SQRT_MACHINE_EPSILON/PetscSqrtReal((PetscReal)eps->ncv);
389: eta = PetscPowReal(PETSC_MACHINE_EPSILON,3.0/4.0)/PetscSqrtReal((PetscReal)eps->ncv);
390: if (anorm < 0.0) {
391: anorm = 1.0;
392: estimate_anorm = PETSC_TRUE;
393: }
394: for (i=0;i<m-k;i++) omega[i] = omega_old[i] = 0.0;
395: for (i=0;i<k;i++) which[i] = PETSC_TRUE;
397: BVSetActiveColumns(eps->V,0,m);
398: for (j=k;j<m;j++) {
399: BVGetColumn(eps->V,j,&vj);
400: BVGetColumn(eps->V,j+1,&vj1);
401: STApply(eps->st,vj,vj1);
402: BVRestoreColumn(eps->V,j,&vj);
403: BVRestoreColumn(eps->V,j+1,&vj1);
404: if (fro) {
405: /* Lanczos step with full reorthogonalization */
406: BVOrthogonalizeColumn(eps->V,j+1,hwork,&norm,breakdown);
407: alpha[j] = PetscRealPart(hwork[j]);
408: } else {
409: /* Lanczos step */
410: which[j] = PETSC_TRUE;
411: if (j-2>=k) which[j-2] = PETSC_FALSE;
412: BVOrthogonalizeSomeColumn(eps->V,j+1,which,hwork,&norm,breakdown);
413: alpha[j] = PetscRealPart(hwork[j]);
414: beta[j] = norm;
416: /* Estimate ||A|| if needed */
417: if (estimate_anorm) {
418: if (j>k) anorm = PetscMax(anorm,PetscAbsReal(alpha[j])+norm+beta[j-1]);
419: else anorm = PetscMax(anorm,PetscAbsReal(alpha[j])+norm);
420: }
422: /* Check if reorthogonalization is needed */
423: reorth = PETSC_FALSE;
424: if (j>k) {
425: update_omega(omega,omega_old,j,alpha,beta-1,eps1,anorm);
426: for (i=0;i<j-k;i++) {
427: if (PetscAbsReal(omega[i]) > delta) reorth = PETSC_TRUE;
428: }
429: }
430: if (reorth || force_reorth) {
431: for (i=0;i<k;i++) which2[i] = PETSC_FALSE;
432: for (i=k;i<=j;i++) which2[i] = PETSC_TRUE;
433: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_PERIODIC) {
434: /* Periodic reorthogonalization */
435: if (force_reorth) force_reorth = PETSC_FALSE;
436: else force_reorth = PETSC_TRUE;
437: for (i=0;i<j-k;i++) omega[i] = eps1;
438: } else {
439: /* Partial reorthogonalization */
440: if (force_reorth) force_reorth = PETSC_FALSE;
441: else {
442: force_reorth = PETSC_TRUE;
443: compute_int(which2+k,omega,j-k,delta,eta);
444: for (i=0;i<j-k;i++) {
445: if (which2[i+k]) omega[i] = eps1;
446: }
447: }
448: }
449: BVOrthogonalizeSomeColumn(eps->V,j+1,which2,hwork,&norm,breakdown);
450: }
451: }
453: if (*breakdown || norm < eps->n*anorm*PETSC_MACHINE_EPSILON) {
454: *M = j+1;
455: break;
456: }
457: if (!fro && norm*delta < anorm*eps1) {
458: fro = PETSC_TRUE;
459: PetscInfo1(eps,"Switching to full reorthogonalization at iteration %D\n",eps->its);
460: }
461: beta[j] = norm;
462: BVScaleColumn(eps->V,j+1,1.0/norm);
463: }
465: if (m>100) {
466: PetscFree5(omega,omega_old,which,which2,hwork);
467: }
468: return(0);
469: }
471: /*
472: EPSBasicLanczos - Computes an m-step Lanczos factorization. The first k
473: columns are assumed to be locked and therefore they are not modified. On
474: exit, the following relation is satisfied:
476: OP * V - V * T = f * e_m^T
478: where the columns of V are the Lanczos vectors, T is a tridiagonal matrix,
479: f is the residual vector and e_m is the m-th vector of the canonical basis.
480: The Lanczos vectors (together with vector f) are B-orthogonal (to working
481: accuracy) if full reorthogonalization is being used, otherwise they are
482: (B-)semi-orthogonal. On exit, beta contains the B-norm of f and the next
483: Lanczos vector can be computed as v_{m+1} = f / beta.
485: This function simply calls another function which depends on the selected
486: reorthogonalization strategy.
487: */
488: static PetscErrorCode EPSBasicLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *m,PetscBool *breakdown,PetscReal anorm)
489: {
490: PetscErrorCode ierr;
491: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
492: PetscScalar *T;
493: PetscInt i,n=*m;
494: PetscReal betam;
495: BVOrthogRefineType orthog_ref;
498: switch (lanczos->reorthog) {
499: case EPS_LANCZOS_REORTHOG_LOCAL:
500: EPSLocalLanczos(eps,alpha,beta,k,m,breakdown);
501: break;
502: case EPS_LANCZOS_REORTHOG_FULL:
503: EPSFullLanczos(eps,alpha,beta,k,m,breakdown);
504: break;
505: case EPS_LANCZOS_REORTHOG_SELECTIVE:
506: EPSSelectiveLanczos(eps,alpha,beta,k,m,breakdown,anorm);
507: break;
508: case EPS_LANCZOS_REORTHOG_PERIODIC:
509: case EPS_LANCZOS_REORTHOG_PARTIAL:
510: EPSPartialLanczos(eps,alpha,beta,k,m,breakdown,anorm);
511: break;
512: case EPS_LANCZOS_REORTHOG_DELAYED:
513: PetscMalloc1(n*n,&T);
514: BVGetOrthogonalization(eps->V,NULL,&orthog_ref,NULL,NULL);
515: if (orthog_ref == BV_ORTHOG_REFINE_NEVER) {
516: EPSDelayedArnoldi1(eps,T,n,k,m,&betam,breakdown);
517: } else {
518: EPSDelayedArnoldi(eps,T,n,k,m,&betam,breakdown);
519: }
520: for (i=k;i<n-1;i++) {
521: alpha[i] = PetscRealPart(T[n*i+i]);
522: beta[i] = PetscRealPart(T[n*i+i+1]);
523: }
524: alpha[n-1] = PetscRealPart(T[n*(n-1)+n-1]);
525: beta[n-1] = betam;
526: PetscFree(T);
527: break;
528: }
529: return(0);
530: }
532: PetscErrorCode EPSSolve_Lanczos(EPS eps)
533: {
534: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
536: PetscInt nconv,i,j,k,l,x,n,*perm,restart,ncv=eps->ncv,r,ld;
537: Vec vi,vj,w;
538: Mat U;
539: PetscScalar *Y,*ritz,stmp;
540: PetscReal *d,*e,*bnd,anorm,beta,norm,rtmp,resnorm;
541: PetscBool breakdown;
542: char *conv,ctmp;
545: DSGetLeadingDimension(eps->ds,&ld);
546: PetscMalloc4(ncv,&ritz,ncv,&bnd,ncv,&perm,ncv,&conv);
548: /* The first Lanczos vector is the normalized initial vector */
549: EPSGetStartVector(eps,0,NULL);
551: anorm = -1.0;
552: nconv = 0;
554: /* Restart loop */
555: while (eps->reason == EPS_CONVERGED_ITERATING) {
556: eps->its++;
558: /* Compute an ncv-step Lanczos factorization */
559: n = PetscMin(nconv+eps->mpd,ncv);
560: DSGetArrayReal(eps->ds,DS_MAT_T,&d);
561: e = d + ld;
562: EPSBasicLanczos(eps,d,e,nconv,&n,&breakdown,anorm);
563: beta = e[n-1];
564: DSRestoreArrayReal(eps->ds,DS_MAT_T,&d);
565: DSSetDimensions(eps->ds,n,0,nconv,0);
566: DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
567: BVSetActiveColumns(eps->V,nconv,n);
569: /* Solve projected problem */
570: DSSolve(eps->ds,ritz,NULL);
571: DSSort(eps->ds,ritz,NULL,NULL,NULL,NULL);
572: DSSynchronize(eps->ds,ritz,NULL);
574: /* Estimate ||A|| */
575: for (i=nconv;i<n;i++)
576: anorm = PetscMax(anorm,PetscAbsReal(PetscRealPart(ritz[i])));
578: /* Compute residual norm estimates as beta*abs(Y(m,:)) + eps*||A|| */
579: DSGetArray(eps->ds,DS_MAT_Q,&Y);
580: for (i=nconv;i<n;i++) {
581: resnorm = beta*PetscAbsScalar(Y[n-1+i*ld]) + PETSC_MACHINE_EPSILON*anorm;
582: (*eps->converged)(eps,ritz[i],eps->eigi[i],resnorm,&bnd[i],eps->convergedctx);
583: if (bnd[i]<eps->tol) conv[i] = 'C';
584: else conv[i] = 'N';
585: }
586: DSRestoreArray(eps->ds,DS_MAT_Q,&Y);
588: /* purge repeated ritz values */
589: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
590: for (i=nconv+1;i<n;i++) {
591: if (conv[i] == 'C' && PetscAbsScalar((ritz[i]-ritz[i-1])/ritz[i]) < eps->tol) conv[i] = 'R';
592: }
593: }
595: /* Compute restart vector */
596: if (breakdown) {
597: PetscInfo2(eps,"Breakdown in Lanczos method (it=%D norm=%g)\n",eps->its,(double)beta);
598: } else {
599: restart = nconv;
600: while (restart<n && conv[restart] != 'N') restart++;
601: if (restart >= n) {
602: breakdown = PETSC_TRUE;
603: } else {
604: for (i=restart+1;i<n;i++) {
605: if (conv[i] == 'N') {
606: SlepcSCCompare(eps->sc,ritz[restart],0.0,ritz[i],0.0,&r);
607: if (r>0) restart = i;
608: }
609: }
610: DSGetArray(eps->ds,DS_MAT_Q,&Y);
611: BVMultColumn(eps->V,1.0,0.0,n,Y+restart*ld+nconv);
612: DSRestoreArray(eps->ds,DS_MAT_Q,&Y);
613: }
614: }
616: /* Count and put converged eigenvalues first */
617: for (i=nconv;i<n;i++) perm[i] = i;
618: for (k=nconv;k<n;k++) {
619: if (conv[perm[k]] != 'C') {
620: j = k + 1;
621: while (j<n && conv[perm[j]] != 'C') j++;
622: if (j>=n) break;
623: l = perm[k]; perm[k] = perm[j]; perm[j] = l;
624: }
625: }
627: /* Sort eigenvectors according to permutation */
628: DSGetArray(eps->ds,DS_MAT_Q,&Y);
629: for (i=nconv;i<k;i++) {
630: x = perm[i];
631: if (x != i) {
632: j = i + 1;
633: while (perm[j] != i) j++;
634: /* swap eigenvalues i and j */
635: stmp = ritz[x]; ritz[x] = ritz[i]; ritz[i] = stmp;
636: rtmp = bnd[x]; bnd[x] = bnd[i]; bnd[i] = rtmp;
637: ctmp = conv[x]; conv[x] = conv[i]; conv[i] = ctmp;
638: perm[j] = x; perm[i] = i;
639: /* swap eigenvectors i and j */
640: for (l=0;l<n;l++) {
641: stmp = Y[l+x*ld]; Y[l+x*ld] = Y[l+i*ld]; Y[l+i*ld] = stmp;
642: }
643: }
644: }
645: DSRestoreArray(eps->ds,DS_MAT_Q,&Y);
647: /* compute converged eigenvectors */
648: DSGetMat(eps->ds,DS_MAT_Q,&U);
649: BVMultInPlace(eps->V,U,nconv,k);
650: MatDestroy(&U);
652: /* purge spurious ritz values */
653: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
654: for (i=nconv;i<k;i++) {
655: BVGetColumn(eps->V,i,&vi);
656: VecNorm(vi,NORM_2,&norm);
657: VecScale(vi,1.0/norm);
658: w = eps->work[0];
659: STApply(eps->st,vi,w);
660: VecAXPY(w,-ritz[i],vi);
661: BVRestoreColumn(eps->V,i,&vi);
662: VecNorm(w,NORM_2,&norm);
663: (*eps->converged)(eps,ritz[i],eps->eigi[i],norm,&bnd[i],eps->convergedctx);
664: if (bnd[i]>=eps->tol) conv[i] = 'S';
665: }
666: for (i=nconv;i<k;i++) {
667: if (conv[i] != 'C') {
668: j = i + 1;
669: while (j<k && conv[j] != 'C') j++;
670: if (j>=k) break;
671: /* swap eigenvalues i and j */
672: stmp = ritz[j]; ritz[j] = ritz[i]; ritz[i] = stmp;
673: rtmp = bnd[j]; bnd[j] = bnd[i]; bnd[i] = rtmp;
674: ctmp = conv[j]; conv[j] = conv[i]; conv[i] = ctmp;
675: /* swap eigenvectors i and j */
676: BVGetColumn(eps->V,i,&vi);
677: BVGetColumn(eps->V,j,&vj);
678: VecSwap(vi,vj);
679: BVRestoreColumn(eps->V,i,&vi);
680: BVRestoreColumn(eps->V,j,&vj);
681: }
682: }
683: k = i;
684: }
686: /* store ritz values and estimated errors */
687: for (i=nconv;i<n;i++) {
688: eps->eigr[i] = ritz[i];
689: eps->errest[i] = bnd[i];
690: }
691: nconv = k;
692: EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,n);
693: (*eps->stopping)(eps,eps->its,eps->max_it,nconv,eps->nev,&eps->reason,eps->stoppingctx);
695: if (eps->reason == EPS_CONVERGED_ITERATING) { /* copy restart vector */
696: BVCopyColumn(eps->V,n,nconv);
697: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL && !breakdown) {
698: /* Reorthonormalize restart vector */
699: BVOrthonormalizeColumn(eps->V,nconv,PETSC_FALSE,NULL,&breakdown);
700: }
701: if (breakdown) {
702: /* Use random vector for restarting */
703: PetscInfo(eps,"Using random vector for restart\n");
704: EPSGetStartVector(eps,nconv,&breakdown);
705: }
706: if (breakdown) { /* give up */
707: eps->reason = EPS_DIVERGED_BREAKDOWN;
708: PetscInfo(eps,"Unable to generate more start vectors\n");
709: }
710: }
711: }
712: eps->nconv = nconv;
714: PetscFree4(ritz,bnd,perm,conv);
715: return(0);
716: }
718: PetscErrorCode EPSSetFromOptions_Lanczos(PetscOptionItems *PetscOptionsObject,EPS eps)
719: {
720: PetscErrorCode ierr;
721: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
722: PetscBool flg;
723: EPSLanczosReorthogType reorthog;
726: PetscOptionsHead(PetscOptionsObject,"EPS Lanczos Options");
728: PetscOptionsEnum("-eps_lanczos_reorthog","Lanczos reorthogonalization","EPSLanczosSetReorthog",EPSLanczosReorthogTypes,(PetscEnum)lanczos->reorthog,(PetscEnum*)&reorthog,&flg);
729: if (flg) { EPSLanczosSetReorthog(eps,reorthog); }
731: PetscOptionsTail();
732: return(0);
733: }
735: static PetscErrorCode EPSLanczosSetReorthog_Lanczos(EPS eps,EPSLanczosReorthogType reorthog)
736: {
737: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
740: switch (reorthog) {
741: case EPS_LANCZOS_REORTHOG_LOCAL:
742: case EPS_LANCZOS_REORTHOG_FULL:
743: case EPS_LANCZOS_REORTHOG_DELAYED:
744: case EPS_LANCZOS_REORTHOG_SELECTIVE:
745: case EPS_LANCZOS_REORTHOG_PERIODIC:
746: case EPS_LANCZOS_REORTHOG_PARTIAL:
747: if (lanczos->reorthog != reorthog) {
748: lanczos->reorthog = reorthog;
749: eps->state = EPS_STATE_INITIAL;
750: }
751: break;
752: default:
753: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
754: }
755: return(0);
756: }
758: /*@
759: EPSLanczosSetReorthog - Sets the type of reorthogonalization used during the Lanczos
760: iteration.
762: Logically Collective on EPS
764: Input Parameters:
765: + eps - the eigenproblem solver context
766: - reorthog - the type of reorthogonalization
768: Options Database Key:
769: . -eps_lanczos_reorthog - Sets the reorthogonalization type (either 'local', 'selective',
770: 'periodic', 'partial', 'full' or 'delayed')
772: Level: advanced
774: .seealso: EPSLanczosGetReorthog(), EPSLanczosReorthogType
775: @*/
776: PetscErrorCode EPSLanczosSetReorthog(EPS eps,EPSLanczosReorthogType reorthog)
777: {
783: PetscTryMethod(eps,"EPSLanczosSetReorthog_C",(EPS,EPSLanczosReorthogType),(eps,reorthog));
784: return(0);
785: }
787: static PetscErrorCode EPSLanczosGetReorthog_Lanczos(EPS eps,EPSLanczosReorthogType *reorthog)
788: {
789: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
792: *reorthog = lanczos->reorthog;
793: return(0);
794: }
796: /*@
797: EPSLanczosGetReorthog - Gets the type of reorthogonalization used during
798: the Lanczos iteration.
800: Not Collective
802: Input Parameter:
803: . eps - the eigenproblem solver context
805: Output Parameter:
806: . reorthog - the type of reorthogonalization
808: Level: advanced
810: .seealso: EPSLanczosSetReorthog(), EPSLanczosReorthogType
811: @*/
812: PetscErrorCode EPSLanczosGetReorthog(EPS eps,EPSLanczosReorthogType *reorthog)
813: {
819: PetscUseMethod(eps,"EPSLanczosGetReorthog_C",(EPS,EPSLanczosReorthogType*),(eps,reorthog));
820: return(0);
821: }
823: PetscErrorCode EPSReset_Lanczos(EPS eps)
824: {
826: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
829: BVDestroy(&lanczos->AV);
830: lanczos->allocsize = 0;
831: return(0);
832: }
834: PetscErrorCode EPSDestroy_Lanczos(EPS eps)
835: {
839: PetscFree(eps->data);
840: PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",NULL);
841: PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",NULL);
842: return(0);
843: }
845: PetscErrorCode EPSView_Lanczos(EPS eps,PetscViewer viewer)
846: {
848: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
849: PetscBool isascii;
852: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
853: if (isascii) {
854: PetscViewerASCIIPrintf(viewer," %s reorthogonalization\n",EPSLanczosReorthogTypes[lanczos->reorthog]);
855: }
856: return(0);
857: }
859: PETSC_EXTERN PetscErrorCode EPSCreate_Lanczos(EPS eps)
860: {
861: EPS_LANCZOS *ctx;
865: PetscNewLog(eps,&ctx);
866: eps->data = (void*)ctx;
868: eps->useds = PETSC_TRUE;
870: eps->ops->solve = EPSSolve_Lanczos;
871: eps->ops->setup = EPSSetUp_Lanczos;
872: eps->ops->setfromoptions = EPSSetFromOptions_Lanczos;
873: eps->ops->destroy = EPSDestroy_Lanczos;
874: eps->ops->reset = EPSReset_Lanczos;
875: eps->ops->view = EPSView_Lanczos;
876: eps->ops->backtransform = EPSBackTransform_Default;
877: eps->ops->computevectors = EPSComputeVectors_Hermitian;
879: PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",EPSLanczosSetReorthog_Lanczos);
880: PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",EPSLanczosGetReorthog_Lanczos);
881: return(0);
882: }