#include <solvernonlinear.hpp>
Non linear solver base interface.
This class provides a uniform interface for nonlinear solvers. This base class is overloaded to provide nonlinear solvers from different packages like PETSC.
Public Types | |
Typedefs | |
| typedef SolverNonLinear< T > | self_type |
|
typedef boost::shared_ptr < SolverNonLinear< T > > | self_ptrtype |
| typedef self_type | solvernonlinear_type |
|
typedef boost::shared_ptr < self_type > | solvernonlinear_ptrtype |
| typedef T | value_type |
| typedef type_traits< T >::real_type | real_type |
|
typedef boost::shared_ptr < Preconditioner< T > > | preconditioner_ptrtype |
|
typedef boost::shared_ptr < Vector< value_type > > | vector_ptrtype |
|
typedef boost::shared_ptr < MatrixSparse< value_type > > | sparse_matrix_ptrtype |
| typedef ublas::matrix< value_type > | dense_matrix_type |
| typedef ublas::vector< value_type > | dense_vector_type |
|
typedef Eigen::Map < Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > > | map_dense_matrix_type |
|
typedef Eigen::Map < Eigen::Matrix< double, Eigen::Dynamic, 1 > > | map_dense_vector_type |
|
typedef boost::function< void(const vector_ptrtype &X, vector_ptrtype &R)> | residual_function_type |
|
typedef boost::function< void(const vector_ptrtype &X, sparse_matrix_ptrtype &J)> | jacobian_function_type |
|
typedef boost::function< void(const vector_ptrtype &X, vector_ptrtype &R, sparse_matrix_ptrtype &J)> | matvec_function_type |
|
typedef boost::function< void(dense_vector_type const &X, dense_vector_type &R)> | dense_residual_function_type |
|
typedef boost::function< void(dense_vector_type const &X, dense_matrix_type &J)> | dense_jacobian_function_type |
|
typedef boost::function< void(dense_vector_type const &X, dense_vector_type &R, dense_matrix_type &J)> | dense_matvec_function_type |
|
typedef boost::function< void(map_dense_vector_type const &X, map_dense_vector_type &R)> | map_dense_residual_function_type |
|
typedef boost::function< void(map_dense_vector_type const &X, map_dense_matrix_type &J)> | map_dense_jacobian_function_type |
|
typedef boost::function< void(map_dense_vector_type const &X, map_dense_vector_type &R, map_dense_matrix_type &J)> | map_dense_matvec_function_type |
Public Member Functions | |
| double | atoleranceKSP () const |
| double | dtoleranceKSP () const |
| size_type | maxitKSP () const |
| double | rtoleranceKSP () const |
| void | setAtoleranceKSP (double tol) |
| void | setDtoleranceKSP (double tol) |
| void | setMaxitKSP (size_type n) |
| void | setRtoleranceKSP (double tol) |
| void | setShowKSPConvergedReason (bool b) |
| void | setShowKSPMonitor (bool b) |
| void | setShowSNESConvergedReason (bool b) |
| void | setShowSNESMonitor (bool b) |
| bool | showKSPConvergedReason () const |
| bool | showKSPMonitor () const |
| bool | showSNESConvergedReason () const |
| bool | showSNESMonitor () const |
Accessors | |
| WorldComm const & | comm () const |
| WorldComm const & | worldComm () const |
| bool | initialized () const |
| virtual void | clear () |
| virtual MatrixStructure | precMatrixStructure () const |
| SolverNonLinearType | getType () const |
| double | getAbsoluteResidualTol () const |
| double | getRelativeResidualTol () const |
| double | getAbsoluteSolutionTol () const |
| uint | getNbItMax () const |
| int | reuseJacobian () const |
| int | reusePreconditioner () const |
| std::string const & | prefix () const |
Mutators | |
| void | setPrefix (std::string const &p) |
| virtual void | setPrecMatrixStructure (MatrixStructure mstruct) |
| void | setType (SolverNonLinearType snl_type) |
| SolverNonLinearType | nlSolverType () const |
| void | setKspSolverType (const SolverType st) |
| SolverType | kspSolverType () const |
| PreconditionerType | preconditionerType () const |
| void | setPreconditionerType (const PreconditionerType pct) |
| void | attachPreconditioner (preconditioner_ptrtype preconditioner) |
| void | setMatSolverPackageType (const MatSolverPackageType mspackt) |
| MatSolverPackageType | matSolverPackageType () const |
| virtual void | setReuse (int jac=1, int prec=1) |
| void | setRelativeResidualTol (double tol) |
| void | setAbsoluteResidualTol (double tol) |
| void | setAbsoluteSolutionTol (double tol) |
| void | setNbItMax (uint n) |
Protected Attributes | |
| double | M_absoluteResidualTol |
| double | M_absoluteSolutionTol |
| double | M_atoleranceKSP |
| double | M_dtoleranceKSP |
| bool | M_is_initialized |
| SolverType | M_kspSolver_type |
| MatSolverPackageType | M_matSolverPackage_type |
| size_type | M_maxitKSP |
| uint | M_nbItMax |
| MatrixStructure | M_prec_matrix_structure |
| preconditioner_ptrtype | M_preconditioner |
| PreconditionerType | M_preconditioner_type |
| std::string | M_prefix |
| double | M_relativeResidualTol |
| int | M_reuse_jac |
| int | M_reuse_prec |
| double | M_rtoleranceKSP |
| bool | M_showKSPConvergedReason |
| bool | M_showKSPMonitor |
| bool | M_showSNESConvergedReason |
| bool | M_showSNESMonitor |
| SolverNonLinearType | M_snl_type |
| WorldComm | M_worldComm |
Constructors, destructor | |
| SolverNonLinear (WorldComm const &worldComm=Environment::worldComm()) | |
| SolverNonLinear (SolverNonLinear const &) | |
| virtual | ~SolverNonLinear () |
| virtual void | init ()=0 |
| static solvernonlinear_ptrtype | build (po::variables_map const &vm, std::string const &prefix="", WorldComm const &worldComm=Environment::worldComm()) |
| static solvernonlinear_ptrtype | build (SolverPackage solver_package, WorldComm const &worldComm=Environment::worldComm()) |
Methods | |
| residual_function_type | residual |
| jacobian_function_type | jacobian |
| matvec_function_type | matvec |
| dense_residual_function_type | dense_residual |
| map_dense_residual_function_type | map_dense_residual |
| dense_jacobian_function_type | dense_jacobian |
| map_dense_jacobian_function_type | map_dense_jacobian |
| dense_matvec_function_type | dense_matvec |
| map_dense_matvec_function_type | map_dense_matvec |
| virtual std::pair< int, real_type > | solve (sparse_matrix_ptrtype &, vector_ptrtype &, vector_ptrtype &, const double, const unsigned int)=0 |
| virtual std::pair< unsigned int, real_type > | solve (dense_matrix_type &, dense_vector_type &, dense_vector_type &, const double, const unsigned int)=0 |
| virtual std::pair< unsigned int, real_type > | solve (map_dense_matrix_type &, map_dense_vector_type &, map_dense_vector_type &, const double, const unsigned int)=0 |
|
inline |
Constructor. Initializes Solver data structures
|
inline |
copy constructor
|
inlinevirtual |
Destructor.
|
inline |
KSP absolute tolerance
References Feel::SolverNonLinear< T >::M_atoleranceKSP.
|
inline |
Attaches a Preconditioner object to be used by the solver
References Feel::SolverNonLinear< T >::M_is_initialized, Feel::SolverNonLinear< T >::M_preconditioner, and Feel::SolverNonLinear< T >::M_preconditioner_type.
|
static |
Builds a NonlinearSolver using the nonlinear solver package specified by the variables_map vm and prefix
|
static |
Builds a NonlinearSolver using the nonlinear solver package specified by solver_package
|
inlinevirtual |
Release all memory and clear data structures.
|
inline |
|
inline |
KSP divergence tolerance
References Feel::SolverNonLinear< T >::M_dtoleranceKSP.
|
pure virtual |
Initialize data structures if not done so already.
|
inline |
References Feel::SolverNonLinear< T >::M_is_initialized.
|
inline |
Returns the type of solver to use.
References Feel::SolverNonLinear< T >::M_kspSolver_type.
|
inline |
Returns the type of preconditioner to use.
References Feel::SolverNonLinear< T >::M_matSolverPackage_type.
|
inline |
KSP maximum number of iterations
References Feel::SolverNonLinear< T >::M_maxitKSP.
|
inline |
Returns the type of solver to use.
References Feel::SolverNonLinear< T >::M_snl_type.
|
inlinevirtual |
|
inline |
Returns the type of preconditioner to use.
References Feel::SolverNonLinear< T >::M_preconditioner, and Feel::SolverNonLinear< T >::M_preconditioner_type.
|
inline |
|
inline |
KSP relative tolerance
References Feel::SolverNonLinear< T >::M_rtoleranceKSP.
|
inline |
Sets the type of solver to use.
References Feel::SolverNonLinear< T >::M_kspSolver_type.
|
inline |
Sets the type of preconditioner to use.
References Feel::SolverNonLinear< T >::M_matSolverPackage_type.
|
inline |
Define the number max of iteration
References Feel::SolverNonLinear< T >::M_nbItMax.
|
inlinevirtual |
References Feel::SolverNonLinear< T >::M_preconditioner.
|
inline |
Sets the type of preconditioner to use.
References Feel::SolverNonLinear< T >::M_preconditioner, and Feel::SolverNonLinear< T >::M_preconditioner_type.
|
inline |
set the prefix of the solver (typically for command line options)
|
inline |
Define values of tolerance for the non linear solver
References Feel::SolverNonLinear< T >::M_relativeResidualTol.
|
inlinevirtual |
set reuse jacobian and/or preconditioner
References Feel::SolverNonLinear< T >::M_reuse_jac, and Feel::SolverNonLinear< T >::M_reuse_prec.
|
inline |
Select type of non linear solver : LINEAR_SEARCH, TRUST_REGION, ...
References Feel::SolverNonLinear< T >::M_snl_type.
|
inline |
show KSP converged reason
|
inline |
show KSP monitor
|
inline |
show SNES converged reason
|
inline |
show SNES monitor
|
pure virtual |
Solves a sparse nonlinear system.
|
pure virtual |
Solves a dense nonlinear system.
|
pure virtual |
Solves a dense nonlinear system ( using eigen).
| dense_jacobian_function_type Feel::SolverNonLinear< T >::dense_jacobian |
Function that computes the Jacobian J(X) of the nonlinear system at the input iterate X.
| dense_matvec_function_type Feel::SolverNonLinear< T >::dense_matvec |
Function that computes either the residual
or the Jacobian
of the nonlinear system at the input iterate
. Note that either R or J could be XSNULL.
| dense_residual_function_type Feel::SolverNonLinear< T >::dense_residual |
Function that computes the residual R(X) of the nonlinear system at the input iterate X.
| jacobian_function_type Feel::SolverNonLinear< T >::jacobian |
Function that computes the Jacobian J(X) of the nonlinear system at the input iterate X.
|
protected |
Absolute tolerances between successive iteration
|
protected |
KSP absolute tolerance
|
protected |
KSP divergence tolerance
|
protected |
Flag indicating if the data structures have been initialized.
|
protected |
Enum stating which type of iterative linear solver to use.
|
protected |
Enum the software that is used to perform the factorization
|
protected |
KSP maximum number of iterations
|
protected |
The maximum numbers of allowable nonlinear iterations
|
protected |
Holds the Preconditioner object to be used for the linear solves.
|
protected |
Enum statitng with type of preconditioner to use.
|
protected |
Two differents tolerances on the residual for the resolution of non linear system
|
protected |
reuse jac level
|
protected |
reuse preconditioner level
|
protected |
KSP relative tolerance
|
protected |
Define the type of non linear solver
| map_dense_jacobian_function_type Feel::SolverNonLinear< T >::map_dense_jacobian |
Function that computes the Jacobian J(X) of the nonlinear system at the input iterate X using eigen.
| map_dense_matvec_function_type Feel::SolverNonLinear< T >::map_dense_matvec |
Function that computes either the residual
or the Jacobian
of the nonlinear system at the input iterate
. Note that either R or J could be XSNULL.
| map_dense_residual_function_type Feel::SolverNonLinear< T >::map_dense_residual |
Function that computes the residual R(X) of the nonlinear system at the input iterate X using eigen.
| matvec_function_type Feel::SolverNonLinear< T >::matvec |
Function that computes either the residual
or the Jacobian
of the nonlinear system at the input iterate
. Note that either R or J could be XSNULL.
| residual_function_type Feel::SolverNonLinear< T >::residual |
Function that computes the residual R(X) of the nonlinear system at the input iterate X.
1.8.5