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