42 #ifndef STOKHOS_SGQUADMODELEVALUATOR_HPP 43 #define STOKHOS_SGQUADMODELEVALUATOR_HPP 45 #include "EpetraExt_ModelEvaluator.h" 47 #include "Teuchos_RCP.hpp" 48 #include "Teuchos_Array.hpp" 68 const Teuchos::RCP<EpetraExt::ModelEvaluator>&
me);
74 Teuchos::RCP<const Epetra_Map>
get_x_map()
const;
77 Teuchos::RCP<const Epetra_Map>
get_f_map()
const;
80 Teuchos::RCP<const Epetra_Map>
get_p_map(
int l)
const;
83 Teuchos::RCP<const Epetra_Map>
get_g_map(
int l)
const;
86 Teuchos::RCP<const Teuchos::Array<std::string> >
90 Teuchos::RCP<const Epetra_Vector>
get_x_init()
const;
93 Teuchos::RCP<const Epetra_Vector>
get_p_init(
int l)
const;
96 Teuchos::RCP<Epetra_Operator>
create_W()
const;
105 void evalModel(
const InArgs& inArgs,
const OutArgs& outArgs)
const;
112 Teuchos::RCP<EpetraExt::ModelEvaluator>
me;
124 Teuchos::RCP<Epetra_Vector>
x_qp;
127 Teuchos::Array< Teuchos::RCP<Epetra_Vector> >
p_qp;
130 Teuchos::RCP<Epetra_Vector>
f_qp;
133 Teuchos::RCP<Epetra_Operator>
W_qp;
136 Teuchos::Array<EpetraExt::ModelEvaluator::Derivative>
dfdp_qp;
139 Teuchos::Array< Teuchos::RCP<Epetra_Vector> >
g_qp;
142 Teuchos::Array<EpetraExt::ModelEvaluator::Derivative>
dgdx_qp;
145 Teuchos::Array<EpetraExt::ModelEvaluator::Derivative>
dgdx_dot_qp;
148 Teuchos::Array< Teuchos::Array<EpetraExt::ModelEvaluator::Derivative> >
dgdp_qp;
154 #endif // STOKHOS_SGQUADMODELEVALUATOR_HPP int num_p
Number of parameter vectors.
int num_g
Number of response vectors.
Teuchos::RCP< const Epetra_Vector > get_x_init() const
Return initial solution.
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
Return parameter vector map.
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Evaluate model on InArgs.
ModelEvaluator adaptor that implements the stochastic Galerkin residual and Jacobian computations usi...
Teuchos::Array< EpetraExt::ModelEvaluator::Derivative > dgdx_dot_qp
Response derivative.
InArgs createInArgs() const
Create InArgs.
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
Return array of parameter names.
OutArgs createOutArgs() const
Create OutArgs.
Teuchos::Array< EpetraExt::ModelEvaluator::Derivative > dgdx_qp
Response derivative.
Teuchos::RCP< const Epetra_Map > get_x_map() const
Return solution vector map.
Teuchos::RCP< const Epetra_Map > get_f_map() const
Return residual vector map.
Teuchos::Array< Teuchos::RCP< Epetra_Vector > > g_qp
Response vectors.
Teuchos::RCP< Epetra_Operator > create_W() const
Create W = alpha*M + beta*J matrix.
Teuchos::Array< Teuchos::RCP< Epetra_Vector > > p_qp
Parameter vectors.
Top-level namespace for Stokhos classes and functions.
Teuchos::Array< EpetraExt::ModelEvaluator::Derivative > dfdp_qp
Residual derivatives.
Teuchos::RCP< Epetra_Vector > x_dot_qp
Time derivative vector.
Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
Return initial parameters.
Teuchos::RCP< Epetra_Vector > f_qp
Residual vector.
Teuchos::RCP< Epetra_Vector > x_qp
Solution vector.
Teuchos::RCP< Epetra_Operator > W_qp
W operator.
Teuchos::Array< Teuchos::Array< EpetraExt::ModelEvaluator::Derivative > > dgdp_qp
Response sensitivities.
Teuchos::RCP< EpetraExt::ModelEvaluator > me
Underlying model evaluator.
SGQuadModelEvaluator(const Teuchos::RCP< EpetraExt::ModelEvaluator > &me)
Teuchos::RCP< const Epetra_Map > get_g_map(int l) const
Return observation vector map.