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, , and the residual, , 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 in the perturbation calculation. | |
Scalar | getLambda () const |
Change the value of in the perturbation calculation. | |
void | setUserDefinedDelta (double delta) |
Change the value of in the perturbation calculation. | |
Scalar | getUserDefinedDelta () const |
Returns the user defined delta, , used for the perturbation. | |
Scalar | getDelta () const |
Returns the last used value of delta 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 > ¶mList) |
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< ParameterList > | getNonconstParameterList () |
RCP< ParameterList > | unsetParameterList () |
RCP< const ParameterList > | getParameterList () 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 > ¶mList) |
RCP< ParameterList > | getMyNonconstParamList () |
RCP< const ParameterList > | getMyParamList () const |
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 in the linear solve iteration sequence. This product can approximated by directional derivatives:
Forward (1st order accurate):
Central (2nd order accurate):
where is the Jacobian, is the vector that the Jacobian is to be applied to, is the solution vector, is the residual evaluation at ,and is a scalar perturbation calculated by a variety of formulas (see MatrixFreeJacobianOperator::E_PerturbationType for references).
NOTE: The base evaluation and its corresponding base solution is assumed to be supplied from an exterior source and that is up to date wrt . 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.
enum NOX::Thyra::MatrixFreeJacobianOperator::E_BaseEvaluationType |
Defined where to get the base objects for the solution, , and the residual, , 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. |
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. |