libMesh
nonlinear_implicit_system.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 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 // Local includes
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"
29 
30 namespace libMesh
31 {
32 
34  const std::string & name_in,
35  const unsigned int number_in) :
36 
37  Parent (es, name_in, number_in),
38  nonlinear_solver (NonlinearSolver<Number>::build(*this)),
39  diff_solver (),
40  _n_nonlinear_iterations (0),
41  _final_nonlinear_residual (1.e20)
42 {
43  // Set default parameters
44  // These were chosen to match the Petsc defaults
45  es.parameters.set<Real> ("linear solver tolerance") = 1e-5;
46  es.parameters.set<Real> ("linear solver minimum tolerance") = 1e-5;
47  es.parameters.set<unsigned int>("linear solver maximum iterations") = 10000;
48 
49  es.parameters.set<unsigned int>("nonlinear solver maximum iterations") = 50;
50  es.parameters.set<unsigned int>("nonlinear solver maximum function evaluations") = 10000;
51 
52  es.parameters.set<Real>("nonlinear solver absolute residual tolerance") = 1e-35;
53  es.parameters.set<Real>("nonlinear solver relative residual tolerance") = 1e-8;
54  es.parameters.set<Real>("nonlinear solver divergence tolerance") = 1e+4;
55  es.parameters.set<Real>("nonlinear solver absolute step tolerance") = 1e-8;
56  es.parameters.set<Real>("nonlinear solver relative step tolerance") = 1e-8;
57 
58  es.parameters.set<bool>("reuse preconditioner") = false;
59  es.parameters.set<unsigned int>("reuse preconditioner maximum linear iterations") = 1;
60 
61  if (this->has_static_condensation())
63 }
64 
65 
66 
68 
69 
70 
72 {
75 }
76 
77 
78 
80 {
81  // clear the nonlinear solver
82  nonlinear_solver->clear();
83 
84  // FIXME - this is necessary for petsc_auto_fieldsplit
85  // nonlinear_solver->init_names(*this);
86 
87  // clear the parent data
88  Parent::clear();
89 }
90 
91 
92 
94 {
95  // re-initialize the nonlinear solver interface
96  nonlinear_solver->clear();
97 
98  // force the solver to get a new preconditioner, in
99  // case reuse was set
100  nonlinear_solver->force_new_preconditioner();
101 
102  // FIXME - this is necessary for petsc_auto_fieldsplit
103  // nonlinear_solver->init_names(*this);
104 
105  if (diff_solver.get())
106  diff_solver->reinit();
107 
108  // initialize parent data
109  Parent::reinit();
110 }
111 
112 
113 
115 {
116  // Get a reference to the EquationSystems
117  const EquationSystems & es =
118  this->get_equation_systems();
119 
120  // Get the user-specified nonlinear solver tolerances
121  const unsigned int maxits = parameters.have_parameter<unsigned int>("nonlinear solver maximum iterations") ?
122  parameters.get<unsigned int>("nonlinear solver maximum iterations") :
123  es.parameters.get<unsigned int>("nonlinear solver maximum iterations");
124 
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");
128 
129  const double abs_resid_tol = parameters.have_parameter<Real>("nonlinear solver absolute residual tolerance") ?
130  double(parameters.get<Real>("nonlinear solver absolute residual tolerance")) :
131  double(es.parameters.get<Real>("nonlinear solver absolute residual tolerance"));
132 
133  const double rel_resid_tol = parameters.have_parameter<Real>("nonlinear solver relative residual tolerance") ?
134  double(parameters.get<Real>("nonlinear solver relative residual tolerance")) :
135  double(es.parameters.get<Real>("nonlinear solver relative residual tolerance"));
136 
137  const double div_tol = parameters.have_parameter<Real>("nonlinear solver divergence tolerance") ?
138  double(parameters.get<Real>("nonlinear solver divergence tolerance")) :
139  double(es.parameters.get<Real>("nonlinear solver divergence tolerance"));
140 
141  const double abs_step_tol = parameters.have_parameter<Real>("nonlinear solver absolute step tolerance") ?
142  double(parameters.get<Real>("nonlinear solver absolute step tolerance")) :
143  double(es.parameters.get<Real>("nonlinear solver absolute step tolerance"));
144 
145  const double rel_step_tol = parameters.have_parameter<Real>("nonlinear solver relative step tolerance")?
146  double(parameters.get<Real>("nonlinear solver relative step tolerance")) :
147  double(es.parameters.get<Real>("nonlinear solver relative step tolerance"));
148 
149  // Get the user-specified linear solver tolerances
150  const auto [maxlinearits, linear_tol] = this->Parent::get_linear_solve_parameters();
151 
152  const double linear_min_tol = parameters.have_parameter<Real>("linear solver minimum tolerance") ?
153  double(parameters.get<Real>("linear solver minimum tolerance")) :
154  double(es.parameters.get<Real>("linear solver minimum tolerance"));
155 
156  const bool reuse_preconditioner = parameters.have_parameter<unsigned int>("reuse preconditioner") ?
157  parameters.get<unsigned int>("reuse preconditioner") :
158  es.parameters.get<bool>("reuse preconditioner");
159  const unsigned int reuse_preconditioner_max_linear_its =
160  parameters.have_parameter<unsigned int>("reuse preconditioner maximum linear iterations") ?
161  parameters.get<unsigned int>("reuse preconditioner maximum linear iterations") :
162  es.parameters.get<unsigned int>("reuse preconditioner maximum linear iterations");
163 
164  // Set all the parameters on the NonlinearSolver
165  nonlinear_solver->max_nonlinear_iterations = maxits;
166  nonlinear_solver->max_function_evaluations = maxfuncs;
167  nonlinear_solver->absolute_residual_tolerance = abs_resid_tol;
168  nonlinear_solver->relative_residual_tolerance = rel_resid_tol;
169  nonlinear_solver->divergence_tolerance = div_tol;
170  nonlinear_solver->absolute_step_tolerance = abs_step_tol;
171  nonlinear_solver->relative_step_tolerance = rel_step_tol;
172  nonlinear_solver->max_linear_iterations = maxlinearits;
173  nonlinear_solver->initial_linear_tolerance = linear_tol;
174  nonlinear_solver->minimum_linear_tolerance = linear_min_tol;
175  nonlinear_solver->set_reuse_preconditioner(reuse_preconditioner);
176  nonlinear_solver->set_reuse_preconditioner_max_linear_its(reuse_preconditioner_max_linear_its);
177 
178  if (diff_solver.get())
179  {
180  diff_solver->max_nonlinear_iterations = maxits;
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;
185  diff_solver->max_linear_iterations = maxlinearits;
186  diff_solver->initial_linear_tolerance = linear_tol;
187  diff_solver->minimum_linear_tolerance = linear_min_tol;
188  }
189 }
190 
191 
192 
194 {
195  // Log how long the nonlinear solve takes.
196  LOG_SCOPE("solve()", "System");
197 
198  this->set_solver_parameters();
199 
200  if (diff_solver.get())
201  {
202  diff_solver->solve();
203 
204  // Store the number of nonlinear iterations required to
205  // solve and the final residual.
206  _n_nonlinear_iterations = diff_solver->total_outer_iterations();
207  _final_nonlinear_residual = 0.; // FIXME - support this!
208  }
209  else
210  {
211  if (this->prefix_with_name())
212  nonlinear_solver->init(this->prefix().c_str());
213  else
214  nonlinear_solver->init();
215 
216  // FIXME - this is necessary for petsc_auto_fieldsplit
217  // nonlinear_solver->init_names(*this);
218 
219  // Solve the nonlinear system.
220  // Store the number of nonlinear iterations required to
221  // solve and the final residual.
223  nonlinear_solver->solve (*matrix, *solution, *rhs,
224  nonlinear_solver->relative_residual_tolerance,
225  nonlinear_solver->max_linear_iterations);
226  }
227 
228  // Update the system after the solve
229  this->update();
230 }
231 
232 
233 
234 std::pair<unsigned int, Real> NonlinearImplicitSystem::get_linear_solve_parameters() const
235 {
236  if (diff_solver.get())
237  return std::make_pair(this->diff_solver->max_linear_iterations,
238  this->diff_solver->relative_residual_tolerance);
239  return std::make_pair(this->nonlinear_solver->max_linear_iterations,
240  this->nonlinear_solver->relative_residual_tolerance);
241 }
242 
243 
244 
245 void NonlinearImplicitSystem::assembly(bool get_residual,
246  bool get_jacobian,
247  bool /*apply_heterogeneous_constraints*/,
248  bool /*apply_no_constraints*/)
249 {
250  // Get current_local_solution in sync
251  this->update();
252 
253  //-----------------------------------------------------------------------------
254  // if the user has provided both function pointers and objects only the pointer
255  // will be used, so catch that as an error
256  libmesh_error_msg_if(nonlinear_solver->jacobian && nonlinear_solver->jacobian_object,
257  "ERROR: cannot specify both a function and object to compute the Jacobian!");
258 
259  libmesh_error_msg_if(nonlinear_solver->residual && nonlinear_solver->residual_object,
260  "ERROR: cannot specify both a function and object to compute the Residual!");
261 
262  libmesh_error_msg_if(nonlinear_solver->matvec && nonlinear_solver->residual_and_jacobian_object,
263  "ERROR: cannot specify both a function and object to compute the combined Residual & Jacobian!");
264 
265 
266  if (get_jacobian)
267  {
268  if (nonlinear_solver->jacobian != nullptr)
269  nonlinear_solver->jacobian (*current_local_solution.get(), *matrix, *this);
270 
271  else if (nonlinear_solver->jacobian_object != nullptr)
272  nonlinear_solver->jacobian_object->jacobian (*current_local_solution.get(), *matrix, *this);
273 
274  else if (nonlinear_solver->matvec != nullptr)
275  nonlinear_solver->matvec (*current_local_solution.get(), get_residual ? rhs : nullptr, matrix, *this);
276 
277  else if (nonlinear_solver->residual_and_jacobian_object != nullptr)
278  nonlinear_solver->residual_and_jacobian_object->residual_and_jacobian (*current_local_solution.get(), get_residual ? rhs : nullptr, matrix, *this);
279 
280  else
281  libmesh_error_msg("Error! Unable to compute residual and/or Jacobian!");
282  }
283 
284  if (get_residual)
285  {
286  if (nonlinear_solver->residual != nullptr)
287  nonlinear_solver->residual (*current_local_solution.get(), *rhs, *this);
288 
289  else if (nonlinear_solver->residual_object != nullptr)
290  nonlinear_solver->residual_object->residual (*current_local_solution.get(), *rhs, *this);
291 
292  else if (nonlinear_solver->matvec != nullptr)
293  {
294  // we might have already grabbed the residual and jacobian together
295  if (!get_jacobian)
296  nonlinear_solver->matvec (*current_local_solution.get(), rhs, nullptr, *this);
297  }
298 
299  else if (nonlinear_solver->residual_and_jacobian_object != nullptr)
300  {
301  // we might have already grabbed the residual and jacobian together
302  if (!get_jacobian)
303  nonlinear_solver->residual_and_jacobian_object->residual_and_jacobian (*current_local_solution.get(), rhs, nullptr, *this);
304  }
305 
306  else
307  libmesh_error_msg("Error! Unable to compute residual and/or Jacobian!");
308  }
309  else
310  libmesh_assert(get_jacobian); // I can't believe you really wanted to assemble *nothing*
311 }
312 
313 
314 
315 
317 {
318  return nonlinear_solver->get_current_nonlinear_iteration_number();
319 }
320 
321 
322 
323 } // namespace libMesh
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
Definition: parameters.h:395
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.
Definition: system.C:454
const EquationSystems & get_equation_systems() const
Definition: system.h:721
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
Definition: system.C:2871
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...
Definition: system.h:1526
void setup_static_condensation_preconditioner(T &solver)
Sets up the static condensation preconditioner for the supplied solver.
const T & get(std::string_view) const
Definition: parameters.h:426
unsigned int _n_nonlinear_iterations
The number of nonlinear iterations required to solve the nonlinear system R(x)=0. ...
std::string prefix() const
Definition: system.h:1938
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
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.
Definition: system.h:1593
libmesh_assert(ctx)
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.
Definition: system.C:510
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 &)
Definition: parameters.h:469
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.
Definition: system.h:1605
virtual std::pair< unsigned int, Real > get_linear_solve_parameters() const override
bool prefix_with_name() const
Definition: system.h:1932