Actual source code: ks-lrep.c
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: "krylovschur"
13: Method: thick-restarted Lanczos for Linear Response eigenvalue problems
15: References:
17: [1] Z. Teng, R.-C. Li, "Convergence analysis of Lanczos-type methods for the
18: linear response eigenvalue problem", J. Comput. Appl. Math. 247, 2013.
20: [2] F. Alvarruiz, B. Mellado-Pinto, J. E. Roman, "Restarted Lanczos methods
21: for the linear response eigenvalue problem", in preparation, 2026.
23: */
24: #include <slepc/private/epsimpl.h>
25: #include "krylovschur.h"
27: static PetscErrorCode Orthog_Teng(Vec x,BV U,BV V,PetscInt j,PetscScalar *h,PetscScalar *c,PetscBool *breakdown)
28: {
29: PetscInt i;
31: PetscFunctionBegin;
32: PetscCall(BVSetActiveColumns(U,0,j));
33: PetscCall(BVSetActiveColumns(V,0,j));
34: /* c = U^* x */
35: PetscCall(BVDotVec(U,x,c));
36: /* x = x-V*c */
37: PetscCall(BVMultVec(V,-1.0,1.0,x,c));
38: /* accumulate orthog coeffs into h */
39: for (i=0;i<2*j;i++) h[i] += c[i];
40: PetscFunctionReturn(PETSC_SUCCESS);
41: }
43: /* Orthogonalize vector x against first j vectors in U and V
44: v is column j-1 of V */
45: static PetscErrorCode OrthogonalizeVector_Teng(Vec x,BV U,BV V,PetscInt j,Vec u,PetscReal *beta,PetscInt k,PetscScalar *h,PetscBool *breakdown)
46: {
47: PetscReal alpha;
48: PetscInt i,l;
50: PetscFunctionBegin;
51: PetscCall(PetscArrayzero(h,2*j));
53: /* Local orthogonalization */
54: l = j==k+1?0:j-2; /* 1st column to orthogonalize against */
55: PetscCall(VecDotRealPart(x,u,&alpha));
56: for (i=l;i<j-1;i++) h[i] = beta[i];
57: h[j-1] = alpha;
58: /* x = x-V(:,l:j-1)*h(l:j-1) */
59: PetscCall(BVSetActiveColumns(V,l,j));
60: PetscCall(BVMultVec(V,-1.0,1.0,x,h+l));
62: /* Full orthogonalization */
63: PetscCall(Orthog_Teng(x,U,V,j,h,h+2*j,breakdown));
64: PetscFunctionReturn(PETSC_SUCCESS);
65: }
67: static PetscErrorCode EPSLREPLanczos_Teng(EPS eps,Mat K,Mat M,BV U,BV V,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscInt *min,PetscBool *breakdown)
68: {
69: PetscInt j,m = *min;
70: Vec u,v,uh,vh;
71: PetscReal beta0;
72: PetscScalar *hwork,lhwork[100],gamma;
74: PetscFunctionBegin;
75: if (4*m > 100) PetscCall(PetscMalloc1(4*m,&hwork));
76: else hwork = lhwork;
78: /* Normalize initial vector */
79: if (k==0) {
80: if (eps->nini==0) PetscCall(BVSetRandomColumn(V,0));
81: PetscCall(BVGetColumn(U,0,&u));
82: PetscCall(BVGetColumn(V,0,&v));
83: PetscCall(MatMult(M,v,u));
84: PetscCall(VecDot(u,v,&gamma));
85: beta0 = PetscSqrtReal(PetscRealPart(gamma)); /* TODO check beta_0=0 */
86: PetscCall(VecScale(u,1.0/beta0));
87: PetscCall(VecScale(v,1.0/beta0));
88: PetscCall(BVRestoreColumn(U,0,&u));
89: PetscCall(BVRestoreColumn(V,0,&v));
90: }
92: for (j=k;j<m;j++) {
93: /* j+1 columns (indices 0 to j) have been computed */
94: PetscCall(BVGetColumn(U,j+1,&uh));
95: PetscCall(BVGetColumn(V,j+1,&vh));
96: PetscCall(BVGetColumn(U,j,&u));
97: PetscCall(MatMult(K,u,vh));
98: PetscCall(OrthogonalizeVector_Teng(vh,U,V,j+1,u,beta,k,hwork,breakdown));
99: alpha[j] = PetscRealPart(hwork[j]);
100: PetscCall(MatMult(M,vh,uh));
101: PetscCall(VecDot(uh,vh,&gamma));
102: beta[j] = PetscSqrtReal(PetscRealPart(gamma)); /* TODO check beta_j=0 */
103: PetscCall(VecScale(uh,1.0/beta[j]));
104: PetscCall(VecScale(vh,1.0/beta[j]));
105: PetscCall(BVRestoreColumn(U,j+1,&uh));
106: PetscCall(BVRestoreColumn(V,j+1,&vh));
107: PetscCall(BVRestoreColumn(U,j,&u));
108: }
109: if (4*m > 100) PetscCall(PetscFree(hwork));
110: PetscFunctionReturn(PETSC_SUCCESS);
111: }
113: static PetscErrorCode EPSComputeVectors_LREP_Teng(EPS eps)
114: {
115: Mat H;
116: Vec u,v,w;
117: BV U,V;
118: IS is[2];
119: PetscInt k;
120: PetscScalar lambda;
121: PetscBool reduced;
123: PetscFunctionBegin;
124: PetscCall(STGetMatrix(eps->st,0,&H));
125: PetscCall(MatNestGetISs(H,is,NULL));
126: PetscCall(SlepcCheckMatLREPReduced(H,&reduced));
127: PetscCall(BVGetSplitRows(eps->V,is[0],is[1],&V,&U));
128: for (k=0;k<eps->nconv;k++) {
129: PetscCall(BVGetColumn(V,k,&v));
130: /* approx eigenvector is [eigr[k]*v; u] */
131: lambda = eps->eigr[k];
132: PetscCall(STBackTransform(eps->st,1,&lambda,&eps->eigi[k]));
133: PetscCall(VecScale(v,lambda));
134: PetscCall(BVRestoreColumn(V,k,&v));
135: }
136: if (!reduced) {
137: /* the eigenvector [v;u] = J*[y;x] where [y;x] is the reduced eigenvector
138: and J = 1/sqrt(2)[I I; I -I], i.e, v=1/sqrt(2)*(y+x) u=1/sqrt(2)*(y-x) */
139: PetscCall(BVCreateVec(V,&w));
140: for (k=0;k<eps->nconv;k++) {
141: PetscCall(BVGetColumn(U,k,&u));
142: PetscCall(BVGetColumn(V,k,&v));
143: PetscCall(VecCopy(u,w));
144: PetscCall(VecCopy(v,u));
145: PetscCall(VecAXPY(u,-1.0,w));
146: PetscCall(VecAXPY(v,1.0,w));
147: PetscCall(BVRestoreColumn(U,k,&u));
148: PetscCall(BVRestoreColumn(V,k,&v));
149: }
150: PetscCall(VecDestroy(&w));
151: }
152: PetscCall(BVRestoreSplitRows(eps->V,is[0],is[1],&V,&U));
153: /* Normalize eigenvectors */
154: PetscCall(BVSetActiveColumns(eps->V,0,eps->nconv));
155: PetscCall(BVNormalize(eps->V,NULL));
156: PetscFunctionReturn(PETSC_SUCCESS);
157: }
159: PetscErrorCode EPSSetUp_KrylovSchur_LREP(EPS eps)
160: {
161: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
162: PetscBool flg;
164: PetscFunctionBegin;
165: PetscCheck((eps->problem_type==EPS_LREP),PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Problem type should be LREP");
166: EPSCheckUnsupportedCondition(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_EXTRACTION | EPS_FEATURE_BALANCE,PETSC_TRUE," with LREP structure");
167: PetscCall(EPSSetDimensions_Default(eps,&eps->nev,&eps->ncv,&eps->mpd));
168: 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");
169: if (eps->max_it==PETSC_DETERMINE) eps->max_it = PetscMax(100,2*eps->n/eps->ncv)*((eps->stop==EPS_STOP_THRESHOLD)?10:1);
171: PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&flg));
172: PetscCheck(flg,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Krylov-Schur LREP only supports shift ST");
173: if (!eps->which) eps->which = EPS_SMALLEST_MAGNITUDE;
175: if (!ctx->keep) ctx->keep = 0.5;
176: PetscCall(STSetStructured(eps->st,PETSC_FALSE));
178: PetscCall(EPSAllocateSolution(eps,1));
179: /* Teng */
180: eps->ops->solve = EPSSolve_KrylovSchur_LREP_Teng;
181: eps->ops->computevectors = EPSComputeVectors_LREP_Teng;
182: PetscCall(DSSetType(eps->ds,DSHEP));
183: PetscCall(DSSetCompact(eps->ds,PETSC_TRUE));
184: PetscCall(DSSetExtraRow(eps->ds,PETSC_TRUE));
185: PetscCall(DSAllocate(eps->ds,eps->ncv+1));
186: PetscFunctionReturn(PETSC_SUCCESS);
187: }
189: PetscErrorCode EPSSolve_KrylovSchur_LREP_Teng(EPS eps)
190: {
191: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
192: PetscInt i,k,l,ld,nv,nconv=0,nevsave,ma,na,Ma,Na;
193: Mat H,Q,K,M,A,B;
194: BV U,V;
195: IS is[2];
196: PetscReal *a,*b,beta;
197: PetscBool reduced,breakdown=PETSC_FALSE;
198: const PetscScalar scal[] = { 1.0, -1.0 };
200: PetscFunctionBegin;
201: PetscCall(DSGetLeadingDimension(eps->ds,&ld));
203: /* Extract matrix blocks */
204: PetscCall(STGetMatrix(eps->st,0,&H));
205: PetscCall(MatNestGetISs(H,is,NULL));
206: PetscCall(SlepcCheckMatLREPReduced(H,&reduced));
207: if (reduced) {
208: PetscCall(MatNestGetSubMat(H,0,1,&K));
209: PetscCall(MatNestGetSubMat(H,1,0,&M));
210: } else {
211: PetscCall(MatNestGetSubMat(H,0,0,&A));
212: PetscCall(MatNestGetSubMat(H,0,1,&B));
213: PetscCall(MatGetSize(A,&Ma,&Na));
214: PetscCall(MatGetLocalSize(A,&ma,&na));
215: /* K = A-B */
216: PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&K));
217: PetscCall(MatSetSizes(K,ma,na,Ma,Na));
218: PetscCall(MatSetType(K,MATCOMPOSITE));
219: PetscCall(MatCompositeAddMat(K,A));
220: PetscCall(MatCompositeAddMat(K,B));
221: PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
222: PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
223: PetscCall(MatCompositeSetScalings(K,scal));
224: /* M = A+B */
225: PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&M));
226: PetscCall(MatSetSizes(M,ma,na,Ma,Na));
227: PetscCall(MatSetType(M,MATCOMPOSITE));
228: PetscCall(MatCompositeAddMat(M,A));
229: PetscCall(MatCompositeAddMat(M,B));
230: PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
231: PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
232: }
234: /* Get the split bases */
235: PetscCall(BVGetSplitRows(eps->V,is[0],is[1],&V,&U));
237: nevsave = eps->nev;
238: eps->nev = (eps->nev+1)/2;
239: l = 0;
241: /* Restart loop */
242: while (eps->reason == EPS_CONVERGED_ITERATING) {
243: eps->its++;
245: /* Compute an nv-step Lanczos factorization */
246: nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
247: PetscCall(DSSetDimensions(eps->ds,nv,eps->nconv,eps->nconv+l));
248: PetscCall(DSGetArrayReal(eps->ds,DS_MAT_T,&a));
249: b = a + ld;
250: PetscCall(EPSLREPLanczos_Teng(eps,K,M,U,V,a,b,eps->nconv+l,&nv,&breakdown));
251: beta = b[nv-1];
252: PetscCall(DSRestoreArrayReal(eps->ds,DS_MAT_T,&a));
253: PetscCall(DSSetDimensions(eps->ds,nv,eps->nconv,eps->nconv+l));
254: PetscCall(DSSetState(eps->ds,l?DS_STATE_RAW:DS_STATE_INTERMEDIATE));
255: PetscCall(BVSetActiveColumns(eps->V,eps->nconv,nv));
257: /* Solve projected problem */
258: PetscCall(DSSolve(eps->ds,eps->eigr,eps->eigi));
259: PetscCall(DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL));
260: PetscCall(DSUpdateExtraRow(eps->ds));
261: PetscCall(DSSynchronize(eps->ds,eps->eigr,eps->eigi));
263: /* Check convergence */
264: for (i=0;i<nv;i++) eps->eigr[i] = PetscSqrtReal(PetscRealPart(eps->eigr[i]));
265: PetscCall(EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta,0.0,1.0,&k));
266: EPSSetCtxThreshold(eps,eps->eigr,eps->eigi,eps->errest,k,nv);
267: PetscCall((*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx));
268: nconv = k;
270: /* Update l */
271: if (eps->reason != EPS_CONVERGED_ITERATING || breakdown || k==nv) l = 0;
272: else l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
273: if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
274: if (l) PetscCall(PetscInfo(eps,"Preparing to restart keeping l=%" PetscInt_FMT " vectors\n",l));
276: if (eps->reason == EPS_CONVERGED_ITERATING) {
277: PetscCheck(!breakdown,PetscObjectComm((PetscObject)eps),PETSC_ERR_CONV_FAILED,"Breakdown in LREP Krylov-Schur (beta=%g)",(double)beta);
278: /* Prepare the Rayleigh quotient for restart */
279: PetscCall(DSTruncate(eps->ds,k+l,PETSC_FALSE));
280: }
281: /* Update the corresponding vectors
282: U(:,idx) = U*Q(:,idx), V(:,idx) = V*Q(:,idx) */
283: PetscCall(DSGetMat(eps->ds,DS_MAT_Q,&Q));
284: PetscCall(BVMultInPlace(U,Q,eps->nconv,k+l));
285: PetscCall(BVMultInPlace(V,Q,eps->nconv,k+l));
286: PetscCall(DSRestoreMat(eps->ds,DS_MAT_Q,&Q));
288: if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
289: PetscCall(BVCopyColumn(eps->V,nv,k+l));
290: if (eps->stop==EPS_STOP_THRESHOLD && nv-k<5) { /* reallocate */
291: eps->ncv = eps->mpd+k;
292: PetscCall(BVRestoreSplitRows(eps->V,is[0],is[1],&V,&U));
293: PetscCall(EPSReallocateSolution(eps,eps->ncv+1));
294: PetscCall(BVGetSplitRows(eps->V,is[0],is[1],&V,&U));
295: for (i=nv;i<eps->ncv;i++) eps->perm[i] = i;
296: PetscCall(DSReallocate(eps->ds,eps->ncv+1));
297: PetscCall(DSGetLeadingDimension(eps->ds,&ld));
298: }
299: }
300: eps->nconv = k;
301: PetscCall(EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nv));
302: }
304: eps->nev = nevsave;
306: PetscCall(DSTruncate(eps->ds,eps->nconv,PETSC_TRUE));
307: PetscCall(BVRestoreSplitRows(eps->V,is[0],is[1],&V,&U));
308: if (!reduced) {
309: PetscCall(MatDestroy(&K));
310: PetscCall(MatDestroy(&M));
311: }
312: PetscFunctionReturn(PETSC_SUCCESS);
313: }