NOX  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Types | Public Member Functions | Protected Attributes | List of all members
NOX::Thyra::MatrixFreeJacobianOperator< Scalar > Class Template Reference

Concrete implementation of a Thyra::LinearOpBase object that approximates a Jacobian operator based on the Jacobian-Free Newton-Krylov method (see Knoll and Keyes Journal of Computational Physics 193 (2004) 357-397 for details). More...

#include <NOX_Thyra_MatrixFreeJacobianOperator.hpp>

Inheritance diagram for NOX::Thyra::MatrixFreeJacobianOperator< Scalar >:
Inheritance graph
[legend]
Collaboration diagram for NOX::Thyra::MatrixFreeJacobianOperator< Scalar >:
Collaboration graph
[legend]

Public Types

enum  E_DifferenceType { Forward =0, Centered =1 }
 Define types for use of the perturbation parameter $ \delta$.
 
enum  E_PerturbationType { SalingerLOCA =0, KelleySalingerPawlowski =1, KnollKeyes =2, UserDefined =3 }
 Defines the algorithm for computing . More...
 
enum  E_BaseEvaluationType { RawThyra, NoxGroup }
 Defined where to get the base objects for the solution, $x$, and the residual, $ f(x) $, and the perturbed residual evaluation from. In many cases the base values have been precomputed and can just be used to eliminate redundant residual evaluations. More...
 

Public Member Functions

 MatrixFreeJacobianOperator (Teuchos::ParameterList &printParams)
 
void setBaseInArgs (const Teuchos::RCP< ::Thyra::ModelEvaluatorBase::InArgs< Scalar > > &base_in_args)
 
void setLambda (double lambda)
 Change the value of $ \lambda $ in the perturbation calculation.
 
Scalar getLambda () const
 Change the value of $ \lambda $ in the perturbation calculation.
 
void setUserDefinedDelta (double delta)
 Change the value of $ \delta $ in the perturbation calculation.
 
Scalar getUserDefinedDelta () const
 Returns the user defined delta, $ \delta $, used for the perturbation.
 
Scalar getDelta () const
 Returns the last used value of delta $ \lambda $ in the perturbation calculation.
 
Setup the base evaluation sources
void setBaseEvaluationToRawThyra (const Teuchos::RCP< const ::Thyra::VectorBase< Scalar > > &x_base, const Teuchos::RCP< const ::Thyra::VectorBase< Scalar > > &f_base, const Teuchos::RCP< ::Thyra::ModelEvaluator< Scalar > > model)
 
void setBaseEvaluationToNOXGroup (const Teuchos::RCP< const NOX::Abstract::Group > &base_group)
 
Derived from Teuchos::ParameterListAcceptor
void setParameterList (const Teuchos::RCP< Teuchos::ParameterList > &paramList)
 
Teuchos::RCP< const
Teuchos::ParameterList
getValidParameters () const
 
Derived from Thyra::LinearOpBase
Teuchos::RCP< const
::Thyra::VectorSpaceBase
< Scalar > > 
range () const
 
Teuchos::RCP< const
::Thyra::VectorSpaceBase
< Scalar > > 
domain () const
 
Teuchos::RCP< const
::Thyra::LinearOpBase< Scalar > > 
clone () const
 
bool opSupportedImpl (::Thyra::EOpTransp M_trans) const
 
void applyImpl (const ::Thyra::EOpTransp M_trans, const ::Thyra::MultiVectorBase< Scalar > &y, const Teuchos::Ptr< ::Thyra::MultiVectorBase< Scalar > > &u, const Scalar alpha, const Scalar beta) const
 
- Public Member Functions inherited from Teuchos::ParameterListAcceptorDefaultBase
RCP< ParameterListgetNonconstParameterList ()
 
RCP< ParameterListunsetParameterList ()
 
RCP< const ParameterListgetParameterList () const
 

Protected Attributes

bool setup_called_
 true if the algorithm has been setup using the setPArameterList() method
 
Teuchos::RCP< const
::Thyra::ModelEvaluator
< Scalar > > 
model_
 User provided interface function.
 
NOX::Utils utils_
 Printing utilities.
 
Teuchos::RCP
< Teuchos::ParameterList
valid_params_
 A list of valid parameters.
 
E_DifferenceType difference_type_
 
E_PerturbationType perturbation_type_
 
E_BaseEvaluationType base_evaluation_type_
 
Scalar lambda_
 Scale factor for eta calculation.
 
Scalar delta_
 Perturbation value to use in the directional derivative.
 
Scalar user_defined_delta_
 Perturbation value to use in the directional derivative.
 
Teuchos::RCP< const
::Thyra::VectorBase< Scalar > > 
x_base_
 
Teuchos::RCP< const
::Thyra::VectorBase< Scalar > > 
f_base_
 
Teuchos::RCP
< ::Thyra::VectorBase< Scalar > > 
x_perturb_
 
Teuchos::RCP
< ::Thyra::VectorBase< Scalar > > 
f_perturb_
 
Teuchos::RCP
< ::Thyra::VectorBase< Scalar > > 
f2_perturb_
 
Teuchos::RCP< const
NOX::Thyra::Group
base_group_
 
Teuchos::RCP
< ::Thyra::ModelEvaluatorBase::InArgs
< Scalar > > 
in_args_
 
Teuchos::RCP
< ::Thyra::ModelEvaluatorBase::OutArgs
< Scalar > > 
out_args_
 

Additional Inherited Members

- Protected Member Functions inherited from Teuchos::ParameterListAcceptorDefaultBase
void setMyParamList (const RCP< ParameterList > &paramList)
 
RCP< ParameterListgetMyNonconstParamList ()
 
RCP< const ParameterListgetMyParamList () const
 

Detailed Description

template<typename Scalar>
class NOX::Thyra::MatrixFreeJacobianOperator< Scalar >

Concrete implementation of a Thyra::LinearOpBase object that approximates a Jacobian operator based on the Jacobian-Free Newton-Krylov method (see Knoll and Keyes Journal of Computational Physics 193 (2004) 357-397 for details).

This operator approximates the Jacobian action on a vector using rediual evaluations. It is used by Newton-Krylov solvers since these methods only require the matrix-vector product $Jy$ in the linear solve iteration sequence. This product can approximated by directional derivatives:

Forward (1st order accurate):

\[ Jy \approx \frac{F(x + \delta y) - F(x)}{\delta} \]

Central (2nd order accurate):

\[ Jy \approx \frac{F(x + \delta y) - F(x - \delta y)}{\delta} \]

where $J$ is the Jacobian, $y$ is the vector that the Jacobian is to be applied to, $x$ is the solution vector, $F(x)$ is the residual evaluation at $ x $,and $\delta$ is a scalar perturbation calculated by a variety of formulas (see MatrixFreeJacobianOperator::E_PerturbationType for references).

NOTE: The base evaluation $f(x)$ and its corresponding base solution $x$ is assumed to be supplied from an exterior source and that $f(x)$ is up to date wrt $x$. Otherwise this algorithm would be inefficient due to extra base evaluations for every Jacobian-vector product. A call must be made to either setBaseEvaluationToRawThyra() or setBaseEvaluationToNOXGroup() to set up the base objects. If a NOX Group is used, the evaluation will be ensured to be correct because of the mechanics for isF(). In the case of raw objects, the user must take care to keep the base values up-to-date.

Member Enumeration Documentation

Defined where to get the base objects for the solution, $x$, and the residual, $ f(x) $, and the perturbed residual evaluation from. In many cases the base values have been precomputed and can just be used to eliminate redundant residual evaluations.

Enumerator
RawThyra 

User defines thyra objects for solution vector, x, residual vector, f, and residual evaluations.

NoxGroup 

Use the nox group registered with this object.

Defines the algorithm for computing .

Enumerator
SalingerLOCA 

Use the analogous algorithm by Salinger in LOCA v1.0 manual SAND2002-0396 p. 28 eqn. 2.43.

KelleySalingerPawlowski 

Use algorithm developed for NOX by Kelley, Salinger and Pawlowski in 2001.

KnollKeyes 

Knoll and Keyes in JCP 193 (2004) 357-397.

UserDefined 

Use a constant value defined by the user.


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