The multi-dimensional Newton method. More...

#include <newtonmethod.hh>

Inheritance diagram for Opm::NewtonMethod< TypeTag >:
Inheritance graph

Public Member Functions

 NewtonMethod (Simulator &simulator)
 
void finishInit ()
 Finialize the construction of the object. More...
 
bool converged () const
 Returns true if the error of the solution is below the tolerance. More...
 
Problem & problem ()
 Returns a reference to the object describing the current physical problem. More...
 
const Problem & problem () const
 Returns a reference to the object describing the current physical problem. More...
 
Model & model ()
 Returns a reference to the numeric model. More...
 
const Model & model () const
 Returns a reference to the numeric model. More...
 
int numIterations () const
 Returns the number of iterations done since the Newton method was invoked. More...
 
void setIterationIndex (int value)
 Set the index of current iteration. More...
 
Scalar tolerance () const
 Return the current tolerance at which the Newton method considers itself to be converged. More...
 
void setTolerance (Scalar value)
 Set the current tolerance at which the Newton method considers itself to be converged. More...
 
bool apply ()
 Run the Newton method. More...
 
Scalar suggestTimeStepSize (Scalar oldDt) const
 Suggest a new time-step size based on the old time-step size. More...
 
std::ostringstream & endIterMsg ()
 Message that should be printed for the user after the end of an iteration. More...
 
void eraseMatrix ()
 Causes the solve() method to discared the structure of the linear system of equations the next time it is called. More...
 
LinearSolverBackend & linearSolver ()
 Returns the linear solver backend object for external use. More...
 
const LinearSolverBackend & linearSolver () const
 Returns the linear solver backend object for external use. More...
 
const TimerprePostProcessTimer () const
 
const TimerlinearizeTimer () const
 
const TimersolveTimer () const
 
const TimerupdateTimer () const
 

Static Public Member Functions

static void registerParameters ()
 Register all run-time parameters for the Newton method. More...
 

Protected Member Functions

bool verbose_ () const
 Returns true if the Newton method ought to be chatty. More...
 
void begin_ (const SolutionVector &)
 Called before the Newton method is applied to an non-linear system of equations. More...
 
void beginIteration_ ()
 Indicates the beginning of a Newton iteration. More...
 
void linearizeDomain_ ()
 Linearize the global non-linear system of equations associated with the spatial domain. More...
 
void linearizeAuxiliaryEquations_ ()
 
void preSolve_ (const SolutionVector &, const GlobalEqVector &currentResidual)
 
void postSolve_ (const SolutionVector &, const GlobalEqVector &, GlobalEqVector &solutionUpdate)
 Update the error of the solution given the previous iteration. More...
 
void update_ (SolutionVector &nextSolution, const SolutionVector &currentSolution, const GlobalEqVector &solutionUpdate, const GlobalEqVector &currentResidual)
 Update the current solution with a delta vector. More...
 
void updateConstraintDof_ (unsigned, PrimaryVariables &nextValue, const Constraints &constraints)
 Update the primary variables for a degree of freedom which is constraint. More...
 
void updatePrimaryVariables_ (unsigned, PrimaryVariables &nextValue, const PrimaryVariables &currentValue, const EqVector &update, const EqVector &)
 Update a single primary variables object. More...
 
void writeConvergence_ (const SolutionVector &currentSolution, const GlobalEqVector &solutionUpdate)
 Write the convergence behaviour of the newton method to disk. More...
 
void endIteration_ (const SolutionVector &, const SolutionVector &)
 Indicates that one Newton iteration was finished. More...
 
bool proceed_ () const
 Returns true iff another Newton iteration should be done. More...
 
void end_ ()
 Indicates that we're done solving the non-linear system of equations. More...
 
void failed_ ()
 Called if the Newton method broke down. More...
 
void succeeded_ ()
 Called if the Newton method was successful. More...
 
int targetIterations_ () const
 
int maxIterations_ () const
 

Static Protected Member Functions

static bool enableConstraints_ ()
 

Protected Attributes

Simulator & simulator_
 
Timer prePostProcessTimer_
 
Timer linearizeTimer_
 
Timer solveTimer_
 
Timer updateTimer_
 
std::ostringstream endIterMsgStream_
 
Scalar error_
 
Scalar lastError_
 
Scalar tolerance_
 
int numIterations_
 
LinearSolverBackend linearSolver_
 
CollectiveCommunication comm_
 
ConvergenceWriter convergenceWriter_
 

Detailed Description

template<class TypeTag>
class Opm::NewtonMethod< TypeTag >

The multi-dimensional Newton method.

This class uses static polymorphism to allow implementations to implement different update/convergence strategies.

Constructor & Destructor Documentation

◆ NewtonMethod()

Member Function Documentation

◆ apply()

◆ begin_()

template<class TypeTag >
void Opm::NewtonMethod< TypeTag >::begin_ ( const SolutionVector &  )
inlineprotected

Called before the Newton method is applied to an non-linear system of equations.

Parameters
uThe initial solution

References Opm::NewtonMethod< TypeTag >::convergenceWriter_, and Opm::NewtonMethod< TypeTag >::numIterations_.

◆ beginIteration_()

◆ converged()

template<class TypeTag >
bool Opm::NewtonMethod< TypeTag >::converged ( ) const
inline

Returns true if the error of the solution is below the tolerance.

References Opm::NewtonMethod< TypeTag >::error_, and Opm::NewtonMethod< TypeTag >::tolerance().

Referenced by Opm::NewtonMethod< TypeTag >::apply(), and Opm::NewtonMethod< TypeTag >::proceed_().

◆ enableConstraints_()

template<class TypeTag >
static bool Opm::NewtonMethod< TypeTag >::enableConstraints_ ( )
inlinestaticprotected

◆ end_()

template<class TypeTag >
void Opm::NewtonMethod< TypeTag >::end_ ( )
inlineprotected

Indicates that we're done solving the non-linear system of equations.

References Opm::NewtonMethod< TypeTag >::convergenceWriter_.

◆ endIteration_()

template<class TypeTag >
void Opm::NewtonMethod< TypeTag >::endIteration_ ( const SolutionVector &  ,
const SolutionVector &   
)
inlineprotected

Indicates that one Newton iteration was finished.

Parameters
nextSolutionThe solution after the current Newton iteration
currentSolutionThe solution at the beginning of the current Newton iteration

References Opm::NewtonMethod< TypeTag >::endIterMsg(), Opm::NewtonMethod< TypeTag >::error_, Opm::NewtonMethod< TypeTag >::numIterations_, Opm::NewtonMethod< TypeTag >::problem(), Opm::NewtonMethod< TypeTag >::simulator_, and Opm::NewtonMethod< TypeTag >::verbose_().

◆ endIterMsg()

template<class TypeTag >
std::ostringstream & Opm::NewtonMethod< TypeTag >::endIterMsg ( )
inline

Message that should be printed for the user after the end of an iteration.

References Opm::NewtonMethod< TypeTag >::endIterMsgStream_.

Referenced by Opm::NewtonMethod< TypeTag >::endIteration_().

◆ eraseMatrix()

template<class TypeTag >
void Opm::NewtonMethod< TypeTag >::eraseMatrix ( )
inline

Causes the solve() method to discared the structure of the linear system of equations the next time it is called.

References Opm::NewtonMethod< TypeTag >::linearSolver_.

◆ failed_()

template<class TypeTag >
void Opm::NewtonMethod< TypeTag >::failed_ ( )
inlineprotected

Called if the Newton method broke down.

This method is called after end_()

References Opm::NewtonMethod< TypeTag >::numIterations_, and Opm::NewtonMethod< TypeTag >::targetIterations_().

◆ finishInit()

template<class TypeTag >
void Opm::NewtonMethod< TypeTag >::finishInit ( )
inline

Finialize the construction of the object.

At this point, it can be assumed that all objects featured by the simulator have been allocated. (But not that they have been fully initialized yet.)

◆ linearizeAuxiliaryEquations_()

template<class TypeTag >
void Opm::NewtonMethod< TypeTag >::linearizeAuxiliaryEquations_ ( )
inlineprotected

◆ linearizeDomain_()

template<class TypeTag >
void Opm::NewtonMethod< TypeTag >::linearizeDomain_ ( )
inlineprotected

Linearize the global non-linear system of equations associated with the spatial domain.

References Opm::NewtonMethod< TypeTag >::model().

◆ linearizeTimer()

template<class TypeTag >
const Timer & Opm::NewtonMethod< TypeTag >::linearizeTimer ( ) const
inline

◆ linearSolver() [1/2]

template<class TypeTag >
LinearSolverBackend & Opm::NewtonMethod< TypeTag >::linearSolver ( )
inline

Returns the linear solver backend object for external use.

References Opm::NewtonMethod< TypeTag >::linearSolver_.

◆ linearSolver() [2/2]

template<class TypeTag >
const LinearSolverBackend & Opm::NewtonMethod< TypeTag >::linearSolver ( ) const
inline

Returns the linear solver backend object for external use.

References Opm::NewtonMethod< TypeTag >::linearSolver_.

◆ maxIterations_()

template<class TypeTag >
int Opm::NewtonMethod< TypeTag >::maxIterations_ ( ) const
inlineprotected

◆ model() [1/2]

◆ model() [2/2]

template<class TypeTag >
const Model & Opm::NewtonMethod< TypeTag >::model ( ) const
inline

Returns a reference to the numeric model.

References Opm::NewtonMethod< TypeTag >::simulator_.

◆ numIterations()

template<class TypeTag >
int Opm::NewtonMethod< TypeTag >::numIterations ( ) const
inline

Returns the number of iterations done since the Newton method was invoked.

References Opm::NewtonMethod< TypeTag >::numIterations_.

Referenced by Opm::NewtonMethod< TypeTag >::proceed_().

◆ postSolve_()

template<class TypeTag >
void Opm::NewtonMethod< TypeTag >::postSolve_ ( const SolutionVector &  ,
const GlobalEqVector &  ,
GlobalEqVector &  solutionUpdate 
)
inlineprotected

Update the error of the solution given the previous iteration.

For our purposes, the error of a solution is defined as the maximum of the weighted residual of a given solution.

Parameters
currentSolutionThe solution at the beginning the current iteration
currentResidualThe residual (i.e., right-hand-side) of the current iteration's solution.
solutionUpdateThe difference between the current and the next solution

References Opm::NewtonMethod< TypeTag >::model(), and Opm::NewtonMethod< TypeTag >::simulator_.

◆ prePostProcessTimer()

template<class TypeTag >
const Timer & Opm::NewtonMethod< TypeTag >::prePostProcessTimer ( ) const
inline

◆ preSolve_()

template<class TypeTag >
void Opm::NewtonMethod< TypeTag >::preSolve_ ( const SolutionVector &  ,
const GlobalEqVector &  currentResidual 
)
inlineprotected

◆ problem() [1/2]

template<class TypeTag >
Problem & Opm::NewtonMethod< TypeTag >::problem ( )
inline

◆ problem() [2/2]

template<class TypeTag >
const Problem & Opm::NewtonMethod< TypeTag >::problem ( ) const
inline

Returns a reference to the object describing the current physical problem.

References Opm::NewtonMethod< TypeTag >::simulator_.

◆ proceed_()

template<class TypeTag >
bool Opm::NewtonMethod< TypeTag >::proceed_ ( ) const
inlineprotected

◆ registerParameters()

template<class TypeTag >
static void Opm::NewtonMethod< TypeTag >::registerParameters ( )
inlinestatic

Register all run-time parameters for the Newton method.

Referenced by Opm::FvBaseDiscretization< TypeTag >::registerParameters().

◆ setIterationIndex()

template<class TypeTag >
void Opm::NewtonMethod< TypeTag >::setIterationIndex ( int  value)
inline

Set the index of current iteration.

Normally this does not need to be called, but if the non-linear solver is implemented externally, it needs to be set in order for the model to do the Right Thing (TM) while linearizing.

References Opm::NewtonMethod< TypeTag >::numIterations_.

◆ setTolerance()

template<class TypeTag >
void Opm::NewtonMethod< TypeTag >::setTolerance ( Scalar  value)
inline

Set the current tolerance at which the Newton method considers itself to be converged.

References Opm::NewtonMethod< TypeTag >::tolerance_.

◆ solveTimer()

template<class TypeTag >
const Timer & Opm::NewtonMethod< TypeTag >::solveTimer ( ) const
inline

◆ succeeded_()

template<class TypeTag >
void Opm::NewtonMethod< TypeTag >::succeeded_ ( )
inlineprotected

Called if the Newton method was successful.

This method is called after end_()

◆ suggestTimeStepSize()

template<class TypeTag >
Scalar Opm::NewtonMethod< TypeTag >::suggestTimeStepSize ( Scalar  oldDt) const
inline

Suggest a new time-step size based on the old time-step size.

The default behavior is to suggest the old time-step size scaled by the ratio between the target iterations and the iterations required to actually solve the last time-step.

References Opm::NewtonMethod< TypeTag >::numIterations_, Opm::NewtonMethod< TypeTag >::problem(), and Opm::NewtonMethod< TypeTag >::targetIterations_().

◆ targetIterations_()

template<class TypeTag >
int Opm::NewtonMethod< TypeTag >::targetIterations_ ( ) const
inlineprotected

◆ tolerance()

template<class TypeTag >
Scalar Opm::NewtonMethod< TypeTag >::tolerance ( ) const
inline

Return the current tolerance at which the Newton method considers itself to be converged.

References Opm::NewtonMethod< TypeTag >::tolerance_.

Referenced by Opm::NewtonMethod< TypeTag >::converged().

◆ update_()

template<class TypeTag >
void Opm::NewtonMethod< TypeTag >::update_ ( SolutionVector &  nextSolution,
const SolutionVector &  currentSolution,
const GlobalEqVector &  solutionUpdate,
const GlobalEqVector &  currentResidual 
)
inlineprotected

Update the current solution with a delta vector.

Different update strategies, such as chopped updates can be implemented by overriding this method. The default behavior is use the standard Newton-Raphson update strategy, i.e.

\[ u^{k+1} = u^k - \Delta u^k \]

Parameters
nextSolutionThe solution vector after the current iteration
currentSolutionThe solution vector after the last iteration
solutionUpdateThe delta vector as calculated by solving the linear system of equations
currentResidualThe residual vector of the current Newton-Raphson iteraton

References Opm::NewtonMethod< TypeTag >::enableConstraints_(), and Opm::NewtonMethod< TypeTag >::model().

Referenced by Opm::FvBaseNewtonMethod< TypeTag >::update_().

◆ updateConstraintDof_()

template<class TypeTag >
void Opm::NewtonMethod< TypeTag >::updateConstraintDof_ ( unsigned  ,
PrimaryVariables &  nextValue,
const Constraints &  constraints 
)
inlineprotected

Update the primary variables for a degree of freedom which is constraint.

◆ updatePrimaryVariables_()

template<class TypeTag >
void Opm::NewtonMethod< TypeTag >::updatePrimaryVariables_ ( unsigned  ,
PrimaryVariables &  nextValue,
const PrimaryVariables &  currentValue,
const EqVector &  update,
const EqVector &   
)
inlineprotected

Update a single primary variables object.

◆ updateTimer()

template<class TypeTag >
const Timer & Opm::NewtonMethod< TypeTag >::updateTimer ( ) const
inline

◆ verbose_()

template<class TypeTag >
bool Opm::NewtonMethod< TypeTag >::verbose_ ( ) const
inlineprotected

Returns true if the Newton method ought to be chatty.

References Opm::NewtonMethod< TypeTag >::comm_.

Referenced by Opm::NewtonMethod< TypeTag >::apply(), and Opm::NewtonMethod< TypeTag >::endIteration_().

◆ writeConvergence_()

template<class TypeTag >
void Opm::NewtonMethod< TypeTag >::writeConvergence_ ( const SolutionVector &  currentSolution,
const GlobalEqVector &  solutionUpdate 
)
inlineprotected

Write the convergence behaviour of the newton method to disk.

This method is called as part of the update proceedure.

References Opm::NewtonMethod< TypeTag >::convergenceWriter_.

Member Data Documentation

◆ comm_

template<class TypeTag >
CollectiveCommunication Opm::NewtonMethod< TypeTag >::comm_
protected

◆ convergenceWriter_

template<class TypeTag >
ConvergenceWriter Opm::NewtonMethod< TypeTag >::convergenceWriter_
protected

◆ endIterMsgStream_

template<class TypeTag >
std::ostringstream Opm::NewtonMethod< TypeTag >::endIterMsgStream_
protected

◆ error_

◆ lastError_

◆ linearizeTimer_

template<class TypeTag >
Timer Opm::NewtonMethod< TypeTag >::linearizeTimer_
protected

◆ linearSolver_

template<class TypeTag >
LinearSolverBackend Opm::NewtonMethod< TypeTag >::linearSolver_
protected

◆ numIterations_

◆ prePostProcessTimer_

template<class TypeTag >
Timer Opm::NewtonMethod< TypeTag >::prePostProcessTimer_
protected

◆ simulator_

◆ solveTimer_

template<class TypeTag >
Timer Opm::NewtonMethod< TypeTag >::solveTimer_
protected

◆ tolerance_

◆ updateTimer_

template<class TypeTag >
Timer Opm::NewtonMethod< TypeTag >::updateTimer_
protected

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