Actual source code: lanczos.c
slepc-3.13.1 2020-04-12
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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 https://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: case EPS_ALL:
54: SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
55: default: ; /* default case to remove warning */
56: }
57: if (!eps->extraction) {
58: EPSSetExtraction(eps,EPS_RITZ);
59: } else if (eps->extraction!=EPS_RITZ) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
60: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
62: EPSAllocateSolution(eps,1);
63: EPS_SetInnerProduct(eps);
64: if (lanczos->reorthog != EPS_LANCZOS_REORTHOG_FULL) {
65: BVGetOrthogonalization(eps->V,NULL,&refine,&eta,&btype);
66: BVSetOrthogonalization(eps->V,BV_ORTHOG_MGS,refine,eta,btype);
67: PetscInfo(eps,"Switching to MGS orthogonalization\n");
68: }
69: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_SELECTIVE) {
70: if (!lanczos->allocsize) {
71: BVDuplicate(eps->V,&lanczos->AV);
72: BVGetSizes(lanczos->AV,NULL,NULL,&lanczos->allocsize);
73: } else { /* make sure V and AV have the same size */
74: BVGetSizes(eps->V,NULL,NULL,&lanczos->allocsize);
75: BVResize(lanczos->AV,lanczos->allocsize,PETSC_FALSE);
76: }
77: }
79: DSSetType(eps->ds,DSHEP);
80: DSSetCompact(eps->ds,PETSC_TRUE);
81: DSAllocate(eps->ds,eps->ncv+1);
82: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
83: EPSSetWorkVecs(eps,1);
84: }
86: if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
87: if (eps->isgeneralized && eps->ishermitian && !eps->ispositive) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method does not work for indefinite problems");
88: return(0);
89: }
91: /*
92: EPSLocalLanczos - Local reorthogonalization.
94: This is the simplest variant. At each Lanczos step, the corresponding Lanczos vector
95: is orthogonalized with respect to the two previous Lanczos vectors, according to
96: the three term Lanczos recurrence. WARNING: This variant does not track the loss of
97: orthogonality that occurs in finite-precision arithmetic and, therefore, the
98: generated vectors are not guaranteed to be (semi-)orthogonal.
99: */
100: static PetscErrorCode EPSLocalLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *M,PetscBool *breakdown)
101: {
103: PetscInt i,j,m = *M;
104: Mat Op;
105: PetscBool *which,lwhich[100];
106: PetscScalar *hwork,lhwork[100];
109: if (m > 100) {
110: PetscMalloc2(m,&which,m,&hwork);
111: } else {
112: which = lwhich;
113: hwork = lhwork;
114: }
115: for (i=0;i<k;i++) which[i] = PETSC_TRUE;
117: BVSetActiveColumns(eps->V,0,m);
118: STGetOperator(eps->st,&Op);
119: for (j=k;j<m;j++) {
120: BVMatMultColumn(eps->V,Op,j);
121: which[j] = PETSC_TRUE;
122: if (j-2>=k) which[j-2] = PETSC_FALSE;
123: BVOrthogonalizeSomeColumn(eps->V,j+1,which,hwork,beta+j,breakdown);
124: alpha[j] = PetscRealPart(hwork[j]);
125: if (*breakdown) {
126: *M = j+1;
127: break;
128: } else {
129: BVScaleColumn(eps->V,j+1,1/beta[j]);
130: }
131: }
132: STRestoreOperator(eps->st,&Op);
133: if (m > 100) {
134: PetscFree2(which,hwork);
135: }
136: return(0);
137: }
139: /*
140: DenseTridiagonal - Solves a real tridiagonal Hermitian Eigenvalue Problem.
142: Input Parameters:
143: + n - dimension of the eigenproblem
144: . D - pointer to the array containing the diagonal elements
145: - E - pointer to the array containing the off-diagonal elements
147: Output Parameters:
148: + w - pointer to the array to store the computed eigenvalues
149: - V - pointer to the array to store the eigenvectors
151: Notes:
152: If V is NULL then the eigenvectors are not computed.
154: This routine use LAPACK routines xSTEVR.
155: */
156: static PetscErrorCode DenseTridiagonal(PetscInt n_,PetscReal *D,PetscReal *E,PetscReal *w,PetscScalar *V)
157: {
159: PetscReal abstol = 0.0,vl,vu,*work;
160: PetscBLASInt il,iu,m,*isuppz,n,lwork,*iwork,liwork,info;
161: const char *jobz;
162: #if defined(PETSC_USE_COMPLEX)
163: PetscInt i,j;
164: PetscReal *VV=NULL;
165: #endif
168: PetscBLASIntCast(n_,&n);
169: PetscBLASIntCast(20*n_,&lwork);
170: PetscBLASIntCast(10*n_,&liwork);
171: if (V) {
172: jobz = "V";
173: #if defined(PETSC_USE_COMPLEX)
174: PetscMalloc1(n*n,&VV);
175: #endif
176: } else jobz = "N";
177: PetscMalloc3(2*n,&isuppz,lwork,&work,liwork,&iwork);
178: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
179: #if defined(PETSC_USE_COMPLEX)
180: PetscStackCallBLAS("LAPACKstevr",LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,VV,&n,isuppz,work,&lwork,iwork,&liwork,&info));
181: #else
182: PetscStackCallBLAS("LAPACKstevr",LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info));
183: #endif
184: PetscFPTrapPop();
185: SlepcCheckLapackInfo("stevr",info);
186: #if defined(PETSC_USE_COMPLEX)
187: if (V) {
188: for (i=0;i<n;i++)
189: for (j=0;j<n;j++)
190: V[i*n+j] = VV[i*n+j];
191: PetscFree(VV);
192: }
193: #endif
194: PetscFree3(isuppz,work,iwork);
195: return(0);
196: }
198: /*
199: EPSSelectiveLanczos - Selective reorthogonalization.
200: */
201: static PetscErrorCode EPSSelectiveLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscReal anorm)
202: {
204: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
205: PetscInt i,j,m = *M,n,nritz=0,nritzo;
206: Vec vj1,av;
207: Mat Op;
208: PetscReal *d,*e,*ritz,norm;
209: PetscScalar *Y,*hwork;
210: PetscBool *which;
213: PetscCalloc6(m+1,&d,m,&e,m,&ritz,m*m,&Y,m,&which,m,&hwork);
214: for (i=0;i<k;i++) which[i] = PETSC_TRUE;
215: STGetOperator(eps->st,&Op);
217: for (j=k;j<m;j++) {
218: BVSetActiveColumns(eps->V,0,m);
220: /* Lanczos step */
221: BVMatMultColumn(eps->V,Op,j);
222: which[j] = PETSC_TRUE;
223: if (j-2>=k) which[j-2] = PETSC_FALSE;
224: BVOrthogonalizeSomeColumn(eps->V,j+1,which,hwork,&norm,breakdown);
225: alpha[j] = PetscRealPart(hwork[j]);
226: beta[j] = norm;
227: if (*breakdown) {
228: *M = j+1;
229: break;
230: }
232: /* Compute eigenvalues and eigenvectors Y of the tridiagonal block */
233: n = j-k+1;
234: for (i=0;i<n;i++) {
235: d[i] = alpha[i+k];
236: e[i] = beta[i+k];
237: }
238: DenseTridiagonal(n,d,e,ritz,Y);
240: /* Estimate ||A|| */
241: for (i=0;i<n;i++)
242: if (PetscAbsReal(ritz[i]) > anorm) anorm = PetscAbsReal(ritz[i]);
244: /* Compute nearly converged Ritz vectors */
245: nritzo = 0;
246: for (i=0;i<n;i++) {
247: if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm) nritzo++;
248: }
249: if (nritzo>nritz) {
250: nritz = 0;
251: for (i=0;i<n;i++) {
252: if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm) {
253: BVSetActiveColumns(eps->V,k,k+n);
254: BVGetColumn(lanczos->AV,nritz,&av);
255: BVMultVec(eps->V,1.0,0.0,av,Y+i*n);
256: BVRestoreColumn(lanczos->AV,nritz,&av);
257: nritz++;
258: }
259: }
260: }
261: if (nritz > 0) {
262: BVGetColumn(eps->V,j+1,&vj1);
263: BVSetActiveColumns(lanczos->AV,0,nritz);
264: BVOrthogonalizeVec(lanczos->AV,vj1,hwork,&norm,breakdown);
265: BVRestoreColumn(eps->V,j+1,&vj1);
266: if (*breakdown) {
267: *M = j+1;
268: break;
269: }
270: }
271: BVScaleColumn(eps->V,j+1,1.0/norm);
272: }
274: STRestoreOperator(eps->st,&Op);
275: PetscFree6(d,e,ritz,Y,which,hwork);
276: return(0);
277: }
279: static void update_omega(PetscReal *omega,PetscReal *omega_old,PetscInt j,PetscReal *alpha,PetscReal *beta,PetscReal eps1,PetscReal anorm)
280: {
281: PetscInt k;
282: PetscReal T,binv;
285: /* Estimate of contribution to roundoff errors from A*v
286: fl(A*v) = A*v + f,
287: where ||f|| \approx eps1*||A||.
288: For a full matrix A, a rule-of-thumb estimate is eps1 = sqrt(n)*eps */
289: T = eps1*anorm;
290: binv = 1.0/beta[j+1];
292: /* Update omega(1) using omega(0)==0 */
293: omega_old[0]= beta[1]*omega[1] + (alpha[0]-alpha[j])*omega[0] - beta[j]*omega_old[0];
294: if (omega_old[0] > 0) omega_old[0] = binv*(omega_old[0] + T);
295: else omega_old[0] = binv*(omega_old[0] - T);
297: /* Update remaining components */
298: for (k=1;k<j-1;k++) {
299: 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];
300: if (omega_old[k] > 0) omega_old[k] = binv*(omega_old[k] + T);
301: else omega_old[k] = binv*(omega_old[k] - T);
302: }
303: omega_old[j-1] = binv*T;
305: /* Swap omega and omega_old */
306: for (k=0;k<j;k++) {
307: omega[k] = omega_old[k];
308: omega_old[k] = omega[k];
309: }
310: omega[j] = eps1;
311: PetscFunctionReturnVoid();
312: }
314: static void compute_int(PetscBool *which,PetscReal *mu,PetscInt j,PetscReal delta,PetscReal eta)
315: {
316: PetscInt i,k,maxpos;
317: PetscReal max;
318: PetscBool found;
321: /* initialize which */
322: found = PETSC_FALSE;
323: maxpos = 0;
324: max = 0.0;
325: for (i=0;i<j;i++) {
326: if (PetscAbsReal(mu[i]) >= delta) {
327: which[i] = PETSC_TRUE;
328: found = PETSC_TRUE;
329: } else which[i] = PETSC_FALSE;
330: if (PetscAbsReal(mu[i]) > max) {
331: maxpos = i;
332: max = PetscAbsReal(mu[i]);
333: }
334: }
335: if (!found) which[maxpos] = PETSC_TRUE;
337: for (i=0;i<j;i++) {
338: if (which[i]) {
339: /* find left interval */
340: for (k=i;k>=0;k--) {
341: if (PetscAbsReal(mu[k])<eta || which[k]) break;
342: else which[k] = PETSC_TRUE;
343: }
344: /* find right interval */
345: for (k=i+1;k<j;k++) {
346: if (PetscAbsReal(mu[k])<eta || which[k]) break;
347: else which[k] = PETSC_TRUE;
348: }
349: }
350: }
351: PetscFunctionReturnVoid();
352: }
354: /*
355: EPSPartialLanczos - Partial reorthogonalization.
356: */
357: static PetscErrorCode EPSPartialLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscReal anorm)
358: {
360: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
361: PetscInt i,j,m = *M;
362: Mat Op;
363: PetscReal norm,*omega,lomega[100],*omega_old,lomega_old[100],eps1,delta,eta;
364: PetscBool *which,lwhich[100],*which2,lwhich2[100];
365: PetscBool reorth = PETSC_FALSE,force_reorth = PETSC_FALSE;
366: PetscBool fro = PETSC_FALSE,estimate_anorm = PETSC_FALSE;
367: PetscScalar *hwork,lhwork[100];
370: if (m>100) {
371: PetscMalloc5(m,&omega,m,&omega_old,m,&which,m,&which2,m,&hwork);
372: } else {
373: omega = lomega;
374: omega_old = lomega_old;
375: which = lwhich;
376: which2 = lwhich2;
377: hwork = lhwork;
378: }
380: eps1 = PetscSqrtReal((PetscReal)eps->n)*PETSC_MACHINE_EPSILON/2;
381: delta = PETSC_SQRT_MACHINE_EPSILON/PetscSqrtReal((PetscReal)eps->ncv);
382: eta = PetscPowReal(PETSC_MACHINE_EPSILON,3.0/4.0)/PetscSqrtReal((PetscReal)eps->ncv);
383: if (anorm < 0.0) {
384: anorm = 1.0;
385: estimate_anorm = PETSC_TRUE;
386: }
387: for (i=0;i<PetscMax(100,m);i++) omega[i] = omega_old[i] = 0.0;
388: for (i=0;i<k;i++) which[i] = PETSC_TRUE;
390: BVSetActiveColumns(eps->V,0,m);
391: STGetOperator(eps->st,&Op);
392: for (j=k;j<m;j++) {
393: BVMatMultColumn(eps->V,Op,j);
394: if (fro) {
395: /* Lanczos step with full reorthogonalization */
396: BVOrthogonalizeColumn(eps->V,j+1,hwork,&norm,breakdown);
397: alpha[j] = PetscRealPart(hwork[j]);
398: } else {
399: /* Lanczos step */
400: which[j] = PETSC_TRUE;
401: if (j-2>=k) which[j-2] = PETSC_FALSE;
402: BVOrthogonalizeSomeColumn(eps->V,j+1,which,hwork,&norm,breakdown);
403: alpha[j] = PetscRealPart(hwork[j]);
404: beta[j] = norm;
406: /* Estimate ||A|| if needed */
407: if (estimate_anorm) {
408: if (j>k) anorm = PetscMax(anorm,PetscAbsReal(alpha[j])+norm+beta[j-1]);
409: else anorm = PetscMax(anorm,PetscAbsReal(alpha[j])+norm);
410: }
412: /* Check if reorthogonalization is needed */
413: reorth = PETSC_FALSE;
414: if (j>k) {
415: update_omega(omega,omega_old,j,alpha,beta-1,eps1,anorm);
416: for (i=0;i<j-k;i++) {
417: if (PetscAbsReal(omega[i]) > delta) reorth = PETSC_TRUE;
418: }
419: }
420: if (reorth || force_reorth) {
421: for (i=0;i<k;i++) which2[i] = PETSC_FALSE;
422: for (i=k;i<=j;i++) which2[i] = PETSC_TRUE;
423: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_PERIODIC) {
424: /* Periodic reorthogonalization */
425: if (force_reorth) force_reorth = PETSC_FALSE;
426: else force_reorth = PETSC_TRUE;
427: for (i=0;i<j-k;i++) omega[i] = eps1;
428: } else {
429: /* Partial reorthogonalization */
430: if (force_reorth) force_reorth = PETSC_FALSE;
431: else {
432: force_reorth = PETSC_TRUE;
433: compute_int(which2+k,omega,j-k,delta,eta);
434: for (i=0;i<j-k;i++) {
435: if (which2[i+k]) omega[i] = eps1;
436: }
437: }
438: }
439: BVOrthogonalizeSomeColumn(eps->V,j+1,which2,hwork,&norm,breakdown);
440: }
441: }
443: if (*breakdown || norm < eps->n*anorm*PETSC_MACHINE_EPSILON) {
444: *M = j+1;
445: break;
446: }
447: if (!fro && norm*delta < anorm*eps1) {
448: fro = PETSC_TRUE;
449: PetscInfo1(eps,"Switching to full reorthogonalization at iteration %D\n",eps->its);
450: }
451: beta[j] = norm;
452: BVScaleColumn(eps->V,j+1,1.0/norm);
453: }
455: STRestoreOperator(eps->st,&Op);
456: if (m>100) {
457: PetscFree5(omega,omega_old,which,which2,hwork);
458: }
459: return(0);
460: }
462: /*
463: EPSBasicLanczos - Computes an m-step Lanczos factorization. The first k
464: columns are assumed to be locked and therefore they are not modified. On
465: exit, the following relation is satisfied:
467: OP * V - V * T = f * e_m^T
469: where the columns of V are the Lanczos vectors, T is a tridiagonal matrix,
470: f is the residual vector and e_m is the m-th vector of the canonical basis.
471: The Lanczos vectors (together with vector f) are B-orthogonal (to working
472: accuracy) if full reorthogonalization is being used, otherwise they are
473: (B-)semi-orthogonal. On exit, beta contains the B-norm of f and the next
474: Lanczos vector can be computed as v_{m+1} = f / beta.
476: This function simply calls another function which depends on the selected
477: reorthogonalization strategy.
478: */
479: static PetscErrorCode EPSBasicLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *m,PetscBool *breakdown,PetscReal anorm)
480: {
481: PetscErrorCode ierr;
482: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
483: PetscScalar *T;
484: PetscInt i,n=*m;
485: PetscReal betam;
486: BVOrthogRefineType orthog_ref;
487: Mat Op;
490: switch (lanczos->reorthog) {
491: case EPS_LANCZOS_REORTHOG_LOCAL:
492: EPSLocalLanczos(eps,alpha,beta,k,m,breakdown);
493: break;
494: case EPS_LANCZOS_REORTHOG_FULL:
495: STGetOperator(eps->st,&Op);
496: BVMatLanczos(eps->V,Op,alpha,beta,k,m,breakdown);
497: STRestoreOperator(eps->st,&Op);
498: break;
499: case EPS_LANCZOS_REORTHOG_SELECTIVE:
500: EPSSelectiveLanczos(eps,alpha,beta,k,m,breakdown,anorm);
501: break;
502: case EPS_LANCZOS_REORTHOG_PERIODIC:
503: case EPS_LANCZOS_REORTHOG_PARTIAL:
504: EPSPartialLanczos(eps,alpha,beta,k,m,breakdown,anorm);
505: break;
506: case EPS_LANCZOS_REORTHOG_DELAYED:
507: PetscMalloc1(n*n,&T);
508: BVGetOrthogonalization(eps->V,NULL,&orthog_ref,NULL,NULL);
509: if (orthog_ref == BV_ORTHOG_REFINE_NEVER) {
510: EPSDelayedArnoldi1(eps,T,n,k,m,&betam,breakdown);
511: } else {
512: EPSDelayedArnoldi(eps,T,n,k,m,&betam,breakdown);
513: }
514: for (i=k;i<n-1;i++) {
515: alpha[i] = PetscRealPart(T[n*i+i]);
516: beta[i] = PetscRealPart(T[n*i+i+1]);
517: }
518: alpha[n-1] = PetscRealPart(T[n*(n-1)+n-1]);
519: beta[n-1] = betam;
520: PetscFree(T);
521: break;
522: }
523: return(0);
524: }
526: PetscErrorCode EPSSolve_Lanczos(EPS eps)
527: {
528: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
530: PetscInt nconv,i,j,k,l,x,n,*perm,restart,ncv=eps->ncv,r,ld;
531: Vec vi,vj,w;
532: Mat U;
533: PetscScalar *Y,*ritz,stmp;
534: PetscReal *d,*e,*bnd,anorm,beta,norm,rtmp,resnorm;
535: PetscBool breakdown;
536: char *conv,ctmp;
539: DSGetLeadingDimension(eps->ds,&ld);
540: PetscMalloc4(ncv,&ritz,ncv,&bnd,ncv,&perm,ncv,&conv);
542: /* The first Lanczos vector is the normalized initial vector */
543: EPSGetStartVector(eps,0,NULL);
545: anorm = -1.0;
546: nconv = 0;
548: /* Restart loop */
549: while (eps->reason == EPS_CONVERGED_ITERATING) {
550: eps->its++;
552: /* Compute an ncv-step Lanczos factorization */
553: n = PetscMin(nconv+eps->mpd,ncv);
554: DSGetArrayReal(eps->ds,DS_MAT_T,&d);
555: e = d + ld;
556: EPSBasicLanczos(eps,d,e,nconv,&n,&breakdown,anorm);
557: beta = e[n-1];
558: DSRestoreArrayReal(eps->ds,DS_MAT_T,&d);
559: DSSetDimensions(eps->ds,n,0,nconv,0);
560: DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
561: BVSetActiveColumns(eps->V,nconv,n);
563: /* Solve projected problem */
564: DSSolve(eps->ds,ritz,NULL);
565: DSSort(eps->ds,ritz,NULL,NULL,NULL,NULL);
566: DSSynchronize(eps->ds,ritz,NULL);
568: /* Estimate ||A|| */
569: for (i=nconv;i<n;i++)
570: anorm = PetscMax(anorm,PetscAbsReal(PetscRealPart(ritz[i])));
572: /* Compute residual norm estimates as beta*abs(Y(m,:)) + eps*||A|| */
573: DSGetArray(eps->ds,DS_MAT_Q,&Y);
574: for (i=nconv;i<n;i++) {
575: resnorm = beta*PetscAbsScalar(Y[n-1+i*ld]) + PETSC_MACHINE_EPSILON*anorm;
576: (*eps->converged)(eps,ritz[i],eps->eigi[i],resnorm,&bnd[i],eps->convergedctx);
577: if (bnd[i]<eps->tol) conv[i] = 'C';
578: else conv[i] = 'N';
579: }
580: DSRestoreArray(eps->ds,DS_MAT_Q,&Y);
582: /* purge repeated ritz values */
583: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
584: for (i=nconv+1;i<n;i++) {
585: if (conv[i] == 'C' && PetscAbsScalar((ritz[i]-ritz[i-1])/ritz[i]) < eps->tol) conv[i] = 'R';
586: }
587: }
589: /* Compute restart vector */
590: if (breakdown) {
591: PetscInfo2(eps,"Breakdown in Lanczos method (it=%D norm=%g)\n",eps->its,(double)beta);
592: } else {
593: restart = nconv;
594: while (restart<n && conv[restart] != 'N') restart++;
595: if (restart >= n) {
596: breakdown = PETSC_TRUE;
597: } else {
598: for (i=restart+1;i<n;i++) {
599: if (conv[i] == 'N') {
600: SlepcSCCompare(eps->sc,ritz[restart],0.0,ritz[i],0.0,&r);
601: if (r>0) restart = i;
602: }
603: }
604: DSGetArray(eps->ds,DS_MAT_Q,&Y);
605: BVMultColumn(eps->V,1.0,0.0,n,Y+restart*ld+nconv);
606: DSRestoreArray(eps->ds,DS_MAT_Q,&Y);
607: }
608: }
610: /* Count and put converged eigenvalues first */
611: for (i=nconv;i<n;i++) perm[i] = i;
612: for (k=nconv;k<n;k++) {
613: if (conv[perm[k]] != 'C') {
614: j = k + 1;
615: while (j<n && conv[perm[j]] != 'C') j++;
616: if (j>=n) break;
617: l = perm[k]; perm[k] = perm[j]; perm[j] = l;
618: }
619: }
621: /* Sort eigenvectors according to permutation */
622: DSGetArray(eps->ds,DS_MAT_Q,&Y);
623: for (i=nconv;i<k;i++) {
624: x = perm[i];
625: if (x != i) {
626: j = i + 1;
627: while (perm[j] != i) j++;
628: /* swap eigenvalues i and j */
629: stmp = ritz[x]; ritz[x] = ritz[i]; ritz[i] = stmp;
630: rtmp = bnd[x]; bnd[x] = bnd[i]; bnd[i] = rtmp;
631: ctmp = conv[x]; conv[x] = conv[i]; conv[i] = ctmp;
632: perm[j] = x; perm[i] = i;
633: /* swap eigenvectors i and j */
634: for (l=0;l<n;l++) {
635: stmp = Y[l+x*ld]; Y[l+x*ld] = Y[l+i*ld]; Y[l+i*ld] = stmp;
636: }
637: }
638: }
639: DSRestoreArray(eps->ds,DS_MAT_Q,&Y);
641: /* compute converged eigenvectors */
642: DSGetMat(eps->ds,DS_MAT_Q,&U);
643: BVMultInPlace(eps->V,U,nconv,k);
644: MatDestroy(&U);
646: /* purge spurious ritz values */
647: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
648: for (i=nconv;i<k;i++) {
649: BVGetColumn(eps->V,i,&vi);
650: VecNorm(vi,NORM_2,&norm);
651: VecScale(vi,1.0/norm);
652: w = eps->work[0];
653: STApply(eps->st,vi,w);
654: VecAXPY(w,-ritz[i],vi);
655: BVRestoreColumn(eps->V,i,&vi);
656: VecNorm(w,NORM_2,&norm);
657: (*eps->converged)(eps,ritz[i],eps->eigi[i],norm,&bnd[i],eps->convergedctx);
658: if (bnd[i]>=eps->tol) conv[i] = 'S';
659: }
660: for (i=nconv;i<k;i++) {
661: if (conv[i] != 'C') {
662: j = i + 1;
663: while (j<k && conv[j] != 'C') j++;
664: if (j>=k) break;
665: /* swap eigenvalues i and j */
666: stmp = ritz[j]; ritz[j] = ritz[i]; ritz[i] = stmp;
667: rtmp = bnd[j]; bnd[j] = bnd[i]; bnd[i] = rtmp;
668: ctmp = conv[j]; conv[j] = conv[i]; conv[i] = ctmp;
669: /* swap eigenvectors i and j */
670: BVGetColumn(eps->V,i,&vi);
671: BVGetColumn(eps->V,j,&vj);
672: VecSwap(vi,vj);
673: BVRestoreColumn(eps->V,i,&vi);
674: BVRestoreColumn(eps->V,j,&vj);
675: }
676: }
677: k = i;
678: }
680: /* store ritz values and estimated errors */
681: for (i=nconv;i<n;i++) {
682: eps->eigr[i] = ritz[i];
683: eps->errest[i] = bnd[i];
684: }
685: nconv = k;
686: EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,n);
687: (*eps->stopping)(eps,eps->its,eps->max_it,nconv,eps->nev,&eps->reason,eps->stoppingctx);
689: if (eps->reason == EPS_CONVERGED_ITERATING) { /* copy restart vector */
690: BVCopyColumn(eps->V,n,nconv);
691: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL && !breakdown) {
692: /* Reorthonormalize restart vector */
693: BVOrthonormalizeColumn(eps->V,nconv,PETSC_FALSE,NULL,&breakdown);
694: }
695: if (breakdown) {
696: /* Use random vector for restarting */
697: PetscInfo(eps,"Using random vector for restart\n");
698: EPSGetStartVector(eps,nconv,&breakdown);
699: }
700: if (breakdown) { /* give up */
701: eps->reason = EPS_DIVERGED_BREAKDOWN;
702: PetscInfo(eps,"Unable to generate more start vectors\n");
703: }
704: }
705: }
706: eps->nconv = nconv;
708: PetscFree4(ritz,bnd,perm,conv);
709: return(0);
710: }
712: PetscErrorCode EPSSetFromOptions_Lanczos(PetscOptionItems *PetscOptionsObject,EPS eps)
713: {
714: PetscErrorCode ierr;
715: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
716: PetscBool flg;
717: EPSLanczosReorthogType reorthog;
720: PetscOptionsHead(PetscOptionsObject,"EPS Lanczos Options");
722: PetscOptionsEnum("-eps_lanczos_reorthog","Lanczos reorthogonalization","EPSLanczosSetReorthog",EPSLanczosReorthogTypes,(PetscEnum)lanczos->reorthog,(PetscEnum*)&reorthog,&flg);
723: if (flg) { EPSLanczosSetReorthog(eps,reorthog); }
725: PetscOptionsTail();
726: return(0);
727: }
729: static PetscErrorCode EPSLanczosSetReorthog_Lanczos(EPS eps,EPSLanczosReorthogType reorthog)
730: {
731: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
734: switch (reorthog) {
735: case EPS_LANCZOS_REORTHOG_LOCAL:
736: case EPS_LANCZOS_REORTHOG_FULL:
737: case EPS_LANCZOS_REORTHOG_DELAYED:
738: case EPS_LANCZOS_REORTHOG_SELECTIVE:
739: case EPS_LANCZOS_REORTHOG_PERIODIC:
740: case EPS_LANCZOS_REORTHOG_PARTIAL:
741: if (lanczos->reorthog != reorthog) {
742: lanczos->reorthog = reorthog;
743: eps->state = EPS_STATE_INITIAL;
744: }
745: break;
746: default:
747: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
748: }
749: return(0);
750: }
752: /*@
753: EPSLanczosSetReorthog - Sets the type of reorthogonalization used during the Lanczos
754: iteration.
756: Logically Collective on eps
758: Input Parameters:
759: + eps - the eigenproblem solver context
760: - reorthog - the type of reorthogonalization
762: Options Database Key:
763: . -eps_lanczos_reorthog - Sets the reorthogonalization type (either 'local', 'selective',
764: 'periodic', 'partial', 'full' or 'delayed')
766: Level: advanced
768: .seealso: EPSLanczosGetReorthog(), EPSLanczosReorthogType
769: @*/
770: PetscErrorCode EPSLanczosSetReorthog(EPS eps,EPSLanczosReorthogType reorthog)
771: {
777: PetscTryMethod(eps,"EPSLanczosSetReorthog_C",(EPS,EPSLanczosReorthogType),(eps,reorthog));
778: return(0);
779: }
781: static PetscErrorCode EPSLanczosGetReorthog_Lanczos(EPS eps,EPSLanczosReorthogType *reorthog)
782: {
783: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
786: *reorthog = lanczos->reorthog;
787: return(0);
788: }
790: /*@
791: EPSLanczosGetReorthog - Gets the type of reorthogonalization used during
792: the Lanczos iteration.
794: Not Collective
796: Input Parameter:
797: . eps - the eigenproblem solver context
799: Output Parameter:
800: . reorthog - the type of reorthogonalization
802: Level: advanced
804: .seealso: EPSLanczosSetReorthog(), EPSLanczosReorthogType
805: @*/
806: PetscErrorCode EPSLanczosGetReorthog(EPS eps,EPSLanczosReorthogType *reorthog)
807: {
813: PetscUseMethod(eps,"EPSLanczosGetReorthog_C",(EPS,EPSLanczosReorthogType*),(eps,reorthog));
814: return(0);
815: }
817: PetscErrorCode EPSReset_Lanczos(EPS eps)
818: {
820: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
823: BVDestroy(&lanczos->AV);
824: lanczos->allocsize = 0;
825: return(0);
826: }
828: PetscErrorCode EPSDestroy_Lanczos(EPS eps)
829: {
833: PetscFree(eps->data);
834: PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",NULL);
835: PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",NULL);
836: return(0);
837: }
839: PetscErrorCode EPSView_Lanczos(EPS eps,PetscViewer viewer)
840: {
842: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
843: PetscBool isascii;
846: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
847: if (isascii) {
848: PetscViewerASCIIPrintf(viewer," %s reorthogonalization\n",EPSLanczosReorthogTypes[lanczos->reorthog]);
849: }
850: return(0);
851: }
853: SLEPC_EXTERN PetscErrorCode EPSCreate_Lanczos(EPS eps)
854: {
855: EPS_LANCZOS *ctx;
859: PetscNewLog(eps,&ctx);
860: eps->data = (void*)ctx;
862: eps->useds = PETSC_TRUE;
864: eps->ops->solve = EPSSolve_Lanczos;
865: eps->ops->setup = EPSSetUp_Lanczos;
866: eps->ops->setfromoptions = EPSSetFromOptions_Lanczos;
867: eps->ops->destroy = EPSDestroy_Lanczos;
868: eps->ops->reset = EPSReset_Lanczos;
869: eps->ops->view = EPSView_Lanczos;
870: eps->ops->backtransform = EPSBackTransform_Default;
871: eps->ops->computevectors = EPSComputeVectors_Hermitian;
873: PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",EPSLanczosSetReorthog_Lanczos);
874: PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",EPSLanczosGetReorthog_Lanczos);
875: return(0);
876: }