Line data Source code
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 : #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" 26 : #include "libmesh/reference_counted_object.h" 27 : #include "libmesh/diff_solver.h" 28 : 29 : // C++ includes 30 : 31 : namespace libMesh 32 : { 33 : 34 : /** 35 : * This class defines a solver which uses the default 36 : * libMesh linear solver in a quasiNewton method to handle a 37 : * DifferentiableSystem 38 : * 39 : * This class is part of the new DifferentiableSystem framework, 40 : * which is still experimental. Users of this framework should 41 : * beware of bugs and future API changes. 42 : * 43 : * \author Roy H. Stogner 44 : * \date 2006 45 : */ 46 740 : class NewtonSolver : public DiffSolver 47 : { 48 : public: 49 : /** 50 : * Constructor. Requires a reference to the system 51 : * to be solved. 52 : */ 53 : explicit 54 : NewtonSolver (sys_type & system); 55 : 56 : /** 57 : * Destructor. 58 : */ 59 : virtual ~NewtonSolver (); 60 : 61 : typedef DiffSolver Parent; 62 : 63 : /** 64 : * The initialization function. This method is used to 65 : * initialize internal data structures before a simulation begins. 66 : */ 67 : virtual void init () override; 68 : 69 : /** 70 : * The reinitialization function. This method is used after 71 : * changes in the mesh. 72 : */ 73 : virtual void reinit () override; 74 : 75 : /** 76 : * This method performs a solve, using an inexact Newton-Krylov 77 : * method with line search. 78 : */ 79 : virtual unsigned int solve () override; 80 : 81 : LinearSolver<Number> & get_linear_solver() 82 : { libmesh_assert(_linear_solver); 83 : return *_linear_solver; 84 : } 85 : 86 : const LinearSolver<Number> & get_linear_solver() const 87 : { libmesh_assert(_linear_solver); 88 : return *_linear_solver; 89 : } 90 : 91 : /** 92 : * If this is set to true, the solver is forced to test the residual 93 : * after each Newton step, and to reduce the length of its steps 94 : * whenever necessary to avoid a residual increase. 95 : * It is currently set to true by default; set it to false to 96 : * avoid unnecessary residual assembly on well-behaved systems. 97 : */ 98 : bool require_residual_reduction; 99 : 100 : /** 101 : * If this is set to true, the solver is forced to test the residual 102 : * after each Newton step, and to reduce the length of its steps 103 : * whenever necessary to avoid an infinite or NaN residual. 104 : * It is currently set to true by default; set it to false to 105 : * avoid unnecessary residual assembly on well-behaved systems. 106 : */ 107 : bool require_finite_residual; 108 : 109 : /** 110 : * If require_residual_reduction is true, the solver may reduce step 111 : * lengths when required. If so, brent_line_search is an option. 112 : * If brent_line_search is set to false, the solver reduces the 113 : * length of its steps by 1/2 iteratively until it finds residual 114 : * reduction. If true, step lengths are first reduced by 1/2 or 115 : * more to find some residual reduction, then Brent's method is used 116 : * to find as much residual reduction as possible. 117 : * 118 : * brent_line_search is currently set to true by default. 119 : */ 120 : bool brent_line_search; 121 : 122 : /** 123 : * If set to true, check for convergence of the linear solve. If no 124 : * convergence is acquired during the linear solve, the nonlinear solve 125 : * fails with DiffSolver::DIVERGED_LINEAR_SOLVER_FAILURE. 126 : * Enabled by default as nonlinear convergence is very difficult, if the 127 : * linear solver is not converged. 128 : */ 129 : bool track_linear_convergence; 130 : 131 : /** 132 : * If the quasi-Newton step length must be reduced to below this 133 : * factor to give a residual reduction, then the Newton solver 134 : * dies with an error message. 135 : * It is currently set to 1e-5 by default. 136 : */ 137 : Real minsteplength; 138 : 139 : /** 140 : * The tolerance for linear solves is kept below this multiplier (which 141 : * defaults to 1e-3) times the norm of the current nonlinear residual 142 : */ 143 : double linear_tolerance_multiplier; 144 : 145 : protected: 146 : 147 : /** 148 : * The \p LinearSolver defines the interface used to 149 : * solve the linear_implicit system. This class handles all the 150 : * details of interfacing with various linear algebra packages 151 : * like PETSc or LASPACK. 152 : */ 153 : std::unique_ptr<LinearSolver<Number>> _linear_solver; 154 : 155 : /** 156 : * This does a line search in the direction opposite \p linear_solution 157 : * to try and minimize the residual of \p newton_iterate. 158 : * \p newton_iterate is moved to the end of the quasiNewton step. 159 : * 160 : * \returns The substep size. 161 : */ 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 : 168 : /** 169 : * This prints output for the convergence criteria based on 170 : * by the given residual and step size. 171 : */ 172 : void print_convergence(unsigned int step_num, 173 : Real current_residual, 174 : Real step_norm, 175 : bool linear_solve_finished); 176 : 177 : /** 178 : * \returns \p true if a convergence criterion has been passed 179 : * by the given residual and step size; false otherwise. 180 : */ 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