OGS 6.2.1-499-g3b941532c.dirty.20191012113459
NumLib::NonlinearSolver< NonlinearSolverTag::Newton > Class Template Referencefinal

Detailed Description

template<>
class NumLib::NonlinearSolver< NonlinearSolverTag::Newton >

Find a solution to a nonlinear equation using the Newton-Raphson method.

Definition at line 83 of file NonlinearSolver.h.

#include <NonlinearSolver.h>

Inheritance diagram for NumLib::NonlinearSolver< NonlinearSolverTag::Newton >:
Collaboration diagram for NumLib::NonlinearSolver< NonlinearSolverTag::Newton >:

Public Types

using System = NonlinearSystem< NonlinearSolverTag::Newton >
 Type of the nonlinear equation system to be solved. More...
 

Public Member Functions

 NonlinearSolver (GlobalLinearSolver &linear_solver, int const maxiter, double const damping=1.0)
 
void setEquationSystem (System &eq, ConvergenceCriterion &conv_crit)
 
void assemble (std::vector< GlobalVector *> const &x, int const process_id) const override
 
NonlinearSolverStatus solve (std::vector< GlobalVector *> &x, std::function< void(int, GlobalVector const &)> const &postIterationCallback, int const process_id) override
 
- Public Member Functions inherited from NumLib::NonlinearSolverBase
virtual ~NonlinearSolverBase ()=default
 

Private Attributes

GlobalLinearSolver & _linear_solver
 
System_equation_system = nullptr
 
ConvergenceCriterion_convergence_criterion = nullptr
 
int const _maxiter
 maximum number of iterations More...
 
double const _damping
 
std::size_t _res_id = 0u
 ID of the residual vector. More...
 
std::size_t _J_id = 0u
 ID of the Jacobian matrix. More...
 
std::size_t _minus_delta_x_id = 0u
 ID of the $ -\Delta x$ vector. More...
 
std::size_t _x_new_id
 ID of the vector storing $ x - (-\Delta x) $. More...
 

Member Typedef Documentation

◆ System

Type of the nonlinear equation system to be solved.

Definition at line 88 of file NonlinearSolver.h.

Constructor & Destructor Documentation

◆ NonlinearSolver()

NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::NonlinearSolver ( GlobalLinearSolver &  linear_solver,
int const  maxiter,
double const  damping = 1.0 
)
inlineexplicit

Constructs a new instance.

Parameters
linear_solverthe linear solver used by this nonlinear solver.
maxiterthe maximum number of iterations used to solve the equation.
dampingA positive damping factor.
See also
_damping

Definition at line 98 of file NonlinearSolver.h.

101  : _linear_solver(linear_solver), _maxiter(maxiter), _damping(damping)
102  {
103  }
int const _maxiter
maximum number of iterations

Member Function Documentation

◆ assemble()

void NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::assemble ( std::vector< GlobalVector *> const &  x,
int const  process_id 
) const
overridevirtual

Only assemble the equation system.

Note
This method is needed to preload CrankNicolson time discretization scheme. It is not used for the general solver steps; in those only the solve() method is needed.
Parameters
xthe state at which the equation system will be assembled.
process_idthe process' id which will be assembled.

Implements NumLib::NonlinearSolverBase.

Definition at line 185 of file NonlinearSolver.cpp.

187 {
188  _equation_system->assemble(x, process_id);
189  // TODO if the equation system would be reset to nullptr after each
190  // assemble() or solve() call, the user would be forced to set the
191  // equation every time and could not forget it.
192 }
virtual void assemble(std::vector< GlobalVector *> const &x, int const process_id)=0

◆ setEquationSystem()

void NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::setEquationSystem ( System eq,
ConvergenceCriterion conv_crit 
)
inline

Set the nonlinear equation system that will be solved. TODO doc

Definition at line 107 of file NonlinearSolver.h.

◆ solve()

NonlinearSolverStatus NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::solve ( std::vector< GlobalVector *> &  x,
std::function< void(int, GlobalVector const &)> const &  postIterationCallback,
int const  process_id 
)
overridevirtual

Assemble and solve the equation system.

Parameters
xin: the initial guess, out: the solution.
postIterationCallbackcalled after each iteration if set.
process_idusually used in staggered schemes.
Return values
trueif the equation system could be solved
falseotherwise

Implements NumLib::NonlinearSolverBase.

Definition at line 194 of file NonlinearSolver.cpp.

References MathLib::LinAlg::axpy(), MathLib::LinAlg::copy(), BaseLib::RunTime::elapsed(), NumLib::FAILURE, NumLib::MatrixProvider::getMatrix(), NumLib::VectorProvider::getVector(), NumLib::GlobalVectorProvider::provider, NumLib::GlobalMatrixProvider::provider, NumLib::MatrixProvider::releaseMatrix(), NumLib::VectorProvider::releaseVector(), NumLib::REPEAT_ITERATION, BaseLib::RunTime::start(), and NumLib::SUCCESS.

198 {
199  namespace LinAlg = MathLib::LinAlg;
200  auto& sys = *_equation_system;
201 
203  _res_id);
204  auto& minus_delta_x =
207  auto& J =
209 
210  bool error_norms_met = false;
211 
212  // TODO be more efficient
213  // init minus_delta_x to the right size
214  LinAlg::copy(*x[process_id], minus_delta_x);
215 
217 
218  int iteration = 1;
219  for (; iteration <= _maxiter;
220  ++iteration, _convergence_criterion->reset())
221  {
222  BaseLib::RunTime timer_dirichlet;
223  double time_dirichlet = 0.0;
224 
225  BaseLib::RunTime time_iteration;
226  time_iteration.start();
227 
228  timer_dirichlet.start();
229  sys.computeKnownSolutions(*x[process_id], process_id);
230  sys.applyKnownSolutions(*x[process_id]);
231  time_dirichlet += timer_dirichlet.elapsed();
232 
233  sys.preIteration(iteration, *x[process_id]);
234 
235  BaseLib::RunTime time_assembly;
236  time_assembly.start();
237  sys.assemble(x, process_id);
238  sys.getResidual(*x[process_id], res);
239  sys.getJacobian(J);
240  INFO("[time] Assembly took %g s.", time_assembly.elapsed());
241 
242  minus_delta_x.setZero();
243 
244  timer_dirichlet.start();
245  sys.applyKnownSolutionsNewton(J, res, minus_delta_x);
246  time_dirichlet += timer_dirichlet.elapsed();
247  INFO("[time] Applying Dirichlet BCs took %g s.", time_dirichlet);
248 
249  if (!sys.isLinear() && _convergence_criterion->hasResidualCheck())
250  {
252  }
253 
254  BaseLib::RunTime time_linear_solver;
255  time_linear_solver.start();
256  bool iteration_succeeded = _linear_solver.solve(J, res, minus_delta_x);
257  INFO("[time] Linear solver took %g s.", time_linear_solver.elapsed());
258 
259  if (!iteration_succeeded)
260  {
261  ERR("Newton: The linear solver failed.");
262  }
263  else
264  {
265  // TODO could be solved in a better way
266  // cf.
267  // http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecWAXPY.html
269  *x[process_id], _x_new_id);
270  LinAlg::axpy(x_new, -_damping, minus_delta_x);
271 
272  if (postIterationCallback)
273  {
274  postIterationCallback(iteration, x_new);
275  }
276 
277  switch(sys.postIteration(x_new))
278  {
280  break;
282  ERR("Newton: The postIteration() hook reported a "
283  "non-recoverable error.");
284  iteration_succeeded = false;
285  break;
287  INFO(
288  "Newton: The postIteration() hook decided that this "
289  "iteration"
290  " has to be repeated.");
291  // TODO introduce some onDestroy hook.
293  continue; // That throws the iteration result away.
294  }
295 
296  LinAlg::copy(x_new, *x[process_id]); // copy new solution to x
298  }
299 
300  if (!iteration_succeeded)
301  {
302  // Don't compute further error norms, but break here.
303  error_norms_met = false;
304  break;
305  }
306 
307  if (sys.isLinear()) {
308  error_norms_met = true;
309  } else {
311  // Note: x contains the new solution!
312  _convergence_criterion->checkDeltaX(minus_delta_x,
313  *x[process_id]);
314  }
315 
316  error_norms_met = _convergence_criterion->isSatisfied();
317  }
318 
319  INFO("[time] Iteration #%u took %g s.", iteration,
320  time_iteration.elapsed());
321 
322  if (error_norms_met)
323  {
324  break;
325  }
326 
327  // Avoid increment of the 'iteration' if the error norms are not met,
328  // but maximum number of iterations is reached.
329  if (iteration >= _maxiter)
330  {
331  break;
332  }
333  }
334 
335  if (iteration > _maxiter)
336  {
337  ERR("Newton: Could not solve the given nonlinear system within %u "
338  "iterations",
339  _maxiter);
340  }
341 
345  minus_delta_x);
346 
347  return {error_norms_met, iteration};
348 }
std::size_t _J_id
ID of the Jacobian matrix.
virtual void releaseMatrix(GlobalMatrix const &A)=0
virtual GlobalVector & getVector()=0
Get an uninitialized vector.
static NUMLIB_EXPORT MatrixProvider & provider
virtual bool hasResidualCheck() const =0
static NUMLIB_EXPORT VectorProvider & provider
std::size_t _x_new_id
ID of the vector storing .
int const _maxiter
maximum number of iterations
virtual void checkDeltaX(GlobalVector const &minus_delta_x, GlobalVector const &x)=0
virtual GlobalMatrix & getMatrix()=0
Get an uninitialized matrix.
Count the running time.
Definition: RunTime.h:31
void start()
Start the timer.
Definition: RunTime.h:35
std::size_t _res_id
ID of the residual vector.
virtual bool hasDeltaXCheck() const =0
virtual void releaseVector(GlobalVector const &x)=0
virtual bool isSatisfied() const
Tell if the convergence criterion is satisfied.
void axpy(MatrixOrVector &y, double const a, MatrixOrVector const &x)
Computes .
Definition: LinAlg.h:58
void copy(MatrixOrVector const &x, MatrixOrVector &y)
Copies x to y.
Definition: LinAlg.h:37
virtual void checkResidual(GlobalVector const &residual)=0
Check if the residual satisfies the convergence criterion.
double elapsed() const
Get the elapsed time after started.
Definition: RunTime.h:51

Member Data Documentation

◆ _convergence_criterion

ConvergenceCriterion* NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_convergence_criterion = nullptr
private

Definition at line 127 of file NonlinearSolver.h.

◆ _damping

double const NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_damping
private

A positive damping factor. The default value 1.0 gives a non-damped Newton method. Common values are in the range 0.5 to 0.7 for somewhat conservative method and seldom become smaller than 0.2 for very conservative approach.

Definition at line 134 of file NonlinearSolver.h.

◆ _equation_system

System* NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_equation_system = nullptr
private

Definition at line 124 of file NonlinearSolver.h.

◆ _J_id

std::size_t NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_J_id = 0u
private

ID of the Jacobian matrix.

Definition at line 137 of file NonlinearSolver.h.

◆ _linear_solver

GlobalLinearSolver& NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_linear_solver
private

Definition at line 123 of file NonlinearSolver.h.

◆ _maxiter

int const NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_maxiter
private

maximum number of iterations

Definition at line 128 of file NonlinearSolver.h.

◆ _minus_delta_x_id

std::size_t NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_minus_delta_x_id = 0u
private

ID of the $ -\Delta x$ vector.

Definition at line 138 of file NonlinearSolver.h.

◆ _res_id

std::size_t NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_res_id = 0u
private

ID of the residual vector.

Definition at line 136 of file NonlinearSolver.h.

◆ _x_new_id

std::size_t NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::_x_new_id
private
Initial value:
=
0u

ID of the vector storing $ x - (-\Delta x) $.

Definition at line 139 of file NonlinearSolver.h.


The documentation for this class was generated from the following files: