Actual source code: lyapii.c
slepc-3.13.0 2020-03-31
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;
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: EPSAllocateSolution(eps,0);
92: EPSSetWorkVecs(eps,3);
93: return(0);
94: }
96: static PetscErrorCode MatMult_EPSLyapIIOperator(Mat M,Vec x,Vec r)
97: {
98: PetscErrorCode ierr;
99: EPS_LYAPII_MATSHELL *matctx;
102: MatShellGetContext(M,(void**)&matctx);
103: MatMult(matctx->S,x,r);
104: BVOrthogonalizeVec(matctx->Q,r,NULL,NULL,NULL);
105: return(0);
106: }
108: static PetscErrorCode MatDestroy_EPSLyapIIOperator(Mat M)
109: {
110: PetscErrorCode ierr;
111: EPS_LYAPII_MATSHELL *matctx;
114: MatShellGetContext(M,(void**)&matctx);
115: MatDestroy(&matctx->S);
116: BVDestroy(&matctx->Q);
117: PetscFree(matctx);
118: return(0);
119: }
121: static PetscErrorCode MatCreateVecs_EPSLyapIIOperator(Mat M,Vec *right,Vec *left)
122: {
123: PetscErrorCode ierr;
124: EPS_LYAPII_MATSHELL *matctx;
127: MatShellGetContext(M,(void**)&matctx);
128: MatCreateVecs(matctx->S,right,left);
129: return(0);
130: }
132: static PetscErrorCode MatMult_EigOperator(Mat M,Vec x,Vec y)
133: {
134: PetscErrorCode ierr;
135: EPS_EIG_MATSHELL *matctx;
136: PetscInt n;
137: #if defined(PETSC_USE_COMPLEX)
138: Mat F;
139: IS perm;
140: #else
141: PetscScalar *S,*Y,*C,zero=0.0,done=1.0,dtwo=2.0;
142: const PetscScalar *X;
143: PetscBLASInt n_;
144: #endif
147: MatShellGetContext(M,(void**)&matctx);
149: #if defined(PETSC_USE_COMPLEX)
150: MatMult(matctx->B,x,matctx->w);
151: MatGetSize(matctx->A,&n,NULL);
152: ISCreateStride(PETSC_COMM_SELF,n,0,1,&perm);
153: MatDuplicate(matctx->A,MAT_COPY_VALUES,&F);
154: MatLUFactor(F,perm,perm,0);
155: MatSolve(F,matctx->w,y);
156: MatDestroy(&F);
157: ISDestroy(&perm);
158: #else
159: VecGetArrayRead(x,&X);
160: VecGetArray(y,&Y);
161: MatDenseGetArray(matctx->S,&S);
163: n = matctx->n;
164: PetscCalloc1(n*n,&C);
165: PetscBLASIntCast(n,&n_);
167: /* C = 2*S*X*S.' */
168: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&dtwo,S,&n_,X,&n_,&zero,Y,&n_));
169: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&n_,&n_,&n_,&done,Y,&n_,S,&n_,&zero,C,&n_));
171: /* Solve S*Y + Y*S' = -C */
172: LMEDenseLyapunov(matctx->lme,n,S,n,C,n,Y,n);
174: PetscFree(C);
175: VecRestoreArrayRead(x,&X);
176: VecRestoreArray(y,&Y);
177: MatDenseRestoreArray(matctx->S,&S);
178: #endif
179: return(0);
180: }
182: static PetscErrorCode MatDestroy_EigOperator(Mat M)
183: {
184: PetscErrorCode ierr;
185: EPS_EIG_MATSHELL *matctx;
188: MatShellGetContext(M,(void**)&matctx);
189: #if defined(PETSC_USE_COMPLEX)
190: MatDestroy(&matctx->A);
191: MatDestroy(&matctx->B);
192: VecDestroy(&matctx->w);
193: #endif
194: PetscFree(matctx);
195: return(0);
196: }
198: /*
199: EV2x2: solve the eigenproblem for a 2x2 matrix M
200: */
201: static PetscErrorCode EV2x2(PetscScalar *M,PetscInt ld,PetscScalar *wr,PetscScalar *wi,PetscScalar *vec)
202: {
204: PetscScalar work[10];
205: PetscBLASInt lwork=10,ld_;
206: #if !defined(PETSC_HAVE_ESSL)
207: PetscBLASInt two=2,info;
208: #if defined(PETSC_USE_COMPLEX)
209: PetscReal rwork[4];
210: #endif
211: #else
212: PetscInt i;
213: PetscBLASInt idummy,io=1;
214: PetscScalar wri[4];
215: #endif
218: PetscBLASIntCast(ld,&ld_);
219: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
220: #if !defined(PETSC_HAVE_ESSL)
221: #if !defined(PETSC_USE_COMPLEX)
222: PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","V",&two,M,&ld_,wr,wi,NULL,&ld_,vec,&ld_,work,&lwork,&info));
223: #else
224: PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","V",&two,M,&ld_,wr,NULL,&ld_,vec,&ld_,work,&lwork,rwork,&info));
225: #endif
226: SlepcCheckLapackInfo("geev",info);
227: #else /* defined(PETSC_HAVE_ESSL) */
228: PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_(&io,M,&ld_,wri,vec,ld_,&idummy,&ld_,work,&lwork));
229: #if !defined(PETSC_USE_COMPLEX)
230: for (i=0;i<2;i++) {
231: wr[i] = wri[2*i];
232: wi[i] = wri[2*i+1];
233: }
234: #else
235: for (i=0;i<2;i++) wr[i] = wri[i];
236: #endif
237: #endif
238: PetscFPTrapPop();
239: return(0);
240: }
242: /*
243: LyapIIBuildRHS: prepare the right-hand side of the Lyapunov equation SY + YS' = -2*S*Z*S'
244: in factored form:
245: if (V) U=sqrt(2)*S*V (uses 1 work vector)
246: else U=sqrt(2)*S*U (uses 2 work vectors)
247: where U,V are assumed to have rk columns.
248: */
249: static PetscErrorCode LyapIIBuildRHS(Mat S,PetscInt rk,Mat U,BV V,Vec *work)
250: {
252: PetscScalar *array,*uu;
253: PetscInt i,nloc;
254: Vec v,u=work[0];
257: MatGetLocalSize(U,&nloc,NULL);
258: for (i=0;i<rk;i++) {
259: MatDenseGetColumn(U,i,&array);
260: if (V) {
261: BVGetColumn(V,i,&v);
262: } else {
263: v = work[1];
264: VecPlaceArray(v,array);
265: }
266: MatMult(S,v,u);
267: if (V) {
268: BVRestoreColumn(V,i,&v);
269: } else {
270: VecResetArray(v);
271: }
272: VecScale(u,PetscSqrtReal(2.0));
273: VecGetArray(u,&uu);
274: PetscMemcpy(array,uu,nloc*sizeof(PetscScalar));
275: VecRestoreArray(u,&uu);
276: MatDenseRestoreColumn(U,&array);
277: }
278: return(0);
279: }
281: /*
282: LyapIIBuildEigenMat: create shell matrix Op=A\B with A = kron(I,S)+kron(S,I), B = -2*kron(S,S)
283: where S is a sequential square dense matrix of order n.
284: v0 is the initial vector, should have the form v0 = w*w' (for instance 1*1')
285: */
286: static PetscErrorCode LyapIIBuildEigenMat(LME lme,Mat S,Mat *Op,Vec *v0)
287: {
288: PetscErrorCode ierr;
289: PetscInt n,m;
290: PetscBool create=PETSC_FALSE;
291: EPS_EIG_MATSHELL *matctx;
292: #if defined(PETSC_USE_COMPLEX)
293: PetscScalar theta,*aa,*bb,*ss;
294: PetscInt i,j,f,c,off,ld;
295: #endif
298: MatGetSize(S,&n,NULL);
299: if (!*Op) create=PETSC_TRUE;
300: else {
301: MatGetSize(*Op,&m,NULL);
302: if (m!=n*n) create=PETSC_TRUE;
303: }
304: if (create) {
305: MatDestroy(Op);
306: VecDestroy(v0);
307: PetscNew(&matctx);
308: #if defined(PETSC_USE_COMPLEX)
309: MatCreateSeqDense(PETSC_COMM_SELF,n*n,n*n,NULL,&matctx->A);
310: MatCreateSeqDense(PETSC_COMM_SELF,n*n,n*n,NULL,&matctx->B);
311: MatCreateVecs(matctx->A,NULL,&matctx->w);
312: #endif
313: MatCreateShell(PETSC_COMM_SELF,n*n,n*n,PETSC_DETERMINE,PETSC_DETERMINE,matctx,Op);
314: MatShellSetOperation(*Op,MATOP_MULT,(void(*)(void))MatMult_EigOperator);
315: MatShellSetOperation(*Op,MATOP_DESTROY,(void(*)(void))MatDestroy_EigOperator);
316: MatCreateVecs(*Op,NULL,v0);
317: } else {
318: MatShellGetContext(*Op,(void**)&matctx);
319: #if defined(PETSC_USE_COMPLEX)
320: MatZeroEntries(matctx->A);
321: #endif
322: }
323: #if defined(PETSC_USE_COMPLEX)
324: MatDenseGetArray(matctx->A,&aa);
325: MatDenseGetArray(matctx->B,&bb);
326: MatDenseGetArray(S,&ss);
327: ld = n*n;
328: for (f=0;f<n;f++) {
329: off = f*n+f*n*ld;
330: for (i=0;i<n;i++) for (j=0;j<n;j++) aa[off+i+j*ld] = ss[i+j*n];
331: for (c=0;c<n;c++) {
332: off = f*n+c*n*ld;
333: theta = ss[f+c*n];
334: for (i=0;i<n;i++) aa[off+i+i*ld] += theta;
335: for (i=0;i<n;i++) for (j=0;j<n;j++) bb[off+i+j*ld] = -2*theta*ss[i+j*n];
336: }
337: }
338: MatDenseRestoreArray(matctx->A,&aa);
339: MatDenseRestoreArray(matctx->B,&bb);
340: MatDenseRestoreArray(S,&ss);
341: #endif
342: matctx->lme = lme;
343: matctx->S = S;
344: matctx->n = n;
345: VecSet(*v0,1.0);
346: return(0);
347: }
349: PetscErrorCode EPSSolve_LyapII(EPS eps)
350: {
351: PetscErrorCode ierr;
352: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
353: PetscInt i,ldds,rk,nloc,mloc,nv,idx,k;
354: Vec v,w,z=eps->work[0],v0=NULL;
355: Mat S,C,Ux[2],Y,Y1,R,U,W,X,Op=NULL;
356: BV V;
357: BVOrthogType type;
358: BVOrthogRefineType refine;
359: PetscScalar *eigr,*eigi,*array,er,ei,*uu,*pV,*s,*xx,*aa,pM[4],vec[4];
360: PetscReal eta;
361: EPS epsrr;
362: PetscReal norm;
363: EPS_LYAPII_MATSHELL *matctx;
366: DSGetLeadingDimension(ctx->ds,&ldds);
368: /* Operator for the Lyapunov equation */
369: PetscNew(&matctx);
370: STGetOperator(eps->st,&matctx->S);
371: MatGetLocalSize(matctx->S,&mloc,&nloc);
372: MatCreateShell(PetscObjectComm((PetscObject)eps),mloc,nloc,PETSC_DETERMINE,PETSC_DETERMINE,matctx,&S);
373: BVDuplicateResize(eps->V,eps->nev+1,&matctx->Q);
374: MatShellSetOperation(S,MATOP_MULT,(void(*)(void))MatMult_EPSLyapIIOperator);
375: MatShellSetOperation(S,MATOP_DESTROY,(void(*)(void))MatDestroy_EPSLyapIIOperator);
376: MatShellSetOperation(S,MATOP_CREATE_VECS,(void(*)(void))MatCreateVecs_EPSLyapIIOperator);
377: LMESetCoefficients(ctx->lme,S,NULL,NULL,NULL);
379: /* Right-hand side */
380: BVDuplicateResize(eps->V,ctx->rkl,&V);
381: BVGetOrthogonalization(V,&type,&refine,&eta,NULL);
382: BVSetOrthogonalization(V,type,refine,eta,BV_ORTHOG_BLOCK_TSQR);
383: MatCreateDense(PetscObjectComm((PetscObject)eps),eps->nloc,PETSC_DECIDE,PETSC_DECIDE,1,NULL,&Ux[0]);
384: MatCreateDense(PetscObjectComm((PetscObject)eps),eps->nloc,PETSC_DECIDE,PETSC_DECIDE,2,NULL,&Ux[1]);
385: nv = ctx->rkl;
386: PetscMalloc3(nv,&s,nv*nv,&eigr,nv*nv,&eigi);
388: /* Initialize first column */
389: EPSGetStartVector(eps,0,NULL);
390: BVGetColumn(eps->V,0,&v);
391: BVInsertVec(V,0,v);
392: BVRestoreColumn(eps->V,0,&v);
393: LyapIIBuildRHS(S,1,Ux[0],V,eps->work);
394: idx = 0;
396: /* EPS for rank reduction */
397: EPSCreate(PETSC_COMM_SELF,&epsrr);
398: EPSSetOptionsPrefix(epsrr,((PetscObject)eps)->prefix);
399: EPSAppendOptionsPrefix(epsrr,"eps_lyapii_");
400: EPSSetDimensions(epsrr,1,PETSC_DEFAULT,PETSC_DEFAULT);
401: EPSSetTolerances(epsrr,PETSC_MACHINE_EPSILON*100,PETSC_DEFAULT);
403: while (eps->reason == EPS_CONVERGED_ITERATING) {
404: eps->its++;
406: /* Matrix for placing the solution of the Lyapunov equation (an alias of V) */
407: BVGetArray(V,&pV);
408: PetscMemzero(pV,nv*eps->nloc*sizeof(PetscScalar));
409: MatCreateDense(PetscObjectComm((PetscObject)eps),eps->nloc,PETSC_DECIDE,PETSC_DECIDE,nv,pV,&Y1);
410: MatCreateLRC(NULL,Y1,NULL,NULL,&Y);
411: MatDestroy(&Y1);
412: LMESetSolution(ctx->lme,Y);
414: /* Solve the Lyapunov equation SY + YS' = -2*S*Z*S' */
415: MatCreateLRC(NULL,Ux[idx],NULL,NULL,&C);
416: LMESetRHS(ctx->lme,C);
417: MatDestroy(&C);
418: LMESolve(ctx->lme);
419: BVRestoreArray(V,&pV);
420: MatDestroy(&Y);
422: /* SVD of the solution: [Q,R]=qr(V); [U,Sigma,~]=svd(R) */
423: DSSetDimensions(ctx->ds,nv,nv,0,0);
424: DSGetMat(ctx->ds,DS_MAT_A,&R);
425: BVSetActiveColumns(V,0,nv);
426: BVOrthogonalize(V,R);
427: DSRestoreMat(ctx->ds,DS_MAT_A,&R);
428: DSSetState(ctx->ds,DS_STATE_RAW);
429: DSSolve(ctx->ds,s,NULL);
431: /* Determine rank */
432: rk = nv;
433: for (i=1;i<nv;i++) if (PetscAbsScalar(s[i]/s[0])<PETSC_SQRT_MACHINE_EPSILON) {rk=i; break;}
434: PetscInfo1(eps,"The computed solution of the Lyapunov equation has rank %D\n",rk);
435: rk = PetscMin(rk,ctx->rkc);
436: DSGetMat(ctx->ds,DS_MAT_U,&U);
437: BVMultInPlace(V,U,0,rk);
438: BVSetActiveColumns(V,0,rk);
439: MatDestroy(&U);
441: /* Rank reduction */
442: DSSetDimensions(ctx->ds,rk,rk,0,0);
443: DSGetMat(ctx->ds,DS_MAT_A,&W);
444: BVMatProject(V,S,V,W);
445: LyapIIBuildEigenMat(ctx->lme,W,&Op,&v0); /* Op=A\B, A=kron(I,S)+kron(S,I), B=-2*kron(S,S) */
446: EPSSetOperators(epsrr,Op,NULL);
447: EPSSetInitialSpace(epsrr,1,&v0);
448: EPSSolve(epsrr);
449: MatDestroy(&W);
450: EPSComputeVectors(epsrr);
451: /* Copy first eigenvector, vec(A)=x */
452: BVGetArray(epsrr->V,&xx);
453: DSGetArray(ctx->ds,DS_MAT_A,&aa);
454: for (i=0;i<rk;i++) {
455: PetscMemcpy(aa+i*ldds,xx+i*rk,rk*sizeof(PetscScalar));
456: }
457: DSRestoreArray(ctx->ds,DS_MAT_A,&aa);
458: BVRestoreArray(epsrr->V,&xx);
459: DSSetState(ctx->ds,DS_STATE_RAW);
460: /* Compute [U,Sigma,~] = svd(A), its rank should be 1 or 2 */
461: DSSolve(ctx->ds,s,NULL);
462: if (PetscAbsScalar(s[1]/s[0])<PETSC_SQRT_MACHINE_EPSILON) rk=1;
463: else rk = 2;
464: PetscInfo1(eps,"The eigenvector has rank %D\n",rk);
465: DSGetMat(ctx->ds,DS_MAT_U,&U);
466: BVMultInPlace(V,U,0,rk);
467: MatDestroy(&U);
469: /* Save V in Ux */
470: idx = (rk==2)?1:0;
471: for (i=0;i<rk;i++) {
472: BVGetColumn(V,i,&v);
473: VecGetArray(v,&uu);
474: MatDenseGetColumn(Ux[idx],i,&array);
475: PetscMemcpy(array,uu,eps->nloc*sizeof(PetscScalar));
476: MatDenseRestoreColumn(Ux[idx],&array);
477: VecRestoreArray(v,&uu);
478: BVRestoreColumn(V,i,&v);
479: }
481: /* Eigenpair approximation */
482: BVGetColumn(V,0,&v);
483: MatMult(S,v,z);
484: VecDot(z,v,pM);
485: BVRestoreColumn(V,0,&v);
486: if (rk>1) {
487: BVGetColumn(V,1,&w);
488: VecDot(z,w,pM+1);
489: MatMult(S,w,z);
490: VecDot(z,w,pM+3);
491: BVGetColumn(V,0,&v);
492: VecDot(z,v,pM+2);
493: BVRestoreColumn(V,0,&v);
494: BVRestoreColumn(V,1,&w);
495: EV2x2(pM,2,eigr,eigi,vec);
496: MatCreateSeqDense(PETSC_COMM_SELF,2,2,vec,&X);
497: BVSetActiveColumns(V,0,rk);
498: BVMultInPlace(V,X,0,rk);
499: MatDestroy(&X);
500: #if !defined(PETSC_USE_COMPLEX)
501: norm = eigr[0]*eigr[0]+eigi[0]*eigi[0];
502: er = eigr[0]/norm; ei = -eigi[0]/norm;
503: #else
504: er =1.0/eigr[0]; ei = 0.0;
505: #endif
506: } else {
507: eigr[0] = pM[0]; eigi[0] = 0.0;
508: er = 1.0/eigr[0]; ei = 0.0;
509: }
510: BVGetColumn(V,0,&v);
511: if (eigi[0]!=0.0) {
512: BVGetColumn(V,1,&w);
513: } else w = NULL;
514: eps->eigr[eps->nconv] = eigr[0]; eps->eigi[eps->nconv] = eigi[0];
515: EPSComputeResidualNorm_Private(eps,PETSC_FALSE,er,ei,v,w,eps->work,&norm);
516: BVRestoreColumn(V,0,&v);
517: if (w) {
518: BVRestoreColumn(V,1,&w);
519: }
520: (*eps->converged)(eps,er,ei,norm,&eps->errest[eps->nconv],eps->convergedctx);
521: k = 0;
522: if (eps->errest[eps->nconv]<eps->tol) {
523: k++;
524: if (rk==2) {
525: #if !defined (PETSC_USE_COMPLEX)
526: eps->eigr[eps->nconv+k] = eigr[0]; eps->eigi[eps->nconv+k] = -eigi[0];
527: #else
528: eps->eigr[eps->nconv+k] = PetscConj(eps->eigr[eps->nconv]);
529: #endif
530: k++;
531: }
532: /* Store converged eigenpairs and vectors for deflation */
533: for (i=0;i<k;i++) {
534: BVGetColumn(V,i,&v);
535: BVInsertVec(eps->V,eps->nconv+i,v);
536: BVInsertVec(matctx->Q,eps->nconv+i,v);
537: BVRestoreColumn(V,i,&v);
538: }
539: eps->nconv += k;
540: if (eps->nconv<eps->nev) {
541: BVSetActiveColumns(matctx->Q,eps->nconv-rk,eps->nconv);
542: BVOrthogonalize(matctx->Q,NULL);
543: idx = 0;
544: BVSetRandomColumn(V,0);
545: BVNormColumn(V,0,NORM_2,&norm);
546: BVScaleColumn(V,0,1.0/norm);
547: LyapIIBuildRHS(S,1,Ux[idx],V,eps->work);
548: }
549: } else {
550: /* Prepare right-hand side */
551: LyapIIBuildRHS(S,rk,Ux[idx],NULL,eps->work);
552: }
553: (*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx);
554: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->nconv+1);
555: }
556: STRestoreOperator(eps->st,&matctx->S);
557: MatDestroy(&S);
558: MatDestroy(&Ux[0]);
559: MatDestroy(&Ux[1]);
560: MatDestroy(&Op);
561: VecDestroy(&v0);
562: BVDestroy(&V);
563: EPSDestroy(&epsrr);
564: PetscFree3(s,eigr,eigi);
565: return(0);
566: }
568: PetscErrorCode EPSSetFromOptions_LyapII(PetscOptionItems *PetscOptionsObject,EPS eps)
569: {
571: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
572: PetscInt k,array[2]={PETSC_DEFAULT,PETSC_DEFAULT};
573: PetscBool flg;
576: PetscOptionsHead(PetscOptionsObject,"EPS Lyapunov Inverse Iteration Options");
578: k = 2;
579: PetscOptionsIntArray("-eps_lyapii_ranks","Ranks for Lyapunov equation (one or two comma-separated integers)","EPSLyapIISetRanks",array,&k,&flg);
580: if (flg) {
581: EPSLyapIISetRanks(eps,array[0],array[1]);
582: }
584: PetscOptionsTail();
586: if (!ctx->lme) { EPSLyapIIGetLME(eps,&ctx->lme); }
587: LMESetFromOptions(ctx->lme);
588: return(0);
589: }
591: static PetscErrorCode EPSLyapIISetRanks_LyapII(EPS eps,PetscInt rkc,PetscInt rkl)
592: {
593: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
596: if (rkc==PETSC_DEFAULT) rkc = 10;
597: if (rkc<2) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The compressed rank %D must be larger than 1",rkc);
598: if (rkl==PETSC_DEFAULT) rkl = 3*rkc;
599: 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);
600: if (rkc != ctx->rkc) {
601: ctx->rkc = rkc;
602: eps->state = EPS_STATE_INITIAL;
603: }
604: if (rkl != ctx->rkl) {
605: ctx->rkl = rkl;
606: eps->state = EPS_STATE_INITIAL;
607: }
608: return(0);
609: }
611: /*@
612: EPSLyapIISetRanks - Set the ranks used in the solution of the Lyapunov equation.
614: Collective on EPS
616: Input Parameters:
617: + eps - the eigenproblem solver context
618: . rkc - the compressed rank
619: - rkl - the Lyapunov rank
621: Options Database Key:
622: . -eps_lyapii_ranks <rkc,rkl> - Sets the rank parameters
624: Notes:
625: Lyapunov inverse iteration needs to solve a large-scale Lyapunov equation
626: at each iteration of the eigensolver. For this, an iterative solver (LME)
627: is used, which requires to prescribe the rank of the solution matrix X. This
628: is the meaning of parameter rkl. Later, this matrix is compressed into
629: another matrix of rank rkc. If not provided, rkl is a small multiple of rkc.
631: Level: intermediate
633: .seealso: EPSLyapIIGetRanks()
634: @*/
635: PetscErrorCode EPSLyapIISetRanks(EPS eps,PetscInt rkc,PetscInt rkl)
636: {
643: PetscTryMethod(eps,"EPSLyapIISetRanks_C",(EPS,PetscInt,PetscInt),(eps,rkc,rkl));
644: return(0);
645: }
647: static PetscErrorCode EPSLyapIIGetRanks_LyapII(EPS eps,PetscInt *rkc,PetscInt *rkl)
648: {
649: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
652: if (rkc) *rkc = ctx->rkc;
653: if (rkl) *rkl = ctx->rkl;
654: return(0);
655: }
657: /*@
658: EPSLyapIIGetRanks - Return the rank values used for the Lyapunov step.
660: Not Collective
662: Input Parameter:
663: . eps - the eigenproblem solver context
665: Output Parameters:
666: + rkc - the compressed rank
667: - rkl - the Lyapunov rank
669: Level: intermediate
671: .seealso: EPSLyapIISetRanks()
672: @*/
673: PetscErrorCode EPSLyapIIGetRanks(EPS eps,PetscInt *rkc,PetscInt *rkl)
674: {
679: PetscUseMethod(eps,"EPSLyapIIGetRanks_C",(EPS,PetscInt*,PetscInt*),(eps,rkc,rkl));
680: return(0);
681: }
683: static PetscErrorCode EPSLyapIISetLME_LyapII(EPS eps,LME lme)
684: {
686: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
689: PetscObjectReference((PetscObject)lme);
690: LMEDestroy(&ctx->lme);
691: ctx->lme = lme;
692: PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->lme);
693: eps->state = EPS_STATE_INITIAL;
694: return(0);
695: }
697: /*@
698: EPSLyapIISetLME - Associate a linear matrix equation solver object (LME) to the
699: eigenvalue solver.
701: Collective on EPS
703: Input Parameters:
704: + eps - the eigenproblem solver context
705: - lme - the linear matrix equation solver object
707: Level: advanced
709: .seealso: EPSLyapIIGetLME()
710: @*/
711: PetscErrorCode EPSLyapIISetLME(EPS eps,LME lme)
712: {
719: PetscTryMethod(eps,"EPSLyapIISetLME_C",(EPS,LME),(eps,lme));
720: return(0);
721: }
723: static PetscErrorCode EPSLyapIIGetLME_LyapII(EPS eps,LME *lme)
724: {
726: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
729: if (!ctx->lme) {
730: LMECreate(PetscObjectComm((PetscObject)eps),&ctx->lme);
731: LMESetOptionsPrefix(ctx->lme,((PetscObject)eps)->prefix);
732: LMEAppendOptionsPrefix(ctx->lme,"eps_lyapii_");
733: PetscObjectIncrementTabLevel((PetscObject)ctx->lme,(PetscObject)eps,1);
734: PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->lme);
735: }
736: *lme = ctx->lme;
737: return(0);
738: }
740: /*@
741: EPSLyapIIGetLME - Retrieve the linear matrix equation solver object (LME)
742: associated with the eigenvalue solver.
744: Not Collective
746: Input Parameter:
747: . eps - the eigenproblem solver context
749: Output Parameter:
750: . lme - the linear matrix equation solver object
752: Level: advanced
754: .seealso: EPSLyapIISetLME()
755: @*/
756: PetscErrorCode EPSLyapIIGetLME(EPS eps,LME *lme)
757: {
763: PetscUseMethod(eps,"EPSLyapIIGetLME_C",(EPS,LME*),(eps,lme));
764: return(0);
765: }
767: PetscErrorCode EPSView_LyapII(EPS eps,PetscViewer viewer)
768: {
770: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
771: PetscBool isascii;
774: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
775: if (isascii) {
776: PetscViewerASCIIPrintf(viewer," ranks: for Lyapunov solver=%D, after compression=%D\n",ctx->rkl,ctx->rkc);
777: if (!ctx->lme) { EPSLyapIIGetLME(eps,&ctx->lme); }
778: PetscViewerASCIIPushTab(viewer);
779: LMEView(ctx->lme,viewer);
780: PetscViewerASCIIPopTab(viewer);
781: }
782: return(0);
783: }
785: PetscErrorCode EPSReset_LyapII(EPS eps)
786: {
788: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
791: if (!ctx->lme) { LMEReset(ctx->lme); }
792: return(0);
793: }
795: PetscErrorCode EPSDestroy_LyapII(EPS eps)
796: {
798: EPS_LYAPII *ctx = (EPS_LYAPII*)eps->data;
801: LMEDestroy(&ctx->lme);
802: DSDestroy(&ctx->ds);
803: PetscFree(eps->data);
804: PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIISetLME_C",NULL);
805: PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIIGetLME_C",NULL);
806: PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIISetRanks_C",NULL);
807: PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIIGetRanks_C",NULL);
808: return(0);
809: }
811: PetscErrorCode EPSSetDefaultST_LyapII(EPS eps)
812: {
816: if (!((PetscObject)eps->st)->type_name) {
817: STSetType(eps->st,STSINVERT);
818: }
819: return(0);
820: }
822: PETSC_EXTERN PetscErrorCode EPSCreate_LyapII(EPS eps)
823: {
824: EPS_LYAPII *ctx;
828: PetscNewLog(eps,&ctx);
829: eps->data = (void*)ctx;
831: eps->ops->solve = EPSSolve_LyapII;
832: eps->ops->setup = EPSSetUp_LyapII;
833: eps->ops->setfromoptions = EPSSetFromOptions_LyapII;
834: eps->ops->reset = EPSReset_LyapII;
835: eps->ops->destroy = EPSDestroy_LyapII;
836: eps->ops->view = EPSView_LyapII;
837: eps->ops->setdefaultst = EPSSetDefaultST_LyapII;
838: eps->ops->backtransform = EPSBackTransform_Default;
840: PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIISetLME_C",EPSLyapIISetLME_LyapII);
841: PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIIGetLME_C",EPSLyapIIGetLME_LyapII);
842: PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIISetRanks_C",EPSLyapIISetRanks_LyapII);
843: PetscObjectComposeFunction((PetscObject)eps,"EPSLyapIIGetRanks_C",EPSLyapIIGetRanks_LyapII);
844: return(0);
845: }