Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
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: "lyapii"
12 :
13 : Method: Lyapunov inverse iteration
14 :
15 : Algorithm:
16 :
17 : Lyapunov inverse iteration using LME solvers
18 :
19 : References:
20 :
21 : [1] H.C. Elman and M. Wu, "Lyapunov inverse iteration for computing a
22 : few rightmost eigenvalues of large generalized eigenvalue problems",
23 : SIAM J. Matrix Anal. Appl. 34(4):1685-1707, 2013.
24 :
25 : [2] K. Meerbergen and A. Spence, "Inverse iteration for purely imaginary
26 : eigenvalues with application to the detection of Hopf bifurcations in
27 : large-scale problems", SIAM J. Matrix Anal. Appl. 31:1982-1999, 2010.
28 : */
29 :
30 : #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
31 : #include <slepcblaslapack.h>
32 :
33 : typedef struct {
34 : LME lme; /* Lyapunov solver */
35 : DS ds; /* used to compute the SVD for compression */
36 : PetscInt rkl; /* prescribed rank for the Lyapunov solver */
37 : PetscInt rkc; /* the compressed rank, cannot be larger than rkl */
38 : } EPS_LYAPII;
39 :
40 : typedef struct {
41 : Mat S; /* the operator matrix, S=A^{-1}*B */
42 : BV Q; /* orthogonal basis of converged eigenvectors */
43 : } EPS_LYAPII_MATSHELL;
44 :
45 : typedef struct {
46 : Mat S; /* the matrix from which the implicit operator is built */
47 : PetscInt n; /* the size of matrix S, the operator is nxn */
48 : LME lme; /* dummy LME object */
49 : #if defined(PETSC_USE_COMPLEX)
50 : Mat A,B,F;
51 : Vec w;
52 : #endif
53 : } EPS_EIG_MATSHELL;
54 :
55 3 : static PetscErrorCode EPSSetUp_LyapII(EPS eps)
56 : {
57 3 : PetscRandom rand;
58 3 : EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
59 :
60 3 : PetscFunctionBegin;
61 3 : EPSCheckSinvert(eps);
62 3 : if (eps->ncv!=PETSC_DEFAULT) {
63 0 : PetscCheck(eps->ncv>=eps->nev+1,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv must be at least nev+1");
64 3 : } else eps->ncv = eps->nev+1;
65 3 : if (eps->mpd!=PETSC_DEFAULT) PetscCall(PetscInfo(eps,"Warning: parameter mpd ignored\n"));
66 3 : if (eps->max_it==PETSC_DEFAULT) eps->max_it = PetscMax(1000*eps->nev,100*eps->n);
67 3 : if (!eps->which) eps->which=EPS_LARGEST_REAL;
68 3 : PetscCheck(eps->which==EPS_LARGEST_REAL,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver supports only largest real eigenvalues");
69 3 : EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_EXTRACTION | EPS_FEATURE_TWOSIDED);
70 :
71 3 : if (!ctx->rkc) ctx->rkc = 10;
72 3 : if (!ctx->rkl) ctx->rkl = 3*ctx->rkc;
73 3 : if (!ctx->lme) PetscCall(EPSLyapIIGetLME(eps,&ctx->lme));
74 3 : PetscCall(LMESetProblemType(ctx->lme,LME_LYAPUNOV));
75 3 : PetscCall(LMESetErrorIfNotConverged(ctx->lme,PETSC_TRUE));
76 :
77 3 : if (!ctx->ds) {
78 3 : PetscCall(DSCreate(PetscObjectComm((PetscObject)eps),&ctx->ds));
79 3 : PetscCall(DSSetType(ctx->ds,DSSVD));
80 : }
81 3 : PetscCall(DSAllocate(ctx->ds,ctx->rkl));
82 :
83 3 : PetscCall(DSSetType(eps->ds,DSNHEP));
84 3 : PetscCall(DSAllocate(eps->ds,eps->ncv));
85 :
86 3 : PetscCall(EPSAllocateSolution(eps,0));
87 3 : PetscCall(BVGetRandomContext(eps->V,&rand)); /* make sure the random context is available when duplicating */
88 3 : PetscCall(EPSSetWorkVecs(eps,3));
89 3 : PetscFunctionReturn(PETSC_SUCCESS);
90 : }
91 :
92 2905 : static PetscErrorCode MatMult_EPSLyapIIOperator(Mat M,Vec x,Vec r)
93 : {
94 2905 : EPS_LYAPII_MATSHELL *matctx;
95 :
96 2905 : PetscFunctionBegin;
97 2905 : PetscCall(MatShellGetContext(M,&matctx));
98 2905 : PetscCall(MatMult(matctx->S,x,r));
99 2905 : PetscCall(BVOrthogonalizeVec(matctx->Q,r,NULL,NULL,NULL));
100 2905 : PetscFunctionReturn(PETSC_SUCCESS);
101 : }
102 :
103 3 : static PetscErrorCode MatDestroy_EPSLyapIIOperator(Mat M)
104 : {
105 3 : EPS_LYAPII_MATSHELL *matctx;
106 :
107 3 : PetscFunctionBegin;
108 3 : PetscCall(MatShellGetContext(M,&matctx));
109 3 : PetscCall(MatDestroy(&matctx->S));
110 3 : PetscCall(PetscFree(matctx));
111 3 : PetscFunctionReturn(PETSC_SUCCESS);
112 : }
113 :
114 619 : static PetscErrorCode MatMult_EigOperator(Mat M,Vec x,Vec y)
115 : {
116 619 : EPS_EIG_MATSHELL *matctx;
117 : #if !defined(PETSC_USE_COMPLEX)
118 619 : PetscInt n,lds;
119 619 : PetscScalar *Y,*C,zero=0.0,done=1.0,dtwo=2.0;
120 619 : const PetscScalar *S,*X;
121 619 : PetscBLASInt n_,lds_;
122 : #endif
123 :
124 619 : PetscFunctionBegin;
125 619 : PetscCall(MatShellGetContext(M,&matctx));
126 :
127 : #if defined(PETSC_USE_COMPLEX)
128 : PetscCall(MatMult(matctx->B,x,matctx->w));
129 : PetscCall(MatSolve(matctx->F,matctx->w,y));
130 : #else
131 619 : PetscCall(VecGetArrayRead(x,&X));
132 619 : PetscCall(VecGetArray(y,&Y));
133 619 : PetscCall(MatDenseGetArrayRead(matctx->S,&S));
134 619 : PetscCall(MatDenseGetLDA(matctx->S,&lds));
135 :
136 619 : n = matctx->n;
137 619 : PetscCall(PetscCalloc1(n*n,&C));
138 619 : PetscCall(PetscBLASIntCast(n,&n_));
139 619 : PetscCall(PetscBLASIntCast(lds,&lds_));
140 :
141 : /* C = 2*S*X*S.' */
142 619 : PetscCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&dtwo,S,&lds_,X,&n_,&zero,Y,&n_));
143 619 : PetscCallBLAS("BLASgemm",BLASgemm_("N","T",&n_,&n_,&n_,&done,Y,&n_,S,&lds_,&zero,C,&n_));
144 :
145 : /* Solve S*Y + Y*S' = -C */
146 619 : PetscCall(LMEDenseLyapunov(matctx->lme,n,(PetscScalar*)S,lds,C,n,Y,n));
147 :
148 619 : PetscCall(PetscFree(C));
149 619 : PetscCall(VecRestoreArrayRead(x,&X));
150 619 : PetscCall(VecRestoreArray(y,&Y));
151 619 : PetscCall(MatDenseRestoreArrayRead(matctx->S,&S));
152 : #endif
153 619 : PetscFunctionReturn(PETSC_SUCCESS);
154 : }
155 :
156 5 : static PetscErrorCode MatDestroy_EigOperator(Mat M)
157 : {
158 5 : EPS_EIG_MATSHELL *matctx;
159 :
160 5 : PetscFunctionBegin;
161 5 : PetscCall(MatShellGetContext(M,&matctx));
162 : #if defined(PETSC_USE_COMPLEX)
163 : PetscCall(MatDestroy(&matctx->A));
164 : PetscCall(MatDestroy(&matctx->B));
165 : PetscCall(MatDestroy(&matctx->F));
166 : PetscCall(VecDestroy(&matctx->w));
167 : #else
168 5 : PetscCall(MatDestroy(&matctx->S));
169 : #endif
170 5 : PetscCall(PetscFree(matctx));
171 5 : PetscFunctionReturn(PETSC_SUCCESS);
172 : }
173 :
174 : /*
175 : EV2x2: solve the eigenproblem for a 2x2 matrix M
176 : */
177 23 : static PetscErrorCode EV2x2(PetscScalar *M,PetscInt ld,PetscScalar *wr,PetscScalar *wi,PetscScalar *vec)
178 : {
179 23 : PetscBLASInt lwork=10,ld_;
180 23 : PetscScalar work[10];
181 23 : PetscBLASInt two=2,info;
182 : #if defined(PETSC_USE_COMPLEX)
183 : PetscReal rwork[6];
184 : #endif
185 :
186 23 : PetscFunctionBegin;
187 23 : PetscCall(PetscBLASIntCast(ld,&ld_));
188 23 : PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
189 : #if !defined(PETSC_USE_COMPLEX)
190 23 : PetscCallBLAS("LAPACKgeev",LAPACKgeev_("N","V",&two,M,&ld_,wr,wi,NULL,&ld_,vec,&ld_,work,&lwork,&info));
191 : #else
192 : PetscCallBLAS("LAPACKgeev",LAPACKgeev_("N","V",&two,M,&ld_,wr,NULL,&ld_,vec,&ld_,work,&lwork,rwork,&info));
193 : #endif
194 23 : SlepcCheckLapackInfo("geev",info);
195 23 : PetscCall(PetscFPTrapPop());
196 23 : PetscFunctionReturn(PETSC_SUCCESS);
197 : }
198 :
199 : /*
200 : LyapIIBuildRHS: prepare the right-hand side of the Lyapunov equation SY + YS' = -2*S*Z*S'
201 : in factored form:
202 : if (V) U=sqrt(2)*S*V (uses 1 work vector)
203 : else U=sqrt(2)*S*U (uses 2 work vectors)
204 : where U,V are assumed to have rk columns.
205 : */
206 26 : static PetscErrorCode LyapIIBuildRHS(Mat S,PetscInt rk,Mat U,BV V,Vec *work)
207 : {
208 26 : PetscScalar *array,*uu;
209 26 : PetscInt i,nloc;
210 26 : Vec v,u=work[0];
211 :
212 26 : PetscFunctionBegin;
213 26 : PetscCall(MatGetLocalSize(U,&nloc,NULL));
214 69 : for (i=0;i<rk;i++) {
215 43 : PetscCall(MatDenseGetColumn(U,i,&array));
216 43 : if (V) PetscCall(BVGetColumn(V,i,&v));
217 : else {
218 36 : v = work[1];
219 36 : PetscCall(VecPlaceArray(v,array));
220 : }
221 43 : PetscCall(MatMult(S,v,u));
222 43 : if (V) PetscCall(BVRestoreColumn(V,i,&v));
223 36 : else PetscCall(VecResetArray(v));
224 43 : PetscCall(VecScale(u,PETSC_SQRT2));
225 43 : PetscCall(VecGetArray(u,&uu));
226 43 : PetscCall(PetscArraycpy(array,uu,nloc));
227 43 : PetscCall(VecRestoreArray(u,&uu));
228 43 : PetscCall(MatDenseRestoreColumn(U,&array));
229 : }
230 26 : PetscFunctionReturn(PETSC_SUCCESS);
231 : }
232 :
233 : /*
234 : LyapIIBuildEigenMat: create shell matrix Op=A\B with A = kron(I,S)+kron(S,I), B = -2*kron(S,S)
235 : where S is a sequential square dense matrix of order n.
236 : v0 is the initial vector, should have the form v0 = w*w' (for instance 1*1')
237 : */
238 26 : static PetscErrorCode LyapIIBuildEigenMat(LME lme,Mat S,Mat *Op,Vec *v0)
239 : {
240 26 : PetscInt n,m;
241 26 : PetscBool create=PETSC_FALSE;
242 26 : EPS_EIG_MATSHELL *matctx;
243 : #if defined(PETSC_USE_COMPLEX)
244 : PetscScalar theta,*aa,*bb;
245 : const PetscScalar *ss;
246 : PetscInt i,j,f,c,off,ld,lds;
247 : IS perm;
248 : #endif
249 :
250 26 : PetscFunctionBegin;
251 26 : PetscCall(MatGetSize(S,&n,NULL));
252 26 : if (!*Op) create=PETSC_TRUE;
253 : else {
254 23 : PetscCall(MatGetSize(*Op,&m,NULL));
255 23 : if (m!=n*n) create=PETSC_TRUE;
256 : }
257 21 : if (create) {
258 5 : PetscCall(MatDestroy(Op));
259 5 : PetscCall(VecDestroy(v0));
260 5 : PetscCall(PetscNew(&matctx));
261 : #if defined(PETSC_USE_COMPLEX)
262 : PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n*n,n*n,NULL,&matctx->A));
263 : PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n*n,n*n,NULL,&matctx->B));
264 : PetscCall(MatCreateVecs(matctx->A,NULL,&matctx->w));
265 : #else
266 5 : PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&matctx->S));
267 : #endif
268 5 : PetscCall(MatCreateShell(PETSC_COMM_SELF,n*n,n*n,PETSC_DETERMINE,PETSC_DETERMINE,matctx,Op));
269 5 : PetscCall(MatShellSetOperation(*Op,MATOP_MULT,(void(*)(void))MatMult_EigOperator));
270 5 : PetscCall(MatShellSetOperation(*Op,MATOP_DESTROY,(void(*)(void))MatDestroy_EigOperator));
271 5 : PetscCall(MatCreateVecs(*Op,NULL,v0));
272 : } else {
273 21 : PetscCall(MatShellGetContext(*Op,&matctx));
274 : #if defined(PETSC_USE_COMPLEX)
275 : PetscCall(MatZeroEntries(matctx->A));
276 : #endif
277 : }
278 : #if defined(PETSC_USE_COMPLEX)
279 : PetscCall(MatDenseGetArray(matctx->A,&aa));
280 : PetscCall(MatDenseGetArray(matctx->B,&bb));
281 : PetscCall(MatDenseGetArrayRead(S,&ss));
282 : PetscCall(MatDenseGetLDA(S,&lds));
283 : ld = n*n;
284 : for (f=0;f<n;f++) {
285 : off = f*n+f*n*ld;
286 : for (i=0;i<n;i++) for (j=0;j<n;j++) aa[off+i+j*ld] = ss[i+j*lds];
287 : for (c=0;c<n;c++) {
288 : off = f*n+c*n*ld;
289 : theta = ss[f+c*lds];
290 : for (i=0;i<n;i++) aa[off+i+i*ld] += theta;
291 : for (i=0;i<n;i++) for (j=0;j<n;j++) bb[off+i+j*ld] = -2*theta*ss[i+j*lds];
292 : }
293 : }
294 : PetscCall(MatDenseRestoreArray(matctx->A,&aa));
295 : PetscCall(MatDenseRestoreArray(matctx->B,&bb));
296 : PetscCall(MatDenseRestoreArrayRead(S,&ss));
297 : PetscCall(ISCreateStride(PETSC_COMM_SELF,n*n,0,1,&perm));
298 : PetscCall(MatDestroy(&matctx->F));
299 : PetscCall(MatDuplicate(matctx->A,MAT_COPY_VALUES,&matctx->F));
300 : PetscCall(MatLUFactor(matctx->F,perm,perm,NULL));
301 : PetscCall(ISDestroy(&perm));
302 : #else
303 26 : PetscCall(MatCopy(S,matctx->S,SAME_NONZERO_PATTERN));
304 : #endif
305 26 : matctx->lme = lme;
306 26 : matctx->n = n;
307 26 : PetscCall(VecSet(*v0,1.0));
308 26 : PetscFunctionReturn(PETSC_SUCCESS);
309 : }
310 :
311 3 : static PetscErrorCode EPSSolve_LyapII(EPS eps)
312 : {
313 3 : EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
314 3 : PetscInt i,ldds,rk,nloc,mloc,nv,idx,k;
315 3 : Vec v,w,z=eps->work[0],v0=NULL;
316 3 : Mat S,C,Ux[2],Y,Y1,R,U,W,X,Op=NULL;
317 3 : BV V;
318 3 : BVOrthogType type;
319 3 : BVOrthogRefineType refine;
320 3 : PetscScalar eigr[2],eigi[2],*array,er,ei,*uu,*s,*xx,*aa,pM[4],vec[4];
321 3 : PetscReal eta;
322 3 : EPS epsrr;
323 3 : PetscReal norm;
324 3 : EPS_LYAPII_MATSHELL *matctx;
325 :
326 3 : PetscFunctionBegin;
327 3 : PetscCall(DSGetLeadingDimension(ctx->ds,&ldds));
328 :
329 : /* Operator for the Lyapunov equation */
330 3 : PetscCall(PetscNew(&matctx));
331 3 : PetscCall(STGetOperator(eps->st,&matctx->S));
332 3 : PetscCall(MatGetLocalSize(matctx->S,&mloc,&nloc));
333 3 : PetscCall(MatCreateShell(PetscObjectComm((PetscObject)eps),mloc,nloc,PETSC_DETERMINE,PETSC_DETERMINE,matctx,&S));
334 3 : matctx->Q = eps->V;
335 3 : PetscCall(MatShellSetOperation(S,MATOP_MULT,(void(*)(void))MatMult_EPSLyapIIOperator));
336 3 : PetscCall(MatShellSetOperation(S,MATOP_DESTROY,(void(*)(void))MatDestroy_EPSLyapIIOperator));
337 3 : PetscCall(LMESetCoefficients(ctx->lme,S,NULL,NULL,NULL));
338 :
339 : /* Right-hand side */
340 3 : PetscCall(BVDuplicateResize(eps->V,ctx->rkl,&V));
341 3 : PetscCall(BVGetOrthogonalization(V,&type,&refine,&eta,NULL));
342 3 : PetscCall(BVSetOrthogonalization(V,type,refine,eta,BV_ORTHOG_BLOCK_TSQR));
343 3 : PetscCall(MatCreateDense(PetscObjectComm((PetscObject)eps),eps->nloc,PETSC_DECIDE,PETSC_DECIDE,1,NULL,&Ux[0]));
344 3 : PetscCall(MatCreateDense(PetscObjectComm((PetscObject)eps),eps->nloc,PETSC_DECIDE,PETSC_DECIDE,2,NULL,&Ux[1]));
345 3 : nv = ctx->rkl;
346 3 : PetscCall(PetscMalloc1(nv,&s));
347 :
348 : /* Initialize first column */
349 3 : PetscCall(EPSGetStartVector(eps,0,NULL));
350 3 : PetscCall(BVGetColumn(eps->V,0,&v));
351 3 : PetscCall(BVInsertVec(V,0,v));
352 3 : PetscCall(BVRestoreColumn(eps->V,0,&v));
353 3 : PetscCall(BVSetActiveColumns(eps->V,0,0)); /* no deflation at the beginning */
354 3 : PetscCall(LyapIIBuildRHS(S,1,Ux[0],V,eps->work));
355 3 : idx = 0;
356 :
357 : /* EPS for rank reduction */
358 3 : PetscCall(EPSCreate(PETSC_COMM_SELF,&epsrr));
359 3 : PetscCall(EPSSetOptionsPrefix(epsrr,((PetscObject)eps)->prefix));
360 3 : PetscCall(EPSAppendOptionsPrefix(epsrr,"eps_lyapii_"));
361 3 : PetscCall(EPSSetDimensions(epsrr,1,PETSC_DEFAULT,PETSC_DEFAULT));
362 3 : PetscCall(EPSSetTolerances(epsrr,PETSC_MACHINE_EPSILON*100,PETSC_DEFAULT));
363 :
364 29 : while (eps->reason == EPS_CONVERGED_ITERATING) {
365 26 : eps->its++;
366 :
367 : /* Matrix for placing the solution of the Lyapunov equation (an alias of V) */
368 26 : PetscCall(BVSetActiveColumns(V,0,nv));
369 26 : PetscCall(BVGetMat(V,&Y1));
370 26 : PetscCall(MatZeroEntries(Y1));
371 26 : PetscCall(MatCreateLRC(NULL,Y1,NULL,NULL,&Y));
372 26 : PetscCall(LMESetSolution(ctx->lme,Y));
373 :
374 : /* Solve the Lyapunov equation SY + YS' = -2*S*Z*S' */
375 26 : PetscCall(MatCreateLRC(NULL,Ux[idx],NULL,NULL,&C));
376 26 : PetscCall(LMESetRHS(ctx->lme,C));
377 26 : PetscCall(MatDestroy(&C));
378 26 : PetscCall(LMESolve(ctx->lme));
379 26 : PetscCall(BVRestoreMat(V,&Y1));
380 26 : PetscCall(MatDestroy(&Y));
381 :
382 : /* SVD of the solution: [Q,R]=qr(V); [U,Sigma,~]=svd(R) */
383 26 : PetscCall(DSSetDimensions(ctx->ds,nv,0,0));
384 26 : PetscCall(DSSVDSetDimensions(ctx->ds,nv));
385 26 : PetscCall(DSGetMat(ctx->ds,DS_MAT_A,&R));
386 26 : PetscCall(BVOrthogonalize(V,R));
387 26 : PetscCall(DSRestoreMat(ctx->ds,DS_MAT_A,&R));
388 26 : PetscCall(DSSetState(ctx->ds,DS_STATE_RAW));
389 26 : PetscCall(DSSolve(ctx->ds,s,NULL));
390 :
391 : /* Determine rank */
392 479 : rk = nv;
393 479 : for (i=1;i<nv;i++) if (PetscAbsScalar(s[i]/s[0])<PETSC_SQRT_MACHINE_EPSILON) {rk=i; break;}
394 26 : PetscCall(PetscInfo(eps,"The computed solution of the Lyapunov equation has rank %" PetscInt_FMT "\n",rk));
395 26 : rk = PetscMin(rk,ctx->rkc);
396 26 : PetscCall(DSGetMat(ctx->ds,DS_MAT_U,&U));
397 26 : PetscCall(BVMultInPlace(V,U,0,rk));
398 26 : PetscCall(DSRestoreMat(ctx->ds,DS_MAT_U,&U));
399 26 : PetscCall(BVSetActiveColumns(V,0,rk));
400 :
401 : /* Rank reduction */
402 26 : PetscCall(DSSetDimensions(ctx->ds,rk,0,0));
403 26 : PetscCall(DSSVDSetDimensions(ctx->ds,rk));
404 26 : PetscCall(DSGetMat(ctx->ds,DS_MAT_A,&W));
405 26 : PetscCall(BVMatProject(V,S,V,W));
406 26 : PetscCall(LyapIIBuildEigenMat(ctx->lme,W,&Op,&v0)); /* Op=A\B, A=kron(I,S)+kron(S,I), B=-2*kron(S,S) */
407 26 : PetscCall(DSRestoreMat(ctx->ds,DS_MAT_A,&W));
408 26 : PetscCall(EPSSetOperators(epsrr,Op,NULL));
409 26 : PetscCall(EPSSetInitialSpace(epsrr,1,&v0));
410 26 : PetscCall(EPSSolve(epsrr));
411 26 : PetscCall(EPSComputeVectors(epsrr));
412 : /* Copy first eigenvector, vec(A)=x */
413 26 : PetscCall(BVGetArray(epsrr->V,&xx));
414 26 : PetscCall(DSGetArray(ctx->ds,DS_MAT_A,&aa));
415 259 : for (i=0;i<rk;i++) PetscCall(PetscArraycpy(aa+i*ldds,xx+i*rk,rk));
416 26 : PetscCall(DSRestoreArray(ctx->ds,DS_MAT_A,&aa));
417 26 : PetscCall(BVRestoreArray(epsrr->V,&xx));
418 26 : PetscCall(DSSetState(ctx->ds,DS_STATE_RAW));
419 : /* Compute [U,Sigma,~] = svd(A), its rank should be 1 or 2 */
420 26 : PetscCall(DSSolve(ctx->ds,s,NULL));
421 26 : if (PetscAbsScalar(s[1]/s[0])<PETSC_SQRT_MACHINE_EPSILON) rk=1;
422 23 : else rk = 2;
423 26 : PetscCall(PetscInfo(eps,"The eigenvector has rank %" PetscInt_FMT "\n",rk));
424 26 : PetscCall(DSGetMat(ctx->ds,DS_MAT_U,&U));
425 26 : PetscCall(BVMultInPlace(V,U,0,rk));
426 26 : PetscCall(DSRestoreMat(ctx->ds,DS_MAT_U,&U));
427 :
428 : /* Save V in Ux */
429 26 : idx = (rk==2)?1:0;
430 75 : for (i=0;i<rk;i++) {
431 49 : PetscCall(BVGetColumn(V,i,&v));
432 49 : PetscCall(VecGetArray(v,&uu));
433 49 : PetscCall(MatDenseGetColumn(Ux[idx],i,&array));
434 49 : PetscCall(PetscArraycpy(array,uu,eps->nloc));
435 49 : PetscCall(MatDenseRestoreColumn(Ux[idx],&array));
436 49 : PetscCall(VecRestoreArray(v,&uu));
437 49 : PetscCall(BVRestoreColumn(V,i,&v));
438 : }
439 :
440 : /* Eigenpair approximation */
441 26 : PetscCall(BVGetColumn(V,0,&v));
442 26 : PetscCall(MatMult(S,v,z));
443 26 : PetscCall(VecDot(z,v,pM));
444 26 : PetscCall(BVRestoreColumn(V,0,&v));
445 26 : if (rk>1) {
446 23 : PetscCall(BVGetColumn(V,1,&w));
447 23 : PetscCall(VecDot(z,w,pM+1));
448 23 : PetscCall(MatMult(S,w,z));
449 23 : PetscCall(VecDot(z,w,pM+3));
450 23 : PetscCall(BVGetColumn(V,0,&v));
451 23 : PetscCall(VecDot(z,v,pM+2));
452 23 : PetscCall(BVRestoreColumn(V,0,&v));
453 23 : PetscCall(BVRestoreColumn(V,1,&w));
454 23 : PetscCall(EV2x2(pM,2,eigr,eigi,vec));
455 23 : PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,2,2,vec,&X));
456 23 : PetscCall(BVSetActiveColumns(V,0,rk));
457 23 : PetscCall(BVMultInPlace(V,X,0,rk));
458 23 : PetscCall(MatDestroy(&X));
459 : #if !defined(PETSC_USE_COMPLEX)
460 23 : norm = eigr[0]*eigr[0]+eigi[0]*eigi[0];
461 23 : er = eigr[0]/norm; ei = -eigi[0]/norm;
462 : #else
463 : er =1.0/eigr[0]; ei = 0.0;
464 : #endif
465 : } else {
466 3 : eigr[0] = pM[0]; eigi[0] = 0.0;
467 3 : er = 1.0/eigr[0]; ei = 0.0;
468 : }
469 26 : PetscCall(BVGetColumn(V,0,&v));
470 26 : if (eigi[0]!=0.0) PetscCall(BVGetColumn(V,1,&w));
471 3 : else w = NULL;
472 26 : eps->eigr[eps->nconv] = eigr[0]; eps->eigi[eps->nconv] = eigi[0];
473 26 : PetscCall(EPSComputeResidualNorm_Private(eps,PETSC_FALSE,er,ei,v,w,eps->work,&norm));
474 26 : PetscCall(BVRestoreColumn(V,0,&v));
475 26 : if (w) PetscCall(BVRestoreColumn(V,1,&w));
476 26 : PetscCall((*eps->converged)(eps,er,ei,norm,&eps->errest[eps->nconv],eps->convergedctx));
477 26 : k = 0;
478 26 : if (eps->errest[eps->nconv]<eps->tol) {
479 7 : k++;
480 7 : if (rk==2) {
481 : #if !defined (PETSC_USE_COMPLEX)
482 6 : eps->eigr[eps->nconv+k] = eigr[0]; eps->eigi[eps->nconv+k] = -eigi[0];
483 : #else
484 : eps->eigr[eps->nconv+k] = PetscConj(eps->eigr[eps->nconv]);
485 : #endif
486 6 : k++;
487 : }
488 : /* Store converged eigenpairs and vectors for deflation */
489 20 : for (i=0;i<k;i++) {
490 13 : PetscCall(BVGetColumn(V,i,&v));
491 13 : PetscCall(BVInsertVec(eps->V,eps->nconv+i,v));
492 13 : PetscCall(BVRestoreColumn(V,i,&v));
493 : }
494 7 : eps->nconv += k;
495 7 : PetscCall(BVSetActiveColumns(eps->V,eps->nconv-rk,eps->nconv));
496 7 : PetscCall(BVOrthogonalize(eps->V,NULL));
497 7 : PetscCall(DSSetDimensions(eps->ds,eps->nconv,0,0));
498 7 : PetscCall(DSGetMat(eps->ds,DS_MAT_A,&W));
499 7 : PetscCall(BVMatProject(eps->V,matctx->S,eps->V,W));
500 7 : PetscCall(DSRestoreMat(eps->ds,DS_MAT_A,&W));
501 7 : if (eps->nconv<eps->nev) {
502 4 : idx = 0;
503 4 : PetscCall(BVSetRandomColumn(V,0));
504 4 : PetscCall(BVNormColumn(V,0,NORM_2,&norm));
505 4 : PetscCall(BVScaleColumn(V,0,1.0/norm));
506 4 : PetscCall(LyapIIBuildRHS(S,1,Ux[idx],V,eps->work));
507 : }
508 : } else {
509 : /* Prepare right-hand side */
510 19 : PetscCall(LyapIIBuildRHS(S,rk,Ux[idx],NULL,eps->work));
511 : }
512 26 : PetscCall((*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx));
513 29 : PetscCall(EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->nconv+1));
514 : }
515 3 : PetscCall(STRestoreOperator(eps->st,&matctx->S));
516 3 : PetscCall(MatDestroy(&S));
517 3 : PetscCall(MatDestroy(&Ux[0]));
518 3 : PetscCall(MatDestroy(&Ux[1]));
519 3 : PetscCall(MatDestroy(&Op));
520 3 : PetscCall(VecDestroy(&v0));
521 3 : PetscCall(BVDestroy(&V));
522 3 : PetscCall(EPSDestroy(&epsrr));
523 3 : PetscCall(PetscFree(s));
524 3 : PetscFunctionReturn(PETSC_SUCCESS);
525 : }
526 :
527 3 : static PetscErrorCode EPSSetFromOptions_LyapII(EPS eps,PetscOptionItems *PetscOptionsObject)
528 : {
529 3 : EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
530 3 : PetscInt k,array[2]={PETSC_DEFAULT,PETSC_DEFAULT};
531 3 : PetscBool flg;
532 :
533 3 : PetscFunctionBegin;
534 3 : PetscOptionsHeadBegin(PetscOptionsObject,"EPS Lyapunov Inverse Iteration Options");
535 :
536 3 : k = 2;
537 3 : PetscCall(PetscOptionsIntArray("-eps_lyapii_ranks","Ranks for Lyapunov equation (one or two comma-separated integers)","EPSLyapIISetRanks",array,&k,&flg));
538 3 : if (flg) PetscCall(EPSLyapIISetRanks(eps,array[0],array[1]));
539 :
540 3 : PetscOptionsHeadEnd();
541 :
542 3 : if (!ctx->lme) PetscCall(EPSLyapIIGetLME(eps,&ctx->lme));
543 3 : PetscCall(LMESetFromOptions(ctx->lme));
544 3 : PetscFunctionReturn(PETSC_SUCCESS);
545 : }
546 :
547 1 : static PetscErrorCode EPSLyapIISetRanks_LyapII(EPS eps,PetscInt rkc,PetscInt rkl)
548 : {
549 1 : EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
550 :
551 1 : PetscFunctionBegin;
552 1 : if (rkc==PETSC_DEFAULT) rkc = 10;
553 1 : PetscCheck(rkc>1,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The compressed rank %" PetscInt_FMT " must be larger than 1",rkc);
554 1 : if (rkl==PETSC_DEFAULT) rkl = 3*rkc;
555 1 : PetscCheck(rkl>=rkc,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The Lyapunov rank %" PetscInt_FMT " cannot be smaller than the compressed rank %" PetscInt_FMT,rkl,rkc);
556 1 : if (rkc != ctx->rkc) {
557 1 : ctx->rkc = rkc;
558 1 : eps->state = EPS_STATE_INITIAL;
559 : }
560 1 : if (rkl != ctx->rkl) {
561 1 : ctx->rkl = rkl;
562 1 : eps->state = EPS_STATE_INITIAL;
563 : }
564 1 : PetscFunctionReturn(PETSC_SUCCESS);
565 : }
566 :
567 : /*@
568 : EPSLyapIISetRanks - Set the ranks used in the solution of the Lyapunov equation.
569 :
570 : Logically Collective
571 :
572 : Input Parameters:
573 : + eps - the eigenproblem solver context
574 : . rkc - the compressed rank
575 : - rkl - the Lyapunov rank
576 :
577 : Options Database Key:
578 : . -eps_lyapii_ranks <rkc,rkl> - Sets the rank parameters
579 :
580 : Notes:
581 : Lyapunov inverse iteration needs to solve a large-scale Lyapunov equation
582 : at each iteration of the eigensolver. For this, an iterative solver (LME)
583 : is used, which requires to prescribe the rank of the solution matrix X. This
584 : is the meaning of parameter rkl. Later, this matrix is compressed into
585 : another matrix of rank rkc. If not provided, rkl is a small multiple of rkc.
586 :
587 : Level: intermediate
588 :
589 : .seealso: EPSLyapIIGetRanks()
590 : @*/
591 1 : PetscErrorCode EPSLyapIISetRanks(EPS eps,PetscInt rkc,PetscInt rkl)
592 : {
593 1 : PetscFunctionBegin;
594 1 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
595 4 : PetscValidLogicalCollectiveInt(eps,rkc,2);
596 4 : PetscValidLogicalCollectiveInt(eps,rkl,3);
597 1 : PetscTryMethod(eps,"EPSLyapIISetRanks_C",(EPS,PetscInt,PetscInt),(eps,rkc,rkl));
598 1 : PetscFunctionReturn(PETSC_SUCCESS);
599 : }
600 :
601 1 : static PetscErrorCode EPSLyapIIGetRanks_LyapII(EPS eps,PetscInt *rkc,PetscInt *rkl)
602 : {
603 1 : EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
604 :
605 1 : PetscFunctionBegin;
606 1 : if (rkc) *rkc = ctx->rkc;
607 1 : if (rkl) *rkl = ctx->rkl;
608 1 : PetscFunctionReturn(PETSC_SUCCESS);
609 : }
610 :
611 : /*@
612 : EPSLyapIIGetRanks - Return the rank values used for the Lyapunov step.
613 :
614 : Not Collective
615 :
616 : Input Parameter:
617 : . eps - the eigenproblem solver context
618 :
619 : Output Parameters:
620 : + rkc - the compressed rank
621 : - rkl - the Lyapunov rank
622 :
623 : Level: intermediate
624 :
625 : .seealso: EPSLyapIISetRanks()
626 : @*/
627 1 : PetscErrorCode EPSLyapIIGetRanks(EPS eps,PetscInt *rkc,PetscInt *rkl)
628 : {
629 1 : PetscFunctionBegin;
630 1 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
631 1 : PetscUseMethod(eps,"EPSLyapIIGetRanks_C",(EPS,PetscInt*,PetscInt*),(eps,rkc,rkl));
632 1 : PetscFunctionReturn(PETSC_SUCCESS);
633 : }
634 :
635 0 : static PetscErrorCode EPSLyapIISetLME_LyapII(EPS eps,LME lme)
636 : {
637 0 : EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
638 :
639 0 : PetscFunctionBegin;
640 0 : PetscCall(PetscObjectReference((PetscObject)lme));
641 0 : PetscCall(LMEDestroy(&ctx->lme));
642 0 : ctx->lme = lme;
643 0 : eps->state = EPS_STATE_INITIAL;
644 0 : PetscFunctionReturn(PETSC_SUCCESS);
645 : }
646 :
647 : /*@
648 : EPSLyapIISetLME - Associate a linear matrix equation solver object (LME) to the
649 : eigenvalue solver.
650 :
651 : Collective
652 :
653 : Input Parameters:
654 : + eps - the eigenproblem solver context
655 : - lme - the linear matrix equation solver object
656 :
657 : Level: advanced
658 :
659 : .seealso: EPSLyapIIGetLME()
660 : @*/
661 0 : PetscErrorCode EPSLyapIISetLME(EPS eps,LME lme)
662 : {
663 0 : PetscFunctionBegin;
664 0 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
665 0 : PetscValidHeaderSpecific(lme,LME_CLASSID,2);
666 0 : PetscCheckSameComm(eps,1,lme,2);
667 0 : PetscTryMethod(eps,"EPSLyapIISetLME_C",(EPS,LME),(eps,lme));
668 0 : PetscFunctionReturn(PETSC_SUCCESS);
669 : }
670 :
671 3 : static PetscErrorCode EPSLyapIIGetLME_LyapII(EPS eps,LME *lme)
672 : {
673 3 : EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
674 :
675 3 : PetscFunctionBegin;
676 3 : if (!ctx->lme) {
677 3 : PetscCall(LMECreate(PetscObjectComm((PetscObject)eps),&ctx->lme));
678 3 : PetscCall(LMESetOptionsPrefix(ctx->lme,((PetscObject)eps)->prefix));
679 3 : PetscCall(LMEAppendOptionsPrefix(ctx->lme,"eps_lyapii_"));
680 3 : PetscCall(PetscObjectIncrementTabLevel((PetscObject)ctx->lme,(PetscObject)eps,1));
681 : }
682 3 : *lme = ctx->lme;
683 3 : PetscFunctionReturn(PETSC_SUCCESS);
684 : }
685 :
686 : /*@
687 : EPSLyapIIGetLME - Retrieve the linear matrix equation solver object (LME)
688 : associated with the eigenvalue solver.
689 :
690 : Not Collective
691 :
692 : Input Parameter:
693 : . eps - the eigenproblem solver context
694 :
695 : Output Parameter:
696 : . lme - the linear matrix equation solver object
697 :
698 : Level: advanced
699 :
700 : .seealso: EPSLyapIISetLME()
701 : @*/
702 3 : PetscErrorCode EPSLyapIIGetLME(EPS eps,LME *lme)
703 : {
704 3 : PetscFunctionBegin;
705 3 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
706 3 : PetscAssertPointer(lme,2);
707 3 : PetscUseMethod(eps,"EPSLyapIIGetLME_C",(EPS,LME*),(eps,lme));
708 3 : PetscFunctionReturn(PETSC_SUCCESS);
709 : }
710 :
711 1 : static PetscErrorCode EPSView_LyapII(EPS eps,PetscViewer viewer)
712 : {
713 1 : EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
714 1 : PetscBool isascii;
715 :
716 1 : PetscFunctionBegin;
717 1 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
718 1 : if (isascii) {
719 1 : PetscCall(PetscViewerASCIIPrintf(viewer," ranks: for Lyapunov solver=%" PetscInt_FMT ", after compression=%" PetscInt_FMT "\n",ctx->rkl,ctx->rkc));
720 1 : if (!ctx->lme) PetscCall(EPSLyapIIGetLME(eps,&ctx->lme));
721 1 : PetscCall(PetscViewerASCIIPushTab(viewer));
722 1 : PetscCall(LMEView(ctx->lme,viewer));
723 1 : PetscCall(PetscViewerASCIIPopTab(viewer));
724 : }
725 1 : PetscFunctionReturn(PETSC_SUCCESS);
726 : }
727 :
728 3 : static PetscErrorCode EPSReset_LyapII(EPS eps)
729 : {
730 3 : EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
731 :
732 3 : PetscFunctionBegin;
733 3 : if (!ctx->lme) PetscCall(LMEReset(ctx->lme));
734 3 : PetscFunctionReturn(PETSC_SUCCESS);
735 : }
736 :
737 3 : static PetscErrorCode EPSDestroy_LyapII(EPS eps)
738 : {
739 3 : EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
740 :
741 3 : PetscFunctionBegin;
742 3 : PetscCall(LMEDestroy(&ctx->lme));
743 3 : PetscCall(DSDestroy(&ctx->ds));
744 3 : PetscCall(PetscFree(eps->data));
745 3 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIISetLME_C",NULL));
746 3 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIIGetLME_C",NULL));
747 3 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIISetRanks_C",NULL));
748 3 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIIGetRanks_C",NULL));
749 3 : PetscFunctionReturn(PETSC_SUCCESS);
750 : }
751 :
752 6 : static PetscErrorCode EPSSetDefaultST_LyapII(EPS eps)
753 : {
754 6 : PetscFunctionBegin;
755 6 : if (!((PetscObject)eps->st)->type_name) PetscCall(STSetType(eps->st,STSINVERT));
756 6 : PetscFunctionReturn(PETSC_SUCCESS);
757 : }
758 :
759 3 : SLEPC_EXTERN PetscErrorCode EPSCreate_LyapII(EPS eps)
760 : {
761 3 : EPS_LYAPII *ctx;
762 :
763 3 : PetscFunctionBegin;
764 3 : PetscCall(PetscNew(&ctx));
765 3 : eps->data = (void*)ctx;
766 :
767 3 : eps->useds = PETSC_TRUE;
768 :
769 3 : eps->ops->solve = EPSSolve_LyapII;
770 3 : eps->ops->setup = EPSSetUp_LyapII;
771 3 : eps->ops->setupsort = EPSSetUpSort_Default;
772 3 : eps->ops->setfromoptions = EPSSetFromOptions_LyapII;
773 3 : eps->ops->reset = EPSReset_LyapII;
774 3 : eps->ops->destroy = EPSDestroy_LyapII;
775 3 : eps->ops->view = EPSView_LyapII;
776 3 : eps->ops->setdefaultst = EPSSetDefaultST_LyapII;
777 3 : eps->ops->backtransform = EPSBackTransform_Default;
778 3 : eps->ops->computevectors = EPSComputeVectors_Schur;
779 :
780 3 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIISetLME_C",EPSLyapIISetLME_LyapII));
781 3 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIIGetLME_C",EPSLyapIIGetLME_LyapII));
782 3 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIISetRanks_C",EPSLyapIISetRanks_LyapII));
783 3 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIIGetRanks_C",EPSLyapIIGetRanks_LyapII));
784 3 : PetscFunctionReturn(PETSC_SUCCESS);
785 : }
|