Actual source code: lanczos.c
slepc-3.22.2 2024-12-02
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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>
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: static PetscErrorCode EPSSetUp_Lanczos(EPS eps)
37: {
38: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
39: BVOrthogRefineType refine;
40: BVOrthogBlockType btype;
41: PetscReal eta;
43: PetscFunctionBegin;
44: EPSCheckHermitianDefinite(eps);
45: EPSCheckNotStructured(eps);
46: PetscCall(EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd));
47: PetscCheck(eps->ncv<=eps->nev+eps->mpd,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv must not be larger than nev+mpd");
48: if (eps->max_it==PETSC_DETERMINE) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
49: if (!eps->which) PetscCall(EPSSetWhichEigenpairs_Default(eps));
50: PetscCheck(eps->which!=EPS_ALL,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support computing all eigenvalues");
51: EPSCheckUnsupported(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_EXTRACTION);
52: EPSCheckIgnored(eps,EPS_FEATURE_BALANCE);
54: PetscCheck(lanczos->reorthog!=(EPSLanczosReorthogType)-1,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"You should explicitly provide the reorthogonalization type, e.g., -eps_lanczos_reorthog local. Note that the EPSLANCZOS solver is *NOT RECOMMENDED* for general use, because it uses explicit restart which typically has slow convergence. The recommended solver is EPSKRYLOVSCHUR (the default), which implements Lanczos with thick restart in the case of symmetric/Hermitian problems");
56: PetscCall(EPSAllocateSolution(eps,1));
57: PetscCall(EPS_SetInnerProduct(eps));
58: if (lanczos->reorthog != EPS_LANCZOS_REORTHOG_FULL) {
59: PetscCall(BVGetOrthogonalization(eps->V,NULL,&refine,&eta,&btype));
60: PetscCall(BVSetOrthogonalization(eps->V,BV_ORTHOG_MGS,refine,eta,btype));
61: PetscCall(PetscInfo(eps,"Switching to MGS orthogonalization\n"));
62: }
63: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_SELECTIVE) {
64: if (!lanczos->allocsize) {
65: PetscCall(BVDuplicate(eps->V,&lanczos->AV));
66: PetscCall(BVGetSizes(lanczos->AV,NULL,NULL,&lanczos->allocsize));
67: } else { /* make sure V and AV have the same size */
68: PetscCall(BVGetSizes(eps->V,NULL,NULL,&lanczos->allocsize));
69: PetscCall(BVResize(lanczos->AV,lanczos->allocsize,PETSC_FALSE));
70: }
71: }
73: PetscCall(DSSetType(eps->ds,DSHEP));
74: PetscCall(DSSetCompact(eps->ds,PETSC_TRUE));
75: PetscCall(DSAllocate(eps->ds,eps->ncv+1));
76: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) PetscCall(EPSSetWorkVecs(eps,1));
77: PetscFunctionReturn(PETSC_SUCCESS);
78: }
80: /*
81: EPSLocalLanczos - Local reorthogonalization.
83: This is the simplest variant. At each Lanczos step, the corresponding Lanczos vector
84: is orthogonalized with respect to the two previous Lanczos vectors, according to
85: the three term Lanczos recurrence. WARNING: This variant does not track the loss of
86: orthogonality that occurs in finite-precision arithmetic and, therefore, the
87: generated vectors are not guaranteed to be (semi-)orthogonal.
88: */
89: static PetscErrorCode EPSLocalLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *M,PetscBool *breakdown)
90: {
91: PetscInt i,j,m = *M;
92: Mat Op;
93: PetscBool *which,lwhich[100];
94: PetscScalar *hwork,lhwork[100];
96: PetscFunctionBegin;
97: if (m > 100) PetscCall(PetscMalloc2(m,&which,m,&hwork));
98: else {
99: which = lwhich;
100: hwork = lhwork;
101: }
102: for (i=0;i<k;i++) which[i] = PETSC_TRUE;
104: PetscCall(BVSetActiveColumns(eps->V,0,m));
105: PetscCall(STGetOperator(eps->st,&Op));
106: for (j=k;j<m;j++) {
107: PetscCall(BVMatMultColumn(eps->V,Op,j));
108: which[j] = PETSC_TRUE;
109: if (j-2>=k) which[j-2] = PETSC_FALSE;
110: PetscCall(BVOrthogonalizeSomeColumn(eps->V,j+1,which,hwork,beta+j,breakdown));
111: alpha[j] = PetscRealPart(hwork[j]);
112: if (PetscUnlikely(*breakdown)) {
113: *M = j+1;
114: break;
115: } else PetscCall(BVScaleColumn(eps->V,j+1,1/beta[j]));
116: }
117: PetscCall(STRestoreOperator(eps->st,&Op));
118: if (m > 100) PetscCall(PetscFree2(which,hwork));
119: PetscFunctionReturn(PETSC_SUCCESS);
120: }
122: /*
123: DenseTridiagonal - Solves a real tridiagonal Hermitian Eigenvalue Problem.
125: Input Parameters:
126: + n - dimension of the eigenproblem
127: . D - pointer to the array containing the diagonal elements
128: - E - pointer to the array containing the off-diagonal elements
130: Output Parameters:
131: + w - pointer to the array to store the computed eigenvalues
132: - V - pointer to the array to store the eigenvectors
134: Notes:
135: If V is NULL then the eigenvectors are not computed.
137: This routine use LAPACK routines xSTEVR.
138: */
139: static PetscErrorCode DenseTridiagonal(PetscInt n_,PetscReal *D,PetscReal *E,PetscReal *w,PetscScalar *V)
140: {
141: PetscReal abstol = 0.0,vl,vu,*work;
142: PetscBLASInt il,iu,m,*isuppz,n,lwork,*iwork,liwork,info;
143: const char *jobz;
144: #if defined(PETSC_USE_COMPLEX)
145: PetscInt i,j;
146: PetscReal *VV=NULL;
147: #endif
149: PetscFunctionBegin;
150: PetscCall(PetscBLASIntCast(n_,&n));
151: PetscCall(PetscBLASIntCast(20*n_,&lwork));
152: PetscCall(PetscBLASIntCast(10*n_,&liwork));
153: if (V) {
154: jobz = "V";
155: #if defined(PETSC_USE_COMPLEX)
156: PetscCall(PetscMalloc1(n*n,&VV));
157: #endif
158: } else jobz = "N";
159: PetscCall(PetscMalloc3(2*n,&isuppz,lwork,&work,liwork,&iwork));
160: PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
161: #if defined(PETSC_USE_COMPLEX)
162: PetscCallBLAS("LAPACKstevr",LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,VV,&n,isuppz,work,&lwork,iwork,&liwork,&info));
163: #else
164: PetscCallBLAS("LAPACKstevr",LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info));
165: #endif
166: PetscCall(PetscFPTrapPop());
167: SlepcCheckLapackInfo("stevr",info);
168: #if defined(PETSC_USE_COMPLEX)
169: if (V) {
170: for (i=0;i<n;i++)
171: for (j=0;j<n;j++)
172: V[i*n+j] = VV[i*n+j];
173: PetscCall(PetscFree(VV));
174: }
175: #endif
176: PetscCall(PetscFree3(isuppz,work,iwork));
177: PetscFunctionReturn(PETSC_SUCCESS);
178: }
180: /*
181: EPSSelectiveLanczos - Selective reorthogonalization.
182: */
183: static PetscErrorCode EPSSelectiveLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscReal anorm)
184: {
185: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
186: PetscInt i,j,m = *M,n,nritz=0,nritzo;
187: Vec vj1,av;
188: Mat Op;
189: PetscReal *d,*e,*ritz,norm;
190: PetscScalar *Y,*hwork;
191: PetscBool *which;
193: PetscFunctionBegin;
194: PetscCall(PetscCalloc6(m+1,&d,m,&e,m,&ritz,m*m,&Y,m,&which,m,&hwork));
195: for (i=0;i<k;i++) which[i] = PETSC_TRUE;
196: PetscCall(STGetOperator(eps->st,&Op));
198: for (j=k;j<m;j++) {
199: PetscCall(BVSetActiveColumns(eps->V,0,m));
201: /* Lanczos step */
202: PetscCall(BVMatMultColumn(eps->V,Op,j));
203: which[j] = PETSC_TRUE;
204: if (j-2>=k) which[j-2] = PETSC_FALSE;
205: PetscCall(BVOrthogonalizeSomeColumn(eps->V,j+1,which,hwork,&norm,breakdown));
206: alpha[j] = PetscRealPart(hwork[j]);
207: beta[j] = norm;
208: if (PetscUnlikely(*breakdown)) {
209: *M = j+1;
210: break;
211: }
213: /* Compute eigenvalues and eigenvectors Y of the tridiagonal block */
214: n = j-k+1;
215: for (i=0;i<n;i++) {
216: d[i] = alpha[i+k];
217: e[i] = beta[i+k];
218: }
219: PetscCall(DenseTridiagonal(n,d,e,ritz,Y));
221: /* Estimate ||A|| */
222: for (i=0;i<n;i++)
223: if (PetscAbsReal(ritz[i]) > anorm) anorm = PetscAbsReal(ritz[i]);
225: /* Compute nearly converged Ritz vectors */
226: nritzo = 0;
227: for (i=0;i<n;i++) {
228: if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm) nritzo++;
229: }
230: if (nritzo>nritz) {
231: nritz = 0;
232: for (i=0;i<n;i++) {
233: if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm) {
234: PetscCall(BVSetActiveColumns(eps->V,k,k+n));
235: PetscCall(BVGetColumn(lanczos->AV,nritz,&av));
236: PetscCall(BVMultVec(eps->V,1.0,0.0,av,Y+i*n));
237: PetscCall(BVRestoreColumn(lanczos->AV,nritz,&av));
238: nritz++;
239: }
240: }
241: }
242: if (nritz > 0) {
243: PetscCall(BVGetColumn(eps->V,j+1,&vj1));
244: PetscCall(BVSetActiveColumns(lanczos->AV,0,nritz));
245: PetscCall(BVOrthogonalizeVec(lanczos->AV,vj1,hwork,&norm,breakdown));
246: PetscCall(BVRestoreColumn(eps->V,j+1,&vj1));
247: if (PetscUnlikely(*breakdown)) {
248: *M = j+1;
249: break;
250: }
251: }
252: PetscCall(BVScaleColumn(eps->V,j+1,1.0/norm));
253: }
255: PetscCall(STRestoreOperator(eps->st,&Op));
256: PetscCall(PetscFree6(d,e,ritz,Y,which,hwork));
257: PetscFunctionReturn(PETSC_SUCCESS);
258: }
260: static void update_omega(PetscReal *omega,PetscReal *omega_old,PetscInt j,PetscReal *alpha,PetscReal *beta,PetscReal eps1,PetscReal anorm)
261: {
262: PetscInt k;
263: PetscReal T,binv;
265: PetscFunctionBegin;
266: /* Estimate of contribution to roundoff errors from A*v
267: fl(A*v) = A*v + f,
268: where ||f|| \approx eps1*||A||.
269: For a full matrix A, a rule-of-thumb estimate is eps1 = sqrt(n)*eps */
270: T = eps1*anorm;
271: binv = 1.0/beta[j+1];
273: /* Update omega(1) using omega(0)==0 */
274: omega_old[0]= beta[1]*omega[1] + (alpha[0]-alpha[j])*omega[0] - beta[j]*omega_old[0];
275: if (omega_old[0] > 0) omega_old[0] = binv*(omega_old[0] + T);
276: else omega_old[0] = binv*(omega_old[0] - T);
278: /* Update remaining components */
279: for (k=1;k<j-1;k++) {
280: 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];
281: if (omega_old[k] > 0) omega_old[k] = binv*(omega_old[k] + T);
282: else omega_old[k] = binv*(omega_old[k] - T);
283: }
284: omega_old[j-1] = binv*T;
286: /* Swap omega and omega_old */
287: for (k=0;k<j;k++) {
288: omega[k] = omega_old[k];
289: omega_old[k] = omega[k];
290: }
291: omega[j] = eps1;
292: PetscFunctionReturnVoid();
293: }
295: static void compute_int(PetscBool *which,PetscReal *mu,PetscInt j,PetscReal delta,PetscReal eta)
296: {
297: PetscInt i,k,maxpos;
298: PetscReal max;
299: PetscBool found;
301: PetscFunctionBegin;
302: /* initialize which */
303: found = PETSC_FALSE;
304: maxpos = 0;
305: max = 0.0;
306: for (i=0;i<j;i++) {
307: if (PetscAbsReal(mu[i]) >= delta) {
308: which[i] = PETSC_TRUE;
309: found = PETSC_TRUE;
310: } else which[i] = PETSC_FALSE;
311: if (PetscAbsReal(mu[i]) > max) {
312: maxpos = i;
313: max = PetscAbsReal(mu[i]);
314: }
315: }
316: if (!found) which[maxpos] = PETSC_TRUE;
318: for (i=0;i<j;i++) {
319: if (which[i]) {
320: /* find left interval */
321: for (k=i;k>=0;k--) {
322: if (PetscAbsReal(mu[k])<eta || which[k]) break;
323: else which[k] = PETSC_TRUE;
324: }
325: /* find right interval */
326: for (k=i+1;k<j;k++) {
327: if (PetscAbsReal(mu[k])<eta || which[k]) break;
328: else which[k] = PETSC_TRUE;
329: }
330: }
331: }
332: PetscFunctionReturnVoid();
333: }
335: /*
336: EPSPartialLanczos - Partial reorthogonalization.
337: */
338: static PetscErrorCode EPSPartialLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscReal anorm)
339: {
340: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
341: PetscInt i,j,m = *M;
342: Mat Op;
343: PetscReal norm,*omega,lomega[100],*omega_old,lomega_old[100],eps1,delta,eta;
344: PetscBool *which,lwhich[100],*which2,lwhich2[100];
345: PetscBool reorth = PETSC_FALSE,force_reorth = PETSC_FALSE;
346: PetscBool fro = PETSC_FALSE,estimate_anorm = PETSC_FALSE;
347: PetscScalar *hwork,lhwork[100];
349: PetscFunctionBegin;
350: if (m>100) PetscCall(PetscMalloc5(m,&omega,m,&omega_old,m,&which,m,&which2,m,&hwork));
351: else {
352: omega = lomega;
353: omega_old = lomega_old;
354: which = lwhich;
355: which2 = lwhich2;
356: hwork = lhwork;
357: }
359: eps1 = PetscSqrtReal((PetscReal)eps->n)*PETSC_MACHINE_EPSILON/2;
360: delta = PETSC_SQRT_MACHINE_EPSILON/PetscSqrtReal((PetscReal)eps->ncv);
361: eta = PetscPowReal(PETSC_MACHINE_EPSILON,3.0/4.0)/PetscSqrtReal((PetscReal)eps->ncv);
362: if (anorm < 0.0) {
363: anorm = 1.0;
364: estimate_anorm = PETSC_TRUE;
365: }
366: for (i=0;i<PetscMax(100,m);i++) omega[i] = omega_old[i] = 0.0;
367: for (i=0;i<k;i++) which[i] = PETSC_TRUE;
369: PetscCall(BVSetActiveColumns(eps->V,0,m));
370: PetscCall(STGetOperator(eps->st,&Op));
371: for (j=k;j<m;j++) {
372: PetscCall(BVMatMultColumn(eps->V,Op,j));
373: if (fro) {
374: /* Lanczos step with full reorthogonalization */
375: PetscCall(BVOrthogonalizeColumn(eps->V,j+1,hwork,&norm,breakdown));
376: alpha[j] = PetscRealPart(hwork[j]);
377: } else {
378: /* Lanczos step */
379: which[j] = PETSC_TRUE;
380: if (j-2>=k) which[j-2] = PETSC_FALSE;
381: PetscCall(BVOrthogonalizeSomeColumn(eps->V,j+1,which,hwork,&norm,breakdown));
382: alpha[j] = PetscRealPart(hwork[j]);
383: beta[j] = norm;
385: /* Estimate ||A|| if needed */
386: if (estimate_anorm) {
387: if (j>k) anorm = PetscMax(anorm,PetscAbsReal(alpha[j])+norm+beta[j-1]);
388: else anorm = PetscMax(anorm,PetscAbsReal(alpha[j])+norm);
389: }
391: /* Check if reorthogonalization is needed */
392: reorth = PETSC_FALSE;
393: if (j>k) {
394: update_omega(omega,omega_old,j,alpha,beta-1,eps1,anorm);
395: for (i=0;i<j-k;i++) {
396: if (PetscAbsReal(omega[i]) > delta) reorth = PETSC_TRUE;
397: }
398: }
399: if (reorth || force_reorth) {
400: for (i=0;i<k;i++) which2[i] = PETSC_FALSE;
401: for (i=k;i<=j;i++) which2[i] = PETSC_TRUE;
402: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_PERIODIC) {
403: /* Periodic reorthogonalization */
404: if (force_reorth) force_reorth = PETSC_FALSE;
405: else force_reorth = PETSC_TRUE;
406: for (i=0;i<j-k;i++) omega[i] = eps1;
407: } else {
408: /* Partial reorthogonalization */
409: if (force_reorth) force_reorth = PETSC_FALSE;
410: else {
411: force_reorth = PETSC_TRUE;
412: compute_int(which2+k,omega,j-k,delta,eta);
413: for (i=0;i<j-k;i++) {
414: if (which2[i+k]) omega[i] = eps1;
415: }
416: }
417: }
418: PetscCall(BVOrthogonalizeSomeColumn(eps->V,j+1,which2,hwork,&norm,breakdown));
419: }
420: }
422: if (PetscUnlikely(*breakdown || norm < eps->n*anorm*PETSC_MACHINE_EPSILON)) {
423: *M = j+1;
424: break;
425: }
426: if (!fro && norm*delta < anorm*eps1) {
427: fro = PETSC_TRUE;
428: PetscCall(PetscInfo(eps,"Switching to full reorthogonalization at iteration %" PetscInt_FMT "\n",eps->its));
429: }
430: beta[j] = norm;
431: PetscCall(BVScaleColumn(eps->V,j+1,1.0/norm));
432: }
434: PetscCall(STRestoreOperator(eps->st,&Op));
435: if (m>100) PetscCall(PetscFree5(omega,omega_old,which,which2,hwork));
436: PetscFunctionReturn(PETSC_SUCCESS);
437: }
439: /*
440: EPSBasicLanczos - Computes an m-step Lanczos factorization. The first k
441: columns are assumed to be locked and therefore they are not modified. On
442: exit, the following relation is satisfied:
444: OP * V - V * T = f * e_m^T
446: where the columns of V are the Lanczos vectors, T is a tridiagonal matrix,
447: f is the residual vector and e_m is the m-th vector of the canonical basis.
448: The Lanczos vectors (together with vector f) are B-orthogonal (to working
449: accuracy) if full reorthogonalization is being used, otherwise they are
450: (B-)semi-orthogonal. On exit, beta contains the B-norm of f and the next
451: Lanczos vector can be computed as v_{m+1} = f / beta.
453: This function simply calls another function which depends on the selected
454: reorthogonalization strategy.
455: */
456: static PetscErrorCode EPSBasicLanczos(EPS eps,PetscInt k,PetscInt *m,PetscReal *betam,PetscBool *breakdown,PetscReal anorm)
457: {
458: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
459: PetscScalar *T;
460: PetscInt i,n=*m,ld;
461: PetscReal *alpha,*beta;
462: BVOrthogRefineType orthog_ref;
463: Mat Op,M;
465: PetscFunctionBegin;
466: PetscCall(DSGetLeadingDimension(eps->ds,&ld));
467: switch (lanczos->reorthog) {
468: case EPS_LANCZOS_REORTHOG_LOCAL:
469: PetscCall(DSGetArrayReal(eps->ds,DS_MAT_T,&alpha));
470: beta = alpha + ld;
471: PetscCall(EPSLocalLanczos(eps,alpha,beta,k,m,breakdown));
472: *betam = beta[*m-1];
473: PetscCall(DSRestoreArrayReal(eps->ds,DS_MAT_T,&alpha));
474: break;
475: case EPS_LANCZOS_REORTHOG_FULL:
476: PetscCall(STGetOperator(eps->st,&Op));
477: PetscCall(DSGetMat(eps->ds,DS_MAT_T,&M));
478: PetscCall(BVMatLanczos(eps->V,Op,M,k,m,betam,breakdown));
479: PetscCall(DSRestoreMat(eps->ds,DS_MAT_T,&M));
480: PetscCall(STRestoreOperator(eps->st,&Op));
481: break;
482: case EPS_LANCZOS_REORTHOG_SELECTIVE:
483: PetscCall(DSGetArrayReal(eps->ds,DS_MAT_T,&alpha));
484: beta = alpha + ld;
485: PetscCall(EPSSelectiveLanczos(eps,alpha,beta,k,m,breakdown,anorm));
486: *betam = beta[*m-1];
487: PetscCall(DSRestoreArrayReal(eps->ds,DS_MAT_T,&alpha));
488: break;
489: case EPS_LANCZOS_REORTHOG_PERIODIC:
490: case EPS_LANCZOS_REORTHOG_PARTIAL:
491: PetscCall(DSGetArrayReal(eps->ds,DS_MAT_T,&alpha));
492: beta = alpha + ld;
493: PetscCall(EPSPartialLanczos(eps,alpha,beta,k,m,breakdown,anorm));
494: *betam = beta[*m-1];
495: PetscCall(DSRestoreArrayReal(eps->ds,DS_MAT_T,&alpha));
496: break;
497: case EPS_LANCZOS_REORTHOG_DELAYED:
498: PetscCall(PetscMalloc1(n*n,&T));
499: PetscCall(BVGetOrthogonalization(eps->V,NULL,&orthog_ref,NULL,NULL));
500: if (orthog_ref == BV_ORTHOG_REFINE_NEVER) PetscCall(EPSDelayedArnoldi1(eps,T,n,k,m,betam,breakdown));
501: else PetscCall(EPSDelayedArnoldi(eps,T,n,k,m,betam,breakdown));
502: n = *m;
503: PetscCall(DSGetArrayReal(eps->ds,DS_MAT_T,&alpha));
504: beta = alpha + ld;
505: for (i=k;i<n-1;i++) {
506: alpha[i] = PetscRealPart(T[n*i+i]);
507: beta[i] = PetscRealPart(T[n*i+i+1]);
508: }
509: alpha[n-1] = PetscRealPart(T[n*(n-1)+n-1]);
510: beta[n-1] = *betam;
511: PetscCall(DSRestoreArrayReal(eps->ds,DS_MAT_T,&alpha));
512: PetscCall(PetscFree(T));
513: break;
514: }
515: PetscFunctionReturn(PETSC_SUCCESS);
516: }
518: static PetscErrorCode EPSSolve_Lanczos(EPS eps)
519: {
520: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
521: PetscInt nconv,i,j,k,l,x,n,*perm,restart,ncv=eps->ncv,r,ld;
522: Vec vi,vj,w;
523: Mat U;
524: PetscScalar *Y,*ritz,stmp;
525: PetscReal *bnd,anorm,beta,norm,rtmp,resnorm;
526: PetscBool breakdown;
527: char *conv,ctmp;
529: PetscFunctionBegin;
530: PetscCall(DSGetLeadingDimension(eps->ds,&ld));
531: PetscCall(PetscMalloc4(ncv,&ritz,ncv,&bnd,ncv,&perm,ncv,&conv));
533: /* The first Lanczos vector is the normalized initial vector */
534: PetscCall(EPSGetStartVector(eps,0,NULL));
536: anorm = -1.0;
537: nconv = 0;
539: /* Restart loop */
540: while (eps->reason == EPS_CONVERGED_ITERATING) {
541: eps->its++;
543: /* Compute an ncv-step Lanczos factorization */
544: n = PetscMin(nconv+eps->mpd,ncv);
545: PetscCall(DSSetDimensions(eps->ds,n,nconv,PETSC_DETERMINE));
546: PetscCall(EPSBasicLanczos(eps,nconv,&n,&beta,&breakdown,anorm));
547: PetscCall(DSSetDimensions(eps->ds,n,nconv,0));
548: PetscCall(DSSetState(eps->ds,DS_STATE_INTERMEDIATE));
549: PetscCall(BVSetActiveColumns(eps->V,nconv,n));
551: /* Solve projected problem */
552: PetscCall(DSSolve(eps->ds,ritz,NULL));
553: PetscCall(DSSort(eps->ds,ritz,NULL,NULL,NULL,NULL));
554: PetscCall(DSSynchronize(eps->ds,ritz,NULL));
556: /* Estimate ||A|| */
557: for (i=nconv;i<n;i++)
558: anorm = PetscMax(anorm,PetscAbsReal(PetscRealPart(ritz[i])));
560: /* Compute residual norm estimates as beta*abs(Y(m,:)) + eps*||A|| */
561: PetscCall(DSGetArray(eps->ds,DS_MAT_Q,&Y));
562: for (i=nconv;i<n;i++) {
563: resnorm = beta*PetscAbsScalar(Y[n-1+i*ld]) + PETSC_MACHINE_EPSILON*anorm;
564: PetscCall((*eps->converged)(eps,ritz[i],eps->eigi[i],resnorm,&bnd[i],eps->convergedctx));
565: if (bnd[i]<eps->tol) conv[i] = 'C';
566: else conv[i] = 'N';
567: }
568: PetscCall(DSRestoreArray(eps->ds,DS_MAT_Q,&Y));
570: /* purge repeated ritz values */
571: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
572: for (i=nconv+1;i<n;i++) {
573: if (conv[i] == 'C' && PetscAbsScalar((ritz[i]-ritz[i-1])/ritz[i]) < eps->tol) conv[i] = 'R';
574: }
575: }
577: /* Compute restart vector */
578: if (breakdown) PetscCall(PetscInfo(eps,"Breakdown in Lanczos method (it=%" PetscInt_FMT " norm=%g)\n",eps->its,(double)beta));
579: else {
580: restart = nconv;
581: while (restart<n && conv[restart] != 'N') restart++;
582: if (restart >= n) {
583: breakdown = PETSC_TRUE;
584: } else {
585: for (i=restart+1;i<n;i++) {
586: if (conv[i] == 'N') {
587: PetscCall(SlepcSCCompare(eps->sc,ritz[restart],0.0,ritz[i],0.0,&r));
588: if (r>0) restart = i;
589: }
590: }
591: PetscCall(DSGetArray(eps->ds,DS_MAT_Q,&Y));
592: PetscCall(BVMultColumn(eps->V,1.0,0.0,n,Y+restart*ld+nconv));
593: PetscCall(DSRestoreArray(eps->ds,DS_MAT_Q,&Y));
594: }
595: }
597: /* Count and put converged eigenvalues first */
598: for (i=nconv;i<n;i++) perm[i] = i;
599: for (k=nconv;k<n;k++) {
600: if (conv[perm[k]] != 'C') {
601: j = k + 1;
602: while (j<n && conv[perm[j]] != 'C') j++;
603: if (j>=n) break;
604: l = perm[k]; perm[k] = perm[j]; perm[j] = l;
605: }
606: }
608: /* Sort eigenvectors according to permutation */
609: PetscCall(DSGetArray(eps->ds,DS_MAT_Q,&Y));
610: for (i=nconv;i<k;i++) {
611: x = perm[i];
612: if (x != i) {
613: j = i + 1;
614: while (perm[j] != i) j++;
615: /* swap eigenvalues i and j */
616: stmp = ritz[x]; ritz[x] = ritz[i]; ritz[i] = stmp;
617: rtmp = bnd[x]; bnd[x] = bnd[i]; bnd[i] = rtmp;
618: ctmp = conv[x]; conv[x] = conv[i]; conv[i] = ctmp;
619: perm[j] = x; perm[i] = i;
620: /* swap eigenvectors i and j */
621: for (l=0;l<n;l++) {
622: stmp = Y[l+x*ld]; Y[l+x*ld] = Y[l+i*ld]; Y[l+i*ld] = stmp;
623: }
624: }
625: }
626: PetscCall(DSRestoreArray(eps->ds,DS_MAT_Q,&Y));
628: /* compute converged eigenvectors */
629: PetscCall(DSGetMat(eps->ds,DS_MAT_Q,&U));
630: PetscCall(BVMultInPlace(eps->V,U,nconv,k));
631: PetscCall(DSRestoreMat(eps->ds,DS_MAT_Q,&U));
633: /* purge spurious ritz values */
634: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL) {
635: for (i=nconv;i<k;i++) {
636: PetscCall(BVGetColumn(eps->V,i,&vi));
637: PetscCall(VecNorm(vi,NORM_2,&norm));
638: PetscCall(VecScale(vi,1.0/norm));
639: w = eps->work[0];
640: PetscCall(STApply(eps->st,vi,w));
641: PetscCall(VecAXPY(w,-ritz[i],vi));
642: PetscCall(BVRestoreColumn(eps->V,i,&vi));
643: PetscCall(VecNorm(w,NORM_2,&norm));
644: PetscCall((*eps->converged)(eps,ritz[i],eps->eigi[i],norm,&bnd[i],eps->convergedctx));
645: if (bnd[i]>=eps->tol) conv[i] = 'S';
646: }
647: for (i=nconv;i<k;i++) {
648: if (conv[i] != 'C') {
649: j = i + 1;
650: while (j<k && conv[j] != 'C') j++;
651: if (j>=k) break;
652: /* swap eigenvalues i and j */
653: stmp = ritz[j]; ritz[j] = ritz[i]; ritz[i] = stmp;
654: rtmp = bnd[j]; bnd[j] = bnd[i]; bnd[i] = rtmp;
655: ctmp = conv[j]; conv[j] = conv[i]; conv[i] = ctmp;
656: /* swap eigenvectors i and j */
657: PetscCall(BVGetColumn(eps->V,i,&vi));
658: PetscCall(BVGetColumn(eps->V,j,&vj));
659: PetscCall(VecSwap(vi,vj));
660: PetscCall(BVRestoreColumn(eps->V,i,&vi));
661: PetscCall(BVRestoreColumn(eps->V,j,&vj));
662: }
663: }
664: k = i;
665: }
667: /* store ritz values and estimated errors */
668: for (i=nconv;i<n;i++) {
669: eps->eigr[i] = ritz[i];
670: eps->errest[i] = bnd[i];
671: }
672: nconv = k;
673: PetscCall(EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,n));
674: PetscCall((*eps->stopping)(eps,eps->its,eps->max_it,nconv,eps->nev,&eps->reason,eps->stoppingctx));
676: if (eps->reason == EPS_CONVERGED_ITERATING) { /* copy restart vector */
677: PetscCall(BVCopyColumn(eps->V,n,nconv));
678: if (lanczos->reorthog == EPS_LANCZOS_REORTHOG_LOCAL && !breakdown) {
679: /* Reorthonormalize restart vector */
680: PetscCall(BVOrthonormalizeColumn(eps->V,nconv,PETSC_FALSE,NULL,&breakdown));
681: }
682: if (breakdown) {
683: /* Use random vector for restarting */
684: PetscCall(PetscInfo(eps,"Using random vector for restart\n"));
685: PetscCall(EPSGetStartVector(eps,nconv,&breakdown));
686: }
687: if (PetscUnlikely(breakdown)) { /* give up */
688: eps->reason = EPS_DIVERGED_BREAKDOWN;
689: PetscCall(PetscInfo(eps,"Unable to generate more start vectors\n"));
690: }
691: }
692: }
693: eps->nconv = nconv;
695: PetscCall(PetscFree4(ritz,bnd,perm,conv));
696: PetscFunctionReturn(PETSC_SUCCESS);
697: }
699: static PetscErrorCode EPSSetFromOptions_Lanczos(EPS eps,PetscOptionItems *PetscOptionsObject)
700: {
701: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
702: PetscBool flg;
703: EPSLanczosReorthogType reorthog=EPS_LANCZOS_REORTHOG_LOCAL,curval;
705: PetscFunctionBegin;
706: PetscOptionsHeadBegin(PetscOptionsObject,"EPS Lanczos Options");
708: curval = (lanczos->reorthog==(EPSLanczosReorthogType)-1)? EPS_LANCZOS_REORTHOG_LOCAL: lanczos->reorthog;
709: PetscCall(PetscOptionsEnum("-eps_lanczos_reorthog","Lanczos reorthogonalization","EPSLanczosSetReorthog",EPSLanczosReorthogTypes,(PetscEnum)curval,(PetscEnum*)&reorthog,&flg));
710: if (flg) PetscCall(EPSLanczosSetReorthog(eps,reorthog));
712: PetscOptionsHeadEnd();
713: PetscFunctionReturn(PETSC_SUCCESS);
714: }
716: static PetscErrorCode EPSLanczosSetReorthog_Lanczos(EPS eps,EPSLanczosReorthogType reorthog)
717: {
718: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
720: PetscFunctionBegin;
721: switch (reorthog) {
722: case EPS_LANCZOS_REORTHOG_LOCAL:
723: case EPS_LANCZOS_REORTHOG_FULL:
724: case EPS_LANCZOS_REORTHOG_DELAYED:
725: case EPS_LANCZOS_REORTHOG_SELECTIVE:
726: case EPS_LANCZOS_REORTHOG_PERIODIC:
727: case EPS_LANCZOS_REORTHOG_PARTIAL:
728: if (lanczos->reorthog != reorthog) {
729: lanczos->reorthog = reorthog;
730: eps->state = EPS_STATE_INITIAL;
731: }
732: break;
733: default:
734: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
735: }
736: PetscFunctionReturn(PETSC_SUCCESS);
737: }
739: /*@
740: EPSLanczosSetReorthog - Sets the type of reorthogonalization used during the Lanczos
741: iteration.
743: Logically Collective
745: Input Parameters:
746: + eps - the eigenproblem solver context
747: - reorthog - the type of reorthogonalization
749: Options Database Key:
750: . -eps_lanczos_reorthog - Sets the reorthogonalization type (either 'local', 'selective',
751: 'periodic', 'partial', 'full' or 'delayed')
753: Level: advanced
755: .seealso: EPSLanczosGetReorthog(), EPSLanczosReorthogType
756: @*/
757: PetscErrorCode EPSLanczosSetReorthog(EPS eps,EPSLanczosReorthogType reorthog)
758: {
759: PetscFunctionBegin;
762: PetscTryMethod(eps,"EPSLanczosSetReorthog_C",(EPS,EPSLanczosReorthogType),(eps,reorthog));
763: PetscFunctionReturn(PETSC_SUCCESS);
764: }
766: static PetscErrorCode EPSLanczosGetReorthog_Lanczos(EPS eps,EPSLanczosReorthogType *reorthog)
767: {
768: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
770: PetscFunctionBegin;
771: *reorthog = lanczos->reorthog;
772: PetscFunctionReturn(PETSC_SUCCESS);
773: }
775: /*@
776: EPSLanczosGetReorthog - Gets the type of reorthogonalization used during
777: the Lanczos iteration.
779: Not Collective
781: Input Parameter:
782: . eps - the eigenproblem solver context
784: Output Parameter:
785: . reorthog - the type of reorthogonalization
787: Level: advanced
789: .seealso: EPSLanczosSetReorthog(), EPSLanczosReorthogType
790: @*/
791: PetscErrorCode EPSLanczosGetReorthog(EPS eps,EPSLanczosReorthogType *reorthog)
792: {
793: PetscFunctionBegin;
795: PetscAssertPointer(reorthog,2);
796: PetscUseMethod(eps,"EPSLanczosGetReorthog_C",(EPS,EPSLanczosReorthogType*),(eps,reorthog));
797: PetscFunctionReturn(PETSC_SUCCESS);
798: }
800: static PetscErrorCode EPSReset_Lanczos(EPS eps)
801: {
802: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
804: PetscFunctionBegin;
805: PetscCall(BVDestroy(&lanczos->AV));
806: lanczos->allocsize = 0;
807: PetscFunctionReturn(PETSC_SUCCESS);
808: }
810: static PetscErrorCode EPSDestroy_Lanczos(EPS eps)
811: {
812: PetscFunctionBegin;
813: PetscCall(PetscFree(eps->data));
814: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",NULL));
815: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",NULL));
816: PetscFunctionReturn(PETSC_SUCCESS);
817: }
819: static PetscErrorCode EPSView_Lanczos(EPS eps,PetscViewer viewer)
820: {
821: EPS_LANCZOS *lanczos = (EPS_LANCZOS*)eps->data;
822: PetscBool isascii;
824: PetscFunctionBegin;
825: PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
826: if (isascii) {
827: if (lanczos->reorthog != (EPSLanczosReorthogType)-1) PetscCall(PetscViewerASCIIPrintf(viewer," %s reorthogonalization\n",EPSLanczosReorthogTypes[lanczos->reorthog]));
828: }
829: PetscFunctionReturn(PETSC_SUCCESS);
830: }
832: SLEPC_EXTERN PetscErrorCode EPSCreate_Lanczos(EPS eps)
833: {
834: EPS_LANCZOS *ctx;
836: PetscFunctionBegin;
837: PetscCall(PetscNew(&ctx));
838: eps->data = (void*)ctx;
839: ctx->reorthog = (EPSLanczosReorthogType)-1;
841: eps->useds = PETSC_TRUE;
843: eps->ops->solve = EPSSolve_Lanczos;
844: eps->ops->setup = EPSSetUp_Lanczos;
845: eps->ops->setupsort = EPSSetUpSort_Default;
846: eps->ops->setfromoptions = EPSSetFromOptions_Lanczos;
847: eps->ops->destroy = EPSDestroy_Lanczos;
848: eps->ops->reset = EPSReset_Lanczos;
849: eps->ops->view = EPSView_Lanczos;
850: eps->ops->backtransform = EPSBackTransform_Default;
851: eps->ops->computevectors = EPSComputeVectors_Hermitian;
853: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",EPSLanczosSetReorthog_Lanczos));
854: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",EPSLanczosGetReorthog_Lanczos));
855: PetscFunctionReturn(PETSC_SUCCESS);
856: }