LCOV - code coverage report
Current view: top level - eps/impls/power - power.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 628 680 92.4 %
Date: 2024-11-23 00:34:26 Functions: 37 40 92.5 %
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: "power"
      12             : 
      13             :    Method: Power Iteration
      14             : 
      15             :    Algorithm:
      16             : 
      17             :        This solver implements the power iteration for finding dominant
      18             :        eigenpairs. It also includes the following well-known methods:
      19             :        - Inverse Iteration: when used in combination with shift-and-invert
      20             :          spectral transformation.
      21             :        - Rayleigh Quotient Iteration (RQI): also with shift-and-invert plus
      22             :          a variable shift.
      23             : 
      24             :        It can also be used for nonlinear inverse iteration on the problem
      25             :        A(x)*x=lambda*B(x)*x, where A and B are not constant but depend on x.
      26             : 
      27             :    References:
      28             : 
      29             :        [1] "Single Vector Iteration Methods in SLEPc", SLEPc Technical Report
      30             :            STR-2, available at https://slepc.upv.es.
      31             : */
      32             : 
      33             : #include <slepc/private/epsimpl.h>                /*I "slepceps.h" I*/
      34             : #include <slepcblaslapack.h>
      35             : /* petsc headers */
      36             : #include <petscdm.h>
      37             : #include <petscsnes.h>
      38             : 
      39             : static PetscErrorCode EPSPowerFormFunction_Update(SNES,Vec,Vec,void*);
      40             : static PetscErrorCode EPSSolve_Power(EPS);
      41             : static PetscErrorCode EPSSolve_TS_Power(EPS);
      42             : 
      43             : typedef struct {
      44             :   EPSPowerShiftType shift_type;
      45             :   PetscBool         nonlinear;
      46             :   PetscBool         update;
      47             :   PetscBool         sign_normalization;
      48             :   SNES              snes;
      49             :   PetscErrorCode    (*formFunctionB)(SNES,Vec,Vec,void*);
      50             :   void              *formFunctionBctx;
      51             :   PetscErrorCode    (*formFunctionA)(SNES,Vec,Vec,void*);
      52             :   void              *formFunctionActx;
      53             :   PetscErrorCode    (*formFunctionAB)(SNES,Vec,Vec,Vec,void*);
      54             :   PetscErrorCode    (*formNorm)(SNES,Vec,PetscReal*,void*);
      55             :   void              *formNormCtx;
      56             :   PetscInt          idx;  /* index of the first nonzero entry in the iteration vector */
      57             :   PetscMPIInt       p;    /* process id of the owner of idx */
      58             :   PetscReal         norm0; /* norm of initial vector */
      59             : } EPS_POWER;
      60             : 
      61          82 : static PetscErrorCode SNESMonitor_PowerUpdate(SNES snes,PetscInt its,PetscReal fnorm,void *ctx)
      62             : {
      63          82 :   EPS            eps;
      64             : 
      65          82 :   PetscFunctionBegin;
      66          82 :   PetscCall(PetscObjectQuery((PetscObject)snes,"eps",(PetscObject *)&eps));
      67          82 :   PetscCheck(eps,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_NULL,"No composed EPS");
      68             :   /* Call EPS monitor on each SNES iteration */
      69          82 :   PetscCall(EPSMonitor(eps,its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev)));
      70          82 :   PetscFunctionReturn(PETSC_SUCCESS);
      71             : }
      72             : 
      73          54 : static PetscErrorCode EPSSetUp_Power(EPS eps)
      74             : {
      75          54 :   EPS_POWER      *power = (EPS_POWER*)eps->data;
      76          54 :   STMatMode      mode;
      77          54 :   Mat            A,B,P;
      78          54 :   Vec            res;
      79          54 :   PetscContainer container;
      80          54 :   PetscErrorCode (*formFunctionA)(SNES,Vec,Vec,void*);
      81          54 :   PetscErrorCode (*formJacobianA)(SNES,Vec,Mat,Mat,void*);
      82          54 :   PetscErrorCode (*formNorm)(SNES,Vec,PetscReal*,void*);
      83          54 :   void           *ctx;
      84             : 
      85          54 :   PetscFunctionBegin;
      86          54 :   EPSCheckNotStructured(eps);
      87          54 :   if (eps->ncv!=PETSC_DETERMINE) {
      88           5 :     PetscCheck(eps->ncv>=eps->nev,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv must be at least nev");
      89          49 :   } else eps->ncv = eps->nev;
      90          54 :   if (eps->mpd!=PETSC_DETERMINE) PetscCall(PetscInfo(eps,"Warning: parameter mpd ignored\n"));
      91          54 :   if (eps->max_it==PETSC_DETERMINE) {
      92             :     /* SNES will directly return the solution for us, and we need to do only one iteration */
      93          32 :     if (power->nonlinear && power->update) eps->max_it = 1;
      94          18 :     else eps->max_it = PetscMax(1000*eps->nev,100*eps->n);
      95             :   }
      96          54 :   if (!eps->which) PetscCall(EPSSetWhichEigenpairs_Default(eps));
      97          54 :   PetscCheck(eps->which==EPS_LARGEST_MAGNITUDE || eps->which==EPS_TARGET_MAGNITUDE,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver supports only largest magnitude or target magnitude eigenvalues");
      98          54 :   if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) {
      99           6 :     PetscCheck(!power->nonlinear,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Variable shifts not allowed in nonlinear problems");
     100           6 :     EPSCheckSinvertCayleyCondition(eps,PETSC_TRUE," (with variable shifts)");
     101           6 :     PetscCall(STGetMatMode(eps->st,&mode));
     102           6 :     PetscCheck(mode!=ST_MATMODE_INPLACE,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"ST matrix mode inplace does not work with variable shifts");
     103             :   }
     104          54 :   EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_CONVERGENCE);
     105          54 :   EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION);
     106          54 :   PetscCall(EPSAllocateSolution(eps,0));
     107          54 :   PetscCall(EPS_SetInnerProduct(eps));
     108             : 
     109          54 :   if (power->nonlinear) {
     110          35 :     PetscCheck(eps->nev==1,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Nonlinear inverse iteration cannot compute more than one eigenvalue");
     111          35 :     PetscCall(EPSSetWorkVecs(eps,3));
     112          35 :     PetscCheck(!power->update || eps->max_it==1,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"More than one iteration is not allowed for Newton eigensolver (SNES)");
     113             : 
     114             :     /* Set up SNES solver */
     115             :     /* If SNES was setup, we need to reset it so that it will be consistent with the current EPS */
     116          35 :     if (power->snes) PetscCall(SNESReset(power->snes));
     117           0 :     else PetscCall(EPSPowerGetSNES(eps,&power->snes));
     118             : 
     119             :     /* For nonlinear solver (Newton), we should scale the initial vector back.
     120             :        The initial vector will be scaled in EPSSetUp. */
     121          35 :     if (eps->IS) PetscCall(VecNorm(eps->IS[0],NORM_2,&power->norm0));
     122             : 
     123          35 :     PetscCall(EPSGetOperators(eps,&A,&B));
     124             : 
     125             :     /* form function */
     126          35 :     PetscCall(PetscObjectQueryFunction((PetscObject)A,"formFunction",&formFunctionA));
     127          35 :     PetscCheck(formFunctionA,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER,"For nonlinear inverse iteration you must compose a callback function 'formFunction' in matrix A");
     128          35 :     PetscCall(PetscObjectQuery((PetscObject)A,"formFunctionCtx",(PetscObject*)&container));
     129          35 :     if (container) PetscCall(PetscContainerGetPointer(container,&ctx));
     130           0 :     else ctx = NULL;
     131          35 :     PetscCall(MatCreateVecs(A,&res,NULL));
     132          35 :     power->formFunctionA = formFunctionA;
     133          35 :     power->formFunctionActx = ctx;
     134          35 :     if (power->update) {
     135          16 :       PetscCall(SNESSetFunction(power->snes,res,EPSPowerFormFunction_Update,ctx));
     136          16 :       PetscCall(PetscObjectQueryFunction((PetscObject)A,"formFunctionAB",&power->formFunctionAB));
     137          16 :       PetscCall(SNESMonitorSet(power->snes,SNESMonitor_PowerUpdate,NULL,NULL));
     138             :     }
     139          19 :     else PetscCall(SNESSetFunction(power->snes,res,formFunctionA,ctx));
     140          35 :     PetscCall(VecDestroy(&res));
     141             : 
     142             :     /* form Jacobian */
     143          35 :     PetscCall(PetscObjectQueryFunction((PetscObject)A,"formJacobian",&formJacobianA));
     144          35 :     PetscCheck(formJacobianA,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER,"For nonlinear inverse iteration you must compose a callback function 'formJacobian' in matrix A");
     145          35 :     PetscCall(PetscObjectQuery((PetscObject)A,"formJacobianCtx",(PetscObject*)&container));
     146          35 :     if (container) PetscCall(PetscContainerGetPointer(container,&ctx));
     147           0 :     else ctx = NULL;
     148             :     /* This allows users to compute a different preconditioning matrix than A.
     149             :      * It is useful when A and B are shell matrices.
     150             :      */
     151          35 :     PetscCall(STGetPreconditionerMat(eps->st,&P));
     152             :     /* If users do not provide a matrix, we simply use A */
     153          35 :     PetscCall(SNESSetJacobian(power->snes,A,P? P:A,formJacobianA,ctx));
     154          35 :     PetscCall(SNESSetFromOptions(power->snes));
     155          35 :     PetscCall(SNESSetUp(power->snes));
     156             : 
     157          35 :     ctx = NULL;
     158          35 :     formNorm = NULL;
     159          35 :     if (B) {
     160          35 :       PetscCall(PetscObjectQueryFunction((PetscObject)B,"formFunction",&power->formFunctionB));
     161          35 :       PetscCall(PetscObjectQuery((PetscObject)B,"formFunctionCtx",(PetscObject*)&container));
     162          35 :       if (power->formFunctionB && container) PetscCall(PetscContainerGetPointer(container,&power->formFunctionBctx));
     163           0 :       else power->formFunctionBctx = NULL;
     164             : 
     165             :       /* form norm */
     166          35 :       PetscCall(PetscObjectQueryFunction((PetscObject)B,"formNorm",&formNorm));
     167          35 :       if (formNorm) {
     168          10 :         PetscCall(PetscObjectQuery((PetscObject)B,"formNormCtx",(PetscObject*)&container));
     169          10 :         if (container) PetscCall(PetscContainerGetPointer(container,&ctx));
     170             :       }
     171             :     }
     172          35 :     power->formNorm = formNorm;
     173          35 :     power->formNormCtx = ctx;
     174             :   } else {
     175          19 :     if (eps->twosided) PetscCall(EPSSetWorkVecs(eps,3));
     176          14 :     else PetscCall(EPSSetWorkVecs(eps,2));
     177          19 :     PetscCall(DSSetType(eps->ds,DSNHEP));
     178          19 :     PetscCall(DSAllocate(eps->ds,eps->nev));
     179             :   }
     180             :   /* dispatch solve method */
     181          54 :   if (eps->twosided) {
     182           5 :     PetscCheck(!power->nonlinear,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Nonlinear inverse iteration does not have two-sided variant");
     183           5 :     PetscCheck(power->shift_type!=EPS_POWER_SHIFT_WILKINSON,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Two-sided variant does not support Wilkinson shifts");
     184           5 :     eps->ops->solve = EPSSolve_TS_Power;
     185          49 :   } else eps->ops->solve = EPSSolve_Power;
     186          54 :   PetscFunctionReturn(PETSC_SUCCESS);
     187             : }
     188             : 
     189             : /*
     190             :    Find the index of the first nonzero entry of x, and the process that owns it.
     191             : */
     192          35 : static PetscErrorCode FirstNonzeroIdx(Vec x,PetscInt *idx,PetscMPIInt *p)
     193             : {
     194          35 :   PetscInt          i,first,last,N;
     195          35 :   PetscLayout       map;
     196          35 :   const PetscScalar *xx;
     197             : 
     198          35 :   PetscFunctionBegin;
     199          35 :   PetscCall(VecGetSize(x,&N));
     200          35 :   PetscCall(VecGetOwnershipRange(x,&first,&last));
     201          35 :   PetscCall(VecGetArrayRead(x,&xx));
     202         385 :   for (i=first;i<last;i++) {
     203         385 :     if (PetscAbsScalar(xx[i-first])>10*PETSC_MACHINE_EPSILON) break;
     204             :   }
     205          35 :   PetscCall(VecRestoreArrayRead(x,&xx));
     206          35 :   if (i==last) i=N;
     207         105 :   PetscCallMPI(MPIU_Allreduce(&i,idx,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)x)));
     208          35 :   PetscCheck(*idx!=N,PetscObjectComm((PetscObject)x),PETSC_ERR_PLIB,"Zero vector found");
     209          35 :   PetscCall(VecGetLayout(x,&map));
     210          35 :   PetscCall(PetscLayoutFindOwner(map,*idx,p));
     211          35 :   PetscFunctionReturn(PETSC_SUCCESS);
     212             : }
     213             : 
     214             : /*
     215             :    Normalize a vector x with respect to a given norm as well as, optionally, the
     216             :    sign of the first nonzero entry (assumed to be idx in process p).
     217             : */
     218       29082 : static PetscErrorCode Normalize(Vec x,PetscReal norm,PetscInt idx,PetscMPIInt p,PetscBool sign_normalization,PetscScalar *sign)
     219             : {
     220       29082 :   PetscScalar       alpha=1.0;
     221       29082 :   PetscInt          first,last;
     222       29082 :   PetscReal         absv;
     223       29082 :   const PetscScalar *xx;
     224             : 
     225       29082 :   PetscFunctionBegin;
     226       29082 :   if (sign_normalization) {
     227       28874 :     PetscCall(VecGetOwnershipRange(x,&first,&last));
     228       28874 :     if (idx>=first && idx<last) {
     229       28874 :       PetscCall(VecGetArrayRead(x,&xx));
     230       28874 :       absv = PetscAbsScalar(xx[idx-first]);
     231       28874 :       if (absv>10*PETSC_MACHINE_EPSILON) alpha = xx[idx-first]/absv;
     232       28874 :       PetscCall(VecRestoreArrayRead(x,&xx));
     233             :     }
     234       57748 :     PetscCallMPI(MPI_Bcast(&alpha,1,MPIU_SCALAR,p,PetscObjectComm((PetscObject)x)));
     235             :   }
     236       29082 :   if (sign) *sign = alpha;
     237       29082 :   alpha *= norm;
     238       29082 :   PetscCall(VecScale(x,1.0/alpha));
     239       29082 :   PetscFunctionReturn(PETSC_SUCCESS);
     240             : }
     241             : 
     242         433 : static PetscErrorCode EPSPowerUpdateFunctionB(EPS eps,Vec x,Vec Bx)
     243             : {
     244         433 :   EPS_POWER      *power = (EPS_POWER*)eps->data;
     245         433 :   Mat            B;
     246             : 
     247         433 :   PetscFunctionBegin;
     248         433 :   PetscCall(STResetMatrixState(eps->st));
     249         433 :   PetscCall(EPSGetOperators(eps,NULL,&B));
     250         433 :   if (B) {
     251         433 :     if (power->formFunctionB) PetscCall((*power->formFunctionB)(power->snes,x,Bx,power->formFunctionBctx));
     252           0 :     else PetscCall(MatMult(B,x,Bx));
     253           0 :   } else PetscCall(VecCopy(x,Bx));
     254         433 :   PetscFunctionReturn(PETSC_SUCCESS);
     255             : }
     256             : 
     257         222 : static PetscErrorCode EPSPowerUpdateFunctionA(EPS eps,Vec x,Vec Ax)
     258             : {
     259         222 :   EPS_POWER      *power = (EPS_POWER*)eps->data;
     260         222 :   Mat            A;
     261             : 
     262         222 :   PetscFunctionBegin;
     263         222 :   PetscCall(STResetMatrixState(eps->st));
     264         222 :   PetscCall(EPSGetOperators(eps,&A,NULL));
     265         222 :   PetscCheck(A,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_NULL,"Matrix A is required for an eigenvalue problem");
     266         222 :   if (power->formFunctionA) PetscCall((*power->formFunctionA)(power->snes,x,Ax,power->formFunctionActx));
     267           0 :   else PetscCall(MatMult(A,x,Ax));
     268         222 :   PetscFunctionReturn(PETSC_SUCCESS);
     269             : }
     270             : 
     271         444 : static PetscErrorCode EPSPowerFormFunction_Update(SNES snes,Vec x,Vec y,void *ctx)
     272             : {
     273         444 :   EPS            eps;
     274         444 :   PetscReal      bx;
     275         444 :   Vec            Bx;
     276         444 :   PetscScalar    sign;
     277         444 :   EPS_POWER      *power;
     278             : 
     279         444 :   PetscFunctionBegin;
     280         444 :   PetscCall(PetscObjectQuery((PetscObject)snes,"eps",(PetscObject *)&eps));
     281         444 :   PetscCheck(eps,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_NULL,"No composed EPS");
     282         444 :   PetscCall(STResetMatrixState(eps->st));
     283         444 :   power = (EPS_POWER*)eps->data;
     284         444 :   Bx = eps->work[2];
     285         444 :   if (power->formFunctionAB) PetscCall((*power->formFunctionAB)(snes,x,y,Bx,ctx));
     286             :   else {
     287         222 :     PetscCall(EPSPowerUpdateFunctionA(eps,x,y));
     288         222 :     PetscCall(EPSPowerUpdateFunctionB(eps,x,Bx));
     289             :   }
     290         444 :   if (power->formNorm) PetscCall((*power->formNorm)(snes,Bx,&bx,power->formNormCtx));
     291         312 :   else PetscCall(VecNorm(Bx,NORM_2,&bx));
     292         444 :   PetscCall(Normalize(Bx,bx,power->idx,power->p,power->sign_normalization,&sign));
     293         444 :   PetscCall(VecAXPY(y,-1.0,Bx));
     294             :   /* Keep tracking eigenvalue update. It would be useful when we want to monitor solver progress via snes monitor. */
     295         444 :   eps->eigr[(eps->nconv < eps->nev)? eps->nconv:(eps->nconv-1)] = 1.0/(bx*sign);
     296         444 :   PetscFunctionReturn(PETSC_SUCCESS);
     297             : }
     298             : 
     299             : /*
     300             :    Use SNES to compute y = A^{-1}*B*x for the nonlinear problem
     301             : */
     302         192 : static PetscErrorCode EPSPowerApply_SNES(EPS eps,Vec x,Vec y)
     303             : {
     304         192 :   EPS_POWER      *power = (EPS_POWER*)eps->data;
     305         192 :   Vec            Bx;
     306             : 
     307         192 :   PetscFunctionBegin;
     308         192 :   PetscCall(VecCopy(x,y));
     309         192 :   if (power->update) PetscCall(SNESSolve(power->snes,NULL,y));
     310             :   else {
     311         176 :     Bx = eps->work[2];
     312         176 :     PetscCall(SNESSolve(power->snes,Bx,y));
     313             :   }
     314         192 :   PetscFunctionReturn(PETSC_SUCCESS);
     315             : }
     316             : 
     317             : /*
     318             :    Use nonlinear inverse power to compute an initial guess
     319             : */
     320          14 : static PetscErrorCode EPSPowerComputeInitialGuess_Update(EPS eps)
     321             : {
     322          14 :   EPS            powereps;
     323          14 :   Mat            A,B,P;
     324          14 :   Vec            v1,v2;
     325          14 :   SNES           snes;
     326          14 :   DM             dm,newdm;
     327          14 :   PetscBool      sign_normalization;
     328             : 
     329          14 :   PetscFunctionBegin;
     330          14 :   PetscCall(EPSCreate(PetscObjectComm((PetscObject)eps),&powereps));
     331          14 :   PetscCall(EPSGetOperators(eps,&A,&B));
     332          14 :   PetscCall(EPSSetType(powereps,EPSPOWER));
     333          14 :   PetscCall(EPSSetOperators(powereps,A,B));
     334          14 :   PetscCall(EPSSetTolerances(powereps,1e-6,4));
     335          14 :   PetscCall(EPSSetOptionsPrefix(powereps,((PetscObject)eps)->prefix));
     336          14 :   PetscCall(EPSAppendOptionsPrefix(powereps,"init_"));
     337          14 :   PetscCall(EPSSetProblemType(powereps,EPS_GNHEP));
     338          14 :   PetscCall(EPSSetWhichEigenpairs(powereps,EPS_TARGET_MAGNITUDE));
     339          14 :   PetscCall(EPSPowerSetNonlinear(powereps,PETSC_TRUE));
     340          14 :   PetscCall(STGetPreconditionerMat(eps->st,&P));
     341             :   /* attach dm to initial solve */
     342          14 :   PetscCall(EPSPowerGetSNES(eps,&snes));
     343          14 :   PetscCall(SNESGetDM(snes,&dm));
     344             :   /* use  dmshell to temporarily store snes context */
     345          14 :   PetscCall(DMCreate(PetscObjectComm((PetscObject)eps),&newdm));
     346          14 :   PetscCall(DMSetType(newdm,DMSHELL));
     347          14 :   PetscCall(DMSetUp(newdm));
     348          14 :   PetscCall(DMCopyDMSNES(dm,newdm));
     349          14 :   PetscCall(EPSPowerGetSNES(powereps,&snes));
     350          14 :   PetscCall(SNESSetDM(snes,dm));
     351          14 :   PetscCall(EPSPowerGetSignNormalization(eps,&sign_normalization));
     352          14 :   PetscCall(EPSPowerSetSignNormalization(powereps,sign_normalization));
     353          14 :   PetscCall(EPSSetFromOptions(powereps));
     354          14 :   if (P) PetscCall(STSetPreconditionerMat(powereps->st,P));
     355          14 :   PetscCall(EPSSolve(powereps));
     356          14 :   PetscCall(BVGetColumn(eps->V,0,&v2));
     357          14 :   PetscCall(BVGetColumn(powereps->V,0,&v1));
     358          14 :   PetscCall(VecCopy(v1,v2));
     359          14 :   PetscCall(BVRestoreColumn(powereps->V,0,&v1));
     360          14 :   PetscCall(BVRestoreColumn(eps->V,0,&v2));
     361          14 :   PetscCall(EPSDestroy(&powereps));
     362             :   /* restore context back to the old nonlinear solver */
     363          14 :   PetscCall(DMCopyDMSNES(newdm,dm));
     364          14 :   PetscCall(DMDestroy(&newdm));
     365          14 :   PetscFunctionReturn(PETSC_SUCCESS);
     366             : }
     367             : 
     368          49 : static PetscErrorCode EPSSolve_Power(EPS eps)
     369             : {
     370          49 :   EPS_POWER           *power = (EPS_POWER*)eps->data;
     371          49 :   PetscInt            k,ld;
     372          49 :   Vec                 v,y,e,Bx;
     373          49 :   Mat                 A;
     374          49 :   KSP                 ksp;
     375          49 :   PetscReal           relerr,norm,norm1,rt1,rt2,cs1;
     376          49 :   PetscScalar         theta,rho,delta,sigma,alpha2,beta1,sn1,*T,sign;
     377          49 :   PetscBool           breakdown;
     378          49 :   KSPConvergedReason  reason;
     379          49 :   SNESConvergedReason snesreason;
     380             : 
     381          49 :   PetscFunctionBegin;
     382          49 :   e = eps->work[0];
     383          49 :   y = eps->work[1];
     384          49 :   if (power->nonlinear) Bx = eps->work[2];
     385             :   else Bx = NULL;
     386             : 
     387          49 :   if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) PetscCall(STGetKSP(eps->st,&ksp));
     388          49 :   if (power->nonlinear) {
     389          35 :     PetscCall(PetscObjectCompose((PetscObject)power->snes,"eps",(PetscObject)eps));
     390             :     /* Compute an initial guess only when users do not provide one */
     391          35 :     if (power->update && !eps->nini) PetscCall(EPSPowerComputeInitialGuess_Update(eps));
     392          14 :   } else PetscCall(DSGetLeadingDimension(eps->ds,&ld));
     393          49 :   if (!power->update) PetscCall(EPSGetStartVector(eps,0,NULL));
     394          49 :   if (power->nonlinear) {
     395          35 :     PetscCall(BVGetColumn(eps->V,0,&v));
     396          35 :     if (eps->nini) {
     397             :       /* We scale the initial vector back if the initial vector was provided by users */
     398           2 :       PetscCall(VecScale(v,power->norm0));
     399             :     }
     400             :     /* Do a couple things:
     401             :      * 1) Determine the first non-zero index for Bx and the proc that owns it. This will be
     402             :      *    used if performing normalization by the sign of the first nonzero element of Bx.
     403             :      * 2) Scale Bx by its norm. For non-update power iterations, Bx (in code terms) is used
     404             :      *    as the RHS argument to SNESSolve. And recall that the formula for generalized
     405             :      *    inverse power iterations in this case is: (Ax)_n = (Bx)_{n-1}/|(Bx)_{n-1}| (in
     406             :      *    math terms)
     407             :      */
     408          35 :     PetscCall(EPSPowerUpdateFunctionB(eps,v,Bx));
     409          35 :     PetscCall(BVRestoreColumn(eps->V,0,&v));
     410          35 :     if (power->formNorm) PetscCall((*power->formNorm)(power->snes,Bx,&norm,power->formNormCtx));
     411          25 :     else PetscCall(VecNorm(Bx,NORM_2,&norm));
     412          35 :     PetscCall(FirstNonzeroIdx(Bx,&power->idx,&power->p));
     413          35 :     PetscCall(Normalize(Bx,norm,power->idx,power->p,power->sign_normalization,NULL));
     414             :   }
     415             : 
     416          49 :   PetscCall(STGetShift(eps->st,&sigma));    /* original shift */
     417          49 :   rho = sigma;
     418             : 
     419       28668 :   while (eps->reason == EPS_CONVERGED_ITERATING) {
     420       28619 :     eps->its++;
     421       28619 :     k = eps->nconv;
     422             : 
     423             :     /* y = OP v */
     424       28619 :     PetscCall(BVGetColumn(eps->V,k,&v));
     425       28619 :     if (power->nonlinear) PetscCall(EPSPowerApply_SNES(eps,v,y));
     426       28427 :     else PetscCall(STApply(eps->st,v,y));
     427       28619 :     PetscCall(BVRestoreColumn(eps->V,k,&v));
     428             : 
     429             :     /* purge previously converged eigenvectors */
     430       28619 :     if (PetscUnlikely(power->nonlinear)) {
     431             :       /* We do not need to call this for Newton eigensolver since eigenvalue is
     432             :        * updated in function evaluations.
     433             :        */
     434         192 :       if (!power->update) {
     435         176 :         PetscCall(EPSPowerUpdateFunctionB(eps,y,Bx));
     436         176 :         if (power->formNorm) PetscCall((*power->formNorm)(power->snes,Bx,&norm,power->formNormCtx));
     437         110 :         else PetscCall(VecNorm(Bx,NORM_2,&norm));
     438         176 :         PetscCall(Normalize(Bx,norm,power->idx,power->p,power->sign_normalization,&sign));
     439             :       }
     440             :     } else {
     441       28427 :       PetscCall(DSGetArray(eps->ds,DS_MAT_A,&T));
     442       28427 :       PetscCall(BVSetActiveColumns(eps->V,0,k));
     443       28427 :       PetscCall(BVOrthogonalizeVec(eps->V,y,T+k*ld,&norm,NULL));
     444             :     }
     445             : 
     446             :     /* theta = (v,y)_B */
     447       28619 :     PetscCall(BVSetActiveColumns(eps->V,k,k+1));
     448       28619 :     PetscCall(BVDotVec(eps->V,y,&theta));
     449       28619 :     if (!power->nonlinear) {
     450       28427 :       T[k+k*ld] = theta;
     451       28427 :       PetscCall(DSRestoreArray(eps->ds,DS_MAT_A,&T));
     452             :     }
     453             : 
     454             :     /* Eigenvalue is already stored in function evaluations.
     455             :      * Assign eigenvalue to theta to make the rest of the code consistent
     456             :      */
     457       28619 :     if (power->update) theta = eps->eigr[eps->nconv];
     458       28603 :     else if (power->nonlinear) theta = 1.0/(norm*sign); /* Eigenvalue: 1/|Bx| */
     459             : 
     460       28619 :     if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) { /* direct & inverse iteration */
     461             : 
     462             :       /* approximate eigenvalue is the Rayleigh quotient */
     463       28584 :       eps->eigr[eps->nconv] = theta;
     464             : 
     465             :       /**
     466             :        * If the Newton method (update, SNES) is used, we do not compute "relerr"
     467             :        * since SNES determines the convergence.
     468             :        */
     469       28584 :       if (PetscUnlikely(power->update)) relerr = 0.;
     470             :       else {
     471             :         /* compute relative error as ||y-theta v||_2/|theta| */
     472       28568 :         PetscCall(VecCopy(y,e));
     473       28568 :         PetscCall(BVGetColumn(eps->V,k,&v));
     474       28568 :         PetscCall(VecAXPY(e,power->nonlinear?-1.0:-theta,v));
     475       28568 :         PetscCall(BVRestoreColumn(eps->V,k,&v));
     476       28568 :         PetscCall(VecNorm(e,NORM_2,&relerr));
     477       28568 :         if (PetscUnlikely(power->nonlinear)) relerr *= PetscAbsScalar(theta);
     478       28392 :         else relerr /= PetscAbsScalar(theta);
     479             :       }
     480             : 
     481             :     } else {  /* RQI */
     482             : 
     483             :       /* delta = ||y||_B */
     484          35 :       delta = norm;
     485             : 
     486             :       /* compute relative error */
     487          35 :       if (rho == 0.0) relerr = PETSC_MAX_REAL;
     488          35 :       else relerr = 1.0 / (norm*PetscAbsScalar(rho));
     489             : 
     490             :       /* approximate eigenvalue is the shift */
     491          35 :       eps->eigr[eps->nconv] = rho;
     492             : 
     493             :       /* compute new shift */
     494          35 :       if (relerr<eps->tol) {
     495           8 :         rho = sigma;  /* if converged, restore original shift */
     496           8 :         PetscCall(STSetShift(eps->st,rho));
     497             :       } else {
     498          27 :         rho = rho + PetscConj(theta)/(delta*delta);  /* Rayleigh quotient R(v) */
     499          27 :         if (power->shift_type == EPS_POWER_SHIFT_WILKINSON) {
     500             :           /* beta1 is the norm of the residual associated with R(v) */
     501          10 :           PetscCall(BVGetColumn(eps->V,k,&v));
     502          10 :           PetscCall(VecAXPY(v,-PetscConj(theta)/(delta*delta),y));
     503          10 :           PetscCall(BVRestoreColumn(eps->V,k,&v));
     504          10 :           PetscCall(BVScaleColumn(eps->V,k,1.0/delta));
     505          10 :           PetscCall(BVNormColumn(eps->V,k,NORM_2,&norm1));
     506          10 :           beta1 = norm1;
     507             : 
     508             :           /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
     509          10 :           PetscCall(STGetMatrix(eps->st,0,&A));
     510          10 :           PetscCall(BVGetColumn(eps->V,k,&v));
     511          10 :           PetscCall(MatMult(A,v,e));
     512          10 :           PetscCall(VecDot(v,e,&alpha2));
     513          10 :           PetscCall(BVRestoreColumn(eps->V,k,&v));
     514          10 :           alpha2 = alpha2 / (beta1 * beta1);
     515             : 
     516             :           /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
     517          10 :           PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
     518          10 :           PetscCallBLAS("LAPACKlaev2",LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1));
     519          10 :           PetscCall(PetscFPTrapPop());
     520          10 :           if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
     521           5 :           else rho = rt2;
     522             :         }
     523             :         /* update operator according to new shift */
     524          27 :         PetscCall(KSPSetErrorIfNotConverged(ksp,PETSC_FALSE));
     525          27 :         PetscCall(STSetShift(eps->st,rho));
     526          27 :         PetscCall(KSPGetConvergedReason(ksp,&reason));
     527          27 :         if (reason) {
     528           0 :           PetscCall(PetscInfo(eps,"Factorization failed, repeat with a perturbed shift\n"));
     529           0 :           rho *= 1+10*PETSC_MACHINE_EPSILON;
     530           0 :           PetscCall(STSetShift(eps->st,rho));
     531           0 :           PetscCall(KSPGetConvergedReason(ksp,&reason));
     532           0 :           PetscCheck(!reason,PetscObjectComm((PetscObject)ksp),PETSC_ERR_CONV_FAILED,"Second factorization failed");
     533             :         }
     534          27 :         PetscCall(KSPSetErrorIfNotConverged(ksp,PETSC_TRUE));
     535             :       }
     536             :     }
     537       28619 :     eps->errest[eps->nconv] = relerr;
     538             : 
     539             :     /* normalize */
     540       28619 :     if (!power->nonlinear) PetscCall(Normalize(y,norm,power->idx,power->p,power->sign_normalization,NULL));
     541       28619 :     PetscCall(BVInsertVec(eps->V,k,y));
     542             : 
     543       28619 :     if (PetscUnlikely(power->update)) {
     544          16 :       PetscCall(SNESGetConvergedReason(power->snes,&snesreason));
     545             :       /* For Newton eigensolver, we are ready to return once SNES converged. */
     546          16 :       if (snesreason>0) eps->nconv = 1;
     547       28603 :     } else if (PetscUnlikely(relerr<eps->tol)) {   /* accept eigenpair */
     548          50 :       eps->nconv = eps->nconv + 1;
     549          50 :       if (eps->nconv<eps->nev) {
     550          31 :         PetscCall(EPSGetStartVector(eps,eps->nconv,&breakdown));
     551          31 :         if (breakdown) {
     552           0 :           eps->reason = EPS_DIVERGED_BREAKDOWN;
     553           0 :           PetscCall(PetscInfo(eps,"Unable to generate more start vectors\n"));
     554             :           break;
     555             :         }
     556             :       }
     557             :     }
     558             :     /* For Newton eigensolver, monitor will be called from SNES monitor */
     559       28619 :     if (!power->update) PetscCall(EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev)));
     560       28619 :     PetscCall((*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx));
     561             : 
     562             :     /**
     563             :      * When a customized stopping test is used, and reason can be set to be converged (EPS_CONVERGED_USER).
     564             :      * In this case, we need to increase eps->nconv to "1" so users can retrieve the solution.
     565             :      */
     566       28619 :     if (PetscUnlikely(power->nonlinear && eps->reason>0)) eps->nconv = 1;
     567             :   }
     568             : 
     569          49 :   if (power->nonlinear) PetscCall(PetscObjectCompose((PetscObject)power->snes,"eps",NULL));
     570             :   else {
     571          14 :     PetscCall(DSSetDimensions(eps->ds,eps->nconv,0,0));
     572          14 :     PetscCall(DSSetState(eps->ds,DS_STATE_RAW));
     573             :   }
     574          49 :   PetscFunctionReturn(PETSC_SUCCESS);
     575             : }
     576             : 
     577           5 : static PetscErrorCode EPSSolve_TS_Power(EPS eps)
     578             : {
     579           5 :   EPS_POWER          *power = (EPS_POWER*)eps->data;
     580           5 :   PetscInt           k,ld;
     581           5 :   Vec                v,w,y,e,z;
     582           5 :   KSP                ksp;
     583           5 :   PetscReal          relerr=1.0,relerrl,delta;
     584           5 :   PetscScalar        theta,rho,alpha,sigma;
     585           5 :   PetscBool          breakdown,breakdownl;
     586           5 :   KSPConvergedReason reason;
     587             : 
     588           5 :   PetscFunctionBegin;
     589           5 :   e = eps->work[0];
     590           5 :   v = eps->work[1];
     591           5 :   w = eps->work[2];
     592             : 
     593           5 :   if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) PetscCall(STGetKSP(eps->st,&ksp));
     594           5 :   PetscCall(DSGetLeadingDimension(eps->ds,&ld));
     595           5 :   PetscCall(EPSGetStartVector(eps,0,NULL));
     596           5 :   PetscCall(EPSGetLeftStartVector(eps,0,NULL));
     597           5 :   PetscCall(BVBiorthonormalizeColumn(eps->V,eps->W,0,NULL));
     598           5 :   PetscCall(BVCopyVec(eps->V,0,v));
     599           5 :   PetscCall(BVCopyVec(eps->W,0,w));
     600           5 :   PetscCall(STGetShift(eps->st,&sigma));    /* original shift */
     601           5 :   rho = sigma;
     602             : 
     603           5 :   while (eps->reason == EPS_CONVERGED_ITERATING) {
     604         421 :     eps->its++;
     605         421 :     k = eps->nconv;
     606             : 
     607             :     /* y = OP v, z = OP' w */
     608         421 :     PetscCall(BVGetColumn(eps->V,k,&y));
     609         421 :     PetscCall(STApply(eps->st,v,y));
     610         421 :     PetscCall(BVRestoreColumn(eps->V,k,&y));
     611         421 :     PetscCall(BVGetColumn(eps->W,k,&z));
     612         421 :     PetscCall(STApplyHermitianTranspose(eps->st,w,z));
     613         421 :     PetscCall(BVRestoreColumn(eps->W,k,&z));
     614             : 
     615             :     /* purge previously converged eigenvectors */
     616         421 :     PetscCall(BVBiorthogonalizeColumn(eps->V,eps->W,k));
     617             : 
     618             :     /* theta = (w,y)_B */
     619         421 :     PetscCall(BVSetActiveColumns(eps->V,k,k+1));
     620         421 :     PetscCall(BVDotVec(eps->V,w,&theta));
     621         421 :     theta = PetscConj(theta);
     622             : 
     623         421 :     if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) { /* direct & inverse iteration */
     624             : 
     625             :       /* approximate eigenvalue is the Rayleigh quotient */
     626         398 :       eps->eigr[eps->nconv] = theta;
     627             : 
     628             :       /* compute relative errors as ||y-theta v||_2/|theta| and ||z-conj(theta) w||_2/|theta|*/
     629         398 :       PetscCall(BVCopyVec(eps->V,k,e));
     630         398 :       PetscCall(VecAXPY(e,-theta,v));
     631         398 :       PetscCall(VecNorm(e,NORM_2,&relerr));
     632         398 :       PetscCall(BVCopyVec(eps->W,k,e));
     633         398 :       PetscCall(VecAXPY(e,-PetscConj(theta),w));
     634         398 :       PetscCall(VecNorm(e,NORM_2,&relerrl));
     635         772 :       relerr = PetscMax(relerr,relerrl)/PetscAbsScalar(theta);
     636             :     }
     637             : 
     638             :     /* normalize */
     639         421 :     PetscCall(BVSetActiveColumns(eps->V,k,k+1));
     640         421 :     PetscCall(BVGetColumn(eps->W,k,&z));
     641         421 :     PetscCall(BVDotVec(eps->V,z,&alpha));
     642         421 :     PetscCall(BVRestoreColumn(eps->W,k,&z));
     643         421 :     delta = PetscSqrtReal(PetscAbsScalar(alpha));
     644         421 :     PetscCheck(delta!=0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_CONV_FAILED,"Breakdown in two-sided Power/RQI");
     645         421 :     PetscCall(BVScaleColumn(eps->V,k,1.0/PetscConj(alpha/delta)));
     646         421 :     PetscCall(BVScaleColumn(eps->W,k,1.0/delta));
     647         421 :     PetscCall(BVCopyVec(eps->V,k,v));
     648         421 :     PetscCall(BVCopyVec(eps->W,k,w));
     649             : 
     650         421 :     if (power->shift_type == EPS_POWER_SHIFT_RAYLEIGH) { /* RQI */
     651             : 
     652             :       /* compute relative error */
     653          23 :       if (rho == 0.0) relerr = PETSC_MAX_REAL;
     654          23 :       else relerr = 1.0 / PetscAbsScalar(delta*rho);
     655             : 
     656             :       /* approximate eigenvalue is the shift */
     657          23 :       eps->eigr[eps->nconv] = rho;
     658             : 
     659             :       /* compute new shift */
     660          23 :       if (relerr<eps->tol) {
     661           5 :         rho = sigma;  /* if converged, restore original shift */
     662           5 :         PetscCall(STSetShift(eps->st,rho));
     663             :       } else {
     664          18 :         rho = rho + PetscConj(theta)/(delta*delta);  /* Rayleigh quotient R(v) */
     665             :         /* update operator according to new shift */
     666          18 :         PetscCall(KSPSetErrorIfNotConverged(ksp,PETSC_FALSE));
     667          18 :         PetscCall(STSetShift(eps->st,rho));
     668          18 :         PetscCall(KSPGetConvergedReason(ksp,&reason));
     669          18 :         if (reason) {
     670           1 :           PetscCall(PetscInfo(eps,"Factorization failed, repeat with a perturbed shift\n"));
     671           1 :           rho *= 1+10*PETSC_MACHINE_EPSILON;
     672           1 :           PetscCall(STSetShift(eps->st,rho));
     673           1 :           PetscCall(KSPGetConvergedReason(ksp,&reason));
     674           1 :           PetscCheck(!reason,PetscObjectComm((PetscObject)ksp),PETSC_ERR_CONV_FAILED,"Second factorization failed");
     675             :         }
     676          18 :         PetscCall(KSPSetErrorIfNotConverged(ksp,PETSC_TRUE));
     677             :       }
     678             :     }
     679         421 :     eps->errest[eps->nconv] = relerr;
     680             : 
     681             :     /* if relerr<tol, accept eigenpair */
     682         421 :     if (relerr<eps->tol) {
     683          11 :       eps->nconv = eps->nconv + 1;
     684          11 :       if (eps->nconv<eps->nev) {
     685           6 :         PetscCall(EPSGetStartVector(eps,eps->nconv,&breakdown));
     686           6 :         PetscCall(EPSGetLeftStartVector(eps,eps->nconv,&breakdownl));
     687           6 :         if (breakdown || breakdownl) {
     688           0 :           eps->reason = EPS_DIVERGED_BREAKDOWN;
     689           0 :           PetscCall(PetscInfo(eps,"Unable to generate more start vectors\n"));
     690             :           break;
     691             :         }
     692           6 :         PetscCall(BVBiorthonormalizeColumn(eps->V,eps->W,eps->nconv,NULL));
     693             :       }
     694             :     }
     695         421 :     PetscCall(EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev)));
     696         426 :     PetscCall((*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx));
     697             :   }
     698             : 
     699           5 :   PetscCall(DSSetDimensions(eps->ds,eps->nconv,0,0));
     700           5 :   PetscCall(DSSetState(eps->ds,DS_STATE_RAW));
     701           5 :   PetscFunctionReturn(PETSC_SUCCESS);
     702             : }
     703             : 
     704       29040 : static PetscErrorCode EPSStopping_Power(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ctx)
     705             : {
     706       29040 :   EPS_POWER      *power = (EPS_POWER*)eps->data;
     707       29040 :   SNESConvergedReason snesreason;
     708             : 
     709       29040 :   PetscFunctionBegin;
     710       29040 :   if (PetscUnlikely(power->update)) {
     711          16 :     PetscCall(SNESGetConvergedReason(power->snes,&snesreason));
     712          16 :     if (snesreason < 0) {
     713           0 :       *reason = EPS_DIVERGED_BREAKDOWN;
     714           0 :       PetscFunctionReturn(PETSC_SUCCESS);
     715             :     }
     716             :   }
     717       29040 :   PetscCall(EPSStoppingBasic(eps,its,max_it,nconv,nev,reason,ctx));
     718       29040 :   PetscFunctionReturn(PETSC_SUCCESS);
     719             : }
     720             : 
     721          19 : static PetscErrorCode EPSBackTransform_Power(EPS eps)
     722             : {
     723          19 :   EPS_POWER      *power = (EPS_POWER*)eps->data;
     724             : 
     725          19 :   PetscFunctionBegin;
     726          19 :   if (power->nonlinear) eps->eigr[0] = 1.0/eps->eigr[0];
     727          19 :   else if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) PetscCall(EPSBackTransform_Default(eps));
     728          19 :   PetscFunctionReturn(PETSC_SUCCESS);
     729             : }
     730             : 
     731          50 : static PetscErrorCode EPSSetFromOptions_Power(EPS eps,PetscOptionItems *PetscOptionsObject)
     732             : {
     733          50 :   EPS_POWER         *power = (EPS_POWER*)eps->data;
     734          50 :   PetscBool         flg,val;
     735          50 :   EPSPowerShiftType shift;
     736             : 
     737          50 :   PetscFunctionBegin;
     738          50 :   PetscOptionsHeadBegin(PetscOptionsObject,"EPS Power Options");
     739             : 
     740          50 :     PetscCall(PetscOptionsEnum("-eps_power_shift_type","Shift type","EPSPowerSetShiftType",EPSPowerShiftTypes,(PetscEnum)power->shift_type,(PetscEnum*)&shift,&flg));
     741          50 :     if (flg) PetscCall(EPSPowerSetShiftType(eps,shift));
     742             : 
     743          50 :     PetscCall(PetscOptionsBool("-eps_power_nonlinear","Use nonlinear inverse iteration","EPSPowerSetNonlinear",power->nonlinear,&val,&flg));
     744          50 :     if (flg) PetscCall(EPSPowerSetNonlinear(eps,val));
     745             : 
     746          50 :     PetscCall(PetscOptionsBool("-eps_power_update","Update residual monolithically","EPSPowerSetUpdate",power->update,&val,&flg));
     747          50 :     if (flg) PetscCall(EPSPowerSetUpdate(eps,val));
     748             : 
     749          50 :     PetscCall(PetscOptionsBool("-eps_power_sign_normalization","Normalize Bx with sign of first nonzero entry","EPSPowerSetSignNormalization",power->sign_normalization,&power->sign_normalization,&flg));
     750             : 
     751          50 :   PetscOptionsHeadEnd();
     752          50 :   PetscFunctionReturn(PETSC_SUCCESS);
     753             : }
     754             : 
     755           9 : static PetscErrorCode EPSPowerSetShiftType_Power(EPS eps,EPSPowerShiftType shift)
     756             : {
     757           9 :   EPS_POWER *power = (EPS_POWER*)eps->data;
     758             : 
     759           9 :   PetscFunctionBegin;
     760           9 :   switch (shift) {
     761           9 :     case EPS_POWER_SHIFT_CONSTANT:
     762             :     case EPS_POWER_SHIFT_RAYLEIGH:
     763             :     case EPS_POWER_SHIFT_WILKINSON:
     764           9 :       if (power->shift_type != shift) {
     765           6 :         power->shift_type = shift;
     766           6 :         eps->state = EPS_STATE_INITIAL;
     767             :       }
     768           9 :       break;
     769           0 :     default:
     770           0 :       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid shift type");
     771             :   }
     772           9 :   PetscFunctionReturn(PETSC_SUCCESS);
     773             : }
     774             : 
     775             : /*@
     776             :    EPSPowerSetShiftType - Sets the type of shifts used during the power
     777             :    iteration. This can be used to emulate the Rayleigh Quotient Iteration
     778             :    (RQI) method.
     779             : 
     780             :    Logically Collective
     781             : 
     782             :    Input Parameters:
     783             : +  eps - the eigenproblem solver context
     784             : -  shift - the type of shift
     785             : 
     786             :    Options Database Key:
     787             : .  -eps_power_shift_type - Sets the shift type (either 'constant' or
     788             :                            'rayleigh' or 'wilkinson')
     789             : 
     790             :    Notes:
     791             :    By default, shifts are constant (EPS_POWER_SHIFT_CONSTANT) and the iteration
     792             :    is the simple power method (or inverse iteration if a shift-and-invert
     793             :    transformation is being used).
     794             : 
     795             :    A variable shift can be specified (EPS_POWER_SHIFT_RAYLEIGH or
     796             :    EPS_POWER_SHIFT_WILKINSON). In this case, the iteration behaves rather like
     797             :    a cubic converging method such as RQI.
     798             : 
     799             :    Level: advanced
     800             : 
     801             : .seealso: EPSPowerGetShiftType(), STSetShift(), EPSPowerShiftType
     802             : @*/
     803           9 : PetscErrorCode EPSPowerSetShiftType(EPS eps,EPSPowerShiftType shift)
     804             : {
     805           9 :   PetscFunctionBegin;
     806           9 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     807          27 :   PetscValidLogicalCollectiveEnum(eps,shift,2);
     808           9 :   PetscTryMethod(eps,"EPSPowerSetShiftType_C",(EPS,EPSPowerShiftType),(eps,shift));
     809           9 :   PetscFunctionReturn(PETSC_SUCCESS);
     810             : }
     811             : 
     812           1 : static PetscErrorCode EPSPowerGetShiftType_Power(EPS eps,EPSPowerShiftType *shift)
     813             : {
     814           1 :   EPS_POWER *power = (EPS_POWER*)eps->data;
     815             : 
     816           1 :   PetscFunctionBegin;
     817           1 :   *shift = power->shift_type;
     818           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     819             : }
     820             : 
     821             : /*@
     822             :    EPSPowerGetShiftType - Gets the type of shifts used during the power
     823             :    iteration.
     824             : 
     825             :    Not Collective
     826             : 
     827             :    Input Parameter:
     828             : .  eps - the eigenproblem solver context
     829             : 
     830             :    Output Parameter:
     831             : .  shift - the type of shift
     832             : 
     833             :    Level: advanced
     834             : 
     835             : .seealso: EPSPowerSetShiftType(), EPSPowerShiftType
     836             : @*/
     837           1 : PetscErrorCode EPSPowerGetShiftType(EPS eps,EPSPowerShiftType *shift)
     838             : {
     839           1 :   PetscFunctionBegin;
     840           1 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     841           1 :   PetscAssertPointer(shift,2);
     842           1 :   PetscUseMethod(eps,"EPSPowerGetShiftType_C",(EPS,EPSPowerShiftType*),(eps,shift));
     843           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     844             : }
     845             : 
     846          33 : static PetscErrorCode EPSPowerSetNonlinear_Power(EPS eps,PetscBool nonlinear)
     847             : {
     848          33 :   EPS_POWER *power = (EPS_POWER*)eps->data;
     849             : 
     850          33 :   PetscFunctionBegin;
     851          33 :   if (power->nonlinear != nonlinear) {
     852          33 :     power->nonlinear = nonlinear;
     853          33 :     eps->useds = PetscNot(nonlinear);
     854          33 :     eps->ops->setupsort = nonlinear? NULL: EPSSetUpSort_Default;
     855          33 :     eps->state = EPS_STATE_INITIAL;
     856             :   }
     857          33 :   PetscFunctionReturn(PETSC_SUCCESS);
     858             : }
     859             : 
     860             : /*@
     861             :    EPSPowerSetNonlinear - Sets a flag to indicate that the problem is nonlinear.
     862             : 
     863             :    Logically Collective
     864             : 
     865             :    Input Parameters:
     866             : +  eps - the eigenproblem solver context
     867             : -  nonlinear - whether the problem is nonlinear or not
     868             : 
     869             :    Options Database Key:
     870             : .  -eps_power_nonlinear - Sets the nonlinear flag
     871             : 
     872             :    Notes:
     873             :    If this flag is set, the solver assumes that the problem is nonlinear,
     874             :    that is, the operators that define the eigenproblem are not constant
     875             :    matrices, but depend on the eigenvector, A(x)*x=lambda*B(x)*x. This is
     876             :    different from the case of nonlinearity with respect to the eigenvalue
     877             :    (use the NEP solver class for this kind of problems).
     878             : 
     879             :    The way in which nonlinear operators are specified is very similar to
     880             :    the case of PETSc's SNES solver. The difference is that the callback
     881             :    functions are provided via composed functions "formFunction" and
     882             :    "formJacobian" in each of the matrix objects passed as arguments of
     883             :    EPSSetOperators(). The application context required for these functions
     884             :    can be attached via a composed PetscContainer.
     885             : 
     886             :    Level: advanced
     887             : 
     888             : .seealso: EPSPowerGetNonlinear(), EPSSetOperators()
     889             : @*/
     890          33 : PetscErrorCode EPSPowerSetNonlinear(EPS eps,PetscBool nonlinear)
     891             : {
     892          33 :   PetscFunctionBegin;
     893          33 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     894          99 :   PetscValidLogicalCollectiveBool(eps,nonlinear,2);
     895          33 :   PetscTryMethod(eps,"EPSPowerSetNonlinear_C",(EPS,PetscBool),(eps,nonlinear));
     896          33 :   PetscFunctionReturn(PETSC_SUCCESS);
     897             : }
     898             : 
     899          19 : static PetscErrorCode EPSPowerGetNonlinear_Power(EPS eps,PetscBool *nonlinear)
     900             : {
     901          19 :   EPS_POWER *power = (EPS_POWER*)eps->data;
     902             : 
     903          19 :   PetscFunctionBegin;
     904          19 :   *nonlinear = power->nonlinear;
     905          19 :   PetscFunctionReturn(PETSC_SUCCESS);
     906             : }
     907             : 
     908             : /*@
     909             :    EPSPowerGetNonlinear - Returns a flag indicating if the problem is nonlinear.
     910             : 
     911             :    Not Collective
     912             : 
     913             :    Input Parameter:
     914             : .  eps - the eigenproblem solver context
     915             : 
     916             :    Output Parameter:
     917             : .  nonlinear - the nonlinear flag
     918             : 
     919             :    Level: advanced
     920             : 
     921             : .seealso: EPSPowerSetUpdate(), EPSPowerSetNonlinear()
     922             : @*/
     923          19 : PetscErrorCode EPSPowerGetNonlinear(EPS eps,PetscBool *nonlinear)
     924             : {
     925          19 :   PetscFunctionBegin;
     926          19 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     927          19 :   PetscAssertPointer(nonlinear,2);
     928          19 :   PetscUseMethod(eps,"EPSPowerGetNonlinear_C",(EPS,PetscBool*),(eps,nonlinear));
     929          19 :   PetscFunctionReturn(PETSC_SUCCESS);
     930             : }
     931             : 
     932          14 : static PetscErrorCode EPSPowerSetUpdate_Power(EPS eps,PetscBool update)
     933             : {
     934          14 :   EPS_POWER *power = (EPS_POWER*)eps->data;
     935             : 
     936          14 :   PetscFunctionBegin;
     937          14 :   PetscCheck(power->nonlinear,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"This option does not make sense for linear problems");
     938          14 :   power->update = update;
     939          14 :   eps->state = EPS_STATE_INITIAL;
     940          14 :   PetscFunctionReturn(PETSC_SUCCESS);
     941             : }
     942             : 
     943             : /*@
     944             :    EPSPowerSetUpdate - Sets a flag to indicate that the residual is updated monolithically
     945             :    for nonlinear problems. This potentially has a better convergence.
     946             : 
     947             :    Logically Collective
     948             : 
     949             :    Input Parameters:
     950             : +  eps - the eigenproblem solver context
     951             : -  update - whether the residual is updated monolithically or not
     952             : 
     953             :    Options Database Key:
     954             : .  -eps_power_update - Sets the update flag
     955             : 
     956             :    Level: advanced
     957             : 
     958             : .seealso: EPSPowerGetUpdate(), EPSPowerGetNonlinear(), EPSSetOperators()
     959             : @*/
     960          14 : PetscErrorCode EPSPowerSetUpdate(EPS eps,PetscBool update)
     961             : {
     962          14 :   PetscFunctionBegin;
     963          14 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     964          42 :   PetscValidLogicalCollectiveBool(eps,update,2);
     965          14 :   PetscTryMethod(eps,"EPSPowerSetUpdate_C",(EPS,PetscBool),(eps,update));
     966          14 :   PetscFunctionReturn(PETSC_SUCCESS);
     967             : }
     968             : 
     969          19 : static PetscErrorCode EPSPowerGetUpdate_Power(EPS eps,PetscBool *update)
     970             : {
     971          19 :   EPS_POWER *power = (EPS_POWER*)eps->data;
     972             : 
     973          19 :   PetscFunctionBegin;
     974          19 :   *update = power->update;
     975          19 :   PetscFunctionReturn(PETSC_SUCCESS);
     976             : }
     977             : 
     978             : /*@
     979             :    EPSPowerGetUpdate - Returns a flag indicating if the residual is updated monolithically
     980             :    for nonlinear problems.
     981             : 
     982             :    Not Collective
     983             : 
     984             :    Input Parameter:
     985             : .  eps - the eigenproblem solver context
     986             : 
     987             :    Output Parameter:
     988             : .  update - the update flag
     989             : 
     990             :    Level: advanced
     991             : 
     992             : .seealso: EPSPowerSetUpdate(), EPSPowerSetNonlinear()
     993             : @*/
     994          19 : PetscErrorCode EPSPowerGetUpdate(EPS eps,PetscBool *update)
     995             : {
     996          19 :   PetscFunctionBegin;
     997          19 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     998          19 :   PetscAssertPointer(update,2);
     999          19 :   PetscUseMethod(eps,"EPSPowerGetUpdate_C",(EPS,PetscBool*),(eps,update));
    1000          19 :   PetscFunctionReturn(PETSC_SUCCESS);
    1001             : }
    1002             : 
    1003          33 : static PetscErrorCode EPSPowerSetSignNormalization_Power(EPS eps,PetscBool sign_normalization)
    1004             : {
    1005          33 :   EPS_POWER *power = (EPS_POWER*)eps->data;
    1006             : 
    1007          33 :   PetscFunctionBegin;
    1008          33 :   power->sign_normalization = sign_normalization;
    1009          33 :   PetscFunctionReturn(PETSC_SUCCESS);
    1010             : }
    1011             : 
    1012             : /*@
    1013             :    EPSPowerSetSignNormalization - Sets a flag to indicate whether the Bx vector
    1014             :    should be normalized by the sign of the first non-zero element in the vector.
    1015             :    E.g., if this is true, the post-normalization value of the first non-zero element
    1016             :    in the vector is guaranteed to be positive.
    1017             : 
    1018             :    Logically Collective
    1019             : 
    1020             :    Input Parameters:
    1021             : +  eps - the eigenproblem solver context
    1022             : -  sign_normalization - whether Bx should be multiplied by the sign of the first non-zero
    1023             :                         element when performing normalization steps
    1024             : 
    1025             :    Options Database Key:
    1026             : .  -eps_power_sign_normalization - Sets the sign normalization flag
    1027             : 
    1028             :    Level: advanced
    1029             : 
    1030             : .seealso: EPSPowerGetSignNormalization()
    1031             : @*/
    1032          33 : PetscErrorCode EPSPowerSetSignNormalization(EPS eps,PetscBool sign_normalization)
    1033             : {
    1034          33 :   PetscFunctionBegin;
    1035          33 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1036          99 :   PetscValidLogicalCollectiveBool(eps,sign_normalization,2);
    1037          33 :   PetscTryMethod(eps,"EPSPowerSetSignNormalization_C",(EPS,PetscBool),(eps,sign_normalization));
    1038          33 :   PetscFunctionReturn(PETSC_SUCCESS);
    1039             : }
    1040             : 
    1041          14 : static PetscErrorCode EPSPowerGetSignNormalization_Power(EPS eps,PetscBool *sign_normalization)
    1042             : {
    1043          14 :   EPS_POWER *power = (EPS_POWER*)eps->data;
    1044             : 
    1045          14 :   PetscFunctionBegin;
    1046          14 :   *sign_normalization = power->sign_normalization;
    1047          14 :   PetscFunctionReturn(PETSC_SUCCESS);
    1048             : }
    1049             : 
    1050             : /*@
    1051             :    EPSPowerGetSignNormalization - Returns a flag indicating whether the Bx vector
    1052             :    is normalized by the sign of the first non-zero element in the vector. E.g.,
    1053             :    if this is true, the post-normalization value of the first non-zero element in
    1054             :    the vector is guaranteed to be positive.
    1055             : 
    1056             :    Not Collective
    1057             : 
    1058             :    Input Parameter:
    1059             : .  eps - the eigenproblem solver context
    1060             : 
    1061             :    Output Parameter:
    1062             : .  sign_normalization - the sign normalization flag
    1063             : 
    1064             :    Level: advanced
    1065             : 
    1066             : .seealso: EPSPowerSetSignNormalization()
    1067             : @*/
    1068          14 : PetscErrorCode EPSPowerGetSignNormalization(EPS eps,PetscBool *sign_normalization)
    1069             : {
    1070          14 :   PetscFunctionBegin;
    1071          14 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1072          14 :   PetscAssertPointer(sign_normalization,2);
    1073          14 :   PetscUseMethod(eps,"EPSPowerGetSignNormalization_C",(EPS,PetscBool*),(eps,sign_normalization));
    1074          14 :   PetscFunctionReturn(PETSC_SUCCESS);
    1075             : }
    1076             : 
    1077           0 : static PetscErrorCode EPSPowerSetSNES_Power(EPS eps,SNES snes)
    1078             : {
    1079           0 :   EPS_POWER      *power = (EPS_POWER*)eps->data;
    1080             : 
    1081           0 :   PetscFunctionBegin;
    1082           0 :   PetscCall(PetscObjectReference((PetscObject)snes));
    1083           0 :   PetscCall(SNESDestroy(&power->snes));
    1084           0 :   power->snes = snes;
    1085           0 :   eps->state = EPS_STATE_INITIAL;
    1086           0 :   PetscFunctionReturn(PETSC_SUCCESS);
    1087             : }
    1088             : 
    1089             : /*@
    1090             :    EPSPowerSetSNES - Associate a nonlinear solver object (SNES) to the
    1091             :    eigenvalue solver (to be used in nonlinear inverse iteration).
    1092             : 
    1093             :    Collective
    1094             : 
    1095             :    Input Parameters:
    1096             : +  eps  - the eigenvalue solver
    1097             : -  snes - the nonlinear solver object
    1098             : 
    1099             :    Level: advanced
    1100             : 
    1101             : .seealso: EPSPowerGetSNES()
    1102             : @*/
    1103           0 : PetscErrorCode EPSPowerSetSNES(EPS eps,SNES snes)
    1104             : {
    1105           0 :   PetscFunctionBegin;
    1106           0 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1107           0 :   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
    1108           0 :   PetscCheckSameComm(eps,1,snes,2);
    1109           0 :   PetscTryMethod(eps,"EPSPowerSetSNES_C",(EPS,SNES),(eps,snes));
    1110           0 :   PetscFunctionReturn(PETSC_SUCCESS);
    1111             : }
    1112             : 
    1113          47 : static PetscErrorCode EPSPowerGetSNES_Power(EPS eps,SNES *snes)
    1114             : {
    1115          47 :   EPS_POWER      *power = (EPS_POWER*)eps->data;
    1116             : 
    1117          47 :   PetscFunctionBegin;
    1118          47 :   if (!power->snes) {
    1119          33 :     PetscCall(SNESCreate(PetscObjectComm((PetscObject)eps),&power->snes));
    1120          33 :     PetscCall(PetscObjectIncrementTabLevel((PetscObject)power->snes,(PetscObject)eps,1));
    1121          33 :     PetscCall(SNESSetOptionsPrefix(power->snes,((PetscObject)eps)->prefix));
    1122          33 :     PetscCall(SNESAppendOptionsPrefix(power->snes,"eps_power_"));
    1123          33 :     PetscCall(PetscObjectSetOptions((PetscObject)power->snes,((PetscObject)eps)->options));
    1124             :   }
    1125          47 :   *snes = power->snes;
    1126          47 :   PetscFunctionReturn(PETSC_SUCCESS);
    1127             : }
    1128             : 
    1129             : /*@
    1130             :    EPSPowerGetSNES - Retrieve the nonlinear solver object (SNES) associated
    1131             :    with the eigenvalue solver.
    1132             : 
    1133             :    Not Collective
    1134             : 
    1135             :    Input Parameter:
    1136             : .  eps - the eigenvalue solver
    1137             : 
    1138             :    Output Parameter:
    1139             : .  snes - the nonlinear solver object
    1140             : 
    1141             :    Level: advanced
    1142             : 
    1143             : .seealso: EPSPowerSetSNES()
    1144             : @*/
    1145          47 : PetscErrorCode EPSPowerGetSNES(EPS eps,SNES *snes)
    1146             : {
    1147          47 :   PetscFunctionBegin;
    1148          47 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
    1149          47 :   PetscAssertPointer(snes,2);
    1150          47 :   PetscUseMethod(eps,"EPSPowerGetSNES_C",(EPS,SNES*),(eps,snes));
    1151          47 :   PetscFunctionReturn(PETSC_SUCCESS);
    1152             : }
    1153             : 
    1154          51 : static PetscErrorCode EPSReset_Power(EPS eps)
    1155             : {
    1156          51 :   EPS_POWER      *power = (EPS_POWER*)eps->data;
    1157             : 
    1158          51 :   PetscFunctionBegin;
    1159          51 :   if (power->snes) PetscCall(SNESReset(power->snes));
    1160          51 :   PetscFunctionReturn(PETSC_SUCCESS);
    1161             : }
    1162             : 
    1163          50 : static PetscErrorCode EPSDestroy_Power(EPS eps)
    1164             : {
    1165          50 :   EPS_POWER      *power = (EPS_POWER*)eps->data;
    1166             : 
    1167          50 :   PetscFunctionBegin;
    1168          50 :   if (power->nonlinear) PetscCall(SNESDestroy(&power->snes));
    1169          50 :   PetscCall(PetscFree(eps->data));
    1170          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",NULL));
    1171          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",NULL));
    1172          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetNonlinear_C",NULL));
    1173          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetNonlinear_C",NULL));
    1174          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetUpdate_C",NULL));
    1175          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetUpdate_C",NULL));
    1176          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSignNormalization_C",NULL));
    1177          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSignNormalization_C",NULL));
    1178          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSNES_C",NULL));
    1179          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSNES_C",NULL));
    1180          50 :   PetscFunctionReturn(PETSC_SUCCESS);
    1181             : }
    1182             : 
    1183           0 : static PetscErrorCode EPSView_Power(EPS eps,PetscViewer viewer)
    1184             : {
    1185           0 :   EPS_POWER      *power = (EPS_POWER*)eps->data;
    1186           0 :   PetscBool      isascii;
    1187             : 
    1188           0 :   PetscFunctionBegin;
    1189           0 :   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
    1190           0 :   if (isascii) {
    1191           0 :     if (power->nonlinear) {
    1192           0 :       if (power->sign_normalization) PetscCall(PetscViewerASCIIPrintf(viewer,"  normalizing Bx by the sign of the first nonzero element\n"));
    1193           0 :       else PetscCall(PetscViewerASCIIPrintf(viewer,"  not normalizing Bx by the sign of the first nonzero element\n"));
    1194           0 :       PetscCall(PetscViewerASCIIPrintf(viewer,"  using nonlinear inverse iteration\n"));
    1195           0 :       if (power->update) PetscCall(PetscViewerASCIIPrintf(viewer,"  updating the residual monolithically\n"));
    1196           0 :       if (!power->snes) PetscCall(EPSPowerGetSNES(eps,&power->snes));
    1197           0 :       PetscCall(PetscViewerASCIIPushTab(viewer));
    1198           0 :       PetscCall(SNESView(power->snes,viewer));
    1199           0 :       PetscCall(PetscViewerASCIIPopTab(viewer));
    1200           0 :     } else PetscCall(PetscViewerASCIIPrintf(viewer,"  %s shifts\n",EPSPowerShiftTypes[power->shift_type]));
    1201             :   }
    1202           0 :   PetscFunctionReturn(PETSC_SUCCESS);
    1203             : }
    1204             : 
    1205          40 : static PetscErrorCode EPSComputeVectors_Power(EPS eps)
    1206             : {
    1207          40 :   EPS_POWER      *power = (EPS_POWER*)eps->data;
    1208             : 
    1209          40 :   PetscFunctionBegin;
    1210          40 :   if (eps->twosided) {
    1211           5 :     PetscCall(EPSComputeVectors_Twosided(eps));
    1212           5 :     PetscCall(BVNormalize(eps->V,NULL));
    1213           5 :     PetscCall(BVNormalize(eps->W,NULL));
    1214          35 :   } else if (!power->nonlinear) PetscCall(EPSComputeVectors_Schur(eps));
    1215          40 :   PetscFunctionReturn(PETSC_SUCCESS);
    1216             : }
    1217             : 
    1218         104 : static PetscErrorCode EPSSetDefaultST_Power(EPS eps)
    1219             : {
    1220         104 :   EPS_POWER      *power = (EPS_POWER*)eps->data;
    1221         104 :   KSP            ksp;
    1222         104 :   PC             pc;
    1223             : 
    1224         104 :   PetscFunctionBegin;
    1225         104 :   if (power->nonlinear) {
    1226          68 :     eps->categ=EPS_CATEGORY_PRECOND;
    1227          68 :     PetscCall(STGetKSP(eps->st,&ksp));
    1228             :     /* Set ST as STPRECOND so it can carry one preconditioning matrix
    1229             :      * It is useful when A and B are shell matrices
    1230             :      */
    1231          68 :     PetscCall(STSetType(eps->st,STPRECOND));
    1232          68 :     PetscCall(KSPGetPC(ksp,&pc));
    1233          68 :     PetscCall(PCSetType(pc,PCNONE));
    1234             :   }
    1235         104 :   PetscFunctionReturn(PETSC_SUCCESS);
    1236             : }
    1237             : 
    1238          50 : SLEPC_EXTERN PetscErrorCode EPSCreate_Power(EPS eps)
    1239             : {
    1240          50 :   EPS_POWER      *ctx;
    1241             : 
    1242          50 :   PetscFunctionBegin;
    1243          50 :   PetscCall(PetscNew(&ctx));
    1244          50 :   eps->data = (void*)ctx;
    1245             : 
    1246          50 :   eps->useds = PETSC_TRUE;
    1247          50 :   eps->categ = EPS_CATEGORY_OTHER;
    1248             : 
    1249          50 :   eps->ops->setup          = EPSSetUp_Power;
    1250          50 :   eps->ops->setupsort      = EPSSetUpSort_Default;
    1251          50 :   eps->ops->setfromoptions = EPSSetFromOptions_Power;
    1252          50 :   eps->ops->reset          = EPSReset_Power;
    1253          50 :   eps->ops->destroy        = EPSDestroy_Power;
    1254          50 :   eps->ops->view           = EPSView_Power;
    1255          50 :   eps->ops->backtransform  = EPSBackTransform_Power;
    1256          50 :   eps->ops->computevectors = EPSComputeVectors_Power;
    1257          50 :   eps->ops->setdefaultst   = EPSSetDefaultST_Power;
    1258          50 :   eps->stopping            = EPSStopping_Power;
    1259          50 :   ctx->sign_normalization  = PETSC_TRUE;
    1260             : 
    1261          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",EPSPowerSetShiftType_Power));
    1262          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",EPSPowerGetShiftType_Power));
    1263          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetNonlinear_C",EPSPowerSetNonlinear_Power));
    1264          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetNonlinear_C",EPSPowerGetNonlinear_Power));
    1265          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetUpdate_C",EPSPowerSetUpdate_Power));
    1266          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetUpdate_C",EPSPowerGetUpdate_Power));
    1267          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSignNormalization_C",EPSPowerSetSignNormalization_Power));
    1268          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSignNormalization_C",EPSPowerGetSignNormalization_Power));
    1269          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSNES_C",EPSPowerSetSNES_Power));
    1270          50 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSNES_C",EPSPowerGetSNES_Power));
    1271          50 :   PetscFunctionReturn(PETSC_SUCCESS);
    1272             : }

Generated by: LCOV version 1.14