MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | List of all members
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. More...
 
 KINSolver (MPI_Comm comm, int strategy, bool oper_grad=true)
 Construct a parallel wrapper to SUNDIALS' KINSOL nonlinear solver. More...
 
virtual ~KINSolver ()
 Destroy the associated KINSOL memory. More...
 
virtual void SetOperator (const Operator &op)
 Set the nonlinear Operator of the system and initialize KINSOL. More...
 
virtual void SetSolver (Solver &solver)
 Set the linear solver for inverting the Jacobian. More...
 
virtual void SetPreconditioner (Solver &solver)
 Equivalent to SetSolver(solver). More...
 
void SetScaledStepTol (double sstol)
 Set KINSOL's scaled step tolerance. More...
 
void SetMaxSetupCalls (int max_calls)
 Set maximum number of nonlinear iterations without a Jacobian update. More...
 
void SetMAA (int maa)
 Set the number of acceleration vectors to use with KIN_FP or KIN_PICARD. More...
 
void SetJFNK (bool use_jfnk)
 Set the Jacobian Free Newton Krylov flag. The default is false. More...
 
void SetLSMaxIter (int m)
 Set the maximum number of linear solver iterations. More...
 
void SetLSMaxRestarts (int m)
 Set the maximum number of linear solver restarts. More...
 
virtual void SetPrintLevel (int print_lvl)
 Set the print level for the KINSetPrintLevel function. More...
 
virtual void SetPrintLevel (PrintLevel)
 This method is not supported and will throw an error. More...
 
virtual void Mult (const Vector &b, Vector &x) const
 Solve the nonlinear system \( F(x) = 0 \). More...
 
void Mult (Vector &x, const Vector &x_scale, const Vector &fx_scale) const
 Solve the nonlinear system \( F(x) = 0 \). More...
 
- Public Member Functions inherited from mfem::NewtonSolver
 NewtonSolver ()
 
 NewtonSolver (MPI_Comm comm_)
 
virtual double ComputeScalingFactor (const Vector &x, const Vector &b) const
 This method can be overloaded in derived classes to implement line search algorithms. More...
 
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. More...
 
void SetAdaptiveLinRtol (const int type=2, const double rtol0=0.5, const double rtol_max=0.9, const double alpha=0.5 *(1.0+sqrt(5.0)), const double gamma=1.0)
 Enable adaptive linear solver relative tolerance algorithm. More...
 
- Public Member Functions inherited from mfem::IterativeSolver
 IterativeSolver ()
 
 IterativeSolver (MPI_Comm comm_)
 
void SetMonitor (IterativeSolverMonitor &m)
 Set the iterative solver monitor. More...
 
MPI_Comm GetComm () const
 Return the associated MPI communicator, or MPI_COMM_NULL if no communicator is set. More...
 
void SetRelTol (double rtol)
 
void SetAbsTol (double atol)
 
void SetMaxIter (int max_it)
 
int GetNumIterations () const
 
bool GetConverged () const
 
double GetFinalNorm () const
 
- Public Member Functions inherited from mfem::Solver
 Solver (int s=0, bool iter_mode=false)
 Initialize a square Solver with size s. More...
 
 Solver (int h, int w, bool iter_mode=false)
 Initialize a Solver with height h and width w. More...
 
- 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. More...
 
 Operator (int s=0)
 Construct a square Operator with given size s (default 0). More...
 
 Operator (int h, int w)
 Construct an Operator with the given height (output size) and width (input size). More...
 
int Height () const
 Get the height (size of output) of the Operator. Synonym with NumRows(). More...
 
int NumRows () const
 Get the number of rows (size of output) of the Operator. Synonym with Height(). More...
 
int Width () const
 Get the width (size of input) of the Operator. Synonym with NumCols(). More...
 
int NumCols () const
 Get the number of columns (size of input) of the Operator. Synonym with Width(). More...
 
virtual MemoryClass GetMemoryClass () const
 Return the MemoryClass preferred by the Operator. More...
 
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. More...
 
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. More...
 
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. More...
 
virtual const OperatorGetProlongation () const
 Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity. More...
 
virtual const OperatorGetRestriction () const
 Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity. More...
 
virtual const OperatorGetOutputProlongation () const
 Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator. NULL means identity. More...
 
virtual const OperatorGetOutputRestrictionTranspose () const
 Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type operators. More...
 
virtual const OperatorGetOutputRestriction () const
 Restriction operator from output vectors for the operator to linear algebra (linear system) vectors. NULL means identity. More...
 
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. More...
 
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. More...
 
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(). More...
 
void FormSystemOperator (const Array< int > &ess_tdof_list, Operator *&A)
 Return in A a parallel (on truedofs) version of this square operator. More...
 
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). More...
 
void FormDiscreteOperator (Operator *&A)
 Return in A a parallel (on truedofs) version of this rectangular operator. More...
 
void PrintMatlab (std::ostream &out, int n, int m=0) const
 Prints operator with input size n and output size m in Matlab format. More...
 
virtual void PrintMatlab (std::ostream &out) const
 Prints operator in Matlab format. More...
 
virtual ~Operator ()
 Virtual destructor. More...
 
Type GetType () const
 Return the type ID of the Operator class. More...
 
- Public Member Functions inherited from mfem::SundialsSolver
void * GetMem () const
 Access the SUNDIALS memory structure. More...
 
int GetFlag () const
 Returns the last flag returned by a call to a SUNDIALS function. More...
 

Protected Member Functions

void SetJFNKSolver (Solver &solver)
 
- Protected Member Functions inherited from mfem::NewtonSolver
void AdaptiveLinRtolPreSolve (const Vector &x, const int it, const double fnorm) const
 Method for the adaptive linear solver rtol invoked before the linear solve. More...
 
void AdaptiveLinRtolPostSolve (const Vector &x, const Vector &b, const int it, const double fnorm) const
 Method for the adaptive linear solver rtol invoked after the linear solve. More...
 
- Protected Member Functions inherited from mfem::IterativeSolver
double Dot (const Vector &x, const Vector &y) const
 
double Norm (const Vector &x) const
 
void Monitor (int it, double norm, const Vector &r, const Vector &x, bool final=false) const
 
PrintLevel FromLegacyPrintLevel (int)
 Convert a legacy print level integer to a PrintLevel object. More...
 
- Protected Member Functions inherited from mfem::Operator
void FormConstrainedSystemOperator (const Array< int > &ess_tdof_list, ConstrainedOperator *&Aout)
 see FormSystemOperator() More...
 
void FormRectangularConstrainedSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, RectangularConstrainedOperator *&Aout)
 see FormRectangularSystemOperator() More...
 
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". More...
 
- 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. More...
 
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 \). More...
 
static int GradientMult (N_Vector v, N_Vector Jv, N_Vector u, booleantype *new_u, void *user_data)
 Wrapper to compute the Jacobian-vector product \( J(u) v = Jv \). More...
 
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 \). More...
 
static int LinSysSolve (SUNLinearSolver LS, SUNMatrix J, N_Vector u, N_Vector b, realtype tol)
 Solve the linear system \( J u = b \). More...
 
static int PrecSetup (N_Vector uu, N_Vector uscale, N_Vector fval, N_Vector fscale, void *user_data)
 Setup the preconditioner. More...
 
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 \). More...
 
- 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. More...
 

Protected Attributes

int global_strategy
 KINSOL solution strategy. More...
 
bool use_oper_grad
 use the Jv prod function More...
 
SundialsNVectory_scale
 
SundialsNVectorf_scale
 scaling vectors More...
 
const Operatorjacobian
 stores oper->GetGradient() More...
 
int maa
 number of acceleration vectors More...
 
bool jfnk = false
 enable JFNK More...
 
Vector wrk
 Work vector needed for the JFNK PC. More...
 
int maxli = 5
 Maximum linear iterations. More...
 
int maxlrs = 0
 Maximum linear solver restarts. More...
 
- Protected Attributes inherited from mfem::NewtonSolver
Vector r
 
Vector c
 
Operatorgrad
 
int lin_rtol_type = 0
 
double lin_rtol0
 
double lin_rtol_max
 
double fnorm_last = 0.0
 
double lnorm_last = 0.0
 
double eta_last = 0.0
 
double gamma
 
double 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. More...
 
double rel_tol
 Relative tolerance. More...
 
double abs_tol
 Absolute tolerance. More...
 
int final_iter
 
bool converged
 
double final_norm
 
int print_level = -1
 (DEPRECATED) Legacy print level definition, which is left for compatibility with custom iterative solvers. More...
 
PrintLevel print_options
 Output behavior for the iterative solver. More...
 
- Protected Attributes inherited from mfem::Operator
int height
 Dimension of the output / number of rows in the matrix. More...
 
int width
 Dimension of the input / number of columns in the matrix. More...
 
- Protected Attributes inherited from mfem::SundialsSolver
void * sundials_mem
 SUNDIALS mem structure. More...
 
int flag
 Last flag returned from a call to SUNDIALS. More...
 
bool reinit
 Flag to signal memory reinitialization is need. More...
 
long saved_global_size
 Global vector length on last initialization. More...
 
SundialsNVectorY
 State vector. More...
 
SUNMatrix A
 
SUNMatrix M
 Mass matrix M. More...
 
SUNLinearSolver LSA
 Linear solver for A. More...
 
SUNLinearSolver LSM
 Linear solver for M. More...
 
SUNNonlinearSolver NLS
 Nonlinear solver. More...
 

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
}
 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. More...
 
- Static Protected Attributes inherited from mfem::SundialsSolver
static constexpr double default_rel_tol = 1e-4
 Default scalar relative tolerance. More...
 
static constexpr double default_abs_tol = 1e-9
 Default scalar absolute tolerance. More...
 

Detailed Description

Interface to the KINSOL library – nonlinear solver methods.

Definition at line 713 of file sundials.hpp.

Constructor & Destructor Documentation

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 1747 of file sundials.cpp.

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 1761 of file sundials.cpp.

mfem::KINSolver::~KINSolver ( )
virtual

Destroy the associated KINSOL memory.

Definition at line 2098 of file sundials.cpp.

Member Function Documentation

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

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

Definition at line 1656 of file sundials.cpp.

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 1679 of file sundials.cpp.

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

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

Definition at line 1696 of file sundials.cpp.

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 1642 of file sundials.cpp.

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

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

Reimplemented from mfem::NewtonSolver.

Definition at line 1985 of file sundials.cpp.

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 2037 of file sundials.cpp.

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 1709 of file sundials.cpp.

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 1727 of file sundials.cpp.

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 826 of file sundials.hpp.

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

Definition at line 1925 of file sundials.cpp.

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 830 of file sundials.hpp.

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 834 of file sundials.hpp.

void mfem::KINSolver::SetMAA ( int  maa)

Set the number of acceleration vectors to use with KIN_FP or KIN_PICARD.

The default is 0. @ note This method must be called before SetOperator() to set the maximum size of the acceleration space. The value of maa can be altered after SetOperator() is called but it can't be higher than initial maximum.

Definition at line 1967 of file sundials.cpp.

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 1961 of file sundials.cpp.

void mfem::KINSolver::SetOperator ( const Operator op)
virtual

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::NewtonSolver.

Definition at line 1776 of file sundials.cpp.

virtual void mfem::KINSolver::SetPreconditioner ( Solver solver)
inlinevirtual

Equivalent to SetSolver(solver).

Reimplemented from mfem::IterativeSolver.

Definition at line 797 of file sundials.hpp.

virtual void mfem::KINSolver::SetPrintLevel ( int  print_lvl)
inlinevirtual

Set the print level for the KINSetPrintLevel function.

Reimplemented from mfem::IterativeSolver.

Definition at line 837 of file sundials.hpp.

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

This method is not supported and will throw an error.

Reimplemented from mfem::IterativeSolver.

Definition at line 1979 of file sundials.cpp.

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 1955 of file sundials.cpp.

void mfem::KINSolver::SetSolver ( Solver solver)
virtual

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 1884 of file sundials.cpp.

Member Data Documentation

SundialsNVector * mfem::KINSolver::f_scale
mutableprotected

scaling vectors

Definition at line 718 of file sundials.hpp.

int mfem::KINSolver::global_strategy
protected

KINSOL solution strategy.

Definition at line 716 of file sundials.hpp.

const Operator* mfem::KINSolver::jacobian
protected

stores oper->GetGradient()

Definition at line 719 of file sundials.hpp.

bool mfem::KINSolver::jfnk = false
protected

enable JFNK

Definition at line 721 of file sundials.hpp.

int mfem::KINSolver::maa
protected

number of acceleration vectors

Definition at line 720 of file sundials.hpp.

int mfem::KINSolver::maxli = 5
protected

Maximum linear iterations.

Definition at line 723 of file sundials.hpp.

int mfem::KINSolver::maxlrs = 0
protected

Maximum linear solver restarts.

Definition at line 724 of file sundials.hpp.

bool mfem::KINSolver::use_oper_grad
protected

use the Jv prod function

Definition at line 717 of file sundials.hpp.

Vector mfem::KINSolver::wrk
protected

Work vector needed for the JFNK PC.

Definition at line 722 of file sundials.hpp.

SundialsNVector* mfem::KINSolver::y_scale
mutableprotected

Definition at line 718 of file sundials.hpp.


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