21 #include "libmesh/nonlinear_implicit_system.h"    22 #include "libmesh/diff_solver.h"    23 #include "libmesh/equation_systems.h"    24 #include "libmesh/libmesh_logging.h"    25 #include "libmesh/nonlinear_solver.h"    26 #include "libmesh/sparse_matrix.h"    27 #include "libmesh/static_condensation.h"    28 #include "libmesh/static_condensation_preconditioner.h"    34                                                   const std::string & name_in,
    35                                                   const unsigned int number_in) :
    37   Parent                    (es, name_in, number_in),
    40   _n_nonlinear_iterations   (0),
    41   _final_nonlinear_residual (1.e20)
    47   es.
parameters.
set<
unsigned int>(
"linear solver maximum iterations") = 10000;
    49   es.
parameters.
set<
unsigned int>(
"nonlinear solver maximum iterations") = 50;
    50   es.
parameters.
set<
unsigned int>(
"nonlinear solver maximum function evaluations") = 10000;
    59   es.
parameters.
set<
unsigned int>(
"reuse preconditioner maximum linear iterations") = 1;
   122     parameters.
get<
unsigned int>(
"nonlinear solver maximum iterations") :
   123     es.
parameters.
get<
unsigned int>(
"nonlinear solver maximum iterations");
   125   const unsigned int maxfuncs = 
parameters.
have_parameter<
unsigned int>(
"nonlinear solver maximum function evaluations") ?
   126     parameters.
get<
unsigned int>(
"nonlinear solver maximum function evaluations") :
   127     es.
parameters.
get<
unsigned int>(
"nonlinear solver maximum function evaluations");
   159   const unsigned int reuse_preconditioner_max_linear_its =
   161     parameters.
get<
unsigned int>(
"reuse preconditioner maximum linear iterations") :
   162       es.
parameters.
get<
unsigned int>(
"reuse preconditioner maximum linear iterations");
   176   nonlinear_solver->set_reuse_preconditioner_max_linear_its(reuse_preconditioner_max_linear_its);
   181       diff_solver->absolute_residual_tolerance = abs_resid_tol;
   182       diff_solver->relative_residual_tolerance = rel_resid_tol;
   183       diff_solver->absolute_step_tolerance = abs_step_tol;
   184       diff_solver->relative_step_tolerance = rel_step_tol;
   186       diff_solver->initial_linear_tolerance = linear_tol;
   187       diff_solver->minimum_linear_tolerance = linear_min_tol;
   196   LOG_SCOPE(
"solve()", 
"System");
   237     return std::make_pair(this->
diff_solver->max_linear_iterations,
   238                           this->diff_solver->relative_residual_tolerance);
   240                         this->nonlinear_solver->relative_residual_tolerance);
   257                        "ERROR: cannot specify both a function and object to compute the Jacobian!");
   260                        "ERROR: cannot specify both a function and object to compute the Residual!");
   263                        "ERROR: cannot specify both a function and object to compute the combined Residual & Jacobian!");
   281         libmesh_error_msg(
"Error! Unable to compute residual and/or Jacobian!");
   307         libmesh_error_msg(
"Error! Unable to compute residual and/or Jacobian!");
 void set_solver_parameters()
Copies system parameters into nonlinear solver parameters. 
This is the EquationSystems class. 
std::unique_ptr< NonlinearSolver< Number > > nonlinear_solver
The NonlinearSolver defines the default interface used to solve the nonlinear_implicit system...
bool have_parameter(std::string_view) const
virtual void create_static_condensation() override
Request that static condensation be performed for this system. 
virtual void create_static_condensation() override
Request that static condensation be performed for this system. 
virtual std::pair< unsigned int, Real > get_linear_solve_parameters() const
virtual void reinit()
Reinitializes degrees of freedom and other required data on the current mesh. 
const EquationSystems & get_equation_systems() const
virtual void clear() override
Clear all the data structures associated with the system. 
NumericVector< Number > * rhs
The system matrix. 
This base class can be inherited from to provide interfaces to nonlinear solvers from different packa...
NonlinearImplicitSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor. 
virtual void clear() override
Clear all the data structures associated with the system. 
bool has_static_condensation() const
The libMesh namespace provides an interface to certain functionality in the library. 
Parameters parameters
Parameters for the system. If a parameter is not provided, it should be retrieved from the EquationSy...
void setup_static_condensation_preconditioner(T &solver)
Sets up the static condensation preconditioner for the supplied solver. 
const T & get(std::string_view) const
unsigned int _n_nonlinear_iterations
The number of nonlinear iterations required to solve the nonlinear system R(x)=0. ...
std::string prefix() const
virtual ~NonlinearImplicitSystem()
Manages consistently variables, degrees of freedom, and coefficient vectors. 
Real _final_nonlinear_residual
The final residual for the nonlinear system R(x) 
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values. 
std::unique_ptr< DiffSolver > diff_solver
The DiffSolver defines an optional interface used to solve the nonlinear_implicit system...
unsigned get_current_nonlinear_iteration_number() const
If called during the solve(), for example by the user-specified residual or Jacobian function...
virtual void solve() override
Assembles & solves the nonlinear system R(x) = 0. 
virtual void update()
Update the local values to reflect the solution on neighboring processors. 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void assembly(bool get_residual, bool get_jacobian, bool apply_heterogeneous_constraints=false, bool apply_no_constraints=false) override
Assembles a residual in rhs and/or a jacobian in matrix, as requested. 
T & set(const std::string &)
SparseMatrix< Number > * matrix
The system matrix. 
virtual void reinit() override
Reinitializes the member data fields associated with the system, so that, e.g., assemble() may be use...
Parameters parameters
Data structure holding arbitrary parameters. 
std::unique_ptr< NumericVector< Number > > current_local_solution
All the values I need to compute my contribution to the simulation at hand. 
virtual std::pair< unsigned int, Real > get_linear_solve_parameters() const override
bool prefix_with_name() const