libMesh
nonlinear_implicit_system.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 // C++ includes
21 
22 // Local includes
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"
29 
30 namespace libMesh
31 {
32 
33 // ------------------------------------------------------------
34 // NonlinearImplicitSystem implementation
36  const std::string & name_in,
37  const unsigned int number_in) :
38 
39  Parent (es, name_in, number_in),
40  nonlinear_solver (NonlinearSolver<Number>::build(*this)),
41  diff_solver (),
42  _n_nonlinear_iterations (0),
43  _final_nonlinear_residual (1.e20)
44 {
45  // Set default parameters
46  // These were chosen to match the Petsc defaults
47  es.parameters.set<Real> ("linear solver tolerance") = 1e-5;
48  es.parameters.set<Real> ("linear solver minimum tolerance") = 1e-5;
49  es.parameters.set<unsigned int>("linear solver maximum iterations") = 10000;
50 
51  es.parameters.set<unsigned int>("nonlinear solver maximum iterations") = 50;
52  es.parameters.set<unsigned int>("nonlinear solver maximum function evaluations") = 10000;
53 
54  es.parameters.set<Real>("nonlinear solver absolute residual tolerance") = 1e-35;
55  es.parameters.set<Real>("nonlinear solver relative residual tolerance") = 1e-8;
56  es.parameters.set<Real>("nonlinear solver divergence tolerance") = 1e+4;
57  es.parameters.set<Real>("nonlinear solver absolute step tolerance") = 1e-8;
58  es.parameters.set<Real>("nonlinear solver relative step tolerance") = 1e-8;
59 }
60 
61 
62 
64 {
65  // Clear data
66  this->clear();
67 }
68 
69 
70 
72 {
73  // clear the nonlinear solver
74  nonlinear_solver->clear();
75 
76  // clear the parent data
77  Parent::clear();
78 }
79 
80 
81 
83 {
84  // re-initialize the nonlinear solver interface
85  nonlinear_solver->clear();
86 
87  if (diff_solver.get())
88  diff_solver->reinit();
89 
90  // initialize parent data
92 }
93 
94 
95 
97 {
98  // Get a reference to the EquationSystems
99  const EquationSystems & es =
100  this->get_equation_systems();
101 
102  // Get the user-specified nonlinear solver tolerances
103  const unsigned int maxits =
104  es.parameters.get<unsigned int>("nonlinear solver maximum iterations");
105 
106  const unsigned int maxfuncs =
107  es.parameters.get<unsigned int>("nonlinear solver maximum function evaluations");
108 
109  const double abs_resid_tol =
110  double(es.parameters.get<Real>("nonlinear solver absolute residual tolerance"));
111 
112  const double rel_resid_tol =
113  double(es.parameters.get<Real>("nonlinear solver relative residual tolerance"));
114 
115  const double div_tol =
116  double(es.parameters.get<Real>("nonlinear solver divergence tolerance"));
117 
118  const double abs_step_tol =
119  double(es.parameters.get<Real>("nonlinear solver absolute step tolerance"));
120 
121  const double rel_step_tol =
122  double(es.parameters.get<Real>("nonlinear solver relative step tolerance"));
123 
124  // Get the user-specified linear solver tolerances
125  const unsigned int maxlinearits =
126  es.parameters.get<unsigned int>("linear solver maximum iterations");
127 
128  const double linear_tol =
129  double(es.parameters.get<Real>("linear solver tolerance"));
130 
131  const double linear_min_tol =
132  double(es.parameters.get<Real>("linear solver minimum tolerance"));
133 
134  // Set all the parameters on the NonlinearSolver
135  nonlinear_solver->max_nonlinear_iterations = maxits;
136  nonlinear_solver->max_function_evaluations = maxfuncs;
137  nonlinear_solver->absolute_residual_tolerance = abs_resid_tol;
138  nonlinear_solver->relative_residual_tolerance = rel_resid_tol;
139  nonlinear_solver->divergence_tolerance = div_tol;
140  nonlinear_solver->absolute_step_tolerance = abs_step_tol;
141  nonlinear_solver->relative_step_tolerance = rel_step_tol;
142  nonlinear_solver->max_linear_iterations = maxlinearits;
143  nonlinear_solver->initial_linear_tolerance = linear_tol;
144  nonlinear_solver->minimum_linear_tolerance = linear_min_tol;
145 
146  if (diff_solver.get())
147  {
148  diff_solver->max_nonlinear_iterations = maxits;
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;
153  diff_solver->max_linear_iterations = maxlinearits;
154  diff_solver->initial_linear_tolerance = linear_tol;
155  diff_solver->minimum_linear_tolerance = linear_min_tol;
156  }
157 }
158 
159 
160 
162 {
163  // Log how long the nonlinear solve takes.
164  START_LOG("solve()", "System");
165 
166  this->set_solver_parameters();
167 
168  if (diff_solver.get())
169  {
170  diff_solver->solve();
171 
172  // Store the number of nonlinear iterations required to
173  // solve and the final residual.
174  _n_nonlinear_iterations = diff_solver->total_outer_iterations();
175  _final_nonlinear_residual = 0.; // FIXME - support this!
176  }
177  else
178  {
179  if (libMesh::on_command_line("--solver-system-names"))
180  nonlinear_solver->init((this->name()+"_").c_str());
181  else
182  nonlinear_solver->init();
183 
184  // Solve the nonlinear system.
185  const std::pair<unsigned int, Real> rval =
186  nonlinear_solver->solve (*matrix, *solution, *rhs,
187  nonlinear_solver->relative_residual_tolerance,
188  nonlinear_solver->max_linear_iterations);
189 
190  // Store the number of nonlinear iterations required to
191  // solve and the final residual.
192  _n_nonlinear_iterations = rval.first;
193  _final_nonlinear_residual = rval.second;
194  }
195 
196  // Stop logging the nonlinear solve
197  STOP_LOG("solve()", "System");
198 
199  // Update the system after the solve
200  this->update();
201 }
202 
203 
204 
205 std::pair<unsigned int, Real> NonlinearImplicitSystem::get_linear_solve_parameters() const
206 {
207  if (diff_solver.get())
208  return std::make_pair(this->diff_solver->max_linear_iterations,
209  this->diff_solver->relative_residual_tolerance);
210  return std::make_pair(this->nonlinear_solver->max_linear_iterations,
211  this->nonlinear_solver->relative_residual_tolerance);
212 }
213 
214 
215 
216 void NonlinearImplicitSystem::assembly(bool get_residual,
217  bool get_jacobian,
218  bool /*apply_heterogeneous_constraints*/,
219  bool /*apply_no_constraints*/)
220 {
221  // Get current_local_solution in sync
222  this->update();
223 
224  //-----------------------------------------------------------------------------
225  // if the user has provided both function pointers and objects only the pointer
226  // will be used, so catch that as an error
227  if (nonlinear_solver->jacobian && nonlinear_solver->jacobian_object)
228  libmesh_error_msg("ERROR: cannot specify both a function and object to compute the Jacobian!");
229 
230  if (nonlinear_solver->residual && nonlinear_solver->residual_object)
231  libmesh_error_msg("ERROR: cannot specify both a function and object to compute the Residual!");
232 
233  if (nonlinear_solver->matvec && nonlinear_solver->residual_and_jacobian_object)
234  libmesh_error_msg("ERROR: cannot specify both a function and object to compute the combined Residual & Jacobian!");
235 
236 
237  if (get_jacobian)
238  {
239  if (nonlinear_solver->jacobian != nullptr)
240  nonlinear_solver->jacobian (*current_local_solution.get(), *matrix, *this);
241 
242  else if (nonlinear_solver->jacobian_object != nullptr)
243  nonlinear_solver->jacobian_object->jacobian (*current_local_solution.get(), *matrix, *this);
244 
245  else if (nonlinear_solver->matvec != nullptr)
246  nonlinear_solver->matvec (*current_local_solution.get(), get_residual ? rhs : nullptr, matrix, *this);
247 
248  else if (nonlinear_solver->residual_and_jacobian_object != nullptr)
249  nonlinear_solver->residual_and_jacobian_object->residual_and_jacobian (*current_local_solution.get(), get_residual ? rhs : nullptr, matrix, *this);
250 
251  else
252  libmesh_error_msg("Error! Unable to compute residual and/or Jacobian!");
253  }
254 
255  if (get_residual)
256  {
257  if (nonlinear_solver->residual != nullptr)
258  nonlinear_solver->residual (*current_local_solution.get(), *rhs, *this);
259 
260  else if (nonlinear_solver->residual_object != nullptr)
261  nonlinear_solver->residual_object->residual (*current_local_solution.get(), *rhs, *this);
262 
263  else if (nonlinear_solver->matvec != nullptr)
264  {
265  // we might have already grabbed the residual and jacobian together
266  if (!get_jacobian)
267  nonlinear_solver->matvec (*current_local_solution.get(), rhs, nullptr, *this);
268  }
269 
270  else if (nonlinear_solver->residual_and_jacobian_object != nullptr)
271  {
272  // we might have already grabbed the residual and jacobian together
273  if (!get_jacobian)
274  nonlinear_solver->residual_and_jacobian_object->residual_and_jacobian (*current_local_solution.get(), rhs, nullptr, *this);
275  }
276 
277  else
278  libmesh_error_msg("Error! Unable to compute residual and/or Jacobian!");
279  }
280  else
281  libmesh_assert(get_jacobian); // I can't believe you really wanted to assemble *nothing*
282 }
283 
284 
285 
286 
288 {
289  return nonlinear_solver->get_current_nonlinear_iteration_number();
290 }
291 
292 
293 
294 } // namespace libMesh
libMesh::System
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:100
libMesh::NonlinearImplicitSystem::~NonlinearImplicitSystem
virtual ~NonlinearImplicitSystem()
Destructor.
Definition: nonlinear_implicit_system.C:63
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::System::get_equation_systems
const EquationSystems & get_equation_systems() const
Definition: system.h:720
libMesh::NonlinearImplicitSystem::_final_nonlinear_residual
Real _final_nonlinear_residual
The final residual for the nonlinear system R(x)
Definition: nonlinear_implicit_system.h:304
libMesh::ExplicitSystem::rhs
NumericVector< Number > * rhs
The system matrix.
Definition: explicit_system.h:114
libMesh::NonlinearSolver
This base class can be inherited from to provide interfaces to nonlinear solvers from different packa...
Definition: nonlinear_solver.h:60
libMesh::NonlinearImplicitSystem::reinit
virtual void reinit() override
Reinitializes the member data fields associated with the system, so that, e.g., assemble() may be use...
Definition: nonlinear_implicit_system.C:82
libMesh::ImplicitSystem::reinit
virtual void reinit() override
Reinitializes the member data fields associated with the system, so that, e.g., assemble() may be use...
Definition: implicit_system.C:149
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::NonlinearImplicitSystem::_n_nonlinear_iterations
unsigned int _n_nonlinear_iterations
The number of nonlinear iterations required to solve the nonlinear system R(x)=0.
Definition: nonlinear_implicit_system.h:299
libMesh::NonlinearImplicitSystem::solve
virtual void solve() override
Assembles & solves the nonlinear system R(x) = 0.
Definition: nonlinear_implicit_system.C:161
libMesh::NonlinearImplicitSystem::diff_solver
std::unique_ptr< DiffSolver > diff_solver
The DiffSolver defines an optional interface used to solve the nonlinear_implicit system.
Definition: nonlinear_implicit_system.h:267
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::NonlinearImplicitSystem::set_solver_parameters
void set_solver_parameters()
Copies system parameters into nonlinear solver parameters.
Definition: nonlinear_implicit_system.C:96
libMesh::NonlinearImplicitSystem::assembly
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.
Definition: nonlinear_implicit_system.C:216
libMesh::ImplicitSystem::clear
virtual void clear() override
Clear all the data structures associated with the system.
Definition: implicit_system.C:66
libMesh::NonlinearImplicitSystem::get_linear_solve_parameters
virtual std::pair< unsigned int, Real > get_linear_solve_parameters() const override
Definition: nonlinear_implicit_system.C:205
libMesh::System::current_local_solution
std::unique_ptr< NumericVector< Number > > current_local_solution
All the values I need to compute my contribution to the simulation at hand.
Definition: system.h:1551
libMesh::ImplicitSystem::matrix
SparseMatrix< Number > * matrix
The system matrix.
Definition: implicit_system.h:393
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::Parameters::set
T & set(const std::string &)
Definition: parameters.h:460
libMesh::NonlinearImplicitSystem::nonlinear_solver
std::unique_ptr< NonlinearSolver< Number > > nonlinear_solver
The NonlinearSolver defines the default interface used to solve the nonlinear_implicit system.
Definition: nonlinear_implicit_system.h:261
libMesh::System::solution
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1539
libMesh::System::name
const std::string & name() const
Definition: system.h:2067
libMesh::on_command_line
bool on_command_line(std::string arg)
Definition: libmesh.C:898
libMesh::NonlinearImplicitSystem::NonlinearImplicitSystem
NonlinearImplicitSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
Definition: nonlinear_implicit_system.C:35
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::NonlinearImplicitSystem::get_current_nonlinear_iteration_number
unsigned get_current_nonlinear_iteration_number() const
If called during the solve(), for example by the user-specified residual or Jacobian function,...
Definition: nonlinear_implicit_system.C:287
libMesh::Parameters::get
const T & get(const std::string &) const
Definition: parameters.h:421
libMesh::NonlinearImplicitSystem::clear
virtual void clear() override
Clear all the data structures associated with the system.
Definition: nonlinear_implicit_system.C:71
libMesh::System::update
virtual void update()
Update the local values to reflect the solution on neighboring processors.
Definition: system.C:408
libMesh::EquationSystems::parameters
Parameters parameters
Data structure holding arbitrary parameters.
Definition: equation_systems.h:557