Actual source code: lyapii.c
slepc-3.13.1 2020-04-12
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: SLEPc eigensolver: "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> /*I "slepceps.h" I*/
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: {
58: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
59: PetscBool istrivial,issinv;
62: if (eps->ncv) {
63: if (eps->ncv<eps->nev+1) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev+1");
64: } else eps->ncv = eps->nev+1;
65: if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
66: if (!ctx->rkc) ctx->rkc = 10;
67: if (!ctx->rkl) ctx->rkl = 3*ctx->rkc;
68: if (!eps->max_it) eps->max_it = PetscMax(1000*eps->nev,100*eps->n);
69: if (!eps->which) eps->which=EPS_LARGEST_REAL;
70: if (eps->which!=EPS_LARGEST_REAL) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
71: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
72: if (eps->balance!=EPS_BALANCE_NONE) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Balancing not supported in this solver");
73: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
74: RGIsTrivial(eps->rg,&istrivial);
75: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
77: PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&issinv);
78: if (!issinv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Must use STSINVERT spectral transformation");
80: if (!ctx->lme) { EPSLyapIIGetLME(eps,&ctx->lme); }
81: LMESetProblemType(ctx->lme,LME_LYAPUNOV);
82: LMESetErrorIfNotConverged(ctx->lme,PETSC_TRUE);
84: if (!ctx->ds) {
85: DSCreate(PetscObjectComm((PetscObject)eps),&ctx->ds);
86: PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->ds);
87: DSSetType(ctx->ds,DSSVD);
88: }
89: DSAllocate(ctx->ds,ctx->rkl);
91: DSSetType(eps->ds,DSNHEP);
92: DSAllocate(eps->ds,eps->ncv);
94: EPSAllocateSolution(eps,0);
95: EPSSetWorkVecs(eps,3);
96: return(0);
97: }
99: static PetscErrorCode MatMult_EPSLyapIIOperator(Mat M,Vec x,Vec r)
100: {
101: PetscErrorCode ierr;
102: EPS_LYAPII_MATSHELL *matctx;
105: MatShellGetContext(M,(void**)&matctx);
106: MatMult(matctx->S,x,r);
107: BVOrthogonalizeVec(matctx->Q,r,NULL,NULL,NULL);
108: return(0);
109: }
111: static PetscErrorCode MatDestroy_EPSLyapIIOperator(Mat M)
112: {
113: PetscErrorCode ierr;
114: EPS_LYAPII_MATSHELL *matctx;
117: MatShellGetContext(M,(void**)&matctx);
118: MatDestroy(&matctx->S);
119: PetscFree(matctx);
120: return(0);
121: }
123: static PetscErrorCode MatMult_EigOperator(Mat M,Vec x,Vec y)
124: {
125: PetscErrorCode ierr;
126: EPS_EIG_MATSHELL *matctx;
127: #if !defined(PETSC_USE_COMPLEX)
128: PetscInt n;
129: PetscScalar *S,*Y,*C,zero=0.0,done=1.0,dtwo=2.0;
130: const PetscScalar *X;
131: PetscBLASInt n_;
132: #endif
135: MatShellGetContext(M,(void**)&matctx);
137: #if defined(PETSC_USE_COMPLEX)
138: MatMult(matctx->B,x,matctx->w);
139: MatSolve(matctx->F,matctx->w,y);
140: #else
141: VecGetArrayRead(x,&X);
142: VecGetArray(y,&Y);
143: MatDenseGetArray(matctx->S,&S);
145: n = matctx->n;
146: PetscCalloc1(n*n,&C);
147: PetscBLASIntCast(n,&n_);
149: /* C = 2*S*X*S.' */
150: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&dtwo,S,&n_,X,&n_,&zero,Y,&n_));
151: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&n_,&n_,&n_,&done,Y,&n_,S,&n_,&zero,C,&n_));
153: /* Solve S*Y + Y*S' = -C */
154: LMEDenseLyapunov(matctx->lme,n,S,n,C,n,Y,n);
156: PetscFree(C);
157: VecRestoreArrayRead(x,&X);
158: VecRestoreArray(y,&Y);
159: MatDenseRestoreArray(matctx->S,&S);
160: #endif
161: return(0);
162: }
164: static PetscErrorCode MatDestroy_EigOperator(Mat M)
165: {
166: PetscErrorCode ierr;
167: EPS_EIG_MATSHELL *matctx;
170: MatShellGetContext(M,(void**)&matctx);
171: #if defined(PETSC_USE_COMPLEX)
172: MatDestroy(&matctx->A);
173: MatDestroy(&matctx->B);
174: MatDestroy(&matctx->F);
175: VecDestroy(&matctx->w);
176: #endif
177: PetscFree(matctx);
178: return(0);
179: }
181: /*
182: EV2x2: solve the eigenproblem for a 2x2 matrix M
183: */
184: static PetscErrorCode EV2x2(PetscScalar *M,PetscInt ld,PetscScalar *wr,PetscScalar *wi,PetscScalar *vec)
185: {
187: PetscScalar work[10];
188: PetscBLASInt lwork=10,ld_;
189: #if !defined(PETSC_HAVE_ESSL)
190: PetscBLASInt two=2,info;
191: #if defined(PETSC_USE_COMPLEX)
192: PetscReal rwork[4];
193: #endif
194: #else
195: PetscInt i;
196: PetscBLASInt idummy,io=1;
197: PetscScalar wri[4];
198: #endif
201: PetscBLASIntCast(ld,&ld_);
202: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
203: #if !defined(PETSC_HAVE_ESSL)
204: #if !defined(PETSC_USE_COMPLEX)
205: PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","V",&two,M,&ld_,wr,wi,NULL,&ld_,vec,&ld_,work,&lwork,&info));
206: #else
207: PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","V",&two,M,&ld_,wr,NULL,&ld_,vec,&ld_,work,&lwork,rwork,&info));
208: #endif
209: SlepcCheckLapackInfo("geev",info);
210: #else /* defined(PETSC_HAVE_ESSL) */
211: PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_(&io,M,&ld_,wri,vec,ld_,&idummy,&ld_,work,&lwork));
212: #if !defined(PETSC_USE_COMPLEX)
213: for (i=0;i<2;i++) {
214: wr[i] = wri[2*i];
215: wi[i] = wri[2*i+1];
216: }
217: #else
218: for (i=0;i<2;i++) wr[i] = wri[i];
219: #endif
220: #endif
221: PetscFPTrapPop();
222: return(0);
223: }
225: /*
226: LyapIIBuildRHS: prepare the right-hand side of the Lyapunov equation SY + YS' = -2*S*Z*S'
227: in factored form:
228: if (V) U=sqrt(2)*S*V (uses 1 work vector)
229: else U=sqrt(2)*S*U (uses 2 work vectors)
230: where U,V are assumed to have rk columns.
231: */
232: static PetscErrorCode LyapIIBuildRHS(Mat S,PetscInt rk,Mat U,BV V,Vec *work)
233: {
235: PetscScalar *array,*uu;
236: PetscInt i,nloc;
237: Vec v,u=work[0];
240: MatGetLocalSize(U,&nloc,NULL);
241: for (i=0;i<rk;i++) {
242: MatDenseGetColumn(U,i,&array);
243: if (V) {
244: BVGetColumn(V,i,&v);
245: } else {
246: v = work[1];
247: VecPlaceArray(v,array);
248: }
249: MatMult(S,v,u);
250: if (V) {
251: BVRestoreColumn(V,i,&v);
252: } else {
253: VecResetArray(v);
254: }
255: VecScale(u,PetscSqrtReal(2.0));
256: VecGetArray(u,&uu);
257: PetscMemcpy(array,uu,nloc*sizeof(PetscScalar));
258: VecRestoreArray(u,&uu);
259: MatDenseRestoreColumn(U,&array);
260: }
261: return(0);
262: }
264: /*
265: LyapIIBuildEigenMat: create shell matrix Op=A\B with A = kron(I,S)+kron(S,I), B = -2*kron(S,S)
266: where S is a sequential square dense matrix of order n.
267: v0 is the initial vector, should have the form v0 = w*w' (for instance 1*1')
268: */
269: static PetscErrorCode LyapIIBuildEigenMat(LME lme,Mat S,Mat *Op,Vec *v0)
270: {
271: PetscErrorCode ierr;
272: PetscInt n,m;
273: PetscBool create=PETSC_FALSE;
274: EPS_EIG_MATSHELL *matctx;
275: #if defined(PETSC_USE_COMPLEX)
276: PetscScalar theta,*aa,*bb,*ss;
277: PetscInt i,j,f,c,off,ld;
278: IS perm;
279: #endif
282: MatGetSize(S,&n,NULL);
283: if (!*Op) create=PETSC_TRUE;
284: else {
285: MatGetSize(*Op,&m,NULL);
286: if (m!=n*n) create=PETSC_TRUE;
287: }
288: if (create) {
289: MatDestroy(Op);
290: VecDestroy(v0);
291: PetscNew(&matctx);
292: #if defined(PETSC_USE_COMPLEX)
293: MatCreateSeqDense(PETSC_COMM_SELF,n*n,n*n,NULL,&matctx->A);
294: MatCreateSeqDense(PETSC_COMM_SELF,n*n,n*n,NULL,&matctx->B);
295: MatCreateVecs(matctx->A,NULL,&matctx->w);
296: #endif
297: MatCreateShell(PETSC_COMM_SELF,n*n,n*n,PETSC_DETERMINE,PETSC_DETERMINE,matctx,Op);
298: MatShellSetOperation(*Op,MATOP_MULT,(void(*)(void))MatMult_EigOperator);
299: MatShellSetOperation(*Op,MATOP_DESTROY,(void(*)(void))MatDestroy_EigOperator);
300: MatCreateVecs(*Op,NULL,v0);
301: } else {
302: MatShellGetContext(*Op,(void**)&matctx);
303: #if defined(PETSC_USE_COMPLEX)
304: MatZeroEntries(matctx->A);
305: #endif
306: }
307: #if defined(PETSC_USE_COMPLEX)
308: MatDenseGetArray(matctx->A,&aa);
309: MatDenseGetArray(matctx->B,&bb);
310: MatDenseGetArray(S,&ss);
311: ld = n*n;
312: for (f=0;f<n;f++) {
313: off = f*n+f*n*ld;
314: for (i=0;i<n;i++) for (j=0;j<n;j++) aa[off+i+j*ld] = ss[i+j*n];
315: for (c=0;c<n;c++) {
316: off = f*n+c*n*ld;
317: theta = ss[f+c*n];
318: for (i=0;i<n;i++) aa[off+i+i*ld] += theta;
319: for (i=0;i<n;i++) for (j=0;j<n;j++) bb[off+i+j*ld] = -2*theta*ss[i+j*n];
320: }
321: }
322: MatDenseRestoreArray(matctx->A,&aa);
323: MatDenseRestoreArray(matctx->B,&bb);
324: MatDenseRestoreArray(S,&ss);
325: ISCreateStride(PETSC_COMM_SELF,n*n,0,1,&perm);
326: MatDestroy(&matctx->F);
327: MatDuplicate(matctx->A,MAT_COPY_VALUES,&matctx->F);
328: MatLUFactor(matctx->F,perm,perm,0);
329: ISDestroy(&perm);
330: #endif
331: matctx->lme = lme;
332: matctx->S = S;
333: matctx->n = n;
334: VecSet(*v0,1.0);
335: return(0);
336: }
338: PetscErrorCode EPSSolve_LyapII(EPS eps)
339: {
340: PetscErrorCode ierr;
341: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
342: PetscInt i,ldds,rk,nloc,mloc,nv,idx,k;
343: Vec v,w,z=eps->work[0],v0=NULL;
344: Mat S,C,Ux[2],Y,Y1,R,U,W,X,Op=NULL;
345: BV V;
346: BVOrthogType type;
347: BVOrthogRefineType refine;
348: PetscScalar eigr[2],eigi[2],*array,er,ei,*uu,*s,*xx,*aa,pM[4],vec[4];
349: PetscReal eta;
350: EPS epsrr;
351: PetscReal norm;
352: EPS_LYAPII_MATSHELL *matctx;
355: DSGetLeadingDimension(ctx->ds,&ldds);
357: /* Operator for the Lyapunov equation */
358: PetscNew(&matctx);
359: STGetOperator(eps->st,&matctx->S);
360: MatGetLocalSize(matctx->S,&mloc,&nloc);
361: MatCreateShell(PetscObjectComm((PetscObject)eps),mloc,nloc,PETSC_DETERMINE,PETSC_DETERMINE,matctx,&S);
362: matctx->Q = eps->V;
363: MatShellSetOperation(S,MATOP_MULT,(void(*)(void))MatMult_EPSLyapIIOperator);
364: MatShellSetOperation(S,MATOP_DESTROY,(void(*)(void))MatDestroy_EPSLyapIIOperator);
365: LMESetCoefficients(ctx->lme,S,NULL,NULL,NULL);
367: /* Right-hand side */
368: BVDuplicateResize(eps->V,ctx->rkl,&V);
369: BVGetOrthogonalization(V,&type,&refine,&eta,NULL);
370: BVSetOrthogonalization(V,type,refine,eta,BV_ORTHOG_BLOCK_TSQR);
371: MatCreateDense(PetscObjectComm((PetscObject)eps),eps->nloc,PETSC_DECIDE,PETSC_DECIDE,1,NULL,&Ux[0]);
372: MatCreateDense(PetscObjectComm((PetscObject)eps),eps->nloc,PETSC_DECIDE,PETSC_DECIDE,2,NULL,&Ux[1]);
373: nv = ctx->rkl;
374: PetscMalloc1(nv,&s);
376: /* Initialize first column */
377: EPSGetStartVector(eps,0,NULL);
378: BVGetColumn(eps->V,0,&v);
379: BVInsertVec(V,0,v);
380: BVRestoreColumn(eps->V,0,&v);
381: BVSetActiveColumns(eps->V,0,0); /* no deflation at the beginning */
382: LyapIIBuildRHS(S,1,Ux[0],V,eps->work);
383: idx = 0;
385: /* EPS for rank reduction */
386: EPSCreate(PETSC_COMM_SELF,&epsrr);
387: EPSSetOptionsPrefix(epsrr,((PetscObject)eps)->prefix);
388: EPSAppendOptionsPrefix(epsrr,"eps_lyapii_");
389: EPSSetDimensions(epsrr,1,PETSC_DEFAULT,PETSC_DEFAULT);
390: EPSSetTolerances(epsrr,PETSC_MACHINE_EPSILON*100,PETSC_DEFAULT);
392: while (eps->reason == EPS_CONVERGED_ITERATING) {
393: eps->its++;
395: /* Matrix for placing the solution of the Lyapunov equation (an alias of V) */
396: BVSetActiveColumns(V,0,nv);
397: BVGetMat(V,&Y1);
398: MatZeroEntries(Y1);
399: MatCreateLRC(NULL,Y1,NULL,NULL,&Y);
400: LMESetSolution(ctx->lme,Y);
402: /* Solve the Lyapunov equation SY + YS' = -2*S*Z*S' */
403: MatCreateLRC(NULL,Ux[idx],NULL,NULL,&C);
404: LMESetRHS(ctx->lme,C);
405: MatDestroy(&C);
406: LMESolve(ctx->lme);
407: BVRestoreMat(V,&Y1);
408: MatDestroy(&Y);
410: /* SVD of the solution: [Q,R]=qr(V); [U,Sigma,~]=svd(R) */
411: DSSetDimensions(ctx->ds,nv,nv,0,0);
412: DSGetMat(ctx->ds,DS_MAT_A,&R);
413: BVOrthogonalize(V,R);
414: DSRestoreMat(ctx->ds,DS_MAT_A,&R);
415: DSSetState(ctx->ds,DS_STATE_RAW);
416: DSSolve(ctx->ds,s,NULL);
418: /* Determine rank */
419: rk = nv;
420: for (i=1;i<nv;i++) if (PetscAbsScalar(s[i]/s[0])<PETSC_SQRT_MACHINE_EPSILON) {rk=i; break;}
421: PetscInfo1(eps,"The computed solution of the Lyapunov equation has rank %D\n",rk);
422: rk = PetscMin(rk,ctx->rkc);
423: DSGetMat(ctx->ds,DS_MAT_U,&U);
424: BVMultInPlace(V,U,0,rk);
425: BVSetActiveColumns(V,0,rk);
426: MatDestroy(&U);
428: /* Rank reduction */
429: DSSetDimensions(ctx->ds,rk,rk,0,0);
430: DSGetMat(ctx->ds,DS_MAT_A,&W);
431: BVMatProject(V,S,V,W);
432: LyapIIBuildEigenMat(ctx->lme,W,&Op,&v0); /* Op=A\B, A=kron(I,S)+kron(S,I), B=-2*kron(S,S) */
433: EPSSetOperators(epsrr,Op,NULL);
434: EPSSetInitialSpace(epsrr,1,&v0);
435: EPSSolve(epsrr);
436: MatDestroy(&W);
437: EPSComputeVectors(epsrr);
438: /* Copy first eigenvector, vec(A)=x */
439: BVGetArray(epsrr->V,&xx);
440: DSGetArray(ctx->ds,DS_MAT_A,&aa);
441: for (i=0;i<rk;i++) {
442: PetscMemcpy(aa+i*ldds,xx+i*rk,rk*sizeof(PetscScalar));
443: }
444: DSRestoreArray(ctx->ds,DS_MAT_A,&aa);
445: BVRestoreArray(epsrr->V,&xx);
446: DSSetState(ctx->ds,DS_STATE_RAW);
447: /* Compute [U,Sigma,~] = svd(A), its rank should be 1 or 2 */
448: DSSolve(ctx->ds,s,NULL);
449: if (PetscAbsScalar(s[1]/s[0])<PETSC_SQRT_MACHINE_EPSILON) rk=1;
450: else rk = 2;
451: PetscInfo1(eps,"The eigenvector has rank %D\n",rk);
452: DSGetMat(ctx->ds,DS_MAT_U,&U);
453: BVMultInPlace(V,U,0,rk);
454: MatDestroy(&U);
456: /* Save V in Ux */
457: idx = (rk==2)?1:0;
458: for (i=0;i<rk;i++) {
459: BVGetColumn(V,i,&v);
460: VecGetArray(v,&uu);
461: MatDenseGetColumn(Ux[idx],i,&array);
462: PetscMemcpy(array,uu,eps->nloc*sizeof(PetscScalar));
463: MatDenseRestoreColumn(Ux[idx],&array);
464: VecRestoreArray(v,&uu);
465: BVRestoreColumn(V,i,&v);
466: }
468: /* Eigenpair approximation */
469: BVGetColumn(V,0,&v);
470: MatMult(S,v,z);
471: VecDot(z,v,pM);
472: BVRestoreColumn(V,0,&v);
473: if (rk>1) {
474: BVGetColumn(V,1,&w);
475: VecDot(z,w,pM+1);
476: MatMult(S,w,z);
477: VecDot(z,w,pM+3);
478: BVGetColumn(V,0,&v);
479: VecDot(z,v,pM+2);
480: BVRestoreColumn(V,0,&v);
481: BVRestoreColumn(V,1,&w);
482: EV2x2(pM,2,eigr,eigi,vec);
483: MatCreateSeqDense(PETSC_COMM_SELF,2,2,vec,&X);
484: BVSetActiveColumns(V,0,rk);
485: BVMultInPlace(V,X,0,rk);
486: MatDestroy(&X);
487: #if !defined(PETSC_USE_COMPLEX)
488: norm = eigr[0]*eigr[0]+eigi[0]*eigi[0];
489: er = eigr[0]/norm; ei = -eigi[0]/norm;
490: #else
491: er =1.0/eigr[0]; ei = 0.0;
492: #endif
493: } else {
494: eigr[0] = pM[0]; eigi[0] = 0.0;
495: er = 1.0/eigr[0]; ei = 0.0;
496: }
497: BVGetColumn(V,0,&v);
498: if (eigi[0]!=0.0) {
499: BVGetColumn(V,1,&w);
500: } else w = NULL;
501: eps->eigr[eps->nconv] = eigr[0]; eps->eigi[eps->nconv] = eigi[0];
502: EPSComputeResidualNorm_Private(eps,PETSC_FALSE,er,ei,v,w,eps->work,&norm);
503: BVRestoreColumn(V,0,&v);
504: if (w) {
505: BVRestoreColumn(V,1,&w);
506: }
507: (*eps->converged)(eps,er,ei,norm,&eps->errest[eps->nconv],eps->convergedctx);
508: k = 0;
509: if (eps->errest[eps->nconv]<eps->tol) {
510: k++;
511: if (rk==2) {
512: #if !defined (PETSC_USE_COMPLEX)
513: eps->eigr[eps->nconv+k] = eigr[0]; eps->eigi[eps->nconv+k] = -eigi[0];
514: #else
515: eps->eigr[eps->nconv+k] = PetscConj(eps->eigr[eps->nconv]);
516: #endif
517: k++;
518: }
519: /* Store converged eigenpairs and vectors for deflation */
520: for (i=0;i<k;i++) {
521: BVGetColumn(V,i,&v);
522: BVInsertVec(eps->V,eps->nconv+i,v);
523: BVRestoreColumn(V,i,&v);
524: }
525: eps->nconv += k;
526: BVSetActiveColumns(eps->V,eps->nconv-rk,eps->nconv);
527: BVOrthogonalize(eps->V,NULL);
528: DSSetDimensions(eps->ds,eps->nconv,0,0,0);
529: DSGetMat(eps->ds,DS_MAT_A,&W);
530: BVMatProject(eps->V,matctx->S,eps->V,W);
531: DSRestoreMat(eps->ds,DS_MAT_A,&W);
532: if (eps->nconv<eps->nev) {
533: idx = 0;
534: BVSetRandomColumn(V,0);
535: BVNormColumn(V,0,NORM_2,&norm);
536: BVScaleColumn(V,0,1.0/norm);
537: LyapIIBuildRHS(S,1,Ux[idx],V,eps->work);
538: }
539: } else {
540: /* Prepare right-hand side */
541: LyapIIBuildRHS(S,rk,Ux[idx],NULL,eps->work);
542: }
543: (*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx);
544: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->nconv+1);
545: }
546: STRestoreOperator(eps->st,&matctx->S);
547: MatDestroy(&S);
548: MatDestroy(&Ux[0]);
549: MatDestroy(&Ux[1]);
550: MatDestroy(&Op);
551: VecDestroy(&v0);
552: BVDestroy(&V);
553: EPSDestroy(&epsrr);
554: PetscFree(s);
555: return(0);
556: }
558: PetscErrorCode EPSSetFromOptions_LyapII(PetscOptionItems *PetscOptionsObject,EPS eps)
559: {
561: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
562: PetscInt k,array[2]={PETSC_DEFAULT,PETSC_DEFAULT};
563: PetscBool flg;
566: PetscOptionsHead(PetscOptionsObject,"EPS Lyapunov Inverse Iteration Options");
568: k = 2;
569: PetscOptionsIntArray("-eps_lyapii_ranks","Ranks for Lyapunov equation (one or two comma-separated integers)","EPSLyapIISetRanks",array,&k,&flg);
570: if (flg) {
571: EPSLyapIISetRanks(eps,array[0],array[1]);
572: }
574: PetscOptionsTail();
576: if (!ctx->lme) { EPSLyapIIGetLME(eps,&ctx->lme); }
577: LMESetFromOptions(ctx->lme);
578: return(0);
579: }
581: static PetscErrorCode EPSLyapIISetRanks_LyapII(EPS eps,PetscInt rkc,PetscInt rkl)
582: {
583: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
586: if (rkc==PETSC_DEFAULT) rkc = 10;
587: if (rkc<2) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The compressed rank %D must be larger than 1",rkc);
588: if (rkl==PETSC_DEFAULT) rkl = 3*rkc;
589: if (rkl<rkc) SETERRQ2(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The Lyapunov rank %D cannot be smaller than the compressed rank %D",rkl,rkc);
590: if (rkc != ctx->rkc) {
591: ctx->rkc = rkc;
592: eps->state = EPS_STATE_INITIAL;
593: }
594: if (rkl != ctx->rkl) {
595: ctx->rkl = rkl;
596: eps->state = EPS_STATE_INITIAL;
597: }
598: return(0);
599: }
601: /*@
602: EPSLyapIISetRanks - Set the ranks used in the solution of the Lyapunov equation.
604: Collective on EPS
606: Input Parameters:
607: + eps - the eigenproblem solver context
608: . rkc - the compressed rank
609: - rkl - the Lyapunov rank
611: Options Database Key:
612: . -eps_lyapii_ranks <rkc,rkl> - Sets the rank parameters
614: Notes:
615: Lyapunov inverse iteration needs to solve a large-scale Lyapunov equation
616: at each iteration of the eigensolver. For this, an iterative solver (LME)
617: is used, which requires to prescribe the rank of the solution matrix X. This
618: is the meaning of parameter rkl. Later, this matrix is compressed into
619: another matrix of rank rkc. If not provided, rkl is a small multiple of rkc.
621: Level: intermediate
623: .seealso: EPSLyapIIGetRanks()
624: @*/
625: PetscErrorCode EPSLyapIISetRanks(EPS eps,PetscInt rkc,PetscInt rkl)
626: {
633: PetscTryMethod(eps,"EPSLyapIISetRanks_C",(EPS,PetscInt,PetscInt),(eps,rkc,rkl));
634: return(0);
635: }
637: static PetscErrorCode EPSLyapIIGetRanks_LyapII(EPS eps,PetscInt *rkc,PetscInt *rkl)
638: {
639: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
642: if (rkc) *rkc = ctx->rkc;
643: if (rkl) *rkl = ctx->rkl;
644: return(0);
645: }
647: /*@
648: EPSLyapIIGetRanks - Return the rank values used for the Lyapunov step.
650: Not Collective
652: Input Parameter:
653: . eps - the eigenproblem solver context
655: Output Parameters:
656: + rkc - the compressed rank
657: - rkl - the Lyapunov rank
659: Level: intermediate
661: .seealso: EPSLyapIISetRanks()
662: @*/
663: PetscErrorCode EPSLyapIIGetRanks(EPS eps,PetscInt *rkc,PetscInt *rkl)
664: {
669: PetscUseMethod(eps,"EPSLyapIIGetRanks_C",(EPS,PetscInt*,PetscInt*),(eps,rkc,rkl));
670: return(0);
671: }
673: static PetscErrorCode EPSLyapIISetLME_LyapII(EPS eps,LME lme)
674: {
676: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
679: PetscObjectReference((PetscObject)lme);
680: LMEDestroy(&ctx->lme);
681: ctx->lme = lme;
682: PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->lme);
683: eps->state = EPS_STATE_INITIAL;
684: return(0);
685: }
687: /*@
688: EPSLyapIISetLME - Associate a linear matrix equation solver object (LME) to the
689: eigenvalue solver.
691: Collective on EPS
693: Input Parameters:
694: + eps - the eigenproblem solver context
695: - lme - the linear matrix equation solver object
697: Level: advanced
699: .seealso: EPSLyapIIGetLME()
700: @*/
701: PetscErrorCode EPSLyapIISetLME(EPS eps,LME lme)
702: {
709: PetscTryMethod(eps,"EPSLyapIISetLME_C",(EPS,LME),(eps,lme));
710: return(0);
711: }
713: static PetscErrorCode EPSLyapIIGetLME_LyapII(EPS eps,LME *lme)
714: {
716: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
719: if (!ctx->lme) {
720: LMECreate(PetscObjectComm((PetscObject)eps),&ctx->lme);
721: LMESetOptionsPrefix(ctx->lme,((PetscObject)eps)->prefix);
722: LMEAppendOptionsPrefix(ctx->lme,"eps_lyapii_");
723: PetscObjectIncrementTabLevel((PetscObject)ctx->lme,(PetscObject)eps,1);
724: PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->lme);
725: }
726: *lme = ctx->lme;
727: return(0);
728: }
730: /*@
731: EPSLyapIIGetLME - Retrieve the linear matrix equation solver object (LME)
732: associated with the eigenvalue solver.
734: Not Collective
736: Input Parameter:
737: . eps - the eigenproblem solver context
739: Output Parameter:
740: . lme - the linear matrix equation solver object
742: Level: advanced
744: .seealso: EPSLyapIISetLME()
745: @*/
746: PetscErrorCode EPSLyapIIGetLME(EPS eps,LME *lme)
747: {
753: PetscUseMethod(eps,"EPSLyapIIGetLME_C",(EPS,LME*),(eps,lme));
754: return(0);
755: }
757: PetscErrorCode EPSView_LyapII(EPS eps,PetscViewer viewer)
758: {
760: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
761: PetscBool isascii;
764: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
765: if (isascii) {
766: PetscViewerASCIIPrintf(viewer," ranks: for Lyapunov solver=%D, after compression=%D\n",ctx->rkl,ctx->rkc);
767: if (!ctx->lme) { EPSLyapIIGetLME(eps,&ctx->lme); }
768: PetscViewerASCIIPushTab(viewer);
769: LMEView(ctx->lme,viewer);
770: PetscViewerASCIIPopTab(viewer);
771: }
772: return(0);
773: }
775: PetscErrorCode EPSReset_LyapII(EPS eps)
776: {
778: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
781: if (!ctx->lme) { LMEReset(ctx->lme); }
782: return(0);
783: }
785: PetscErrorCode EPSDestroy_LyapII(EPS eps)
786: {
788: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
791: LMEDestroy(&ctx->lme);
792: DSDestroy(&ctx->ds);
793: PetscFree(eps->data);
794: PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIISetLME_C",NULL);
795: PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIIGetLME_C",NULL);
796: PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIISetRanks_C",NULL);
797: PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIIGetRanks_C",NULL);
798: return(0);
799: }
801: PetscErrorCode EPSSetDefaultST_LyapII(EPS eps)
802: {
806: if (!((PetscObject)eps->st)->type_name) {
807: STSetType(eps->st,STSINVERT);
808: }
809: return(0);
810: }
812: PETSC_EXTERN PetscErrorCode EPSCreate_LyapII(EPS eps)
813: {
814: EPS_LYAPII *ctx;
818: PetscNewLog(eps,&ctx);
819: eps->data = (void*)ctx;
821: eps->useds = PETSC_TRUE;
823: eps->ops->solve = EPSSolve_LyapII;
824: eps->ops->setup = EPSSetUp_LyapII;
825: eps->ops->setfromoptions = EPSSetFromOptions_LyapII;
826: eps->ops->reset = EPSReset_LyapII;
827: eps->ops->destroy = EPSDestroy_LyapII;
828: eps->ops->view = EPSView_LyapII;
829: eps->ops->setdefaultst = EPSSetDefaultST_LyapII;
830: eps->ops->backtransform = EPSBackTransform_Default;
831: eps->ops->computevectors = EPSComputeVectors_Schur;
833: PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIISetLME_C",EPSLyapIISetLME_LyapII);
834: PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIIGetLME_C",EPSLyapIIGetLME_LyapII);
835: PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIISetRanks_C",EPSLyapIISetRanks_LyapII);
836: PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIIGetRanks_C",EPSLyapIIGetRanks_LyapII);
837: return(0);
838: }