Thyra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Related Functions | List of all members
Thyra::NonlinearSolverBase< Scalar > Class Template Referenceabstract

Base class for all nonlinear equation solvers. More...

#include <Thyra_NonlinearSolverBase.hpp>

Inheritance diagram for Thyra::NonlinearSolverBase< Scalar >:
Inheritance graph
[legend]

Related Functions

(Note that these are not member functions.)

template<class Scalar >
const SolveStatus< Scalar > solve (NonlinearSolverBase< Scalar > &nonlinearSolver, VectorBase< Scalar > *x, const SolveCriteria< Scalar > *solveCriteria=NULL, VectorBase< Scalar > *delta=NULL)
 

Pure virtual functions that must be overridden in subclasses

virtual void setModel (const RCP< const ModelEvaluator< Scalar > > &model)=0
 Set the model that defines the nonlinear equations. More...
 
virtual RCP< const
ModelEvaluator< Scalar > > 
getModel () const =0
 Get the model that defines the nonlinear equations. More...
 
virtual SolveStatus< Scalar > solve (VectorBase< Scalar > *x, const SolveCriteria< Scalar > *solveCriteria=NULL, VectorBase< Scalar > *delta=NULL)=0
 Solve a set of nonlinear equations from a given starting point. More...
 

Virtual functions with default implementation

virtual bool supportsCloning () const
 Return if this solver object supports cloning or not. More...
 
virtual RCP
< NonlinearSolverBase< Scalar > > 
cloneNonlinearSolver () const
 Clone the solver algorithm if supported. More...
 
virtual RCP< const VectorBase
< Scalar > > 
get_current_x () const
 Return the current value of the solution x as computed in the last solve() operation if supported. More...
 
virtual bool is_W_current () const
 Returns true if *get_W() is current with respect to *get_current_x(). More...
 
virtual RCP
< LinearOpWithSolveBase
< Scalar > > 
get_nonconst_W (const bool forceUpToDate=false)
 Get a nonconst RCP to the Jacobian if available. More...
 
virtual RCP< const
LinearOpWithSolveBase< Scalar > > 
get_W () const
 Get a const RCP to the Jacobian if available. More...
 
virtual void set_W_is_current (bool W_is_current)
 Set if *get_W() is current with respect to *get_current_x(). More...
 

Detailed Description

template<class Scalar>
class Thyra::NonlinearSolverBase< Scalar >

Base class for all nonlinear equation solvers.

Warning! This interface is highly experimental and general developers should not even consider using it in any way if there is any expectation of code stability!

ToDo: Finish documentation.

ToDo:

Definition at line 75 of file Thyra_NonlinearSolverBase.hpp.

Member Function Documentation

template<class Scalar>
virtual void Thyra::NonlinearSolverBase< Scalar >::setModel ( const RCP< const ModelEvaluator< Scalar > > &  model)
pure virtual

Set the model that defines the nonlinear equations.

After the model is set, only the residual f can change between solves and not the structure of the Jacobian W. If a more significant change to *model occurs, then this function must be called again to reset the model and reinitialize.

Implemented in Thyra::DampenedNewtonNonlinearSolver< Scalar >, and Thyra::LinearNonlinearSolver< Scalar >.

template<class Scalar>
virtual RCP<const ModelEvaluator<Scalar> > Thyra::NonlinearSolverBase< Scalar >::getModel ( ) const
pure virtual

Get the model that defines the nonlinear equations.

Implemented in Thyra::DampenedNewtonNonlinearSolver< Scalar >, and Thyra::LinearNonlinearSolver< Scalar >.

template<class Scalar>
virtual SolveStatus<Scalar> Thyra::NonlinearSolverBase< Scalar >::solve ( VectorBase< Scalar > *  x,
const SolveCriteria< Scalar > *  solveCriteria = NULL,
VectorBase< Scalar > *  delta = NULL 
)
pure virtual

Solve a set of nonlinear equations from a given starting point.

Parameters
x[in/out] On input, *x contains the guess for the solution of f(x)=0 as defined by *this->getModel() On output, *x will contain the best estimate of the solution.
solveCriteria[in] If solveCriteria!=NULL then *solveCriteria will contain the criteria for the nonlinear solve.
Returns
The returned SolveStatus object gives the status of the returned solution *x.

Preconditions:

Implemented in Thyra::DampenedNewtonNonlinearSolver< Scalar >, and Thyra::LinearNonlinearSolver< Scalar >.

template<class Scalar >
bool Thyra::NonlinearSolverBase< Scalar >::supportsCloning ( ) const
virtual

Return if this solver object supports cloning or not.

The default implementation returns false.

Definition at line 240 of file Thyra_NonlinearSolverBase.hpp.

template<class Scalar >
RCP< NonlinearSolverBase< Scalar > > Thyra::NonlinearSolverBase< Scalar >::cloneNonlinearSolver ( ) const
virtual

Clone the solver algorithm if supported.

Postconditions:

Note that cloning a nonlinear solver in this case does not imply that the Jacobian state will be copied as well, shallow or deep. Instead, here cloning means to just clone the solver algorithm and it will do a showllow of the model as well if a model is set. Since the model is stateless, this is okay. Therefore, do not assume that the state of *returnValue is exactly the same as the state of *this. You have been warned!

The default implementation returns Teuchos::null which is consistent with the default implementation of supportsCloning(). If this function is overridden in a base class to support cloning, then supportsCloning() must be overridden to return true.

Definition at line 247 of file Thyra_NonlinearSolverBase.hpp.

template<class Scalar >
RCP< const VectorBase< Scalar > > Thyra::NonlinearSolverBase< Scalar >::get_current_x ( ) const
virtual

Return the current value of the solution x as computed in the last solve() operation if supported.

The default implementation returns return.get()==NULL.

Reimplemented in Thyra::DampenedNewtonNonlinearSolver< Scalar >.

Definition at line 254 of file Thyra_NonlinearSolverBase.hpp.

template<class Scalar >
bool Thyra::NonlinearSolverBase< Scalar >::is_W_current ( ) const
virtual

Returns true if *get_W() is current with respect to *get_current_x().

The default implementation returns false.

Reimplemented in Thyra::DampenedNewtonNonlinearSolver< Scalar >.

Definition at line 260 of file Thyra_NonlinearSolverBase.hpp.

template<class Scalar >
RCP< LinearOpWithSolveBase< Scalar > > Thyra::NonlinearSolverBase< Scalar >::get_nonconst_W ( const bool  forceUpToDate = false)
virtual

Get a nonconst RCP to the Jacobian if available.

Parameters
forceUpToDate[in] If true, then the underlying W will be forced to be up to date w.r.t the current x if a Jacobian exists.

Postconditions:

Through this the RCP returned from this function, a client can change the W object held internally. If the object gets changed the client should call set_W_is_current(false).

The default implementation returns return.get()==NULL.

Reimplemented in Thyra::DampenedNewtonNonlinearSolver< Scalar >, and Thyra::LinearNonlinearSolver< Scalar >.

Definition at line 267 of file Thyra_NonlinearSolverBase.hpp.

template<class Scalar >
RCP< const LinearOpWithSolveBase< Scalar > > Thyra::NonlinearSolverBase< Scalar >::get_W ( ) const
virtual

Get a const RCP to the Jacobian if available.

Through this interface the client should not change the object W.

The default implementation returns return.get()==NULL.

Reimplemented in Thyra::DampenedNewtonNonlinearSolver< Scalar >, and Thyra::LinearNonlinearSolver< Scalar >.

Definition at line 274 of file Thyra_NonlinearSolverBase.hpp.

template<class Scalar >
void Thyra::NonlinearSolverBase< Scalar >::set_W_is_current ( bool  W_is_current)
virtual

Set if *get_W() is current with respect to *get_current_x().

Preconditions: this->get_W().get()!=NULL

The default implementation throwns an exception.

Reimplemented in Thyra::DampenedNewtonNonlinearSolver< Scalar >.

Definition at line 280 of file Thyra_NonlinearSolverBase.hpp.

Friends And Related Function Documentation

template<class Scalar >
const SolveStatus< Scalar > solve ( NonlinearSolverBase< Scalar > &  nonlinearSolver,
VectorBase< Scalar > *  x,
const SolveCriteria< Scalar > *  solveCriteria = NULL,
VectorBase< Scalar > *  delta = NULL 
)
related

Definition at line 225 of file Thyra_NonlinearSolverBase.hpp.


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