NOX  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
NOX::Thyra::Group Class Reference

A concrete implementation of the NOX::Abstract::Group using Thyra. More...

#include <NOX_Thyra_Group.H>

Inheritance diagram for NOX::Thyra::Group:
Inheritance graph
[legend]
Collaboration diagram for NOX::Thyra::Group:
Collaboration graph
[legend]

Public Member Functions

 Group (const NOX::Thyra::Vector &initialGuess, const Teuchos::RCP< const ::Thyra::ModelEvaluator< double > > &model, const Teuchos::RCP< const ::Thyra::VectorBase< double > > &weightVector=Teuchos::null, const Teuchos::RCP< const ::Thyra::VectorBase< double > > &rightWeightVector=Teuchos::null, const Teuchos::RCP<::Thyra::VectorBase< double > > &inv_rightWeightVector=Teuchos::null, const bool rightScalingFirst=false)
 The default constructor that uses the linear solver from the ModelEvaluator. More...
 
 Group (const NOX::Thyra::Vector &initialGuess, const Teuchos::RCP< const ::Thyra::ModelEvaluator< double > > &model, const Teuchos::RCP< ::Thyra::LinearOpBase< double > > &linearOp, const Teuchos::RCP< const ::Thyra::LinearOpWithSolveFactoryBase< double > > &lowsFactory, const Teuchos::RCP< ::Thyra::PreconditionerBase< double > > &precOp, const Teuchos::RCP< ::Thyra::PreconditionerFactoryBase< double > > &precFactory, const Teuchos::RCP< const ::Thyra::VectorBase< double > > &weightVector=Teuchos::null, const Teuchos::RCP< const ::Thyra::VectorBase< double > > &rightWeightVector=Teuchos::null, const Teuchos::RCP<::Thyra::VectorBase< double > > &inv_rightWeightVector=Teuchos::null, const bool rightScalingFirst=false, const bool updatePreconditioner=true, const bool jacobianIsEvaluated=false)
 Power user constructor that takes explicit linear solver objects to handle different combinations. More...
 
 Group (const NOX::Thyra::Group &source, NOX::CopyType type=DeepCopy)
 Copy constructor.
 
 ~Group ()
 Destructor.
 
NOX::Abstract::Groupoperator= (const NOX::Abstract::Group &source)
 Copies the source group into this group. More...
 
NOX::Abstract::Groupoperator= (const NOX::Thyra::Group &source)
 
Teuchos::RCP< const
::Thyra::VectorBase< double > > 
get_current_x () const
 
Teuchos::RCP
< ::Thyra::LinearOpBase
< double > > 
getNonconstJacobianOperator ()
 
Teuchos::RCP< const
::Thyra::LinearOpBase< double > > 
getJacobianOperator () const
 
Teuchos::RCP< const
::Thyra::LinearOpBase< double > > 
getScaledJacobianOperator () const
 
void unscaleJacobianOperator () const
 
Teuchos::RCP
< ::Thyra::LinearOpWithSolveBase
< double > > 
getNonconstJacobian ()
 
Teuchos::RCP< const
::Thyra::LinearOpWithSolveBase
< double > > 
getJacobian () const
 
Teuchos::RCP
< ::Thyra::PreconditionerBase
< double > > 
getNonconstPreconditioner ()
 
Teuchos::RCP< const
::Thyra::PreconditionerBase
< double > > 
getPreconditioner () const
 
void setJacobianOperator (const Teuchos::RCP<::Thyra::LinearOpBase< double >> &op)
 Dangerous power user function for LOCA Householder bordered algorithm.
 
void setPreconditionerMatrix (const Teuchos::RCP< const ::Thyra::DefaultLinearOpSource< double >> &op)
 Dangerous power user function for LOCA Householder bordered algorithm. This is the Matrix M that is used to initialize a stratimikos preconditioner. NOTE: this sets the losb_ object used to update prec_!
 
virtual void logLastLinearSolveStats (NOX::SolverStats &stats) const
 Adds statistics from last linear solve to the SovlerStats object.
 
virtual Teuchos::RCP
< NOX::Abstract::Group
clone (NOX::CopyType type=NOX::DeepCopy) const
 Create a new Group of the same derived type as this one by cloning this one, and return a ref count pointer to the new group. More...
 
void print () const
 Print out the group.
 
Teuchos::RCP< const
::Thyra::ModelEvaluator
< double > > 
getModel () const
 
void enablePseudoTransientTerms (const Teuchos::RCP< const ::Thyra::VectorBase< double >> &x_dot, const double alpha, const double beta, const double t)
 Set the transient terms on the Group and use them in the underlying evalModelImpl() calls. More...
 
void disablePseudoTransientTerms ()
 Disable the pseudo traansient terms in the underlying evalModel() calls. Sets x_dot, alpha, beta and t back to steady state values.
 
bool usingPseudoTransientTerms () const
 Check for whether the pseudo transient support is enabled for residual and Jacobian evaluations.
 
void setBasePoint (const ::Thyra::ModelEvaluatorBase::InArgs< double > &base_point_params)
 
void unsetBasePoint ()
 Unset the base point parameters so that they are not used internally.
 
bool usingBasePoint () const
 Returns true if a base point has been set.
 
"Compute" functions.
void setX (const NOX::Abstract::Vector &y)
 Set the solution vector x to y. More...
 
void setX (const NOX::Thyra::Vector &y)
 See above.
 
void computeX (const NOX::Abstract::Group &grp, const NOX::Abstract::Vector &d, double step)
 Compute x = grp.x + step * d. More...
 
void computeX (const NOX::Thyra::Group &grp, const NOX::Thyra::Vector &d, double step)
 See above.
 
NOX::Abstract::Group::ReturnType computeF ()
 Compute and store F(x). More...
 
NOX::Abstract::Group::ReturnType computeJacobian ()
 Compute and store Jacobian. More...
 
NOX::Abstract::Group::ReturnType computeGradient ()
 Compute and store gradient. More...
 
NOX::Abstract::Group::ReturnType computeNewton (Teuchos::ParameterList &params)
 Compute the Newton direction, using parameters for the linear solve. More...
 
Jacobian operations.

Operations using the Jacobian matrix. These may not be defined in matrix-free scenarios.

NOX::Abstract::Group::ReturnType applyJacobian (const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const
 Applies Jacobian to the given input vector and puts the answer in the result. More...
 
NOX::Abstract::Group::ReturnType applyJacobian (const NOX::Thyra::Vector &input, NOX::Thyra::Vector &result) const
 
NOX::Abstract::Group::ReturnType applyJacobianMultiVector (const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const
 applyJacobian for multiple right-hand sides More...
 
NOX::Abstract::Group::ReturnType applyJacobianTranspose (const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const
 Applies Jacobian-Transpose to the given input vector and puts the answer in the result. More...
 
NOX::Abstract::Group::ReturnType applyJacobianTranspose (const NOX::Thyra::Vector &input, NOX::Thyra::Vector &result) const
 
NOX::Abstract::Group::ReturnType applyJacobianTransposeMultiVector (const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const
 applyJacobianTranspose for multiple right-hand sides More...
 
NOX::Abstract::Group::ReturnType applyJacobianInverse (Teuchos::ParameterList &params, const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const
 Applies the inverse of the Jacobian matrix to the given input vector and puts the answer in result. More...
 
NOX::Abstract::Group::ReturnType applyJacobianInverse (Teuchos::ParameterList &params, const NOX::Thyra::Vector &input, NOX::Thyra::Vector &result) const
 
NOX::Abstract::Group::ReturnType applyJacobianInverseMultiVector (Teuchos::ParameterList &params, const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const
 applyJacobianInverse for multiple right-hand sides More...
 
NOX::Abstract::Group::ReturnType applyRightPreconditioning (bool useTranspose, Teuchos::ParameterList &params, const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const
 Apply right preconditiong to the given input vector. More...
 
"Is" functions

Checks to see if various objects have been computed. Returns true if the corresponding "compute" function has been called since the last update to the solution vector (via instantiation or computeX).

bool isF () const
 Return true if F is valid.
 
bool isJacobian () const
 Return true if the Jacobian is valid. More...
 
bool isGradient () const
 Return true if the gradient is valid. More...
 
bool isNewton () const
 Return true if the Newton direction is valid. More...
 
"Get" functions

Note that these function do not check whether or not the vectors are valid. Must use the "Is" functions for that purpose.

const NOX::Abstract::VectorgetX () const
 Return solution vector.
 
const NOX::Abstract::VectorgetScaledX () const
 
const NOX::Abstract::VectorgetF () const
 Return F(x)
 
double getNormF () const
 Return 2-norm of F(x). More...
 
const NOX::Abstract::VectorgetGradient () const
 Return gradient.
 
const NOX::Abstract::VectorgetNewton () const
 Return Newton direction.
 
Teuchos::RCP< const
NOX::Abstract::Vector
getXPtr () const
 Return RCP to solution vector.
 
Teuchos::RCP< const
NOX::Abstract::Vector
getFPtr () const
 Return RCP to F(x)
 
Teuchos::RCP< const
NOX::Abstract::Vector
getGradientPtr () const
 Return RCP to gradient.
 
Teuchos::RCP< const
NOX::Abstract::Vector
getNewtonPtr () const
 Return RCP to Newton direction.
 
- Public Member Functions inherited from NOX::Abstract::Group
 Group ()
 Constructor. More...
 
virtual Teuchos::RCP
< NOX::Abstract::Group
getNestedGroup ()
 Return an internally stored group from this group. More...
 
virtual Teuchos::RCP< const
NOX::Abstract::Group
getNestedGroup () const
 Return an internally stored group from this group. More...
 
virtual
NOX::Abstract::Group::ReturnType 
applyRightPreconditioningMultiVector (bool useTranspose, Teuchos::ParameterList &params, const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const
 applyRightPreconditioning for multiple right-hand sides More...
 
virtual
NOX::Abstract::Group::ReturnType 
getNormLastLinearSolveResidual (double &residual) const
 Return the norm of the last linear solve residual as the result of either a call to computeNewton() or applyJacobianInverse(). More...
 

Protected Member Functions

void resetIsValidFlags ()
 resets the isValid flags to false
 
NOX::Abstract::Group::ReturnType applyJacobianInverseMultiVector (Teuchos::ParameterList &p, const ::Thyra::MultiVectorBase< double > &input,::Thyra::MultiVectorBase< double > &result) const
 Apply Jacobian inverse using Thyra objects.
 
::Thyra::ESolveMeasureNormType getThyraNormType (const std::string &name) const
 
void updateLOWS () const
 Finalizes LOWS to be a valid solver for the Jacobian.
 
void scaleResidualAndJacobian () const
 
void unscaleResidualAndJacobian () const
 
void computeScaledSolution ()
 

Protected Attributes

Teuchos::RCP< const
::Thyra::ModelEvaluator
< double > > 
model_
 Problem interface.
 
Teuchos::RCP< NOX::Thyra::Vectorx_vec_
 Solution vector.
 
Teuchos::RCP< NOX::Thyra::Vectorf_vec_
 Residual vector.
 
Teuchos::RCP< NOX::Thyra::Vectornewton_vec_
 Newton direction vector.
 
Teuchos::RCP< NOX::Thyra::Vectorgradient_vec_
 Gradient direction vector.
 
Teuchos::RCP
< NOX::SharedObject
< ::Thyra::LinearOpWithSolveBase
< double >, NOX::Thyra::Group > > 
shared_jacobian_
 Shared Jacobian operator with solve.
 
Teuchos::RCP
< ::Thyra::LinearOpBase
< double > > 
lop_
 Jacobian operator.
 
Teuchos::RCP< const
::Thyra::LinearOpWithSolveFactoryBase
< double > > 
lows_factory_
 Thyra LOWS factory for building Jacobians.
 
Teuchos::RCP< const
::Thyra::DefaultLinearOpSource
< double > > 
losb_
 Source base needed to create preconditioner.
 
Teuchos::RCP
< ::Thyra::PreconditionerBase
< double > > 
prec_
 Preconditioner for Jacobian.
 
Teuchos::RCP
< ::Thyra::PreconditionerFactoryBase
< double > > 
prec_factory_
 Preconditioner factory.
 
Teuchos::RCP< const
::Thyra::VectorBase< double > > 
weight_vec_
 Optional wieghting vector for function scaling. NOX assumes that this vector can be updated in between nonlinear iterations. More...
 
Teuchos::RCP< const
::Thyra::VectorBase< double > > 
right_weight_vec_
 Optional wieghting vector for solution (right) scaling.
 
Teuchos::RCP
< ::Thyra::VectorBase< double > > 
inv_weight_vec_
 Inverse of weight vector used to unscale function (left) scaling. NOX assumes that this vector can be updated in between nonlinear iterations.
 
Teuchos::RCP
< ::Thyra::VectorBase< double > > 
inv_right_weight_vec_
 Inverse of weight vector used to unscale solution (right) scaling.
 
Teuchos::RCP< NOX::Thyra::Vectorscaled_x_vec_
 Scaled solution vector scaled by the.
 
bool rightScalingFirst_
 Do right scaling before left scaling?
 
bool updatePreconditioner_
 If set to true, the preconditioner matrix values will be automatically updated via precFactory or ModelEvalautor. If set to false, the user must manually handle updating the preconditioner.
 
NOX::Abstract::Group::ReturnType last_linear_solve_status_
 The status of the last linear solve performed.
 
int last_linear_solve_num_iters_
 Number of iterations for last linear solve performed.
 
double last_linear_solve_achieved_tol_
 The tolerance achieved by the last linear solver.
 
bool use_pseudo_transient_terms_
 
Teuchos::RCP< const
::Thyra::VectorBase< double > > 
x_dot_
 
double alpha_
 
double beta_
 
double t_
 
bool use_base_point_
 
::Thyra::ModelEvaluatorBase::InArgs
< double > 
base_point_
 
IsValid flags

True if the current solution is up-to-date with respect to the currect solution vector.

bool is_valid_f_
 
bool is_valid_jacobian_
 
bool is_valid_newton_dir_
 
bool is_valid_gradient_dir_
 
bool is_valid_lows_
 

Additional Inherited Members

- Public Types inherited from NOX::Abstract::Group
enum  ReturnType {
  Ok, NotDefined, BadDependency, NotConverged,
  Failed
}
 The computation of, say, the Newton direction in computeNewton() may fail in many different ways, so we have included a variety of return codes to describe the failures. Of course, we also have a code for success. More...
 

Detailed Description

A concrete implementation of the NOX::Abstract::Group using Thyra.

NOTE: This Group supports row sum scaling of the function (residual and Jacobian). This is enabled by setting a weight vector on the initial guess vector in the Group constructor. The residual and Jacobian must be scaled before and then unscaled after calls to construct the preconditioner and solve the linear system. This follows the nox epetra group. This design should be changed in a future nox refactor, but requires significant changes to the Group object.

Constructor & Destructor Documentation

NOX::Thyra::Group::Group ( const NOX::Thyra::Vector initialGuess,
const Teuchos::RCP< const ::Thyra::ModelEvaluator< double > > &  model,
const Teuchos::RCP< const ::Thyra::VectorBase< double > > &  weightVector = Teuchos::null,
const Teuchos::RCP< const ::Thyra::VectorBase< double > > &  rightWeightVector = Teuchos::null,
const Teuchos::RCP<::Thyra::VectorBase< double > > &  inv_rightWeightVector = Teuchos::null,
const bool  rightScalingFirst = false 
)

The default constructor that uses the linear solver from the ModelEvaluator.

Most users should use this constructor. It is meant to be used in conjunction with a stratimikos linear solver that is built as part of the input model evaluator. For finer grained control over the use of the preconditioner and for Jacobian-Free Newton-Krylov cases, the power user constructor should be used.

Parameters
[in]initialGuessInitial guess for the solution vector
[in]modelModelEvaluator
[in]weightVectorOptional diagonal weighting vector for the model.
[in]rightWeightVectorOptional solution vector weighting
[in]inv_rightWeightVectorOptional inverse solution vector weighting
[in]rightScalingFirstOptional bool to select if right scaling should be applied before left scaling

References NOX::DeepCopy, f_vec_, gradient_vec_, inv_right_weight_vec_, lop_, losb_, lows_factory_, newton_vec_, nonnull(), Teuchos::nonnull(), prec_, prec_factory_, Teuchos::RCP< T >::ptr(), Teuchos::rcp(), resetIsValidFlags(), right_weight_vec_, rightScalingFirst_, scaled_x_vec_, NOX::ShapeCopy, shared_jacobian_, TEUCHOS_ASSERT, weight_vec_, and x_vec_.

NOX::Thyra::Group::Group ( const NOX::Thyra::Vector initialGuess,
const Teuchos::RCP< const ::Thyra::ModelEvaluator< double > > &  model,
const Teuchos::RCP< ::Thyra::LinearOpBase< double > > &  linearOp,
const Teuchos::RCP< const ::Thyra::LinearOpWithSolveFactoryBase< double > > &  lowsFactory,
const Teuchos::RCP< ::Thyra::PreconditionerBase< double > > &  precOp,
const Teuchos::RCP< ::Thyra::PreconditionerFactoryBase< double > > &  precFactory,
const Teuchos::RCP< const ::Thyra::VectorBase< double > > &  weightVector = Teuchos::null,
const Teuchos::RCP< const ::Thyra::VectorBase< double > > &  rightWeightVector = Teuchos::null,
const Teuchos::RCP<::Thyra::VectorBase< double > > &  inv_rightWeightVector = Teuchos::null,
const bool  rightScalingFirst = false,
const bool  updatePreconditioner = true,
const bool  jacobianIsEvaluated = false 
)

Power user constructor that takes explicit linear solver objects to handle different combinations.

This class allows the user to set user-defined linear operators and preconditioners (and corresponding factories). The user can set the linear_op to be a Jacobian-Free Newton Krylov operator (use the class NOX::Thyra::MatrixFreeJacobianOperator).

Parameters
[in]initialGuess(Required) Initial guess for the solution vector
[in]model(Required) ModelEvaluator
[in]linearOp(Optional) Forward operator for the Jacobian. Must be non-null for Newton-based solvers.
[in]lowsFactory(Optional) Factory for building and updating linear solver.
[in]precOp(Optional) Preconditioner operator. If set to Teuchos::null and a non-null prec_factory exists, the prec_op will be constructed using the preconditioner factory.
[in]precFactory(Optional) Factory for updating the precondiitoner. If set to Teuchos::null and there is a non-null prec_op, then the preconditioner will be updated using the model evaluator as long as the ModelEvaluator::outArgs supports W_prec.
[in]weightVector(Optional) diagonal weighting vector for the model.
[in]rightWeightVectorOptional solution vector weighting
[in]inv_rightWeightVectorOptional inverse solution vector weighting
[in]rightScalingFirstOptional bool to select if right scaling should be applied before left scaling
[in]updatePreconditionerOptional bool to select if the Group should auotmatically update the preconditioner matrix values between Newton iterations
[in]jacobianIsEvaluatedOptional bool, if true this means that the input Jacobian operator (linearOp) has been evaluated externally and is consistent with the initialGuess. In this case, the isJacobian() flag is initialized to true.

References NOX::DeepCopy, f_vec_, gradient_vec_, inv_right_weight_vec_, Teuchos::is_null(), is_null(), lop_, losb_, lows_factory_, newton_vec_, nonnull(), prec_, prec_factory_, Teuchos::RCP< T >::ptr(), Teuchos::rcp(), resetIsValidFlags(), right_weight_vec_, rightScalingFirst_, scaled_x_vec_, NOX::ShapeCopy, shared_jacobian_, TEUCHOS_TEST_FOR_EXCEPTION, weight_vec_, and x_vec_.

Member Function Documentation

NOX::Abstract::Group::ReturnType NOX::Thyra::Group::applyJacobian ( const NOX::Abstract::Vector input,
NOX::Abstract::Vector result 
) const
virtual

Applies Jacobian to the given input vector and puts the answer in the result.

Computes

\[ v = J u, \]

where $ J$ is the Jacobian, $ u$ is the input vector, and $ v$ is the result vector.

Returns

Reimplemented from NOX::Abstract::Group.

References Teuchos::dyn_cast().

NOX::Abstract::Group::ReturnType NOX::Thyra::Group::applyJacobianInverse ( Teuchos::ParameterList params,
const NOX::Abstract::Vector input,
NOX::Abstract::Vector result 
) const
virtual

Applies the inverse of the Jacobian matrix to the given input vector and puts the answer in result.

Computes

\[ v = J^{-1} u, \]

where $ J$ is the Jacobian, $ u$ is the input vector, and $ v$ is the result vector.

The "Tolerance" parameter specifies that the solution should be such that

\[ \frac{\| J v - u \|_2}{\max \{ 1, \|u\|_2\} } < \mbox{Tolerance} \]

Returns

The parameter "Tolerance" may be added/modified in the list of parameters - this is the ideal solution tolerance for an iterative linear solve.

Reimplemented from NOX::Abstract::Group.

References Teuchos::dyn_cast().

NOX::Abstract::Group::ReturnType NOX::Thyra::Group::applyJacobianInverseMultiVector ( Teuchos::ParameterList params,
const NOX::Abstract::MultiVector input,
NOX::Abstract::MultiVector result 
) const
virtual

applyJacobianInverse for multiple right-hand sides

The default implementation here calls applyJacobianInverse() for each right hand side serially but should be overloaded if a block solver is available.

Reimplemented from NOX::Abstract::Group.

References Teuchos::dyn_cast(), and NOX::Thyra::MultiVector::getThyraMultiVector().

NOX::Abstract::Group::ReturnType NOX::Thyra::Group::applyJacobianMultiVector ( const NOX::Abstract::MultiVector input,
NOX::Abstract::MultiVector result 
) const
virtual

applyJacobian for multiple right-hand sides

The default implementation here calls applyJacobian() for each right hand side serially but should be overloaded if a block method is available.

Reimplemented from NOX::Abstract::Group.

References Teuchos::dyn_cast(), NOX::Thyra::MultiVector::getThyraMultiVector(), nonnull(), NOX::Abstract::Group::Ok, Teuchos::RCP< T >::ptr(), and TEUCHOS_TEST_FOR_EXCEPTION.

NOX::Abstract::Group::ReturnType NOX::Thyra::Group::applyJacobianTranspose ( const NOX::Abstract::Vector input,
NOX::Abstract::Vector result 
) const
virtual

Applies Jacobian-Transpose to the given input vector and puts the answer in the result.

Computes

\[ v = J^T u, \]

where $ J$ is the Jacobian, $ u$ is the input vector, and $ v$ is the result vector.

Returns

Reimplemented from NOX::Abstract::Group.

References Teuchos::dyn_cast().

NOX::Abstract::Group::ReturnType NOX::Thyra::Group::applyJacobianTransposeMultiVector ( const NOX::Abstract::MultiVector input,
NOX::Abstract::MultiVector result 
) const
virtual

applyJacobianTranspose for multiple right-hand sides

The default implementation here calls applyJacobianTranspose() for each right hand side serially but should be overloaded if a block method is available.

Reimplemented from NOX::Abstract::Group.

References Teuchos::dyn_cast(), NOX::Abstract::Group::Failed, NOX::Thyra::MultiVector::getThyraMultiVector(), nonnull(), NOX::Abstract::Group::Ok, Teuchos::RCP< T >::ptr(), and TEUCHOS_TEST_FOR_EXCEPTION.

NOX::Abstract::Group::ReturnType NOX::Thyra::Group::applyRightPreconditioning ( bool  useTranspose,
Teuchos::ParameterList params,
const NOX::Abstract::Vector input,
NOX::Abstract::Vector result 
) const
virtual

Apply right preconditiong to the given input vector.

Let $ M$ be a right preconditioner for the Jacobian $ J$; in other words, $ M$ is a matrix such that

\[ JM \approx I. \]

Compute

\[ u = M^{-1} v, \]

where $ u$ is the input vector and $ v$ is the result vector.

If useTranspose is true, then the transpose of the preconditioner is applied:

\[ u = {M^{-1}}^T v, \]

The transpose preconditioner is currently only required for Tensor methods.

The "Tolerance" parameter specifies that the solution should be such that

\[ \frac{\| M v - u \|_2}{\max \{ 1, \|u\|_2\} } < \mbox{Tolerance} \]

Returns

The parameters are from the "Linear %Solver" sublist of the "Direction" sublist that is passed to solver during construction.

Reimplemented from NOX::Abstract::Group.

References NOX::Thyra::Vector::getThyraRCPVector(), is_null(), nonnull(), NOX::Abstract::Group::Ok, and Teuchos::RCP< T >::ptr().

Teuchos::RCP< NOX::Abstract::Group > NOX::Thyra::Group::clone ( NOX::CopyType  type = NOX::DeepCopy) const
virtual

Create a new Group of the same derived type as this one by cloning this one, and return a ref count pointer to the new group.

If type is NOX::DeepCopy, then we need to create an exact replica of "this". Otherwise, if type is NOX::ShapeCopy, we need only replicate the shape of "this" (only the memory is allocated, the values are not copied into the vectors and Jacobian). Returns NULL if clone is not supported.

Note
Any shared data should have its ownership transfered to this group from the source for a NOX::DeepCopy.

Implements NOX::Abstract::Group.

Reimplemented in LOCA::Thyra::Group, and LOCA::Thyra::GroupWrapper.

References Teuchos::rcp().

NOX::Abstract::Group::ReturnType NOX::Thyra::Group::computeF ( )
virtual

Compute and store F(x).

Note
It's generally useful to also compute and store the 2-norm of F(x) at this point for later access by the getNormF() function.
Returns

Implements NOX::Abstract::Group.

References NOX::Abstract::Group::Failed, and NOX::Abstract::Group::Ok.

NOX::Abstract::Group::ReturnType NOX::Thyra::Group::computeGradient ( )
virtual

Compute and store gradient.

We can pose the nonlinear equation problem $ F(x) = 0$ as an optimization problem as follows:

\[ \min f(x) \equiv \frac{1}{2} \|F(x)\|_2^2. \]

In that case, the gradient (of $ f$) is defined as

\[ g \equiv J^T F. \]

Returns

Reimplemented from NOX::Abstract::Group.

References NOX::Abstract::Group::Failed, nonnull(), and TEUCHOS_TEST_FOR_EXCEPTION.

NOX::Abstract::Group::ReturnType NOX::Thyra::Group::computeJacobian ( )
virtual

Compute and store Jacobian.

Recall that

\[ F(x) = \left[ \begin{array}{c} F_1(x) \\ F_2(x) \\ \vdots \\ F_n(x) \\ \end{array} \right]. \]

The Jacobian is denoted by $ J$ and defined by

\[ J_{ij} = \frac{\partial F_i}{\partial x_j} (x). \]

Note
If this is a shared object, this group should take ownership of the Jacobian before it computes it.
Returns

Reimplemented from NOX::Abstract::Group.

References NOX::Abstract::Group::Failed, nonnull(), and NOX::Abstract::Group::Ok.

NOX::Abstract::Group::ReturnType NOX::Thyra::Group::computeNewton ( Teuchos::ParameterList params)
virtual

Compute the Newton direction, using parameters for the linear solve.

The Newton direction is the solution, s, of

\[ J s = -F. \]

The parameters are from the "Linear %Solver" sublist of the "Direction" sublist that is passed to solver during construction.

The "Tolerance" parameter may be added/modified in the sublist of "Linear Solver" parameters that is passed into this function. The solution should be such that

\[ \frac{\| J s - (-F) \|_2}{\max \{ 1, \|F\|_2\} } < \mbox{Tolerance} \]

Returns

Reimplemented from NOX::Abstract::Group.

void NOX::Thyra::Group::computeX ( const NOX::Abstract::Group grp,
const NOX::Abstract::Vector d,
double  step 
)
virtual

Compute x = grp.x + step * d.

Let $ x$ denote this group's solution vector. Let $\hat x$ denote the result of grp.getX(). Then set

\[ x = \hat x + \mbox{step} \; d. \]

Note
This should invalidate the function value, Jacobian, gradient, and Newton direction.
Throw an error if the copy fails.
Returns
Reference to this object

Implements NOX::Abstract::Group.

void NOX::Thyra::Group::enablePseudoTransientTerms ( const Teuchos::RCP< const ::Thyra::VectorBase< double >> &  x_dot,
const double  alpha,
const double  beta,
const double  t 
)

Set the transient terms on the Group and use them in the underlying evalModelImpl() calls.

Parameters
x_dotTime derivative term. Can be set to null.
alphaModel evaluator transient timer derivative multiplier
betaModel evaluator transient Jacobian multiplier
tCurrent time value
double NOX::Thyra::Group::getNormF ( ) const
virtual

Return 2-norm of F(x).

In other words,

\[ \sqrt{\sum_{i=1}^n F_i^2} \]

Implements NOX::Abstract::Group.

const NOX::Abstract::Vector & NOX::Thyra::Group::getScaledX ( ) const
virtual

If right scaling vector exist, return a right scaled vector.

Note
Default to getX

Reimplemented from NOX::Abstract::Group.

References nonnull().

bool NOX::Thyra::Group::isGradient ( ) const
virtual

Return true if the gradient is valid.

Note
Default implementation in NOX::Abstract::Group returns false.

Reimplemented from NOX::Abstract::Group.

bool NOX::Thyra::Group::isJacobian ( ) const
virtual

Return true if the Jacobian is valid.

Note
Default implementation in NOX::Abstract::Group returns false.

Reimplemented from NOX::Abstract::Group.

References nonnull().

Referenced by Group().

bool NOX::Thyra::Group::isNewton ( ) const
virtual

Return true if the Newton direction is valid.

Note
Default implementation in NOX::Abstract::Group returns false.

Reimplemented from NOX::Abstract::Group.

NOX::Abstract::Group & NOX::Thyra::Group::operator= ( const NOX::Abstract::Group source)
virtual

Copies the source group into this group.

Note
Any shared data owned by the source should have its ownership transfered to this group. This may result in a secret modification to the source object.

Implements NOX::Abstract::Group.

Referenced by LOCA::Thyra::Group::operator=().

void NOX::Thyra::Group::setBasePoint ( const ::Thyra::ModelEvaluatorBase::InArgs< double > &  base_point_params)

Set default parameters to be used with inArgs. This is a hack to support a bad design in PIRO. PIRO wraps nox as a model evaluator. The PIRO inArgs suppots all inargs, but nox only knows about nonlinear solver realted inputs. It knows nothing about parameters. PIRO should have a model evaluator wrapper that populates the extra inArgs as needed to override nominal values. Instead, we have to hard code storage of all in arg parameters in nox groups. This is only used in the Thyra::NOXNonlinearSolver and only when wrapped within a PIRO model evaluator that calls the setBasePoint on Thyra::NOXNonlinearSolver. This function should eventually be deprecated in favor of addign a wrapper ME to the piro nox solver class.

void NOX::Thyra::Group::setX ( const NOX::Abstract::Vector y)
virtual

Set the solution vector x to y.

Note
This should invalidate the function value, Jacobian, gradient, and Newton direction.
Throw an error if the copy fails.
Returns
Reference to this object

Implements NOX::Abstract::Group.

Referenced by Thyra::NOXNonlinearSolver::solve().

Member Data Documentation

Teuchos::RCP<const ::Thyra::VectorBase<double> > NOX::Thyra::Group::weight_vec_
protected

Optional wieghting vector for function scaling. NOX assumes that this vector can be updated in between nonlinear iterations.

This is pulled out of the initial guess vector

Referenced by Group().


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