LCOV - code coverage report
Current view: top level - eps/impls/lyapii - lyapii.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 423 443 95.5 %
Date: 2024-11-21 00:40:22 Functions: 21 23 91.3 %
Legend: Lines: hit not hit

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

Generated by: LCOV version 1.14