MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::KINSolver Class Reference

Interface to the KINSOL library – nonlinear solver methods. More...

#include <sundials.hpp>

Inheritance diagram for mfem::KINSolver:
[legend]
Collaboration diagram for mfem::KINSolver:
[legend]

Public Member Functions

 KINSolver (int strategy, bool oper_grad=true)
 Construct a serial wrapper to SUNDIALS' KINSOL nonlinear solver.
 
 KINSolver (MPI_Comm comm, int strategy, bool oper_grad=true)
 Construct a parallel wrapper to SUNDIALS' KINSOL nonlinear solver.
 
virtual ~KINSolver ()
 Destroy the associated KINSOL memory.
 
void SetOperator (const Operator &op) override
 Set the nonlinear Operator of the system and initialize KINSOL.
 
void SetSolver (Solver &solver) override
 Set the linear solver for inverting the Jacobian.
 
void SetPreconditioner (Solver &solver) override
 Equivalent to SetSolver(solver).
 
void SetScaledStepTol (double sstol)
 Set KINSOL's scaled step tolerance.
 
void SetMaxSetupCalls (int max_calls)
 Set maximum number of nonlinear iterations without a Jacobian update.
 
void EnableAndersonAcc (int n, int orth=KIN_ORTH_MGS, int delay=0, double damping=1.0)
 Enable Anderson Acceleration for KIN_FP or KIN_PICARD.
 
void SetDamping (double damping)
 
void SetJFNK (bool use_jfnk)
 Set the Jacobian Free Newton Krylov flag. The default is false.
 
void SetLSMaxIter (int m)
 Set the maximum number of linear solver iterations.
 
void SetLSMaxRestarts (int m)
 Set the maximum number of linear solver restarts.
 
void SetPrintLevel (int print_lvl) override
 Set the print level for the KINSetPrintLevel function.
 
void SetPrintLevel (PrintLevel) override
 This method is not supported and will throw an error.
 
void Mult (const Vector &b, Vector &x) const override
 Solve the nonlinear system \( F(x) = 0 \).
 
void Mult (Vector &x, const Vector &x_scale, const Vector &fx_scale) const
 Solve the nonlinear system \( F(x) = 0 \).
 
- Public Member Functions inherited from mfem::NewtonSolver
 NewtonSolver ()
 
 NewtonSolver (MPI_Comm comm_)
 
void SetOperator (const Operator &op) override
 Also calls SetOperator for the preconditioner if there is one.
 
void Mult (const Vector &b, Vector &x) const override
 Solve the nonlinear system with right-hand side b.
 
virtual real_t ComputeScalingFactor (const Vector &x, const Vector &b) const
 This method can be overloaded in derived classes to implement line search algorithms.
 
virtual void ProcessNewState (const Vector &x) const
 This method can be overloaded in derived classes to perform computations that need knowledge of the newest Newton state.
 
void SetAdaptiveLinRtol (const int type=2, const real_t rtol0=0.5, const real_t rtol_max=0.9, const real_t alpha=0.5 *(1.0+sqrt(5.0)), const real_t gamma=1.0)
 Enable adaptive linear solver relative tolerance algorithm.
 
- Public Member Functions inherited from mfem::IterativeSolver
 IterativeSolver ()
 
 IterativeSolver (MPI_Comm comm_)
 
void SetMonitor (IterativeSolverMonitor &m)
 Set the iterative solver monitor.
 
MPI_Comm GetComm () const
 Return the associated MPI communicator, or MPI_COMM_NULL if no communicator is set.
 
void SetRelTol (real_t rtol)
 
void SetAbsTol (real_t atol)
 
void SetMaxIter (int max_it)
 
int GetNumIterations () const
 Returns the number of iterations taken during the last call to Mult()
 
bool GetConverged () const
 Returns true if the last call to Mult() converged successfully.
 
real_t GetInitialNorm () const
 Returns the initial residual norm from the last call to Mult().
 
real_t GetFinalNorm () const
 Returns the final residual norm after termination of the solver during the last call to Mult().
 
real_t GetFinalRelNorm () const
 Returns the final residual norm after termination of the solver during the last call to Mult(), divided by the initial residual norm. Returns -1 if one of these norms is left undefined by the solver.
 
- Public Member Functions inherited from mfem::Solver
 Solver (int s=0, bool iter_mode=false)
 Initialize a square Solver with size s.
 
 Solver (int h, int w, bool iter_mode=false)
 Initialize a Solver with height h and width w.
 
- Public Member Functions inherited from mfem::Operator
void InitTVectors (const Operator *Po, const Operator *Ri, const Operator *Pi, Vector &x, Vector &b, Vector &X, Vector &B) const
 Initializes memory for true vectors of linear system.
 
 Operator (int s=0)
 Construct a square Operator with given size s (default 0).
 
 Operator (int h, int w)
 Construct an Operator with the given height (output size) and width (input size).
 
int Height () const
 Get the height (size of output) of the Operator. Synonym with NumRows().
 
int NumRows () const
 Get the number of rows (size of output) of the Operator. Synonym with Height().
 
int Width () const
 Get the width (size of input) of the Operator. Synonym with NumCols().
 
int NumCols () const
 Get the number of columns (size of input) of the Operator. Synonym with Width().
 
virtual MemoryClass GetMemoryClass () const
 Return the MemoryClass preferred by the Operator.
 
virtual void MultTranspose (const Vector &x, Vector &y) const
 Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an error.
 
virtual void AddMult (const Vector &x, Vector &y, const real_t a=1.0) const
 Operator application: y+=A(x) (default) or y+=a*A(x).
 
virtual void AddMultTranspose (const Vector &x, Vector &y, const real_t a=1.0) const
 Operator transpose application: y+=A^t(x) (default) or y+=a*A^t(x).
 
virtual void ArrayMult (const Array< const Vector * > &X, Array< Vector * > &Y) const
 Operator application on a matrix: Y=A(X).
 
virtual void ArrayMultTranspose (const Array< const Vector * > &X, Array< Vector * > &Y) const
 Action of the transpose operator on a matrix: Y=A^t(X).
 
virtual void ArrayAddMult (const Array< const Vector * > &X, Array< Vector * > &Y, const real_t a=1.0) const
 Operator application on a matrix: Y+=A(X) (default) or Y+=a*A(X).
 
virtual void ArrayAddMultTranspose (const Array< const Vector * > &X, Array< Vector * > &Y, const real_t a=1.0) const
 Operator transpose application on a matrix: Y+=A^t(X) (default) or Y+=a*A^t(X).
 
virtual OperatorGetGradient (const Vector &x) const
 Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate an error.
 
virtual void AssembleDiagonal (Vector &diag) const
 Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operators. In some cases, only an approximation of the diagonal is computed.
 
virtual const OperatorGetProlongation () const
 Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity.
 
virtual const OperatorGetRestriction () const
 Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity.
 
virtual const OperatorGetOutputProlongation () const
 Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator. NULL means identity.
 
virtual const OperatorGetOutputRestrictionTranspose () const
 Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type operators.
 
virtual const OperatorGetOutputRestriction () const
 Restriction operator from output vectors for the operator to linear algebra (linear system) vectors. NULL means identity.
 
void FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0)
 Form a constrained linear system using a matrix-free approach.
 
void FormRectangularLinearSystem (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B)
 Form a column-constrained linear system using a matrix-free approach.
 
virtual void RecoverFEMSolution (const Vector &X, const Vector &b, Vector &x)
 Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear system obtained from Operator::FormLinearSystem() or Operator::FormRectangularLinearSystem().
 
void FormSystemOperator (const Array< int > &ess_tdof_list, Operator *&A)
 Return in A a parallel (on truedofs) version of this square operator.
 
void FormRectangularSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Operator *&A)
 Return in A a parallel (on truedofs) version of this rectangular operator (including constraints).
 
void FormDiscreteOperator (Operator *&A)
 Return in A a parallel (on truedofs) version of this rectangular operator.
 
void PrintMatlab (std::ostream &out, int n, int m=0) const
 Prints operator with input size n and output size m in Matlab format.
 
virtual void PrintMatlab (std::ostream &out) const
 Prints operator in Matlab format.
 
virtual ~Operator ()
 Virtual destructor.
 
Type GetType () const
 Return the type ID of the Operator class.
 
- Public Member Functions inherited from mfem::SundialsSolver
void * GetMem () const
 Access the SUNDIALS memory structure.
 
int GetFlag () const
 Returns the last flag returned by a call to a SUNDIALS function.
 

Protected Member Functions

void SetJFNKSolver (Solver &solver)
 
- Protected Member Functions inherited from mfem::NewtonSolver
void AdaptiveLinRtolPreSolve (const Vector &x, const int it, const real_t fnorm) const
 Method for the adaptive linear solver rtol invoked before the linear solve.
 
void AdaptiveLinRtolPostSolve (const Vector &x, const Vector &b, const int it, const real_t fnorm) const
 Method for the adaptive linear solver rtol invoked after the linear solve.
 
- Protected Member Functions inherited from mfem::IterativeSolver
virtual real_t Dot (const Vector &x, const Vector &y) const
 Return the standard (l2, i.e., Euclidean) inner product of x and y.
 
real_t Norm (const Vector &x) const
 Return the inner product norm of x, using the inner product defined by Dot()
 
bool Monitor (int it, real_t norm, const Vector &r, const Vector &x, bool final=false) const
 Monitor both the residual r and the solution x.
 
PrintLevel FromLegacyPrintLevel (int)
 Convert a legacy print level integer to a PrintLevel object.
 
- Protected Member Functions inherited from mfem::Operator
void FormConstrainedSystemOperator (const Array< int > &ess_tdof_list, ConstrainedOperator *&Aout)
 see FormSystemOperator()
 
void FormRectangularConstrainedSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, RectangularConstrainedOperator *&Aout)
 see FormRectangularSystemOperator()
 
OperatorSetupRAP (const Operator *Pi, const Operator *Po)
 Returns RAP Operator of this, using input/output Prolongation matrices Pi corresponds to "P", Po corresponds to "Rt".
 
- Protected Member Functions inherited from mfem::SundialsSolver
bool Parallel () const
 
bool Parallel () const
 
 SundialsSolver ()
 Protected constructor: objects of this type should be constructed only as part of a derived class.
 
void AllocateEmptyNVector (N_Vector &y)
 
void AllocateEmptyNVector (N_Vector &y, MPI_Comm comm)
 

Static Protected Member Functions

static int Mult (const N_Vector u, N_Vector fu, void *user_data)
 Wrapper to compute the nonlinear residual \( F(u) = 0 \).
 
static int GradientMult (N_Vector v, N_Vector Jv, N_Vector u, sunbooleantype *new_u, void *user_data)
 Wrapper to compute the Jacobian-vector product \( J(u) v = Jv \).
 
static int LinSysSetup (N_Vector u, N_Vector fu, SUNMatrix J, void *user_data, N_Vector tmp1, N_Vector tmp2)
 Setup the linear system \( J u = b \).
 
static int LinSysSolve (SUNLinearSolver LS, SUNMatrix J, N_Vector u, N_Vector b, sunrealtype tol)
 Solve the linear system \( J u = b \).
 
static int PrecSetup (N_Vector uu, N_Vector uscale, N_Vector fval, N_Vector fscale, void *user_data)
 Setup the preconditioner.
 
static int PrecSolve (N_Vector uu, N_Vector uscale, N_Vector fval, N_Vector fscale, N_Vector vv, void *user_data)
 Solve the preconditioner equation \( Pz = v \).
 
- Static Protected Member Functions inherited from mfem::IterativeSolver
static int GuessLegacyPrintLevel (PrintLevel)
 Use some heuristics to guess a legacy print level corresponding to the given PrintLevel.
 

Protected Attributes

int global_strategy
 KINSOL solution strategy.
 
bool use_oper_grad
 use the Jv prod function
 
SundialsNVectory_scale
 
SundialsNVectorf_scale
 scaling vectors
 
const Operatorjacobian
 stores oper->GetGradient()
 
int aa_n = 0
 number of acceleration vectors
 
int aa_delay
 Anderson Acceleration delay.
 
double aa_damping
 Anderson Acceleration damping.
 
int aa_orth
 Anderson Acceleration orthogonalization routine.
 
double fp_damping = 1.0
 Fixed Point or Picard damping parameter.
 
bool jfnk = false
 enable JFNK
 
Vector wrk
 Work vector needed for the JFNK PC.
 
int maxli = 5
 Maximum linear iterations.
 
int maxlrs = 0
 Maximum linear solver restarts.
 
- Protected Attributes inherited from mfem::NewtonSolver
Vector r
 
Vector c
 
Operatorgrad
 
int lin_rtol_type = 0
 
real_t lin_rtol0
 
real_t lin_rtol_max
 
real_t fnorm_last = 0.0
 
real_t lnorm_last = 0.0
 
real_t eta_last = 0.0
 
real_t gamma
 
real_t alpha
 
- Protected Attributes inherited from mfem::IterativeSolver
const Operatoroper
 
Solverprec
 
IterativeSolverMonitormonitor = nullptr
 
int max_iter
 Limit for the number of iterations the solver is allowed to do.
 
real_t rel_tol
 Relative tolerance.
 
real_t abs_tol
 Absolute tolerance.
 
int final_iter = -1
 
bool converged = false
 
real_t initial_norm = -1.0
 
real_t final_norm = -1.0
 
int print_level = -1
 (DEPRECATED) Legacy print level definition, which is left for compatibility with custom iterative solvers.
 
PrintLevel print_options
 Output behavior for the iterative solver.
 
- Protected Attributes inherited from mfem::Operator
int height
 Dimension of the output / number of rows in the matrix.
 
int width
 Dimension of the input / number of columns in the matrix.
 
- Protected Attributes inherited from mfem::SundialsSolver
void * sundials_mem
 SUNDIALS mem structure.
 
int flag
 Last flag returned from a call to SUNDIALS.
 
bool reinit
 Flag to signal memory reinitialization is need.
 
long saved_global_size
 Global vector length on last initialization.
 
SundialsNVectorY
 State vector.
 
SUNMatrix A
 
SUNMatrix M
 Mass matrix M.
 
SUNLinearSolver LSA
 Linear solver for A.
 
SUNLinearSolver LSM
 Linear solver for M.
 
SUNNonlinearSolver NLS
 Nonlinear solver.
 

Additional Inherited Members

- Public Types inherited from mfem::Operator
enum  DiagonalPolicy { DIAG_ZERO , DIAG_ONE , DIAG_KEEP }
 Defines operator diagonal policy upon elimination of rows and/or columns. More...
 
enum  Type {
  ANY_TYPE , MFEM_SPARSEMAT , Hypre_ParCSR , PETSC_MATAIJ ,
  PETSC_MATIS , PETSC_MATSHELL , PETSC_MATNEST , PETSC_MATHYPRE ,
  PETSC_MATGENERIC , Complex_Operator , MFEM_ComplexSparseMat , Complex_Hypre_ParCSR ,
  Complex_DenseMat , MFEM_Block_Matrix , MFEM_Block_Operator
}
 Enumeration defining IDs for some classes derived from Operator. More...
 
- Public Attributes inherited from mfem::Solver
bool iterative_mode
 If true, use the second argument of Mult() as an initial guess.
 
- Static Protected Attributes inherited from mfem::SundialsSolver
static constexpr double default_rel_tol = 1e-4
 Default scalar relative tolerance.
 
static constexpr double default_abs_tol = 1e-9
 Default scalar absolute tolerance.
 

Detailed Description

Interface to the KINSOL library – nonlinear solver methods.

Definition at line 886 of file sundials.hpp.

Constructor & Destructor Documentation

◆ KINSolver() [1/2]

mfem::KINSolver::KINSolver ( int strategy,
bool oper_grad = true )

Construct a serial wrapper to SUNDIALS' KINSOL nonlinear solver.

Parameters
[in]strategySpecifies the nonlinear solver strategy: KIN_NONE / KIN_LINESEARCH / KIN_PICARD / KIN_FP.
[in]oper_gradSpecifies whether the solver should use its Operator's GetGradient() method to compute the Jacobian of the system.

Definition at line 1996 of file sundials.cpp.

◆ KINSolver() [2/2]

mfem::KINSolver::KINSolver ( MPI_Comm comm,
int strategy,
bool oper_grad = true )

Construct a parallel wrapper to SUNDIALS' KINSOL nonlinear solver.

Parameters
[in]commThe MPI communicator used to partition the system.
[in]strategySpecifies the nonlinear solver strategy: KIN_NONE / KIN_LINESEARCH / KIN_PICARD / KIN_FP.
[in]oper_gradSpecifies whether the solver should use its Operator's GetGradient() method to compute the Jacobian of the system.

Definition at line 2014 of file sundials.cpp.

◆ ~KINSolver()

mfem::KINSolver::~KINSolver ( )
virtual

Destroy the associated KINSOL memory.

Definition at line 2410 of file sundials.cpp.

Member Function Documentation

◆ EnableAndersonAcc()

void mfem::KINSolver::EnableAndersonAcc ( int n,
int orth = KIN_ORTH_MGS,
int delay = 0,
double damping = 1.0 )

Enable Anderson Acceleration for KIN_FP or KIN_PICARD.

Note
Has to be called once before SetOperator() in order to set up the maximum subspace size. Subsequent calls need n less or equal to the initial subspace size.
Parameters
[in]nAnderson Acceleration subspace size
[in]orthAnderson Acceleration orthogonalization routine
[in]delayAnderson Acceleration delay
[in]dampingAnderson Acceleration damping parameter valid from 0 < d <= 1.0. Default is 1.0 (no damping)

Definition at line 2239 of file sundials.cpp.

◆ GradientMult()

int mfem::KINSolver::GradientMult ( N_Vector v,
N_Vector Jv,
N_Vector u,
sunbooleantype * new_u,
void * user_data )
staticprotected

Wrapper to compute the Jacobian-vector product \( J(u) v = Jv \).

Definition at line 1905 of file sundials.cpp.

◆ LinSysSetup()

int mfem::KINSolver::LinSysSetup ( N_Vector u,
N_Vector fu,
SUNMatrix J,
void * user_data,
N_Vector tmp1,
N_Vector tmp2 )
staticprotected

Setup the linear system \( J u = b \).

Definition at line 1928 of file sundials.cpp.

◆ LinSysSolve()

int mfem::KINSolver::LinSysSolve ( SUNLinearSolver LS,
SUNMatrix J,
N_Vector u,
N_Vector b,
sunrealtype tol )
staticprotected

Solve the linear system \( J u = b \).

Definition at line 1945 of file sundials.cpp.

◆ Mult() [1/3]

int mfem::KINSolver::Mult ( const N_Vector u,
N_Vector fu,
void * user_data )
staticprotected

Wrapper to compute the nonlinear residual \( F(u) = 0 \).

Definition at line 1891 of file sundials.cpp.

◆ Mult() [2/3]

void mfem::KINSolver::Mult ( const Vector & b,
Vector & x ) const
overridevirtual

Solve the nonlinear system \( F(x) = 0 \).

This method computes the x_scale and fx_scale vectors and calls the other Mult(Vector&, Vector&, Vector&) const method. The x_scale vector is a vector of ones and values of fx_scale are determined by comparing the chosen relative and functional norm (i.e. absolute) tolerances.

Parameters
[in]bNot used, KINSOL always assumes zero RHS
[in,out]xOn input, initial guess, if iterative_mode = true, otherwise the initial guess is zero; on output, the solution

Implements mfem::Operator.

Definition at line 2294 of file sundials.cpp.

◆ Mult() [3/3]

void mfem::KINSolver::Mult ( Vector & x,
const Vector & x_scale,
const Vector & fx_scale ) const

Solve the nonlinear system \( F(x) = 0 \).

Calls KINSol() to solve the nonlinear system. Before calling KINSol(), this functions uses the data members inherited from class IterativeSolver to set corresponding KINSOL options.

Parameters
[in,out]xOn input, initial guess, if iterative_mode = true, otherwise the initial guess is zero; on output, the solution
[in]x_scaleElements of a diagonal scaling matrix D, s.t. D*x has all elements roughly the same when x is close to a solution
[in]fx_scaleElements of a diagonal scaling matrix E, s.t. D*F(x) has all elements roughly the same when x is not too close to a solution

Definition at line 2346 of file sundials.cpp.

◆ PrecSetup()

int mfem::KINSolver::PrecSetup ( N_Vector uu,
N_Vector uscale,
N_Vector fval,
N_Vector fscale,
void * user_data )
staticprotected

Setup the preconditioner.

Definition at line 1958 of file sundials.cpp.

◆ PrecSolve()

int mfem::KINSolver::PrecSolve ( N_Vector uu,
N_Vector uscale,
N_Vector fval,
N_Vector fscale,
N_Vector vv,
void * user_data )
staticprotected

Solve the preconditioner equation \( Pz = v \).

Definition at line 1976 of file sundials.cpp.

◆ SetDamping()

void mfem::KINSolver::SetDamping ( double damping)

Specifies the value of the damping parameter in the fixed point or Picard iteration. param[in] damping fixed point iteration or Picard damping parameter

Definition at line 2278 of file sundials.cpp.

◆ SetJFNK()

void mfem::KINSolver::SetJFNK ( bool use_jfnk)
inline

Set the Jacobian Free Newton Krylov flag. The default is false.

This flag indicates to use JFNK as the linear solver for KINSOL. This means that the Solver object set in SetSolver() or SetPreconditioner() is used as a preconditioner for an FGMRES algorithm provided by SpFGMR from SUNDIALS. Furthermore, all Jacobian-vector products in the outer Krylov method are approximated by a difference quotient and the relative tolerance for the outer Krylov method is adaptive. See the KINSOL User Manual for details.

Definition at line 1012 of file sundials.hpp.

◆ SetJFNKSolver()

void mfem::KINSolver::SetJFNKSolver ( Solver & solver)
protected

Definition at line 2196 of file sundials.cpp.

◆ SetLSMaxIter()

void mfem::KINSolver::SetLSMaxIter ( int m)
inline

Set the maximum number of linear solver iterations.

Note
Only valid in combination with JFNK

Definition at line 1016 of file sundials.hpp.

◆ SetLSMaxRestarts()

void mfem::KINSolver::SetLSMaxRestarts ( int m)
inline

Set the maximum number of linear solver restarts.

Note
Only valid in combination with JFNK

Definition at line 1020 of file sundials.hpp.

◆ SetMaxSetupCalls()

void mfem::KINSolver::SetMaxSetupCalls ( int max_calls)

Set maximum number of nonlinear iterations without a Jacobian update.

The default is 10.

Note
This method must be called after SetOperator().

Definition at line 2233 of file sundials.cpp.

◆ SetOperator()

void mfem::KINSolver::SetOperator ( const Operator & op)
overridevirtual

Set the nonlinear Operator of the system and initialize KINSOL.

Note
If this method is called a second time with a different problem size, then non-default KINSOL-specific options will be lost and will need to be set again.

Reimplemented from mfem::IterativeSolver.

Definition at line 2033 of file sundials.cpp.

◆ SetPreconditioner()

void mfem::KINSolver::SetPreconditioner ( Solver & solver)
inlineoverridevirtual

Equivalent to SetSolver(solver).

Reimplemented from mfem::IterativeSolver.

Definition at line 974 of file sundials.hpp.

◆ SetPrintLevel() [1/2]

void mfem::KINSolver::SetPrintLevel ( int print_lvl)
inlineoverridevirtual

Set the print level for the KINSetPrintLevel function.

Reimplemented from mfem::IterativeSolver.

Definition at line 1023 of file sundials.hpp.

◆ SetPrintLevel() [2/2]

void mfem::KINSolver::SetPrintLevel ( PrintLevel )
overridevirtual

This method is not supported and will throw an error.

Reimplemented from mfem::IterativeSolver.

Definition at line 2288 of file sundials.cpp.

◆ SetScaledStepTol()

void mfem::KINSolver::SetScaledStepTol ( double sstol)

Set KINSOL's scaled step tolerance.

The default tolerance is \( U^\frac{2}{3} \) , where U = machine unit round-off.

Note
This method must be called after SetOperator().

Definition at line 2227 of file sundials.cpp.

◆ SetSolver()

void mfem::KINSolver::SetSolver ( Solver & solver)
overridevirtual

Set the linear solver for inverting the Jacobian.

Note
This function assumes that Operator::GetGradient(const Vector &) is implemented by the Operator specified by SetOperator(const Operator &).

This method must be called after SetOperator().

Reimplemented from mfem::NewtonSolver.

Definition at line 2155 of file sundials.cpp.

Member Data Documentation

◆ aa_damping

double mfem::KINSolver::aa_damping
protected

Anderson Acceleration damping.

Definition at line 895 of file sundials.hpp.

◆ aa_delay

int mfem::KINSolver::aa_delay
protected

Anderson Acceleration delay.

Definition at line 894 of file sundials.hpp.

◆ aa_n

int mfem::KINSolver::aa_n = 0
protected

number of acceleration vectors

Definition at line 893 of file sundials.hpp.

◆ aa_orth

int mfem::KINSolver::aa_orth
protected

Anderson Acceleration orthogonalization routine.

Definition at line 896 of file sundials.hpp.

◆ f_scale

SundialsNVector * mfem::KINSolver::f_scale
protected

scaling vectors

Definition at line 891 of file sundials.hpp.

◆ fp_damping

double mfem::KINSolver::fp_damping = 1.0
protected

Fixed Point or Picard damping parameter.

Definition at line 897 of file sundials.hpp.

◆ global_strategy

int mfem::KINSolver::global_strategy
protected

KINSOL solution strategy.

Definition at line 889 of file sundials.hpp.

◆ jacobian

const Operator* mfem::KINSolver::jacobian
protected

stores oper->GetGradient()

Definition at line 892 of file sundials.hpp.

◆ jfnk

bool mfem::KINSolver::jfnk = false
protected

enable JFNK

Definition at line 898 of file sundials.hpp.

◆ maxli

int mfem::KINSolver::maxli = 5
protected

Maximum linear iterations.

Definition at line 900 of file sundials.hpp.

◆ maxlrs

int mfem::KINSolver::maxlrs = 0
protected

Maximum linear solver restarts.

Definition at line 901 of file sundials.hpp.

◆ use_oper_grad

bool mfem::KINSolver::use_oper_grad
protected

use the Jv prod function

Definition at line 890 of file sundials.hpp.

◆ wrk

Vector mfem::KINSolver::wrk
protected

Work vector needed for the JFNK PC.

Definition at line 899 of file sundials.hpp.

◆ y_scale

SundialsNVector* mfem::KINSolver::y_scale
mutableprotected

Definition at line 891 of file sundials.hpp.


The documentation for this class was generated from the following files: