Go to the documentation of this file.
   23 #include "libmesh/nonlinear_implicit_system.h" 
   24 #include "libmesh/diff_solver.h" 
   25 #include "libmesh/equation_systems.h" 
   26 #include "libmesh/libmesh_logging.h" 
   27 #include "libmesh/nonlinear_solver.h" 
   28 #include "libmesh/sparse_matrix.h" 
   36                                                   const std::string & name_in,
 
   37                                                   const unsigned int number_in) :
 
   39   Parent                    (es, name_in, number_in),
 
   42   _n_nonlinear_iterations   (0),
 
   43   _final_nonlinear_residual (1.e20)
 
   49   es.
parameters.
set<
unsigned int>(
"linear solver maximum iterations") = 10000;
 
   51   es.
parameters.
set<
unsigned int>(
"nonlinear solver maximum iterations") = 50;
 
   52   es.
parameters.
set<
unsigned int>(
"nonlinear solver maximum function evaluations") = 10000;
 
  103   const unsigned int maxits =
 
  104     es.
parameters.
get<
unsigned int>(
"nonlinear solver maximum iterations");
 
  106   const unsigned int maxfuncs =
 
  107     es.
parameters.
get<
unsigned int>(
"nonlinear solver maximum function evaluations");
 
  109   const double abs_resid_tol =
 
  112   const double rel_resid_tol =
 
  115   const double div_tol =
 
  118   const double abs_step_tol =
 
  121   const double rel_step_tol =
 
  125   const unsigned int maxlinearits =
 
  126     es.
parameters.
get<
unsigned int>(
"linear solver maximum iterations");
 
  128   const double linear_tol =
 
  131   const double linear_min_tol =
 
  149       diff_solver->absolute_residual_tolerance = abs_resid_tol;
 
  150       diff_solver->relative_residual_tolerance = rel_resid_tol;
 
  151       diff_solver->absolute_step_tolerance = abs_step_tol;
 
  152       diff_solver->relative_step_tolerance = rel_step_tol;
 
  154       diff_solver->initial_linear_tolerance = linear_tol;
 
  155       diff_solver->minimum_linear_tolerance = linear_min_tol;
 
  164   START_LOG(
"solve()", 
"System");
 
  185       const std::pair<unsigned int, Real> rval =
 
  197   STOP_LOG(
"solve()", 
"System");
 
  208     return std::make_pair(this->
diff_solver->max_linear_iterations,
 
  209                           this->diff_solver->relative_residual_tolerance);
 
  211                         this->nonlinear_solver->relative_residual_tolerance);
 
  228     libmesh_error_msg(
"ERROR: cannot specify both a function and object to compute the Jacobian!");
 
  231     libmesh_error_msg(
"ERROR: cannot specify both a function and object to compute the Residual!");
 
  234     libmesh_error_msg(
"ERROR: cannot specify both a function and object to compute the combined Residual & Jacobian!");
 
  252         libmesh_error_msg(
"Error! Unable to compute residual and/or Jacobian!");
 
  278         libmesh_error_msg(
"Error! Unable to compute residual and/or Jacobian!");
 
  
Manages consistently variables, degrees of freedom, and coefficient vectors.
 
virtual ~NonlinearImplicitSystem()
Destructor.
 
const EquationSystems & get_equation_systems() const
 
Real _final_nonlinear_residual
The final residual for the nonlinear system R(x)
 
NumericVector< Number > * rhs
The system matrix.
 
This base class can be inherited from to provide interfaces to nonlinear solvers from different packa...
 
virtual void reinit() override
Reinitializes the member data fields associated with the system, so that, e.g., assemble() may be use...
 
virtual void reinit() override
Reinitializes the member data fields associated with the system, so that, e.g., assemble() may be use...
 
The libMesh namespace provides an interface to certain functionality in the library.
 
unsigned int _n_nonlinear_iterations
The number of nonlinear iterations required to solve the nonlinear system R(x)=0.
 
virtual void solve() override
Assembles & solves the nonlinear system R(x) = 0.
 
std::unique_ptr< DiffSolver > diff_solver
The DiffSolver defines an optional interface used to solve the nonlinear_implicit system.
 
void set_solver_parameters()
Copies system parameters into nonlinear solver parameters.
 
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.
 
virtual void clear() override
Clear all the data structures associated with the system.
 
virtual std::pair< unsigned int, Real > get_linear_solve_parameters() const override
 
std::unique_ptr< NumericVector< Number > > current_local_solution
All the values I need to compute my contribution to the simulation at hand.
 
SparseMatrix< Number > * matrix
The system matrix.
 
This is the EquationSystems class.
 
T & set(const std::string &)
 
std::unique_ptr< NonlinearSolver< Number > > nonlinear_solver
The NonlinearSolver defines the default interface used to solve the nonlinear_implicit system.
 
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
 
const std::string & name() const
 
bool on_command_line(std::string arg)
 
NonlinearImplicitSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
unsigned get_current_nonlinear_iteration_number() const
If called during the solve(), for example by the user-specified residual or Jacobian function,...
 
const T & get(const std::string &) const
 
virtual void clear() override
Clear all the data structures associated with the system.
 
virtual void update()
Update the local values to reflect the solution on neighboring processors.
 
Parameters parameters
Data structure holding arbitrary parameters.