44#include "Teuchos_Assert.hpp"
51 const Teuchos::RCP<EpetraExt::ModelEvaluator>& me_,
52 const Teuchos::Array<int>& mp_p_index_map_,
53 const Teuchos::Array<int>& mp_g_index_map_,
54 const Teuchos::Array< Teuchos::RCP<const Epetra_Map> >& base_g_maps_) :
64 InArgs me_inargs =
me->createInArgs();
65 OutArgs me_outargs =
me->createOutArgs();
67 num_g = me_outargs.Ng();
69 TEUCHOS_TEST_FOR_EXCEPTION(
72 <<
"Error! Stokhos::MPInverseModelEvaluator::MPInverseModelEvaluator():"
73 <<
" Base response map array size does not match size of index array!");
78Teuchos::RCP<const Epetra_Map>
85Teuchos::RCP<const Epetra_Map>
92Teuchos::RCP<const Epetra_Map>
96 TEUCHOS_TEST_FOR_EXCEPTION(
97 l >=
num_p || l < 0, std::logic_error,
98 std::endl <<
"Error! Stokhos::MPInverseModelEvaluator::get_p_map():"
99 <<
" Invalid parameter index l = " << l << std::endl);
100 return me->get_p_map(l);
103Teuchos::RCP<const Epetra_Map>
107 TEUCHOS_TEST_FOR_EXCEPTION(
108 l >=
num_g || l < 0, std::logic_error,
109 std::endl <<
"Error! Stokhos::MPInverseModelEvaluator::get_g_map():"
110 <<
" Invalid response index l = " << l << std::endl);
114Teuchos::RCP<const Teuchos::Array<std::string> >
118 TEUCHOS_TEST_FOR_EXCEPTION(
119 l >=
num_p || l < 0, std::logic_error,
120 std::endl <<
"Error! Stokhos::MPInverseModelEvaluator::get_p_names():"
121 <<
" Invalid parameter index l = " << l << std::endl);
122 return me->get_p_names(l);
125Teuchos::RCP<const Epetra_Vector>
129 TEUCHOS_TEST_FOR_EXCEPTION(
130 l >=
num_p || l < 0, std::logic_error,
131 std::endl <<
"Error! Stokhos::MPInverseModelEvaluator::get_p_init():"
132 <<
" Invalid parameter index l = " << l << std::endl);
133 return me->get_p_init(l);
136EpetraExt::ModelEvaluator::InArgs
140 InArgs me_inargs =
me->createInArgs();
142 inArgs.setModelEvalDescription(this->description());
143 inArgs.set_Np(
num_p);
150EpetraExt::ModelEvaluator::OutArgs
153 OutArgsSetup outArgs;
154 OutArgs me_outargs =
me->createOutArgs();
156 outArgs.setModelEvalDescription(this->description());
163 me_outargs.supports(OUT_ARG_DgDp,i,
j));
171 const OutArgs& outArgs)
const
174 InArgs me_inargs =
me->createInArgs();
177 for (
int i=0; i<
num_p; i++)
178 me_inargs.set_p(i, inArgs.get_p(i));
183 if (p_mp != Teuchos::null) {
184 me_inargs.set_p(i+
num_p, p_mp->getBlockVector());
189 OutArgs me_outargs =
me->createOutArgs();
197 mp_vector_t g_mp = outArgs.get_g_mp(ii);
198 if (g_mp != Teuchos::null) {
199 me_outargs.set_g(i, Teuchos::rcp_dynamic_cast<Epetra_Vector>(g_mp->getBlockVector()));
204 if (!outArgs.supports(OUT_ARG_DgDp_mp, ii,
j).none()) {
205 MPDerivative dgdp_mp = outArgs.get_DgDp_mp(ii,
j);
206 Teuchos::RCP<Stokhos::ProductEpetraMultiVector> dgdp_mp_mv =
207 dgdp_mp.getMultiVector();
208 Teuchos::RCP<Epetra_Operator> dgdp_mp_op =
209 dgdp_mp.getLinearOp();
210 if (dgdp_mp_mv != Teuchos::null) {
212 i,
j, Derivative(dgdp_mp_mv->getBlockMultiVector(),
213 dgdp_mp.getMultiVectorOrientation()));
215 else if (dgdp_mp_op != Teuchos::null) {
216 me_outargs.set_DgDp(i,
j, Derivative(dgdp_mp_op));
224 me->evalModel(me_inargs, me_outargs);
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Evaluate model on InArgs.
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > base_g_maps
Base maps of block g vectors.
OutArgs createOutArgs() const
Create OutArgs.
Teuchos::RCP< const Epetra_Map > get_f_map() const
Return residual vector map.
int num_p
Number of parameters.
Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
Return initial parameters.
Teuchos::Array< int > mp_p_index_map
Mapping between multipoint block parameters and mp parameters.
int num_p_mp
Number of multi-point parameter vectors.
Teuchos::Array< int > mp_g_index_map
Mapping between stochastic block responses and sg responses.
MPInverseModelEvaluator(const Teuchos::RCP< EpetraExt::ModelEvaluator > &me, const Teuchos::Array< int > &mp_p_index_map, const Teuchos::Array< int > &mp_g_index_map, const Teuchos::Array< Teuchos::RCP< const Epetra_Map > > &base_g_maps)
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
Return parameter vector map.
int num_g
Number of responses.
Teuchos::RCP< EpetraExt::ModelEvaluator > me
Underlying model evaluator.
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
Return array of parameter names.
InArgs createInArgs() const
Create InArgs.
Teuchos::RCP< const Epetra_Map > get_g_map(int l) const
Return response map.
Teuchos::RCP< const Epetra_Map > get_x_map() const
Return solution vector map.
int num_g_mp
Number of multi-point response vectors.