Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | List of all members
panzer::ModelEvaluator< Scalar > Class Template Reference

#include <Panzer_ModelEvaluator.hpp>

Inherits StateFuncModelEvaluatorBase< Scalar >.

Classes

struct  ParameterObject
 
struct  ResponseObject
 

Public Member Functions

template<typename EvalT >
void disableEvaluationType ()
 
void buildVolumeFieldManagers (const bool value)
 
void buildBCFieldManagers (const bool value)
 
void setupModel (const Teuchos::RCP< panzer::WorksetContainer > &wc, const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &physicsBlocks, const std::vector< panzer::BC > &bcs, const panzer::EquationSetFactory &eqset_factory, const panzer::BCStrategyFactory &bc_factory, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &volume_cm_factory, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &bc_cm_factory, const Teuchos::ParameterList &closure_models, const Teuchos::ParameterList &user_data, bool writeGraph=false, const std::string &graphPrefix="", const Teuchos::ParameterList &me_params=Teuchos::ParameterList())
 
int addParameter (const std::string &name, const Scalar &initial)
 
int addParameter (const Teuchos::Array< std::string > &names, const Teuchos::Array< Scalar > &initialValues)
 
int addDistributedParameter (const std::string &name, const Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > &vs, const Teuchos::RCP< GlobalEvaluationData > &ged, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &initial, const Teuchos::RCP< const GlobalIndexer > &ugi=Teuchos::null)
 
void addNonParameterGlobalEvaluationData (const std::string &name, const Teuchos::RCP< GlobalEvaluationData > &ged)
 
int addFlexibleResponse (const std::string &responseName, const std::vector< WorksetDescriptor > &wkst_desc, const Teuchos::RCP< ResponseMESupportBuilderBase > &builder)
 
template<typename ResponseEvaluatorFactory_BuilderT >
int addResponse (const std::string &responseName, const std::vector< WorksetDescriptor > &wkst_desc, const ResponseEvaluatorFactory_BuilderT &builder)
 
void buildResponses (const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &physicsBlocks, const panzer::EquationSetFactory &eqset_factory, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &cm_factory, const Teuchos::ParameterList &closure_models, const Teuchos::ParameterList &user_data, const bool write_graphviz_file=false, const std::string &graphviz_file_prefix="")
 
void buildResponses (const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &physicsBlocks, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &cm_factory, const Teuchos::ParameterList &closure_models, const Teuchos::ParameterList &user_data, const bool write_graphviz_file=false, const std::string &graphviz_file_prefix="")
 
void buildDistroParamDfDp_RL (const Teuchos::RCP< panzer::WorksetContainer > &wc, const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &physicsBlocks, const std::vector< panzer::BC > &bcs, const panzer::EquationSetFactory &eqset_factory, const panzer::BCStrategyFactory &bc_factory, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &cm_factory, const Teuchos::ParameterList &closure_models, const Teuchos::ParameterList &user_data, const bool write_graphviz_file=false, const std::string &graphviz_file_prefix="")
 
void buildDistroParamDgDp_RL (const Teuchos::RCP< panzer::WorksetContainer > &wc, const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &physicsBlocks, const std::vector< panzer::BC > &bcs, const panzer::EquationSetFactory &eqset_factory, const panzer::BCStrategyFactory &bc_factory, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &cm_factory, const Teuchos::ParameterList &closure_models, const Teuchos::ParameterList &user_data, const bool write_graphviz_file=false, const std::string &graphviz_file_prefix="")
 
void setOneTimeDirichletBeta (const Scalar &beta) const
 
void applyDirichletBCs (const Teuchos::RCP< Thyra::VectorBase< Scalar > > &x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &f) const
 
void setupAssemblyInArgs (const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, panzer::AssemblyEngineInArgs &ae_inargs) const
 
panzer::AssemblyEngine_TemplateManager
< panzer::Traits
getAssemblyEngineTemplateManager () const
 return a copy of the model evaluators template manager, this is shallow class so pass by value More...
 
Teuchos::RCP
< panzer::ResponseLibrary
< panzer::Traits > > 
getResponseLibrary () const
 
int getXTangentVectorIndex (const int index) const
 
int getXDotTangentVectorIndex (const int index) const
 
Teuchos::RCP< const
Thyra::VectorBase< Scalar > > 
get_parameter_vector (int index) const
 
void evalModel_D2gDx2 (int rIndex, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &delta_x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &D2gDx2) const
 
void evalModel_D2gDp2 (int rIndex, int pIndex, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &delta_x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &D2gDp2) const
 
void evalModel_D2gDpDx (int rIndex, int pIndex, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &delta_x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &D2gDpDx) const
 
void evalModel_D2gDxDp (int rIndex, int pIndex, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &delta_p, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &D2gDxDp) const
 
void evalModel_D2fDx2 (const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &delta_x, const Teuchos::RCP< Thyra::LinearOpBase< Scalar > > &D2fDx2) const
 
void evalModel_D2fDp2 (int pIndex, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &delta_x, const Teuchos::RCP< Thyra::LinearOpBase< Scalar > > &D2fDp2) const
 
void evalModel_D2fDpDx (int pIndex, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &delta_x, const Teuchos::RCP< Thyra::LinearOpBase< Scalar > > &D2fDpDx) const
 
void evalModel_D2fDxDp (int pIndex, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &delta_p, const Teuchos::RCP< Thyra::LinearOpBase< Scalar > > &D2fDxDp) const
 

Protected Member Functions

virtual void evalModelImpl_basic (const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
 Evaluate a simple model, meaning a residual and a jacobian, no fancy stochastic galerkin or multipoint. More...
 
virtual void evalModelImpl_basic_g (const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
 Construct a simple response dicatated by this set of out args. More...
 
virtual void evalModelImpl_basic_dgdx (const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
 
virtual void evalModelImpl_basic_dgdp_scalar (const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
 
virtual void evalModelImpl_basic_dgdp_distro (const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
 
virtual void evalModelImpl_basic_dfdp_scalar (const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
 
virtual void evalModelImpl_basic_dfdp_scalar_fd (const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
 
virtual void evalModelImpl_basic_dfdp_distro (const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
 
bool required_basic_g (const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
 Does this set of out args require a simple response? More...
 
bool required_basic_dgdx (const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
 Are their required responses in the out args? DgDx. More...
 
bool required_basic_dgdp_scalar (const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
 Are their required responses in the out args? DgDp. More...
 
bool required_basic_dgdp_distro (const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
 Are their required responses in the out args? DgDp. More...
 
bool required_basic_dfdp_scalar (const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
 Are derivatives of the residual with respect to the scalar parameters in the out args? DfDp. More...
 
bool required_basic_dfdp_distro (const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
 Are derivatives of the residual with respect to the distributed parameters in the out args? DfDp. More...
 
void initializeNominalValues () const
 Initialize the nominal values with good starting conditions. More...
 
void setParameters (const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs) const
 
void resetParameters () const
 

Private Member Functions

Teuchos::RCP< ParameterObjectcreateScalarParameter (const Teuchos::Array< std::string > &names, const Teuchos::Array< Scalar > &in_values) const
 
Teuchos::RCP< ParameterObjectcreateDistributedParameter (const std::string &key, const Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > &vs, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &initial, const Teuchos::RCP< const GlobalIndexer > &ugi) const
 

Private Attributes

double t_init_
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< Scalar > > 
x_space_
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< Scalar > > 
f_space_
 
Thyra::ModelEvaluatorBase::InArgs
< Scalar > 
prototypeInArgs_
 
Thyra::ModelEvaluatorBase::OutArgs
< Scalar > 
prototypeOutArgs_
 
Thyra::ModelEvaluatorBase::InArgs
< Scalar > 
nominalValues_
 
panzer::AssemblyEngine_TemplateManager
< panzer::Traits
ae_tm_
 
std::vector< Teuchos::RCP
< ParameterObject > > 
parameters_
 
std::vector< Teuchos::RCP
< Thyra::VectorSpaceBase
< double > > > 
tangent_space_
 
int num_me_parameters_
 
bool do_fd_dfdp_
 
double fd_perturb_size_
 
bool require_in_args_refresh_
 
bool require_out_args_refresh_
 
Teuchos::RCP
< panzer::ResponseLibrary
< panzer::Traits > > 
responseLibrary_
 
std::vector< Teuchos::RCP
< ResponseObject > > 
responses_
 
Teuchos::RCP< panzer::GlobalDataglobal_data_
 
bool build_transient_support_
 
Teuchos::RCP< const
panzer::LinearObjFactory
< panzer::Traits > > 
lof_
 
Teuchos::RCP
< panzer::LinearObjContainer
ghostedContainer_
 
Teuchos::RCP
< ReadOnlyVector_GlobalEvaluationData
xContainer_
 
Teuchos::RCP
< ReadOnlyVector_GlobalEvaluationData
xdotContainer_
 
Teuchos::RCP< const
Thyra::LinearOpWithSolveFactoryBase
< Scalar > > 
solverFactory_
 
GlobalEvaluationDataContainer nonParamGlobalEvaluationData_
 
GlobalEvaluationDataContainer distrParamGlobalEvaluationData_
 
bool oneTimeDirichletBeta_on_
 
Scalar oneTimeDirichletBeta_
 
bool build_volume_field_managers_
 
bool build_bc_field_managers_
 
std::vector< bool > active_evaluation_types_
 
unsigned long long write_matrix_count_
 

Constructors/Initializers/Accessors

 ModelEvaluator (const Teuchos::RCP< panzer::FieldManagerBuilder > &fmb, const Teuchos::RCP< panzer::ResponseLibrary< panzer::Traits > > &rLibrary, const Teuchos::RCP< const panzer::LinearObjFactory< panzer::Traits > > &lof, const std::vector< Teuchos::RCP< Teuchos::Array< std::string > > > &p_names, const std::vector< Teuchos::RCP< Teuchos::Array< double > > > &p_values, const Teuchos::RCP< const Thyra::LinearOpWithSolveFactoryBase< Scalar > > &solverFactory, const Teuchos::RCP< panzer::GlobalData > &global_data, bool build_transient_support, double t_init)
 
 ModelEvaluator (const Teuchos::RCP< const panzer::LinearObjFactory< panzer::Traits > > &lof, const Teuchos::RCP< const Thyra::LinearOpWithSolveFactoryBase< Scalar > > &solverFactory, const Teuchos::RCP< panzer::GlobalData > &global_data, bool build_transient_support, double t_init)
 
 ModelEvaluator ()
 

Public functions overridden from ModelEvaulator.

Teuchos::RCP< const
Thyra::VectorSpaceBase< Scalar > > 
get_x_space () const override
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< Scalar > > 
get_f_space () const override
 
Teuchos::RCP< const
Teuchos::Array< std::string > > 
get_p_names (int i) const override
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< Scalar > > 
get_p_space (int i) const override
 
Teuchos::ArrayView< const
std::string > 
get_g_names (int i) const override
 
const std::string & get_g_name (int i) const
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< Scalar > > 
get_g_space (int i) const override
 
Teuchos::RCP
< Thyra::LinearOpBase< Scalar > > 
create_W_op () const override
 
Teuchos::RCP< const
Thyra::LinearOpWithSolveFactoryBase
< Scalar > > 
get_W_factory () const override
 
Teuchos::RCP
< Thyra::LinearOpBase< Scalar > > 
create_DfDp_op (int i) const override
 
Thyra::ModelEvaluatorBase::InArgs
< Scalar > 
createInArgs () const override
 
Thyra::ModelEvaluatorBase::InArgs
< Scalar > 
getNominalValues () const override
 

Private functions overridden from ModelEvaulatorDefaultBase.

Thyra::ModelEvaluatorBase::OutArgs
< Scalar > 
createOutArgsImpl () const override
 
virtual void evalModelImpl (const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const override
 

Detailed Description

template<typename Scalar>
class panzer::ModelEvaluator< Scalar >

Definition at line 41 of file Panzer_ModelEvaluator.hpp.

Constructor & Destructor Documentation

template<typename Scalar >
panzer::ModelEvaluator< Scalar >::ModelEvaluator ( const Teuchos::RCP< panzer::FieldManagerBuilder > &  fmb,
const Teuchos::RCP< panzer::ResponseLibrary< panzer::Traits > > &  rLibrary,
const Teuchos::RCP< const panzer::LinearObjFactory< panzer::Traits > > &  lof,
const std::vector< Teuchos::RCP< Teuchos::Array< std::string > > > &  p_names,
const std::vector< Teuchos::RCP< Teuchos::Array< double > > > &  p_values,
const Teuchos::RCP< const Thyra::LinearOpWithSolveFactoryBase< Scalar > > &  solverFactory,
const Teuchos::RCP< panzer::GlobalData > &  global_data,
bool  build_transient_support,
double  t_init 
)

Definition at line 54 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
panzer::ModelEvaluator< Scalar >::ModelEvaluator ( const Teuchos::RCP< const panzer::LinearObjFactory< panzer::Traits > > &  lof,
const Teuchos::RCP< const Thyra::LinearOpWithSolveFactoryBase< Scalar > > &  solverFactory,
const Teuchos::RCP< panzer::GlobalData > &  global_data,
bool  build_transient_support,
double  t_init 
)

Definition at line 116 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
panzer::ModelEvaluator< Scalar >::ModelEvaluator ( )

Definition at line 162 of file Panzer_ModelEvaluator_impl.hpp.

Member Function Documentation

template<typename Scalar >
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > panzer::ModelEvaluator< Scalar >::get_x_space ( ) const
override

Definition at line 171 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > panzer::ModelEvaluator< Scalar >::get_f_space ( ) const
override

Definition at line 179 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
Teuchos::RCP< const Teuchos::Array< std::string > > panzer::ModelEvaluator< Scalar >::get_p_names ( int  i) const
override

Definition at line 186 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > panzer::ModelEvaluator< Scalar >::get_p_space ( int  i) const
override

Definition at line 215 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
Teuchos::ArrayView< const std::string > panzer::ModelEvaluator< Scalar >::get_g_names ( int  i) const
override

Definition at line 232 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
const std::string & panzer::ModelEvaluator< Scalar >::get_g_name ( int  i) const

Definition at line 242 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > panzer::ModelEvaluator< Scalar >::get_g_space ( int  i) const
override

Definition at line 252 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
Teuchos::RCP< Thyra::LinearOpBase< Scalar > > panzer::ModelEvaluator< Scalar >::create_W_op ( ) const
override

Definition at line 699 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
Teuchos::RCP< const Thyra::LinearOpWithSolveFactoryBase< Scalar > > panzer::ModelEvaluator< Scalar >::get_W_factory ( ) const
override

Definition at line 711 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
Teuchos::RCP< Thyra::LinearOpBase< Scalar > > panzer::ModelEvaluator< Scalar >::create_DfDp_op ( int  i) const
override

Definition at line 719 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
Thyra::ModelEvaluatorBase::InArgs< Scalar > panzer::ModelEvaluator< Scalar >::createInArgs ( ) const
override

Definition at line 262 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
Thyra::ModelEvaluatorBase::InArgs< Scalar > panzer::ModelEvaluator< Scalar >::getNominalValues ( ) const
override

Definition at line 269 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
template<typename EvalT >
void panzer::ModelEvaluator< Scalar >::disableEvaluationType ( )
inline

Disable an evaluation type from AssemblyEngine_TemplateManager, FieldManagerBuilder and PhysicsBlock objects. This will prevent the allocation of unused resources.

Definition at line 115 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::buildVolumeFieldManagers ( const bool  value)

If set to false, disables building volume field managers to save memory if not needed. Must be called BEFORE setupModel() is called. Defaults to true.

Definition at line 363 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::buildBCFieldManagers ( const bool  value)

If set to false, disables building bc field managers to save memory if not needed. Must be called BEFORE setupModel() is called. Defaults to true.

Definition at line 370 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::setupModel ( const Teuchos::RCP< panzer::WorksetContainer > &  wc,
const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &  physicsBlocks,
const std::vector< panzer::BC > &  bcs,
const panzer::EquationSetFactory eqset_factory,
const panzer::BCStrategyFactory bc_factory,
const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &  volume_cm_factory,
const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &  bc_cm_factory,
const Teuchos::ParameterList closure_models,
const Teuchos::ParameterList user_data,
bool  writeGraph = false,
const std::string &  graphPrefix = "",
const Teuchos::ParameterList me_params = Teuchos::ParameterList() 
)

Definition at line 377 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
int panzer::ModelEvaluator< Scalar >::addParameter ( const std::string &  name,
const Scalar &  initial 
)

Add a simple (i.e. nondistributed) parameter to the model evaluator.

Note that these parameters will automatically use the parameter library passed into the model evaluator object through the GlobalData.

Parameters
[in]nameName of the parameter
[in]initialInitial (default) condition for this parameter
Returns
The index associated with this parameter for accessing it through the ModelEvaluator interface.
Note
The implementation for this is a call to the Array version of addParameter

Definition at line 770 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
int panzer::ModelEvaluator< Scalar >::addParameter ( const Teuchos::Array< std::string > &  names,
const Teuchos::Array< Scalar > &  initialValues 
)

Add a simple (i.e. nondistributed) parameter to the model evaluator.

Note that these parameters will automatically use the parameter library passed into the model evaluator object through the GlobalData.

Parameters
[in]namesNames of the parameter
[in]initialValuesInitial values for the parameters
Returns
The index associated with this parameter for accessing it through the ModelEvaluator interface.

Definition at line 783 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
int panzer::ModelEvaluator< Scalar >::addDistributedParameter ( const std::string &  name,
const Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > &  vs,
const Teuchos::RCP< GlobalEvaluationData > &  ged,
const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &  initial,
const Teuchos::RCP< const GlobalIndexer > &  ugi = Teuchos::null 
)

Add a distributed parameter to the model evaluator

Distributed parameters are special in that they most likely will require a global to ghost call before being used in the evaluator. This function registers the parameter and any needed machinery to perform the global to ghost call.

Parameters
[in]nameName of the distributed parameter
[in]vsVector space that this corresponds to
[in]gedGlobal evaluation data object that handles ghosting
[in]initialInitial value to use for this parameter (defaults in the equation set)
[in]ugiUnique global indexer used for this parameter. Useful in constructing derivatives.
Returns
The index associated with this parameter for accessing it through the ModelEvaluator interface.

Definition at line 819 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::addNonParameterGlobalEvaluationData ( const std::string &  name,
const Teuchos::RCP< GlobalEvaluationData > &  ged 
)

Add a global evaluation data object that will be filled as a side effect when evalModel is called. This is useful for building things like auxiliary operators used in block preconditioning. This will not be used as a parameter (or response) to the model evaluator.

Parameters
[in]nameName to associate with global evaluation data object
[in]gedPointer to a global evaluation data object

Definition at line 840 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
int panzer::ModelEvaluator< Scalar >::addFlexibleResponse ( const std::string &  responseName,
const std::vector< WorksetDescriptor > &  wkst_desc,
const Teuchos::RCP< ResponseMESupportBuilderBase > &  builder 
)

Add a response specified by a list of WorksetDescriptor objects. The specifics of the response are specified by the response factory builder. This version supports computing derivatives with respect to both the state ('x') and control ('p') variables and is thus ``flexible''.

NOTE: Response factories must use a response of type ResponseMESupportBase. This is how the model evaluator parses and puts responses in the right location. If this condition is violated the evalModel call will fail. Furthermore, this method cannot be called after buildRespones has been called.

Parameters
[in]responseNameName of the response to be added.
[in]wkst_descA vector of descriptors describing the types of elements that make up the response.
[in]builderBuilder that builds the correct response object.
Returns
The index associated with this response for accessing it through the ModelEvaluator interface.

Definition at line 848 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
template<typename ResponseEvaluatorFactory_BuilderT >
int panzer::ModelEvaluator< Scalar >::addResponse ( const std::string &  responseName,
const std::vector< WorksetDescriptor > &  wkst_desc,
const ResponseEvaluatorFactory_BuilderT &  builder 
)

Add a response specified by a list of WorksetDescriptor objects. The specifics of the response are specified by the response factory builder.

NOTE: Response factories must use a response of type ResponseMESupportBase. This is how the model evaluator parses and puts responses in the right location. If this condition is violated the evalModel call will fail. Furthermore, this method cannot be called after buildRespones has been called.

Parameters
[in]responseNameName of the response to be added.
[in]wkst_descA vector of descriptors describing the types of elements that make up the response.
[in]builderBuilder that builds the correct response object.
Returns
The index associated with this response for accessing it through the ModelEvaluator interface.

Definition at line 688 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::buildResponses ( const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &  physicsBlocks,
const panzer::EquationSetFactory eqset_factory,
const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &  cm_factory,
const Teuchos::ParameterList closure_models,
const Teuchos::ParameterList user_data,
const bool  write_graphviz_file = false,
const std::string &  graphviz_file_prefix = "" 
)
inline

Build all the responses set on the model evaluator. Once this method is called no other responses can be added. An exception is thrown if they are.

Definition at line 248 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::buildResponses ( const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &  physicsBlocks,
const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &  cm_factory,
const Teuchos::ParameterList closure_models,
const Teuchos::ParameterList user_data,
const bool  write_graphviz_file = false,
const std::string &  graphviz_file_prefix = "" 
)
inline

Build all the responses set on the model evaluator. Once this method is called no other responses can be added. An exception is thrown if they are.

Definition at line 272 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::buildDistroParamDfDp_RL ( const Teuchos::RCP< panzer::WorksetContainer > &  wc,
const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &  physicsBlocks,
const std::vector< panzer::BC > &  bcs,
const panzer::EquationSetFactory eqset_factory,
const panzer::BCStrategyFactory bc_factory,
const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &  cm_factory,
const Teuchos::ParameterList closure_models,
const Teuchos::ParameterList user_data,
const bool  write_graphviz_file = false,
const std::string &  graphviz_file_prefix = "" 
)

This method builds the response libraries that build the dfdp sensitivities for the distributed parameters if requested. Note that in general the user is expected to call this through setupModel and not call it directly.

Definition at line 2283 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::buildDistroParamDgDp_RL ( const Teuchos::RCP< panzer::WorksetContainer > &  wc,
const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &  physicsBlocks,
const std::vector< panzer::BC > &  bcs,
const panzer::EquationSetFactory eqset_factory,
const panzer::BCStrategyFactory bc_factory,
const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &  cm_factory,
const Teuchos::ParameterList closure_models,
const Teuchos::ParameterList user_data,
const bool  write_graphviz_file = false,
const std::string &  graphviz_file_prefix = "" 
)

This method builds the response libraries that build the dgdp sensitivities for the distributed parameters if requested. This only applies to "flexible" responses. Note that in general the user is expected to call this through setupModel and not call it directly.

Definition at line 2334 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::setOneTimeDirichletBeta ( const Scalar &  beta) const

This function is intended for experts only, it allows for a beta to be set for the dirichlet conditions only. This allows the dirichlet condition to be propagated to the mass matrix. The reason it is one time only is that it breaks encapsulation, and should be only used if absolutely neccessary.

Parameters
[in]betaValue of beta to use.

Definition at line 2402 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::applyDirichletBCs ( const Teuchos::RCP< Thyra::VectorBase< Scalar > > &  x,
const Teuchos::RCP< Thyra::VectorBase< Scalar > > &  f 
) const

Apply the dirichlet boundary conditions to the vector "f" using the "x" values as the current solution.

Definition at line 866 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::setupAssemblyInArgs ( const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
panzer::AssemblyEngineInArgs ae_inargs 
) const

Setup all the assembly input arguments required by "inArgs".

Parameters
[in]inArgsModel evalutor input arguments
in/out]ae_inArgs Assembly engine input arguments.

Definition at line 454 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
panzer::AssemblyEngine_TemplateManager<panzer::Traits> panzer::ModelEvaluator< Scalar >::getAssemblyEngineTemplateManager ( ) const
inline

return a copy of the model evaluators template manager, this is shallow class so pass by value

Returns
The AssemblyEngine template manager

Definition at line 355 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
Teuchos::RCP<panzer::ResponseLibrary<panzer::Traits> > panzer::ModelEvaluator< Scalar >::getResponseLibrary ( ) const
inline

Definition at line 359 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
int panzer::ModelEvaluator< Scalar >::getXTangentVectorIndex ( const int  index) const
inline

Returns the x tangent vector index for a given parameter index

Definition at line 364 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
int panzer::ModelEvaluator< Scalar >::getXDotTangentVectorIndex ( const int  index) const
inline

Returns the xdot tangent vector index for a given parameter index

Definition at line 375 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
Teuchos::RCP<const Thyra::VectorBase<Scalar> > panzer::ModelEvaluator< Scalar >::get_parameter_vector ( int  index) const
inline

Initializes the given vector with current values of the parameters

Definition at line 386 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::evalModel_D2gDx2 ( int  rIndex,
const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &  delta_x,
const Teuchos::RCP< Thyra::VectorBase< Scalar > > &  D2gDx2 
) const

Compute second (x) derivative of the response in the direction delta_x.

Parameters
[in]rIndexResponse to differentiate
[in]inArgsInput arguments that sets the state
[in]delta_xDirection to take the derivative with respect to.
[out]D2gDx2Result vector allocated by get_x_space().

Definition at line 939 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::evalModel_D2gDp2 ( int  rIndex,
int  pIndex,
const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &  delta_x,
const Teuchos::RCP< Thyra::VectorBase< Scalar > > &  D2gDp2 
) const

Compute second (p) derivative of the response in the direction delta_p.

Parameters
[in]rIndexResponse to differentiate
[in]pIndexParameter to differentiate with respect to
[in]inArgsInput arguments that sets the state
[in]delta_pDirection to take the derivative with respect to.
[out]D2gDp2Result vector allocated by get_p_space(pIndex).

Definition at line 1036 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::evalModel_D2gDpDx ( int  rIndex,
int  pIndex,
const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &  delta_x,
const Teuchos::RCP< Thyra::VectorBase< Scalar > > &  D2gDpDx 
) const

Compute second (p) derivative of the response in the direction delta_x.

Parameters
[in]rIndexResponse to differentiate
[in]pIndexParameter to differentiate with respect to
[in]inArgsInput arguments that sets the state
[in]delta_xDirection to take the derivative with respect to.
[out]D2gDpDxResult vector allocated by get_x_space().

Definition at line 1090 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::evalModel_D2gDxDp ( int  rIndex,
int  pIndex,
const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &  delta_p,
const Teuchos::RCP< Thyra::VectorBase< Scalar > > &  D2gDxDp 
) const

Compute second (p) derivative of the response in the direction delta_x.

Parameters
[in]rIndexResponse to differentiate
[in]pIndexParameter to differentiate with respect to
[in]inArgsInput arguments that sets the state
[in]delta_pDirection to take the derivative with respect to.
[out]D2gDxDpResult vector allocated by get_x_space().

Definition at line 986 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::evalModel_D2fDx2 ( const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &  delta_x,
const Teuchos::RCP< Thyra::LinearOpBase< Scalar > > &  D2fDx2 
) const

Compute second (x) derivative of the residual in the direction delta_x.

Parameters
[in]rIndexResponse to differentiate
[in]inArgsInput arguments that sets the state
[in]delta_xDirection to take the derivative with respect to.
[out]D2fDx2Result vector allocated by get_x_space().

Definition at line 1144 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::evalModel_D2fDp2 ( int  pIndex,
const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &  delta_x,
const Teuchos::RCP< Thyra::LinearOpBase< Scalar > > &  D2fDp2 
) const

Compute second (p) derivative of the residual in the direction delta_p.

Parameters
[in]rIndexResponse to differentiate
[in]pIndexParameter to differentiate with respect to
[in]inArgsInput arguments that sets the state
[in]delta_pDirection to take the derivative with respect to.
[out]D2fDp2Result vector allocated by get_p_space(pIndex).

Definition at line 1395 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::evalModel_D2fDpDx ( int  pIndex,
const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &  delta_x,
const Teuchos::RCP< Thyra::LinearOpBase< Scalar > > &  D2fDpDx 
) const

Compute second (p) derivative of the residual in the direction delta_x.

Parameters
[in]rIndexResponse to differentiate
[in]pIndexParameter to differentiate with respect to
[in]inArgsInput arguments that sets the state
[in]delta_xDirection to take the derivative with respect to.
[out]D2fDpDxResult vector allocated by get_x_space().

Definition at line 1344 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::evalModel_D2fDxDp ( int  pIndex,
const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &  delta_p,
const Teuchos::RCP< Thyra::LinearOpBase< Scalar > > &  D2fDxDp 
) const

Compute second (p) derivative of the residual in the direction delta_p.

Parameters
[in]rIndexResponse to differentiate
[in]pIndexParameter to differentiate with respect to
[in]inArgsInput arguments that sets the state
[in]delta_pDirection to take the derivative with respect to.
[out]D2fDxDpResult vector allocated by get_x_space().

Definition at line 1242 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
Thyra::ModelEvaluatorBase::OutArgs< Scalar > panzer::ModelEvaluator< Scalar >::createOutArgsImpl ( ) const
overrideprotected

Definition at line 640 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::evalModelImpl ( const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &  outArgs 
) const
overrideprotectedvirtual

Definition at line 1446 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::evalModelImpl_basic ( const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &  outArgs 
) const
protectedvirtual

Evaluate a simple model, meaning a residual and a jacobian, no fancy stochastic galerkin or multipoint.

Definition at line 1480 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::evalModelImpl_basic_g ( const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &  outArgs 
) const
protectedvirtual

Construct a simple response dicatated by this set of out args.

Definition at line 1647 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::evalModelImpl_basic_dgdx ( const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &  outArgs 
) const
protectedvirtual

handles evaluation of responses dgdx

Note
This method should (basically) be a no-op if required_basic_dgdx(outArgs)==false. However, for efficiency this is not checked.

Definition at line 1684 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::evalModelImpl_basic_dgdp_scalar ( const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &  outArgs 
) const
protectedvirtual

handles evaluation of responses dgdp (scalar) defined as dg/dx * dx/dp + dg/dp

Note
This method should (basically) be a no-op if required_basic_dgdp_scalar(outArgs)==false. However, for efficiency this is not checked.

Definition at line 1730 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::evalModelImpl_basic_dgdp_distro ( const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &  outArgs 
) const
protectedvirtual

handles evaluation of responses dgdp (distributed)

Note
This method should (basically) be a no-op if required_basic_dgdx_distro(outArgs)==false. However, for efficiency this is not checked.

Definition at line 1815 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::evalModelImpl_basic_dfdp_scalar ( const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &  outArgs 
) const
protectedvirtual

handles evaluation of dfdp (tangent) defined as df/dx * dx/dp + df/dp

Note
This method should (basically) be a no-op if required_basic_dfdp_scalar(outArgs)==false. However, for efficiency this is not checked.

Definition at line 1874 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::evalModelImpl_basic_dfdp_scalar_fd ( const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &  outArgs 
) const
protectedvirtual

handles evaluation of dfdp (tangent) defined as df/dx * dx/dp + df/dp using finite-differences

Note
This method should (basically) be a no-op if required_basic_dfdp_scalar(outArgs)==false. However, for efficiency this is not checked.

Definition at line 1990 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::evalModelImpl_basic_dfdp_distro ( const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs,
const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &  outArgs 
) const
protectedvirtual

handles evaluation of dfdp

Note
This method should (basically) be a no-op if required_basic_dfdp_distro(outArgs)==false. However, for efficiency this is not checked.

Definition at line 2094 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
bool panzer::ModelEvaluator< Scalar >::required_basic_g ( const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &  outArgs) const
protected

Does this set of out args require a simple response?

Definition at line 2149 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
bool panzer::ModelEvaluator< Scalar >::required_basic_dgdx ( const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &  outArgs) const
protected

Are their required responses in the out args? DgDx.

Definition at line 2161 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
bool panzer::ModelEvaluator< Scalar >::required_basic_dgdp_scalar ( const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &  outArgs) const
protected

Are their required responses in the out args? DgDp.

Definition at line 2181 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
bool panzer::ModelEvaluator< Scalar >::required_basic_dgdp_distro ( const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &  outArgs) const
protected

Are their required responses in the out args? DgDp.

Definition at line 2207 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
bool panzer::ModelEvaluator< Scalar >::required_basic_dfdp_scalar ( const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &  outArgs) const
protected

Are derivatives of the residual with respect to the scalar parameters in the out args? DfDp.

Definition at line 2233 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
bool panzer::ModelEvaluator< Scalar >::required_basic_dfdp_distro ( const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &  outArgs) const
protected

Are derivatives of the residual with respect to the distributed parameters in the out args? DfDp.

Definition at line 2258 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::initializeNominalValues ( ) const
protected

Initialize the nominal values with good starting conditions.

Definition at line 307 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::setParameters ( const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs) const
protected

Definition at line 2485 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
void panzer::ModelEvaluator< Scalar >::resetParameters ( ) const
protected

Definition at line 2507 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
Teuchos::RCP< typename panzer::ModelEvaluator< Scalar >::ParameterObject > panzer::ModelEvaluator< Scalar >::createScalarParameter ( const Teuchos::Array< std::string > &  names,
const Teuchos::Array< Scalar > &  in_values 
) const
private

Definition at line 2411 of file Panzer_ModelEvaluator_impl.hpp.

template<typename Scalar >
Teuchos::RCP< typename panzer::ModelEvaluator< Scalar >::ParameterObject > panzer::ModelEvaluator< Scalar >::createDistributedParameter ( const std::string &  key,
const Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > &  vs,
const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &  initial,
const Teuchos::RCP< const GlobalIndexer > &  ugi 
) const
private

Definition at line 2461 of file Panzer_ModelEvaluator_impl.hpp.

Member Data Documentation

template<typename Scalar >
double panzer::ModelEvaluator< Scalar >::t_init_
private

Definition at line 634 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> > panzer::ModelEvaluator< Scalar >::x_space_
private

Definition at line 636 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> > panzer::ModelEvaluator< Scalar >::f_space_
private

Definition at line 637 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
Thyra::ModelEvaluatorBase::InArgs<Scalar> panzer::ModelEvaluator< Scalar >::prototypeInArgs_
mutableprivate

Definition at line 639 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
Thyra::ModelEvaluatorBase::OutArgs<Scalar> panzer::ModelEvaluator< Scalar >::prototypeOutArgs_
mutableprivate

Definition at line 640 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
Thyra::ModelEvaluatorBase::InArgs<Scalar> panzer::ModelEvaluator< Scalar >::nominalValues_
mutableprivate

Definition at line 642 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
panzer::AssemblyEngine_TemplateManager<panzer::Traits> panzer::ModelEvaluator< Scalar >::ae_tm_
mutableprivate

Definition at line 644 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
std::vector<Teuchos::RCP<ParameterObject> > panzer::ModelEvaluator< Scalar >::parameters_
private

Definition at line 646 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
std::vector<Teuchos::RCP<Thyra::VectorSpaceBase<double> > > panzer::ModelEvaluator< Scalar >::tangent_space_
private

Definition at line 647 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
int panzer::ModelEvaluator< Scalar >::num_me_parameters_
private

Definition at line 648 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
bool panzer::ModelEvaluator< Scalar >::do_fd_dfdp_
private

Definition at line 649 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
double panzer::ModelEvaluator< Scalar >::fd_perturb_size_
private

Definition at line 650 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
bool panzer::ModelEvaluator< Scalar >::require_in_args_refresh_
mutableprivate

Definition at line 652 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
bool panzer::ModelEvaluator< Scalar >::require_out_args_refresh_
mutableprivate

Definition at line 653 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
Teuchos::RCP<panzer::ResponseLibrary<panzer::Traits> > panzer::ModelEvaluator< Scalar >::responseLibrary_
mutableprivate

Definition at line 656 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
std::vector<Teuchos::RCP<ResponseObject> > panzer::ModelEvaluator< Scalar >::responses_
private

Definition at line 657 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
Teuchos::RCP<panzer::GlobalData> panzer::ModelEvaluator< Scalar >::global_data_
private

Definition at line 659 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
bool panzer::ModelEvaluator< Scalar >::build_transient_support_
private

Definition at line 660 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
Teuchos::RCP<const panzer::LinearObjFactory<panzer::Traits> > panzer::ModelEvaluator< Scalar >::lof_
private

Definition at line 663 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
Teuchos::RCP<panzer::LinearObjContainer> panzer::ModelEvaluator< Scalar >::ghostedContainer_
mutableprivate

Definition at line 664 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
Teuchos::RCP<ReadOnlyVector_GlobalEvaluationData> panzer::ModelEvaluator< Scalar >::xContainer_
mutableprivate

Definition at line 665 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
Teuchos::RCP<ReadOnlyVector_GlobalEvaluationData> panzer::ModelEvaluator< Scalar >::xdotContainer_
mutableprivate

Definition at line 666 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
Teuchos::RCP<const Thyra::LinearOpWithSolveFactoryBase<Scalar> > panzer::ModelEvaluator< Scalar >::solverFactory_
private

Definition at line 669 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
GlobalEvaluationDataContainer panzer::ModelEvaluator< Scalar >::nonParamGlobalEvaluationData_
private

Definition at line 671 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
GlobalEvaluationDataContainer panzer::ModelEvaluator< Scalar >::distrParamGlobalEvaluationData_
private

Definition at line 672 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
bool panzer::ModelEvaluator< Scalar >::oneTimeDirichletBeta_on_
mutableprivate

Definition at line 674 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
Scalar panzer::ModelEvaluator< Scalar >::oneTimeDirichletBeta_
mutableprivate

Definition at line 675 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
bool panzer::ModelEvaluator< Scalar >::build_volume_field_managers_
private

Definition at line 677 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
bool panzer::ModelEvaluator< Scalar >::build_bc_field_managers_
private

Definition at line 678 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
std::vector<bool> panzer::ModelEvaluator< Scalar >::active_evaluation_types_
private

Definition at line 679 of file Panzer_ModelEvaluator.hpp.

template<typename Scalar >
unsigned long long panzer::ModelEvaluator< Scalar >::write_matrix_count_
mutableprivate

Definition at line 681 of file Panzer_ModelEvaluator.hpp.


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