|
NOX Development
|
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>


Public Types | |
| enum | E_DifferenceType { Forward =0 , Centered =1 } |
| Define types for use of the perturbation parameter | |
| 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, | |
Public Member Functions | |
| MatrixFreeJacobianOperator (Teuchos::ParameterList &printParams) | |
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) |
| void | setBaseInArgs (const Teuchos::RCP< ::Thyra::ModelEvaluatorBase::InArgs< Scalar > > &base_in_args) |
Derived from Teuchos::ParameterListAcceptor | |
| void | setParameterList (const Teuchos::RCP< Teuchos::ParameterList > ¶mList) |
| Teuchos::RCP< const Teuchos::ParameterList > | getValidParameters () const |
Derived from Thyra::LinearOpBase | |
| 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_ |
| 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 |
| void | setLambda (double lambda) |
| Change the value of | |
| Scalar | getLambda () const |
| Change the value of | |
| void | setUserDefinedDelta (double delta) |
| Change the value of | |
| Scalar | getUserDefinedDelta () const |
| Returns the user defined delta, | |
| Scalar | getDelta () const |
| Returns the last used value of delta | |
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 
Forward (1st order accurate):
![\[ Jy \approx \frac{F(x + \delta y) - F(x)}{\delta} \]](form_691.png)
Central (2nd order accurate):
![\[ Jy \approx \frac{F(x + \delta y) - F(x - \delta y)}{\delta} \]](form_692.png)
where 





NOTE: The base evaluation 



| enum NOX::Thyra::MatrixFreeJacobianOperator::E_BaseEvaluationType |
Defined where to get the base objects for the solution, 

| 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. |
| enum NOX::Thyra::MatrixFreeJacobianOperator::E_PerturbationType |
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. |