newton_solver.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 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 #ifndef LIBMESH_NEWTON_SOLVER_H
21 #define LIBMESH_NEWTON_SOLVER_H
22 
23 // Local includes
24 #include "libmesh/libmesh_common.h"
25 #include "libmesh/linear_solver.h"
27 #include "libmesh/diff_solver.h"
28 
29 // C++ includes
30 
31 namespace libMesh
32 {
33 
46 class NewtonSolver : public DiffSolver
47 {
48 public:
53  explicit
55 
59  virtual ~NewtonSolver ();
60 
61  typedef DiffSolver Parent;
62 
67  virtual void init () override;
68 
73  virtual void reinit () override;
74 
79  virtual unsigned int solve () override;
80 
82  { libmesh_assert(_linear_solver);
83  return *_linear_solver;
84  }
85 
87  { libmesh_assert(_linear_solver);
88  return *_linear_solver;
89  }
90 
99 
108 
121 
130 
138 
144 
145 protected:
146 
153  std::unique_ptr<LinearSolver<Number>> _linear_solver;
154 
162  Real line_search(Real tol,
163  Real last_residual,
164  Real & current_residual,
165  NumericVector<Number> & newton_iterate,
166  const NumericVector<Number> & linear_solution);
167 
172  void print_convergence(unsigned int step_num,
173  Real current_residual,
174  Real step_norm,
175  bool linear_solve_finished);
176 
181  bool test_convergence(Real current_residual,
182  Real step_norm,
183  bool linear_solve_finished);
184 };
185 
186 
187 } // namespace libMesh
188 
189 
190 #endif // LIBMESH_NEWTON_SOLVER_H
virtual void init() override
NewtonSolver(sys_type &system)
LinearSolver< Number > & get_linear_solver()
Definition: newton_solver.h:81
std::unique_ptr< LinearSolver< Number > > _linear_solver
virtual unsigned int solve() override
const sys_type & system() const
Definition: diff_solver.h:138
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void reinit() override
bool test_convergence(Real current_residual, Real step_norm, bool linear_solve_finished)
const LinearSolver< Number > & get_linear_solver() const
Definition: newton_solver.h:86
Real line_search(Real tol, Real last_residual, Real &current_residual, NumericVector< Number > &newton_iterate, const NumericVector< Number > &linear_solution)
Definition: newton_solver.C:38
void print_convergence(unsigned int step_num, Real current_residual, Real step_norm, bool linear_solve_finished)
Manages consistently variables, degrees of freedom, coefficient vectors, and matrices for implicit sy...