Thyra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | List of all members
Thyra::ModelEvaluatorBase::OutArgs< Scalar > Class Template Reference

Concrete aggregate class for all output arguments computable by a ModelEvaluator subclass object. More...

#include <Thyra_ModelEvaluatorBase_decl.hpp>

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

Public Member Functions

 OutArgs ()
 
int Np () const
 Return the number of parameter subvectors p(l) supported (Np >= 0). More...
 
int Ng () const
 Return the number of axillary response functions g(j)(...) supported (Ng >= 0). More...
 
bool supports (EOutArgsMembers arg) const
 Determine if an input argument is supported or not. More...
 
const DerivativeSupportsupports (EOutArgsDfDp arg, int l) const
 Determine if DfDp(l) is supported or not, where 0 <= l && l < Np(). More...
 
const DerivativeSupportsupports (EOutArgsDgDx_dot arg, int j) const
 Determine if DgDx_dot(j) is supported or not, 0 <= j && j < Ng(). More...
 
const DerivativeSupportsupports (EOutArgsDgDx arg, int j) const
 Determine if DgDx(j) is supported or not, 0 <= j && j < Ng(). More...
 
const DerivativeSupportsupports (EOutArgsDgDp arg, int j, int l) const
 Determine if DgDp(j,l) is supported or not, 0 <= j && j < Ng() and 0 <= l && l < Np(). More...
 
bool supports (EOutArgs_hess_vec_prod_f_xx arg) const
 Determine if hess_vec_prod_f_xx is supported or not. More...
 
bool supports (EOutArgs_hess_vec_prod_f_xp arg, int l) const
 Determine if hess_vec_prod_f_xp(l) is supported or not, 0 <= l && l < Np(). More...
 
bool supports (EOutArgs_hess_vec_prod_f_px arg, int l) const
 Determine if hess_vec_prod_f_px(l) is supported or not, 0 <= l && l < Np(). More...
 
bool supports (EOutArgs_hess_vec_prod_f_pp arg, int l1, int l2) const
 Determine if hess_vec_prod_f_pp(l1,l2) is supported or not, 0 <= l1 && l1 < Np() and 0 <= l2 && l2 < Np(). More...
 
bool supports (EOutArgs_hess_vec_prod_g_xx arg, int j) const
 Determine if hess_vec_prod_g_xx(j,l) is supported or not, 0 <= j && j < Ng(). More...
 
bool supports (EOutArgs_hess_vec_prod_g_xp arg, int j, int l) const
 Determine if hess_vec_prod_g_xp(j,l) is supported or not, 0 <= j && j < Ng() and 0 <= l && l < Np(). More...
 
bool supports (EOutArgs_hess_vec_prod_g_px arg, int j, int l) const
 Determine if hess_vec_prod_g_px(j,l) is supported or not, 0 <= j && j < Ng() and 0 <= l && l < Np(). More...
 
bool supports (EOutArgs_hess_vec_prod_g_pp arg, int j, int l1, int l2) const
 Determine if hess_vec_prod_g_pp(j,l1,l2) is supported or not, 0 <= j && j < Ng(), 0 <= l1 && l1 < Np(), and 0 <= l2 && l2 < Np(). More...
 
bool supports (EOutArgs_hess_f_xx arg) const
 Determine if hess_f_xx is supported or not. More...
 
bool supports (EOutArgs_hess_f_xp arg, int l) const
 Determine if hess_f_xp(l) is supported or not, 0 <= l && l < Np(). More...
 
bool supports (EOutArgs_hess_f_pp arg, int l1, int l2) const
 Determine if hess_f_pp(l1,l2) is supported or not, 0 <= l1 && l1 < Np() and 0 <= l2 && l2 < Np(). More...
 
bool supports (EOutArgs_hess_g_xx arg, int j) const
 Determine if hess_g_xx(j,l) is supported or not, 0 <= j && j < Ng(). More...
 
bool supports (EOutArgs_hess_g_xp arg, int j, int l) const
 Determine if hess_g_xp(j,l) is supported or not, 0 <= j && j < Ng() and 0 <= l && l < Np(). More...
 
bool supports (EOutArgs_hess_g_pp arg, int j, int l1, int l2) const
 Determine if hess_g_pp(j,l1,l2) is supported or not, 0 <= j && j < Ng(), 0 <= l1 && l1 < Np(), and 0 <= l2 && l2 < Np(). More...
 
bool supports (EOutArgs_H_xx arg) const
 Determine if H_xx is supported or not. More...
 
bool supports (EOutArgs_H_xp arg, int l) const
 Determine if H_xp(l) is supported or not, 0 <= l && l < Np(). More...
 
bool supports (EOutArgs_H_pp arg, int l1, int l2) const
 Determine if H_pp(l1,l2) is supported or not, 0 <= l1 && l1 < Np() and 0 <= l2 && l2 < Np(). More...
 
void set_f (const Evaluation< VectorBase< Scalar > > &f)
 Precondition: supports(OUT_ARG_f)==true. More...
 
Evaluation< VectorBase< Scalar > > get_f () const
 Precondition: supports(OUT_ARG_f)==true. More...
 
void set_g (int j, const Evaluation< VectorBase< Scalar > > &g_j)
 Precondition: supports(OUT_ARG_g)==true. More...
 
Evaluation< VectorBase< Scalar > > get_g (int j) const
 Precondition: supports(OUT_ARG_g)==true.. More...
 
void set_W (const RCP< LinearOpWithSolveBase< Scalar > > &W)
 Precondition: supports(OUT_ARG_W)==true. More...
 
RCP< LinearOpWithSolveBase
< Scalar > > 
get_W () const
 Precondition: supports(OUT_ARG_W)==true. More...
 
template<typename ObjectType >
bool supports () const
 Determines if an extended output argument of type ObjectType is supported. More...
 
template<typename ObjectType >
void set (const RCP< const ObjectType > &uo)
 Set an extended output argument of type ObjectType in OutArgs. Precondition: supports()==true. More...
 
template<typename ObjectType >
RCP< const ObjectType > get () const
 Get an extended output argument of type ObjectType from OutArgs. Precondition: supports()==true. More...
 
void set_f_mp (const RCP< Stokhos::ProductEpetraVector > &f_mp)
 Precondition: supports(OUT_ARG_f_mp)==true. More...
 
RCP< Stokhos::ProductEpetraVector > get_f_mp () const
 Precondition: supports(OUT_ARG_f_mp)==true. More...
 
void set_g_mp (int j, const RCP< Stokhos::ProductEpetraVector > &g_mp_j)
 Precondition: supports(OUT_ARG_g_mp)==true. More...
 
RCP< Stokhos::ProductEpetraVector > get_g_mp (int j) const
 Precondition: supports(OUT_ARG_g_mp)==true.. More...
 
void set_W_mp (const RCP< Stokhos::ProductEpetraOperator > &W_mp)
 Precondition: supports(OUT_ARG_W_mp)==true. More...
 
RCP
< Stokhos::ProductEpetraOperator > 
get_W_mp () const
 Precondition: supports(OUT_ARG_W_mp)==true. More...
 
void set_W_op (const RCP< LinearOpBase< Scalar > > &W_op)
 Precondition: supports(OUT_ARG_W_op)==true. More...
 
RCP< LinearOpBase< Scalar > > get_W_op () const
 Precondition: supports(OUT_ARG_W_op)==true. More...
 
void set_W_prec (const RCP< PreconditionerBase< Scalar > > &W_prec)
 Precondition: supports(OUT_ARG_W_op)==true. More...
 
RCP< PreconditionerBase< Scalar > > get_W_prec () const
 Precondition: supports(OUT_ARG_W_op)==true. More...
 
DerivativeProperties get_W_properties () const
 Return the known properties of W (precondition: supports(OUT_ARG_f)==true). More...
 
void set_DfDp (int l, const Derivative< Scalar > &DfDp_l)
 Precondition: supports(OUT_ARG_DfDp,l)==true. More...
 
Derivative< Scalar > get_DfDp (int l) const
 Precondition: supports(OUT_ARG_DfDp,l)==true. More...
 
DerivativeProperties get_DfDp_properties (int l) const
 Return the know properties of DfDp(l) (precondition: supports(OUT_ARG_DfDp,l)==true). More...
 
void set_DgDx_dot (int j, const Derivative< Scalar > &DgDx_dot_j)
 Precondition: supports(OUT_ARG_DgDx_dot,j)==true. More...
 
Derivative< Scalar > get_DgDx_dot (int j) const
 Precondition: supports(OUT_ARG_DgDx_dot,j)==true. More...
 
DerivativeProperties get_DgDx_dot_properties (int j) const
 Return the know properties of DgDx_dot(j) (precondition: supports(OUT_ARG_DgDx_dot,j)==true). More...
 
void set_DgDx (int j, const Derivative< Scalar > &DgDx_j)
 Precondition: supports(OUT_ARG_DgDx,j)==true. More...
 
Derivative< Scalar > get_DgDx (int j) const
 Precondition: supports(OUT_ARG_DgDx,j)==true. More...
 
DerivativeProperties get_DgDx_properties (int j) const
 Return the know properties of DgDx(j) (precondition: supports(OUT_ARG_DgDx,j)==true). More...
 
void set_DgDp (int j, int l, const Derivative< Scalar > &DgDp_j_l)
 Precondition: supports(OUT_ARG_DgDp,j,l)==true. More...
 
Derivative< Scalar > get_DgDp (int j, int l) const
 Precondition: supports(OUT_ARG_DgDp,j,l)==true. More...
 
DerivativeProperties get_DgDp_properties (int j, int l) const
 Return the know properties of DgDp(j,l) (precondition: supports(OUT_ARG_DgDp,j,l)==true). More...
 
void set_hess_vec_prod_f_xx (const RCP< MultiVectorBase< Scalar > > &hess_vec_prod_f_xx)
 Precondition: supports(OUT_ARG_hess_vec_prod_f_xx)==true. More...
 
void set_hess_vec_prod_f_xp (int l, const RCP< MultiVectorBase< Scalar > > &hess_vec_prod_f_xp_l)
 Precondition: supports(OUT_ARG_hess_vec_prod_f_xp,l)==true. More...
 
void set_hess_vec_prod_f_px (int l, const RCP< MultiVectorBase< Scalar > > &hess_vec_prod_f_px_l)
 Precondition: supports(OUT_ARG_hess_vec_prod_f_px,l)==true. More...
 
void set_hess_vec_prod_f_pp (int l1, int l2, const RCP< MultiVectorBase< Scalar > > &hess_vec_prod_f_pp_l1_l2)
 Precondition: supports(OUT_ARG_hess_vec_prod_f_pp,l1,l2)==true. More...
 
RCP< MultiVectorBase< Scalar > > get_hess_vec_prod_f_xx () const
 Precondition: supports(OUT_ARG_hess_vec_prod_f_xx)==true. More...
 
RCP< MultiVectorBase< Scalar > > get_hess_vec_prod_f_xp (int l) const
 Precondition: supports(OUT_ARG_hess_vec_prod_f_xp,l)==true. More...
 
RCP< MultiVectorBase< Scalar > > get_hess_vec_prod_f_px (int l) const
 Precondition: supports(OUT_ARG_hess_vec_prod_f_px,l)==true. More...
 
RCP< MultiVectorBase< Scalar > > get_hess_vec_prod_f_pp (int l1, int l2) const
 Precondition: supports(OUT_ARG_hess_vec_prod_f_pp,l1,l2)==true. More...
 
void set_hess_vec_prod_g_xx (int j, const RCP< MultiVectorBase< Scalar > > &hess_vec_prod_g_xx_j)
 Precondition: supports(OUT_ARG_hess_vec_prod_g_xx,j)==true. More...
 
void set_hess_vec_prod_g_xp (int j, int l, const RCP< MultiVectorBase< Scalar > > &hess_vec_prod_g_xp_j_l)
 Precondition: supports(OUT_ARG_hess_vec_prod_g_xp,j,l)==true. More...
 
void set_hess_vec_prod_g_px (int j, int l, const RCP< MultiVectorBase< Scalar > > &hess_vec_prod_g_px_j_l)
 Precondition: supports(OUT_ARG_hess_vec_prod_g_px,j,l)==true. More...
 
void set_hess_vec_prod_g_pp (int j, int l1, int l2, const RCP< MultiVectorBase< Scalar > > &hess_vec_prod_g_pp_j_l1_l2)
 Precondition: supports(OUT_ARG_hess_vec_prod_g_pp,j,l1,l2)==true. More...
 
RCP< MultiVectorBase< Scalar > > get_hess_vec_prod_g_xx (int j) const
 Precondition: supports(OUT_ARG_hess_vec_prod_g_xx,j)==true. More...
 
RCP< MultiVectorBase< Scalar > > get_hess_vec_prod_g_xp (int j, int l) const
 Precondition: supports(OUT_ARG_hess_vec_prod_g_xp,j,l)==true. More...
 
RCP< MultiVectorBase< Scalar > > get_hess_vec_prod_g_px (int j, int l) const
 Precondition: supports(OUT_ARG_hess_vec_prod_g_px,j,l)==true. More...
 
RCP< MultiVectorBase< Scalar > > get_hess_vec_prod_g_pp (int j, int l1, int l2) const
 Precondition: supports(OUT_ARG_hess_vec_prod_g_pp,j,l1,l2)==true. More...
 
void set_hess_f_xx (const RCP< LinearOpBase< Scalar > > &hess_f_xx)
 Precondition: supports(OUT_ARG_hess_f_xx)==true. More...
 
void set_hess_f_xp (int l, const RCP< LinearOpBase< Scalar > > &hess_f_xp_l)
 Precondition: supports(OUT_ARG_hess_f_xp,l)==true. More...
 
void set_hess_f_pp (int l1, int l2, const RCP< LinearOpBase< Scalar > > &hess_f_pp_l1_l2)
 Precondition: supports(OUT_ARG_hess_f_pp,l1,l2)==true. More...
 
void set_hess_g_xx (int j, const RCP< LinearOpBase< Scalar > > &hess_g_xx_j)
 Precondition: supports(OUT_ARG_hess_g_xx,j)==true. More...
 
void set_hess_g_xp (int j, int l, const RCP< LinearOpBase< Scalar > > &hess_g_xp_j_l)
 Precondition: supports(OUT_ARG_hess_g_xp,j,l)==true. More...
 
void set_hess_g_pp (int j, int l1, int l2, const RCP< LinearOpBase< Scalar > > &hess_g_pp_j_l1_l2)
 Precondition: supports(OUT_ARG_hess_g_pp,j,l1,l2)==true. More...
 
void set_H_xx (const RCP< LinearOpBase< Scalar > > &H_xx)
 Precondition: supports(OUT_ARG_H_xx)==true. More...
 
void set_H_xp (int l, const RCP< LinearOpBase< Scalar > > &H_xp_l)
 Precondition: supports(OUT_ARG_H_xp,l)==true. More...
 
void set_H_pp (int l1, int l2, const RCP< LinearOpBase< Scalar > > &H_pp_l1_l2)
 Precondition: supports(OUT_ARG_H_pp,l1,l2)==true. More...
 
RCP< LinearOpBase< Scalar > > get_hess_f_xx () const
 Precondition: supports(OUT_ARG_hess_f_xx)==true. More...
 
RCP< LinearOpBase< Scalar > > get_hess_f_xp (int l) const
 Precondition: supports(OUT_ARG_hess_f_xp,l)==true. More...
 
RCP< LinearOpBase< Scalar > > get_hess_f_pp (int l1, int l2) const
 Precondition: supports(OUT_ARG_hess_f_pp,l1,l2)==true. More...
 
RCP< LinearOpBase< Scalar > > get_hess_g_xx (int j) const
 Precondition: supports(OUT_ARG_hess_g_xx,j)==true. More...
 
RCP< LinearOpBase< Scalar > > get_hess_g_xp (int j, int l) const
 Precondition: supports(OUT_ARG_hess_g_xp,j,l)==true. More...
 
RCP< LinearOpBase< Scalar > > get_hess_g_pp (int j, int l1, int l2) const
 Precondition: supports(OUT_ARG_hess_g_pp,j,l1,l2)==true. More...
 
RCP< LinearOpBase< Scalar > > get_H_xx () const
 Precondition: supports(OUT_ARG_H_xx)==true. More...
 
RCP< LinearOpBase< Scalar > > get_H_xp (int l) const
 Precondition: supports(OUT_ARG_H_xp,l)==true. More...
 
RCP< LinearOpBase< Scalar > > get_H_pp (int l1, int l2) const
 Precondition: supports(OUT_ARG_H_pp,l1,l2)==true. More...
 
void setArgs (const OutArgs< Scalar > &outArgs, bool ignoreUnsupported=false)
 Set all arguments fron outArgs into *this. More...
 
void setFailed () const
 Set that the evaluation as a whole failed. More...
 
bool isFailed () const
 Return if the evaluation failed or not. More...
 
bool isEmpty () const
 
void assertSameSupport (const OutArgs< Scalar > &outArgs) const
 Assert that two OutArgs objects have the same support. More...
 
std::string modelEvalDescription () const
 
std::string description () const
 
void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
 Create a more detailed description along about this object and the ModelEvaluator that created it. More...
 

Protected Member Functions

void _setModelEvalDescription (const std::string &modelEvalDescription)
 
void _set_Np_Ng (int Np, int Ng)
 
void _setSupports (EOutArgsMembers arg, bool supports)
 
template<typename ObjectType >
void _setSupports (const bool supports)
 
void _setSupports (EOutArgsDfDp arg, int l, const DerivativeSupport &)
 
void _setSupports (EOutArgsDgDx_dot arg, int j, const DerivativeSupport &)
 
void _setSupports (EOutArgsDgDx arg, int j, const DerivativeSupport &)
 
void _setSupports (EOutArgsDgDp arg, int j, int l, const DerivativeSupport &)
 
void _setSupports (EOutArgs_hess_vec_prod_f_xx arg, bool supports)
 
void _setSupports (EOutArgs_hess_vec_prod_f_xp arg, int l, bool supports)
 
void _setSupports (EOutArgs_hess_vec_prod_f_px arg, int l, bool supports)
 
void _setSupports (EOutArgs_hess_vec_prod_f_pp arg, int l1, int l2, bool supports)
 
void _setSupports (EOutArgs_hess_vec_prod_g_xx arg, int j, bool supports)
 
void _setSupports (EOutArgs_hess_vec_prod_g_xp arg, int j, int l, bool supports)
 
void _setSupports (EOutArgs_hess_vec_prod_g_px arg, int j, int l, bool supports)
 
void _setSupports (EOutArgs_hess_vec_prod_g_pp arg, int j, int l1, int l2, bool supports)
 
void _setSupports (EOutArgs_hess_f_xx arg, bool supports)
 
void _setSupports (EOutArgs_hess_f_xp arg, int l, bool supports)
 
void _setSupports (EOutArgs_hess_f_pp arg, int l1, int l2, bool supports)
 
void _setSupports (EOutArgs_hess_g_xx arg, int j, bool supports)
 
void _setSupports (EOutArgs_hess_g_xp arg, int j, int l, bool supports)
 
void _setSupports (EOutArgs_hess_g_pp arg, int j, int l1, int l2, bool supports)
 
void _setSupports (EOutArgs_H_xx arg, bool supports)
 
void _setSupports (EOutArgs_H_xp arg, int l, bool supports)
 
void _setSupports (EOutArgs_H_pp arg, int l1, int l2, bool supports)
 
void _set_W_properties (const DerivativeProperties &properties)
 
void _set_DfDp_properties (int l, const DerivativeProperties &properties)
 
void _set_DgDx_dot_properties (int j, const DerivativeProperties &properties)
 
void _set_DgDx_properties (int j, const DerivativeProperties &properties)
 
void _set_DgDp_properties (int j, int l, const DerivativeProperties &properties)
 
void _setSupports (const OutArgs< Scalar > &inputOutArgs)
 
void _setUnsupportsAndRelated (EInArgsMembers arg)
 
void _setUnsupportsAndRelated (EOutArgsMembers arg)
 
void _setHessianSupports (const bool supports)
 

Detailed Description

template<class Scalar>
class Thyra::ModelEvaluatorBase::OutArgs< Scalar >

Concrete aggregate class for all output arguments computable by a ModelEvaluator subclass object.

Note that const OutArgs object means that a client can not change what output objects are being pointed to but they can still change the states of the contained objects. This is slight variation on the concept of logical const-ness in C++ but it is totally consistent with the vary nature of this class.

In addition to storing the output objects themselves, this class also allows the storage if the properties of some of the objects as well. Therefore, objects of this type are used to communicate a lot of different information about the output functions and derivatives supported by a model. It tells clients what functions and derivatives are supported and what the know properties are.

A client can not directly set what input arguments are supported or not supported. Only a subclass of ModelEvaluator can do that (through the OutArgsSetup subclass).

Extended types: The methods for supports(), set() and get() that are templated on the ObjectType are used to support what are called extended types. This functionality allows developers to inject new objects into the ModelEvaluator evaluations. This can be used for expermenting with new capabilities that require adding new/special arguments to the outArgs. This also can be used to reduce the clutter and complexity of the model evaluator outArgs object. For example, users could create their own outArgs for a stochastic computation:

struct StochasticOutArgs {
StochLib::MultiVector<Scalar> stochastic_f;
...
};

This object could then be used in a model evaluator without changing the base classes in Thyra:

auto outArgs = model->createOutArgs();
assert(outArgs.template supports<StochasticOutArgs<Scalar>>());
RCP<StochasticOutArgs> stochOutArgs = rcp(new StochasticOutArgs<Scalar>);
stochOutArgs->stochastic_f = ... // create and set objects
...
outArgs.set(stochOutArgs);

Definition at line 900 of file Thyra_ModelEvaluatorBase_decl.hpp.

Constructor & Destructor Documentation

template<class Scalar >
Thyra::ModelEvaluatorBase::OutArgs< Scalar >::OutArgs ( )

Definition at line 860 of file Thyra_ModelEvaluatorBase_def.hpp.

Member Function Documentation

template<class Scalar >
int Thyra::ModelEvaluatorBase::OutArgs< Scalar >::Np ( ) const

Return the number of parameter subvectors p(l) supported (Np >= 0).

Definition at line 870 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
int Thyra::ModelEvaluatorBase::OutArgs< Scalar >::Ng ( ) const

Return the number of axillary response functions g(j)(...) supported (Ng >= 0).

Definition at line 875 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgsMembers  arg) const

Determine if an input argument is supported or not.

Definition at line 880 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
const ModelEvaluatorBase::DerivativeSupport & Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgsDfDp  arg,
int  l 
) const

Determine if DfDp(l) is supported or not, where 0 <= l && l < Np().

Definition at line 895 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
const ModelEvaluatorBase::DerivativeSupport & Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgsDgDx_dot  arg,
int  j 
) const

Determine if DgDx_dot(j) is supported or not, 0 <= j && j < Ng().

Definition at line 906 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
const ModelEvaluatorBase::DerivativeSupport & Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgsDgDx  arg,
int  j 
) const

Determine if DgDx(j) is supported or not, 0 <= j && j < Ng().

Definition at line 917 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
const ModelEvaluatorBase::DerivativeSupport & Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgsDgDp  arg,
int  j,
int  l 
) const

Determine if DgDp(j,l) is supported or not, 0 <= j && j < Ng() and 0 <= l && l < Np().

Definition at line 928 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgs_hess_vec_prod_f_xx  arg) const

Determine if hess_vec_prod_f_xx is supported or not.

Definition at line 939 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgs_hess_vec_prod_f_xp  arg,
int  l 
) const

Determine if hess_vec_prod_f_xp(l) is supported or not, 0 <= l && l < Np().

Definition at line 947 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgs_hess_vec_prod_f_px  arg,
int  l 
) const

Determine if hess_vec_prod_f_px(l) is supported or not, 0 <= l && l < Np().

Definition at line 956 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgs_hess_vec_prod_f_pp  arg,
int  l1,
int  l2 
) const

Determine if hess_vec_prod_f_pp(l1,l2) is supported or not, 0 <= l1 && l1 < Np() and 0 <= l2 && l2 < Np().

Definition at line 965 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgs_hess_vec_prod_g_xx  arg,
int  j 
) const

Determine if hess_vec_prod_g_xx(j,l) is supported or not, 0 <= j && j < Ng().

Definition at line 975 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgs_hess_vec_prod_g_xp  arg,
int  j,
int  l 
) const

Determine if hess_vec_prod_g_xp(j,l) is supported or not, 0 <= j && j < Ng() and 0 <= l && l < Np().

Definition at line 984 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgs_hess_vec_prod_g_px  arg,
int  j,
int  l 
) const

Determine if hess_vec_prod_g_px(j,l) is supported or not, 0 <= j && j < Ng() and 0 <= l && l < Np().

Definition at line 994 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgs_hess_vec_prod_g_pp  arg,
int  j,
int  l1,
int  l2 
) const

Determine if hess_vec_prod_g_pp(j,l1,l2) is supported or not, 0 <= j && j < Ng(), 0 <= l1 && l1 < Np(), and 0 <= l2 && l2 < Np().

Definition at line 1004 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgs_hess_f_xx  arg) const

Determine if hess_f_xx is supported or not.

Definition at line 1015 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgs_hess_f_xp  arg,
int  l 
) const

Determine if hess_f_xp(l) is supported or not, 0 <= l && l < Np().

Definition at line 1023 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgs_hess_f_pp  arg,
int  l1,
int  l2 
) const

Determine if hess_f_pp(l1,l2) is supported or not, 0 <= l1 && l1 < Np() and 0 <= l2 && l2 < Np().

Definition at line 1032 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgs_hess_g_xx  arg,
int  j 
) const

Determine if hess_g_xx(j,l) is supported or not, 0 <= j && j < Ng().

Definition at line 1069 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgs_hess_g_xp  arg,
int  j,
int  l 
) const

Determine if hess_g_xp(j,l) is supported or not, 0 <= j && j < Ng() and 0 <= l && l < Np().

Definition at line 1078 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgs_hess_g_pp  arg,
int  j,
int  l1,
int  l2 
) const

Determine if hess_g_pp(j,l1,l2) is supported or not, 0 <= j && j < Ng(), 0 <= l1 && l1 < Np(), and 0 <= l2 && l2 < Np().

Definition at line 1088 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgs_H_xx  arg) const

Determine if H_xx is supported or not.

Definition at line 1042 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgs_H_xp  arg,
int  l 
) const

Determine if H_xp(l) is supported or not, 0 <= l && l < Np().

Definition at line 1050 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( EOutArgs_H_pp  arg,
int  l1,
int  l2 
) const

Determine if H_pp(l1,l2) is supported or not, 0 <= l1 && l1 < Np() and 0 <= l2 && l2 < Np().

Definition at line 1059 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_f ( const Evaluation< VectorBase< Scalar > > &  f)

Precondition: supports(OUT_ARG_f)==true.

Definition at line 1156 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
ModelEvaluatorBase::Evaluation< VectorBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_f ( ) const

Precondition: supports(OUT_ARG_f)==true.

Definition at line 1167 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_g ( int  j,
const Evaluation< VectorBase< Scalar > > &  g_j 
)

Precondition: supports(OUT_ARG_g)==true.

Definition at line 1175 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
ModelEvaluatorBase::Evaluation< VectorBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_g ( int  j) const

Precondition: supports(OUT_ARG_g)==true..

Definition at line 1186 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_W ( const RCP< LinearOpWithSolveBase< Scalar > > &  W)

Precondition: supports(OUT_ARG_W)==true.

Definition at line 1232 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< LinearOpWithSolveBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_W ( ) const

Precondition: supports(OUT_ARG_W)==true.

Definition at line 1243 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
template<typename ObjectType >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::supports ( ) const

Determines if an extended output argument of type ObjectType is supported.

Definition at line 1715 of file Thyra_ModelEvaluatorBase_decl.hpp.

template<class Scalar >
template<typename ObjectType >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set ( const RCP< const ObjectType > &  uo)

Set an extended output argument of type ObjectType in OutArgs. Precondition: supports()==true.

Definition at line 1728 of file Thyra_ModelEvaluatorBase_decl.hpp.

template<class Scalar>
template<typename ObjectType >
RCP<const ObjectType> Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get ( ) const

Get an extended output argument of type ObjectType from OutArgs. Precondition: supports()==true.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_f_mp ( const RCP< Stokhos::ProductEpetraVector > &  f_mp)

Precondition: supports(OUT_ARG_f_mp)==true.

Definition at line 1194 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< Stokhos::ProductEpetraVector > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_f_mp ( ) const

Precondition: supports(OUT_ARG_f_mp)==true.

Definition at line 1205 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_g_mp ( int  j,
const RCP< Stokhos::ProductEpetraVector > &  g_mp_j 
)

Precondition: supports(OUT_ARG_g_mp)==true.

Definition at line 1213 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< Stokhos::ProductEpetraVector > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_g_mp ( int  j) const

Precondition: supports(OUT_ARG_g_mp)==true..

Definition at line 1224 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_W_mp ( const RCP< Stokhos::ProductEpetraOperator > &  W_mp)

Precondition: supports(OUT_ARG_W_mp)==true.

Definition at line 1251 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< Stokhos::ProductEpetraOperator > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_W_mp ( ) const

Precondition: supports(OUT_ARG_W_mp)==true.

Definition at line 1262 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_W_op ( const RCP< LinearOpBase< Scalar > > &  W_op)

Precondition: supports(OUT_ARG_W_op)==true.

Definition at line 1270 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< LinearOpBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_W_op ( ) const

Precondition: supports(OUT_ARG_W_op)==true.

Definition at line 1281 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_W_prec ( const RCP< PreconditionerBase< Scalar > > &  W_prec)

Precondition: supports(OUT_ARG_W_op)==true.

Definition at line 1289 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< PreconditionerBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_W_prec ( ) const

Precondition: supports(OUT_ARG_W_op)==true.

Definition at line 1300 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
ModelEvaluatorBase::DerivativeProperties Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_W_properties ( ) const

Return the known properties of W (precondition: supports(OUT_ARG_f)==true).

Definition at line 1309 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_DfDp ( int  l,
const Derivative< Scalar > &  DfDp_l 
)

Precondition: supports(OUT_ARG_DfDp,l)==true.

Definition at line 1317 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
ModelEvaluatorBase::Derivative< Scalar > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_DfDp ( int  l) const

Precondition: supports(OUT_ARG_DfDp,l)==true.

Definition at line 1328 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
ModelEvaluatorBase::DerivativeProperties Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_DfDp_properties ( int  l) const

Return the know properties of DfDp(l) (precondition: supports(OUT_ARG_DfDp,l)==true).

Definition at line 1337 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_DgDx_dot ( int  j,
const Derivative< Scalar > &  DgDx_dot_j 
)

Precondition: supports(OUT_ARG_DgDx_dot,j)==true.

Definition at line 1373 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
ModelEvaluatorBase::Derivative< Scalar > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_DgDx_dot ( int  j) const

Precondition: supports(OUT_ARG_DgDx_dot,j)==true.

Definition at line 1384 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
ModelEvaluatorBase::DerivativeProperties Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_DgDx_dot_properties ( int  j) const

Return the know properties of DgDx_dot(j) (precondition: supports(OUT_ARG_DgDx_dot,j)==true).

Definition at line 1393 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_DgDx ( int  j,
const Derivative< Scalar > &  DgDx_j 
)

Precondition: supports(OUT_ARG_DgDx,j)==true.

Definition at line 1429 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
ModelEvaluatorBase::Derivative< Scalar > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_DgDx ( int  j) const

Precondition: supports(OUT_ARG_DgDx,j)==true.

Definition at line 1440 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
ModelEvaluatorBase::DerivativeProperties Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_DgDx_properties ( int  j) const

Return the know properties of DgDx(j) (precondition: supports(OUT_ARG_DgDx,j)==true).

Definition at line 1449 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_DgDp ( int  j,
int  l,
const Derivative< Scalar > &  DgDp_j_l 
)

Precondition: supports(OUT_ARG_DgDp,j,l)==true.

Definition at line 1485 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
ModelEvaluatorBase::Derivative< Scalar > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_DgDp ( int  j,
int  l 
) const

Precondition: supports(OUT_ARG_DgDp,j,l)==true.

Definition at line 1496 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
ModelEvaluatorBase::DerivativeProperties Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_DgDp_properties ( int  j,
int  l 
) const

Return the know properties of DgDp(j,l) (precondition: supports(OUT_ARG_DgDp,j,l)==true).

Definition at line 1505 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_hess_vec_prod_f_xx ( const RCP< MultiVectorBase< Scalar > > &  hess_vec_prod_f_xx)

Precondition: supports(OUT_ARG_hess_vec_prod_f_xx)==true.

Definition at line 1541 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_hess_vec_prod_f_xp ( int  l,
const RCP< MultiVectorBase< Scalar > > &  hess_vec_prod_f_xp_l 
)

Precondition: supports(OUT_ARG_hess_vec_prod_f_xp,l)==true.

Definition at line 1550 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_hess_vec_prod_f_px ( int  l,
const RCP< MultiVectorBase< Scalar > > &  hess_vec_prod_f_px_l 
)

Precondition: supports(OUT_ARG_hess_vec_prod_f_px,l)==true.

Definition at line 1559 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_hess_vec_prod_f_pp ( int  l1,
int  l2,
const RCP< MultiVectorBase< Scalar > > &  hess_vec_prod_f_pp_l1_l2 
)

Precondition: supports(OUT_ARG_hess_vec_prod_f_pp,l1,l2)==true.

Definition at line 1568 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< MultiVectorBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_hess_vec_prod_f_xx ( ) const

Precondition: supports(OUT_ARG_hess_vec_prod_f_xx)==true.

Definition at line 1578 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< MultiVectorBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_hess_vec_prod_f_xp ( int  l) const

Precondition: supports(OUT_ARG_hess_vec_prod_f_xp,l)==true.

Definition at line 1586 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< MultiVectorBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_hess_vec_prod_f_px ( int  l) const

Precondition: supports(OUT_ARG_hess_vec_prod_f_px,l)==true.

Definition at line 1594 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< MultiVectorBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_hess_vec_prod_f_pp ( int  l1,
int  l2 
) const

Precondition: supports(OUT_ARG_hess_vec_prod_f_pp,l1,l2)==true.

Definition at line 1602 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_hess_vec_prod_g_xx ( int  j,
const RCP< MultiVectorBase< Scalar > > &  hess_vec_prod_g_xx_j 
)

Precondition: supports(OUT_ARG_hess_vec_prod_g_xx,j)==true.

Definition at line 1610 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_hess_vec_prod_g_xp ( int  j,
int  l,
const RCP< MultiVectorBase< Scalar > > &  hess_vec_prod_g_xp_j_l 
)

Precondition: supports(OUT_ARG_hess_vec_prod_g_xp,j,l)==true.

Definition at line 1619 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_hess_vec_prod_g_px ( int  j,
int  l,
const RCP< MultiVectorBase< Scalar > > &  hess_vec_prod_g_px_j_l 
)

Precondition: supports(OUT_ARG_hess_vec_prod_g_px,j,l)==true.

Definition at line 1628 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_hess_vec_prod_g_pp ( int  j,
int  l1,
int  l2,
const RCP< MultiVectorBase< Scalar > > &  hess_vec_prod_g_pp_j_l1_l2 
)

Precondition: supports(OUT_ARG_hess_vec_prod_g_pp,j,l1,l2)==true.

Definition at line 1637 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< MultiVectorBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_hess_vec_prod_g_xx ( int  j) const

Precondition: supports(OUT_ARG_hess_vec_prod_g_xx,j)==true.

Definition at line 1647 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< MultiVectorBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_hess_vec_prod_g_xp ( int  j,
int  l 
) const

Precondition: supports(OUT_ARG_hess_vec_prod_g_xp,j,l)==true.

Definition at line 1655 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< MultiVectorBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_hess_vec_prod_g_px ( int  j,
int  l 
) const

Precondition: supports(OUT_ARG_hess_vec_prod_g_px,j,l)==true.

Definition at line 1663 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< MultiVectorBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_hess_vec_prod_g_pp ( int  j,
int  l1,
int  l2 
) const

Precondition: supports(OUT_ARG_hess_vec_prod_g_pp,j,l1,l2)==true.

Definition at line 1671 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_hess_f_xx ( const RCP< LinearOpBase< Scalar > > &  hess_f_xx)

Precondition: supports(OUT_ARG_hess_f_xx)==true.

Definition at line 1678 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_hess_f_xp ( int  l,
const RCP< LinearOpBase< Scalar > > &  hess_f_xp_l 
)

Precondition: supports(OUT_ARG_hess_f_xp,l)==true.

Definition at line 1687 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_hess_f_pp ( int  l1,
int  l2,
const RCP< LinearOpBase< Scalar > > &  hess_f_pp_l1_l2 
)

Precondition: supports(OUT_ARG_hess_f_pp,l1,l2)==true.

Definition at line 1696 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_hess_g_xx ( int  j,
const RCP< LinearOpBase< Scalar > > &  hess_g_xx_j 
)

Precondition: supports(OUT_ARG_hess_g_xx,j)==true.

Definition at line 1729 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_hess_g_xp ( int  j,
int  l,
const RCP< LinearOpBase< Scalar > > &  hess_g_xp_j_l 
)

Precondition: supports(OUT_ARG_hess_g_xp,j,l)==true.

Definition at line 1738 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_hess_g_pp ( int  j,
int  l1,
int  l2,
const RCP< LinearOpBase< Scalar > > &  hess_g_pp_j_l1_l2 
)

Precondition: supports(OUT_ARG_hess_g_pp,j,l1,l2)==true.

Definition at line 1747 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_H_xx ( const RCP< LinearOpBase< Scalar > > &  H_xx)

Precondition: supports(OUT_ARG_H_xx)==true.

Definition at line 1780 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_H_xp ( int  l,
const RCP< LinearOpBase< Scalar > > &  H_xp_l 
)

Precondition: supports(OUT_ARG_H_xp,l)==true.

Definition at line 1789 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::set_H_pp ( int  l1,
int  l2,
const RCP< LinearOpBase< Scalar > > &  H_pp_l1_l2 
)

Precondition: supports(OUT_ARG_H_pp,l1,l2)==true.

Definition at line 1798 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< LinearOpBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_hess_f_xx ( ) const

Precondition: supports(OUT_ARG_hess_f_xx)==true.

Definition at line 1706 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< LinearOpBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_hess_f_xp ( int  l) const

Precondition: supports(OUT_ARG_hess_f_xp,l)==true.

Definition at line 1714 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< LinearOpBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_hess_f_pp ( int  l1,
int  l2 
) const

Precondition: supports(OUT_ARG_hess_f_pp,l1,l2)==true.

Definition at line 1722 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< LinearOpBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_hess_g_xx ( int  j) const

Precondition: supports(OUT_ARG_hess_g_xx,j)==true.

Definition at line 1757 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< LinearOpBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_hess_g_xp ( int  j,
int  l 
) const

Precondition: supports(OUT_ARG_hess_g_xp,j,l)==true.

Definition at line 1765 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< LinearOpBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_hess_g_pp ( int  j,
int  l1,
int  l2 
) const

Precondition: supports(OUT_ARG_hess_g_pp,j,l1,l2)==true.

Definition at line 1773 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< LinearOpBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_H_xx ( ) const

Precondition: supports(OUT_ARG_H_xx)==true.

Definition at line 1808 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< LinearOpBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_H_xp ( int  l) const

Precondition: supports(OUT_ARG_H_xp,l)==true.

Definition at line 1816 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
RCP< LinearOpBase< Scalar > > Thyra::ModelEvaluatorBase::OutArgs< Scalar >::get_H_pp ( int  l1,
int  l2 
) const

Precondition: supports(OUT_ARG_H_pp,l1,l2)==true.

Definition at line 1824 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::setArgs ( const OutArgs< Scalar > &  outArgs,
bool  ignoreUnsupported = false 
)

Set all arguments fron outArgs into *this.

If ignoreUnsupported==true, then arguments in outArgs that are not supported in *this will be ignored. Othereise, if an unsupported argument is set, then an exception will be thrown.

Definition at line 1855 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::setFailed ( ) const

Set that the evaluation as a whole failed.

Note that this function is declared as const even through it technically changes the state of *this object. This was done so that this property could be set by a ModelEvaluator subclass in evalModel() which takes a const OutArgs object. This is consistent with the behavior of the rest of a const OutArgs object in that a client is allowed to change the state of objects through a const OutArgs object, they just can't change what objects are pointed to.

Definition at line 2096 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::isFailed ( ) const

Return if the evaluation failed or not.

If the evaluation failed, no assumptions should be made at all about the state of the output objects.

Definition at line 2112 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
bool Thyra::ModelEvaluatorBase::OutArgs< Scalar >::isEmpty ( ) const

Definition at line 2119 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::assertSameSupport ( const OutArgs< Scalar > &  outArgs) const

Assert that two OutArgs objects have the same support.

Definition at line 2152 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
std::string Thyra::ModelEvaluatorBase::OutArgs< Scalar >::modelEvalDescription ( ) const

Definition at line 2226 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
std::string Thyra::ModelEvaluatorBase::OutArgs< Scalar >::description ( ) const
virtual

Reimplemented from Teuchos::Describable.

Definition at line 2233 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::describe ( Teuchos::FancyOStream out,
const Teuchos::EVerbosityLevel  verbLevel 
) const
virtual

Create a more detailed description along about this object and the ModelEvaluator that created it.

Reimplemented from Teuchos::Describable.

Definition at line 2249 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setModelEvalDescription ( const std::string &  modelEvalDescription)
protected

Definition at line 2356 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_set_Np_Ng ( int  Np,
int  Ng 
)
protected

Definition at line 2364 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgsMembers  arg,
bool  supports 
)
protected

Definition at line 2454 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
template<class ObjectType >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( const bool  supports)
protected

Definition at line 1764 of file Thyra_ModelEvaluatorBase_decl.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgsDfDp  arg,
int  l,
const DerivativeSupport supports_in 
)
protected

Definition at line 2467 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgsDgDx_dot  arg,
int  j,
const DerivativeSupport supports_in 
)
protected

Definition at line 2478 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgsDgDx  arg,
int  j,
const DerivativeSupport supports_in 
)
protected

Definition at line 2488 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgsDgDp  arg,
int  j,
int  l,
const DerivativeSupport supports_in 
)
protected

Definition at line 2534 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgs_hess_vec_prod_f_xx  arg,
bool  supports 
)
protected

Definition at line 2498 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgs_hess_vec_prod_f_xp  arg,
int  l,
bool  supports 
)
protected

Definition at line 2506 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgs_hess_vec_prod_f_px  arg,
int  l,
bool  supports 
)
protected

Definition at line 2515 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgs_hess_vec_prod_f_pp  arg,
int  l1,
int  l2,
bool  supports 
)
protected

Definition at line 2524 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgs_hess_vec_prod_g_xx  arg,
int  j,
bool  supports 
)
protected

Definition at line 2545 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgs_hess_vec_prod_g_xp  arg,
int  j,
int  l,
bool  supports 
)
protected

Definition at line 2554 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgs_hess_vec_prod_g_px  arg,
int  j,
int  l,
bool  supports 
)
protected

Definition at line 2564 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgs_hess_vec_prod_g_pp  arg,
int  j,
int  l1,
int  l2,
bool  supports 
)
protected

Definition at line 2574 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgs_hess_f_xx  arg,
bool  supports 
)
protected

Definition at line 2586 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgs_hess_f_xp  arg,
int  l,
bool  supports 
)
protected

Definition at line 2594 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgs_hess_f_pp  arg,
int  l1,
int  l2,
bool  supports 
)
protected

Definition at line 2603 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgs_hess_g_xx  arg,
int  j,
bool  supports 
)
protected

Definition at line 2613 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgs_hess_g_xp  arg,
int  j,
int  l,
bool  supports 
)
protected

Definition at line 2622 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgs_hess_g_pp  arg,
int  j,
int  l1,
int  l2,
bool  supports 
)
protected

Definition at line 2632 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgs_H_xx  arg,
bool  supports 
)
protected

Definition at line 2643 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgs_H_xp  arg,
int  l,
bool  supports 
)
protected

Definition at line 2651 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( EOutArgs_H_pp  arg,
int  l1,
int  l2,
bool  supports 
)
protected

Definition at line 2660 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_set_W_properties ( const DerivativeProperties properties)
protected

Definition at line 2724 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_set_DfDp_properties ( int  l,
const DerivativeProperties properties 
)
protected

Definition at line 2733 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_set_DgDx_dot_properties ( int  j,
const DerivativeProperties properties 
)
protected

Definition at line 2743 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_set_DgDx_properties ( int  j,
const DerivativeProperties properties 
)
protected

Definition at line 2753 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_set_DgDp_properties ( int  j,
int  l,
const DerivativeProperties properties 
)
protected

Definition at line 2763 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar>
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setSupports ( const OutArgs< Scalar > &  inputOutArgs)
protected

Definition at line 2813 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setUnsupportsAndRelated ( EInArgsMembers  arg)
protected

Definition at line 2919 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setUnsupportsAndRelated ( EOutArgsMembers  arg)
protected

Definition at line 2950 of file Thyra_ModelEvaluatorBase_def.hpp.

template<class Scalar >
void Thyra::ModelEvaluatorBase::OutArgs< Scalar >::_setHessianSupports ( const bool  supports)
protected

Definition at line 2983 of file Thyra_ModelEvaluatorBase_def.hpp.


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