MFEM
v4.1.0
Finite element discretization library
|
Interface to the KINSOL library – nonlinear solver methods. More...
#include <sundials.hpp>
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... | |
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... | |
Public Member Functions inherited from mfem::IterativeSolver | |
IterativeSolver () | |
IterativeSolver (MPI_Comm _comm) | |
void | SetRelTol (double rtol) |
void | SetAbsTol (double atol) |
void | SetMaxIter (int max_it) |
void | SetPrintLevel (int print_lvl) |
int | GetNumIterations () const |
int | 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 Operator & | GetGradient (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 const Operator * | GetProlongation () const |
Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity. More... | |
virtual const Operator * | GetRestriction () const |
Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity. More... | |
virtual const Operator * | GetOutputProlongation () const |
Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator. NULL means identity. More... | |
virtual const Operator * | GetOutputRestriction () 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=0, int m=0) const |
Prints operator with input size n and output size m 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 retured by a call to a SUNDIALS function. More... | |
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... | |
Protected Attributes | |
int | global_strategy |
KINSOL solution strategy. More... | |
bool | use_oper_grad |
use the Jv prod function More... | |
N_Vector | y_scale |
N_Vector | f_scale |
scaling vectors More... | |
const Operator * | jacobian |
stores oper->GetGradient() More... | |
int | maa |
number of acceleration vectors More... | |
Protected Attributes inherited from mfem::NewtonSolver | |
Vector | r |
Vector | c |
Protected Attributes inherited from mfem::IterativeSolver | |
const Operator * | oper |
Solver * | prec |
int | max_iter |
int | print_level |
double | rel_tol |
double | abs_tol |
int | final_iter |
int | converged |
double | final_norm |
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... | |
N_Vector | y |
State vector. More... | |
SUNMatrix | A |
Linear system A = I - gamma J, M - gamma J, or J. More... | |
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 | 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... | |
Protected Member Functions inherited from mfem::IterativeSolver | |
double | Dot (const Vector &x, const Vector &y) const |
double | Norm (const Vector &x) const |
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... | |
Operator * | SetupRAP (const Operator *Pi, const Operator *Po) |
Returns RAP Operator of this, taking in input/output Prolongation matrices. 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... | |
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... | |
Interface to the KINSOL library – nonlinear solver methods.
Definition at line 370 of file sundials.hpp.
mfem::KINSolver::KINSolver | ( | int | strategy, |
bool | oper_grad = true |
||
) |
Construct a serial wrapper to SUNDIALS' KINSOL nonlinear solver.
[in] | strategy | Specifies the nonlinear solver strategy: KIN_NONE / KIN_LINESEARCH / KIN_PICARD / KIN_FP. |
[in] | oper_grad | Specifies whether the solver should use its Operator's GetGradient() method to compute the Jacobian of the system. |
Definition at line 1004 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.
[in] | comm | The MPI communicator used to partition the system. |
[in] | strategy | Specifies the nonlinear solver strategy: KIN_NONE / KIN_LINESEARCH / KIN_PICARD / KIN_FP. |
[in] | oper_grad | Specifies whether the solver should use its Operator's GetGradient() method to compute the Jacobian of the system. |
Definition at line 1020 of file sundials.cpp.
|
virtual |
Destroy the associated KINSOL memory.
Definition at line 1349 of file sundials.cpp.
|
staticprotected |
Wrapper to compute the Jacobian-vector product \( J(u) v = Jv \).
Definition at line 951 of file sundials.cpp.
|
staticprotected |
Setup the linear system \( J u = b \).
Definition at line 974 of file sundials.cpp.
|
staticprotected |
Solve the linear system \( J u = b \).
Definition at line 991 of file sundials.cpp.
|
staticprotected |
Wrapper to compute the nonlinear residual \( F(u) = 0 \).
Definition at line 937 of file sundials.cpp.
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.
[in] | b | Not used, KINSOL always assumes zero RHS |
[in,out] | x | On input, initial guess, if iterative_mode = true, otherwise the initial guess is zero; on output, the solution |
Reimplemented from mfem::NewtonSolver.
Definition at line 1244 of file sundials.cpp.
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.
[in,out] | x | On input, initial guess, if iterative_mode = true, otherwise the initial guess is zero; on output, the solution |
[in] | x_scale | Elements 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_scale | Elements 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 1295 of file sundials.cpp.
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 1231 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.
Definition at line 1225 of file sundials.cpp.
|
virtual |
Set the nonlinear Operator of the system and initialize KINSOL.
Reimplemented from mfem::NewtonSolver.
Definition at line 1052 of file sundials.cpp.
|
inlinevirtual |
Equivalent to SetSolver(solver).
Reimplemented from mfem::IterativeSolver.
Definition at line 433 of file sundials.hpp.
void mfem::KINSolver::SetScaledStepTol | ( | double | sstol | ) |
Set KINSOL's scaled step tolerance.
The default tolerance is \( U^\frac{2}{3} \) , where U = machine unit roundoff.
Definition at line 1219 of file sundials.cpp.
|
virtual |
Set the linear solver for inverting the Jacobian.
This method must be called after SetOperator().
Reimplemented from mfem::NewtonSolver.
Definition at line 1185 of file sundials.cpp.
|
mutableprotected |
scaling vectors
Definition at line 375 of file sundials.hpp.
|
protected |
KINSOL solution strategy.
Definition at line 373 of file sundials.hpp.
|
protected |
stores oper->GetGradient()
Definition at line 376 of file sundials.hpp.
|
protected |
number of acceleration vectors
Definition at line 377 of file sundials.hpp.
|
protected |
use the Jv prod function
Definition at line 374 of file sundials.hpp.
|
mutableprotected |
Definition at line 375 of file sundials.hpp.