11 #include "Thyra_LinearOpWithSolveFactoryHelpers.hpp"
14 #include "Thyra_DetachedMultiVectorView.hpp"
15 #include "Thyra_ModelEvaluatorDelegatorBase.hpp"
16 #include "EpetraExt_ModelEvaluatorScalingTools.h"
17 #include "Epetra_RowMatrix.h"
28 const std::string StateFunctionScaling_name =
"State Function Scaling";
34 stateFunctionScalingValidator;
35 const std::string StateFunctionScaling_default =
"None";
41 const EpetraExt::ModelEvaluator::OutArgs &epetraOutArgs
45 const RCP<Epetra_Operator>
46 eW = epetraOutArgs.get_W();
47 const RCP<Epetra_RowMatrix>
48 ermW = Teuchos::rcp_dynamic_cast<Epetra_RowMatrix>(eW,
false);
50 is_null(ermW), std::logic_error,
51 "Thyra::EpetraModelEvaluator::evalModel(...): Error, if\n"
52 "scaling is turned on, the underlying Epetra_Operator created\n"
53 "an initialized by the underlying epetra model evaluator\n"
54 "\"" << epetraOutArgs.modelEvalDescription() <<
"\"\n"
55 "must support the Epetra_RowMatrix interface through a dynamic cast.\n"
65 const EpetraExt::ModelEvaluator &epetraModel
70 eW = epetraModel.create_W();
73 "Error, the call to create_W() returned null on the "
74 "EpetraExt::ModelEvaluator object "
75 "\"" << epetraModel.description() <<
"\". This may mean that "
76 "the underlying model does not support more than one copy of "
92 :nominalValuesAndBoundsAreUpdated_(false), stateFunctionScaling_(STATE_FUNC_SCALING_NONE),
93 currentInArgsOutArgs_(false), finalPointWasSolved_(false)
99 const RCP<LinearOpWithSolveFactoryBase<double> > &W_factory
101 :nominalValuesAndBoundsAreUpdated_(false), stateFunctionScaling_(STATE_FUNC_SCALING_NONE),
102 currentInArgsOutArgs_(false), finalPointWasSolved_(false)
110 const RCP<LinearOpWithSolveFactoryBase<double> > &W_factory
127 EpetraExt::ModelEvaluator::InArgs inArgs =
epetraModel_->createInArgs();
130 for(
int l = 0; l < implicit_cast<int>(
p_space_.
size()); ++l ) {
135 is_null(p_map_l), std::logic_error,
136 "Error, the the map p["<<l<<
"] for the model \""
137 <<epetraModel->description()<<
"\" can not be null!");
144 EpetraExt::ModelEvaluator::OutArgs outArgs =
epetraModel_->createOutArgs();
147 for(
int j = 0; j < implicit_cast<int>(
g_space_.
size()); ++j ) {
172 const ModelEvaluatorBase::InArgs<double>& nominalValues
185 typedef ModelEvaluatorBase MEB;
226 RCP<LinearOpWithSolveFactoryBase<double> > *W_factory
240 const ModelEvaluatorBase::InArgs<double>&
258 std::ostringstream oss;
259 oss <<
"Thyra::EpetraModelEvaluator{";
260 oss <<
"epetraModel=";
265 oss <<
",W_factory=";
287 *
paramList_, StateFunctionScaling_name, StateFunctionScaling_default
291 Teuchos::readVerboseObjectSublist(&*
paramList_,
this);
294 #endif // TEUCHOS_DEBUG
326 using Teuchos::tuple;
327 using Teuchos::rcp_implicit_cast;
333 stateFunctionScalingValidator =
rcp(
334 new StringToIntegralParameterEntryValidator<EStateFunctionScaling>(
340 "Do not scale the state function f(...) in this class.",
342 "Scale the state function f(...) and all its derivatives\n"
343 "using the row sum scaling from the initial Jacobian\n"
344 "W=d(f)/d(x). Note, this only works with Epetra_CrsMatrix\n"
347 tuple<EStateFunctionScaling>(
351 StateFunctionScaling_name
354 pl->
set(StateFunctionScaling_name,StateFunctionScaling_default,
355 "Determines if and how the state function f(...) and all of its\n"
356 "derivatives are scaled. The scaling is done explicitly so there should\n"
357 "be no impact on the meaning of inner products or tolerances for\n"
359 rcp_implicit_cast<const PEV>(stateFunctionScalingValidator)
361 Teuchos::setupVerboseObjectSublist(&*pl);
435 ModelEvaluatorBase::InArgs<double>
443 ModelEvaluatorBase::InArgs<double>
451 ModelEvaluatorBase::InArgs<double>
489 const ModelEvaluatorBase::InArgs<double> &finalPoint,
534 ModelEvaluatorBase::OutArgs<double>
545 const ModelEvaluatorBase::InArgs<double>& inArgs_in,
546 const ModelEvaluatorBase::OutArgs<double>& outArgs
551 using Teuchos::rcp_const_cast;
552 using Teuchos::rcp_dynamic_cast;
555 typedef EpetraExt::ModelEvaluator EME;
573 inArgs.setArgs(inArgs_in);
579 if (inArgs.supports(Thyra::ModelEvaluatorBase::IN_ARG_x_dot)) {
580 if (
is_null(inArgs_in.get_x_dot()))
585 typedef double Scalar;
586 THYRA_MODEL_EVALUATOR_DECORATOR_EVAL_MODEL_GEN_BEGIN(
591 const bool firstTimeStateFuncScaling
598 VOTSLOWSF W_factory_outputTempState(
W_factory_,out,verbLevel);
612 *out <<
"\nSetting-up/creating input arguments ...\n";
616 EME::InArgs epetraScaledInArgs =
epetraModel_->createInArgs();
621 EME::InArgs epetraInArgs =
epetraModel_->createInArgs();
622 EpetraExt::unscaleModelVars(
636 *out <<
"\nSetting-up/creating output arguments ...\n";
641 EME::OutArgs epetraUnscaledOutArgs =
epetraModel_->createOutArgs();
653 &epetraUnscaledOutArgs,
661 if (firstTimeStateFuncScaling) {
668 <<
"\nTime to setup OutArgs = "
676 *out <<
"\nEvaluating the Epetra output functions ...\n";
679 epetraModel_->evalModel(epetraInArgs,epetraUnscaledOutArgs);
684 <<
"\nTime to evaluate Epetra output functions = "
696 *out <<
"\nCompute scale factors if needed ...\n";
699 if (firstTimeStateFuncScaling) {
706 <<
"\nTime to compute scale factors = "
714 *out <<
"\nScale the output objects ...\n";
717 EME::OutArgs epetraOutArgs =
epetraModel_->createOutArgs();
718 bool allFuncsWhereScaled =
false;
719 EpetraExt::scaleModelFuncs(
721 &epetraOutArgs, &allFuncsWhereScaled,
727 !allFuncsWhereScaled, std::logic_error,
728 "Error, we can not currently handle epetra output objects that could not be"
729 " scaled. Special code will have to be added to handle this (i.e. using"
730 " implicit diagonal and multiplied linear operators to implicitly do"
737 <<
"\nTime to scale the output objects = "
746 *out <<
"\nFinish processing and wrapping the output objects ...\n";
750 epetraOutArgs, W_op, efwdW, eW,
757 <<
"\nTime to finish processing and wrapping the output objects = "
764 THYRA_MODEL_EVALUATOR_DECORATOR_EVAL_MODEL_END();
773 const EpetraExt::ModelEvaluator::InArgs &epetraInArgs,
774 ModelEvaluatorBase::InArgs<double> *inArgs
779 typedef ModelEvaluatorBase MEB;
783 if(inArgs->supports(MEB::IN_ARG_x)) {
787 if(inArgs->supports(MEB::IN_ARG_x_dot)) {
791 if(inArgs->supports(MEB::IN_ARG_x_mp)) {
792 inArgs->set_x_mp( epetraInArgs.get_x_mp() );
795 if(inArgs->supports(MEB::IN_ARG_x_dot_mp)) {
796 inArgs->set_x_dot_mp( epetraInArgs.get_x_dot_mp() );
799 const int l_Np = inArgs->Np();
800 for(
int l = 0; l < l_Np; ++l ) {
803 for(
int l = 0; l < l_Np; ++l ) {
804 if(inArgs->supports(MEB::IN_ARG_p_mp, l))
805 inArgs->set_p_mp( l, epetraInArgs.get_p_mp(l) );
808 if(inArgs->supports(MEB::IN_ARG_t)) {
809 inArgs->set_t(epetraInArgs.get_t());
816 const ModelEvaluatorBase::InArgs<double> &inArgs,
817 EpetraExt::ModelEvaluator::InArgs *epetraInArgs
822 using Teuchos::rcp_const_cast;
823 #ifdef HAVE_THYRA_ME_POLYNOMIAL
825 #endif // HAVE_THYRA_ME_POLYNOMIAL
831 if( inArgs.supports(IN_ARG_x_dot) && (x_dot = inArgs.get_x_dot()).
get() ) {
833 epetraInArgs->set_x_dot(e_x_dot);
836 if( inArgs.supports(IN_ARG_x_dot_mp) && (x_dot_mp = inArgs.get_x_dot_mp()).
get() ) {
837 epetraInArgs->set_x_dot_mp(x_dot_mp);
841 if( inArgs.supports(IN_ARG_x) && (x = inArgs.get_x()).
get() ) {
843 epetraInArgs->set_x(e_x);
846 if( inArgs.supports(IN_ARG_x_mp) && (x_mp = inArgs.get_x_mp()).
get() ) {
847 epetraInArgs->set_x_mp(x_mp);
851 for(
int l = 0; l < inArgs.Np(); ++l ) {
852 p_l = inArgs.get_p(l);
856 for(
int l = 0; l < inArgs.Np(); ++l ) {
857 if (inArgs.supports(IN_ARG_p_mp,l)) {
858 p_mp_l = inArgs.get_p_mp(l);
859 if(p_mp_l.
get()) epetraInArgs->set_p_mp(l,p_mp_l);
863 #ifdef HAVE_THYRA_ME_POLYNOMIAL
868 inArgs.supports(IN_ARG_x_dot_poly)
869 && (x_dot_poly = inArgs.get_x_dot_poly()).
get()
873 rcp(
new Polynomial<Epetra_Vector>(x_dot_poly->degree()));
874 for (
unsigned int i=0; i<=x_dot_poly->degree(); i++) {
875 epetra_ptr = rcp_const_cast<Epetra_Vector>(
877 epetra_x_dot_poly->setCoefficientPtr(i,epetra_ptr);
879 epetraInArgs->set_x_dot_poly(epetra_x_dot_poly);
884 inArgs.supports(IN_ARG_x_poly)
885 && (x_poly = inArgs.get_x_poly()).
get()
889 rcp(
new Polynomial<Epetra_Vector>(x_poly->degree()));
890 for (
unsigned int i=0; i<=x_poly->degree(); i++) {
891 epetra_ptr = rcp_const_cast<Epetra_Vector>(
893 epetra_x_poly->setCoefficientPtr(i,epetra_ptr);
895 epetraInArgs->set_x_poly(epetra_x_poly);
898 #endif // HAVE_THYRA_ME_POLYNOMIAL
900 if( inArgs.supports(IN_ARG_t) )
901 epetraInArgs->set_t(inArgs.get_t());
903 if( inArgs.supports(IN_ARG_alpha) )
904 epetraInArgs->set_alpha(inArgs.get_alpha());
906 if( inArgs.supports(IN_ARG_beta) )
907 epetraInArgs->set_beta(inArgs.get_beta());
909 if( inArgs.supports(IN_ARG_step_size) )
910 epetraInArgs->set_step_size(inArgs.get_step_size());
912 if( inArgs.supports(IN_ARG_stage_number) )
913 epetraInArgs->set_stage_number(inArgs.get_stage_number());
918 const ModelEvaluatorBase::OutArgs<double> &outArgs,
919 EpetraExt::ModelEvaluator::OutArgs *epetraUnscaledOutArgs_inout,
920 RCP<LinearOpBase<double> > *W_op_out,
927 using Teuchos::rcp_const_cast;
928 using Teuchos::rcp_dynamic_cast;
943 EpetraExt::ModelEvaluator::OutArgs &epetraUnscaledOutArgs = *epetraUnscaledOutArgs_inout;
951 if( outArgs.supports(OUT_ARG_f) && (f = outArgs.get_f()).
get() )
954 if( outArgs.supports(OUT_ARG_f_mp) && (f_mp = outArgs.get_f_mp()).
get() )
955 epetraUnscaledOutArgs.set_f_mp(f_mp);
961 for(
int j = 0; j < outArgs.Ng(); ++j ) {
962 g_j = outArgs.get_g(j);
966 for(
int j = 0; j < outArgs.Ng(); ++j ) {
967 if (outArgs.supports(OUT_ARG_g_mp,j)) {
968 g_mp_j = outArgs.get_g_mp(j);
969 if(g_mp_j.
get()) epetraUnscaledOutArgs.set_g_mp(j,g_mp_j);
977 if( outArgs.supports(OUT_ARG_W_mp) && (W_mp = outArgs.get_W_mp()).
get() )
978 epetraUnscaledOutArgs.set_W_mp(W_mp);
984 if (outArgs.supports(OUT_ARG_W_op) &&
nonnull(W_op = outArgs.get_W_op())) {
996 eW = efwdW->epetra_op();
997 epetraUnscaledOutArgs.set_W(eW);
1006 Derivative<double> DfDp_l;
1007 for(
int l = 0; l < outArgs.Np(); ++l ) {
1008 if( !outArgs.supports(OUT_ARG_DfDp,l).none()
1009 && !(DfDp_l = outArgs.get_DfDp(l)).isEmpty() )
1014 MPDerivative DfDp_mp_l;
1015 for(
int l = 0; l < outArgs.Np(); ++l ) {
1016 if( !outArgs.supports(OUT_ARG_DfDp_mp,l).none()
1017 && !(DfDp_mp_l = outArgs.get_DfDp_mp(l)).isEmpty() )
1026 Derivative<double> DgDx_dot_j;
1027 for(
int j = 0; j < outArgs.Ng(); ++j ) {
1028 if( !outArgs.supports(OUT_ARG_DgDx_dot,j).none()
1029 && !(DgDx_dot_j = outArgs.get_DgDx_dot(j)).isEmpty() )
1034 MPDerivative DgDx_dot_mp_j;
1035 for(
int j = 0; j < outArgs.Ng(); ++j ) {
1036 if( !outArgs.supports(OUT_ARG_DgDx_dot_mp,j).none()
1037 && !(DgDx_dot_mp_j = outArgs.get_DgDx_dot_mp(j)).isEmpty() )
1046 Derivative<double> DgDx_j;
1047 for(
int j = 0; j < outArgs.Ng(); ++j ) {
1048 if( !outArgs.supports(OUT_ARG_DgDx,j).none()
1049 && !(DgDx_j = outArgs.get_DgDx(j)).isEmpty() )
1054 MPDerivative DgDx_mp_j;
1055 for(
int j = 0; j < outArgs.Ng(); ++j ) {
1056 if( !outArgs.supports(OUT_ARG_DgDx_mp,j).none()
1057 && !(DgDx_mp_j = outArgs.get_DgDx_mp(j)).isEmpty() )
1066 DerivativeSupport DgDp_j_l_support;
1067 Derivative<double> DgDp_j_l;
1068 for (
int j = 0; j < outArgs.Ng(); ++j ) {
1069 for (
int l = 0; l < outArgs.Np(); ++l ) {
1070 if (!(DgDp_j_l_support = outArgs.supports(OUT_ARG_DgDp,j,l)).none()
1071 && !(DgDp_j_l = outArgs.get_DgDp(j,l)).isEmpty() )
1077 DerivativeSupport DgDp_mp_j_l_support;
1078 MPDerivative DgDp_mp_j_l;
1079 for (
int j = 0; j < outArgs.Ng(); ++j ) {
1080 for (
int l = 0; l < outArgs.Np(); ++l ) {
1081 if (!(DgDp_mp_j_l_support = outArgs.supports(OUT_ARG_DgDp_mp,j,l)).none()
1082 && !(DgDp_mp_j_l = outArgs.get_DgDp_mp(j,l)).isEmpty() )
1090 #ifdef HAVE_THYRA_ME_POLYNOMIAL
1094 if( outArgs.supports(OUT_ARG_f_poly) && (f_poly = outArgs.get_f_poly()).
get() )
1098 for (
unsigned int i=0; i<=f_poly->degree(); i++) {
1101 f_poly->getCoefficient(i)));
1102 epetra_f_poly->setCoefficientPtr(i,epetra_ptr);
1104 epetraUnscaledOutArgs.set_f_poly(epetra_f_poly);
1107 #endif // HAVE_THYRA_ME_POLYNOMIAL
1113 EpetraExt::ModelEvaluator::InArgs *epetraInArgs_inout,
1114 EpetraExt::ModelEvaluator::OutArgs *epetraUnscaledOutArgs_inout,
1120 typedef EpetraExt::ModelEvaluator EME;
1122 #ifdef TEUCHOS_DEBUG
1127 EpetraExt::ModelEvaluator::InArgs
1128 &epetraInArgs = *epetraInArgs_inout;
1129 EpetraExt::ModelEvaluator::OutArgs
1130 &epetraUnscaledOutArgs = *epetraUnscaledOutArgs_inout;
1136 epetraUnscaledOutArgs.supports(EME::OUT_ARG_f)
1138 epetraUnscaledOutArgs.funcOrDerivesAreSet(EME::OUT_ARG_f)
1142 epetraUnscaledOutArgs.supports(EME::OUT_ARG_W)
1144 is_null(epetraUnscaledOutArgs.get_W())
1159 <<
"\nCreating a temporary Epetra W to compute scale factors"
1160 <<
" for f(...) ...\n";
1161 epetraUnscaledOutArgs.set_W(create_and_assert_W(*
epetraModel_));
1162 if( epetraInArgs.supports(EME::IN_ARG_beta) )
1163 epetraInArgs.set_beta(1.0);
1164 if( epetraInArgs.supports(EME::IN_ARG_alpha) )
1165 epetraInArgs.set_alpha(0.0);
1166 if( epetraInArgs.supports(EME::IN_ARG_step_size) )
1167 epetraInArgs.set_step_size(0.0);
1168 if( epetraInArgs.supports(EME::IN_ARG_stage_number) )
1169 epetraInArgs.set_stage_number(1.0);
1176 const EpetraExt::ModelEvaluator::OutArgs &epetraUnscaledOutArgs,
1184 using Teuchos::rcp_const_cast;
1195 ermW = get_Epetra_RowMatrix(epetraUnscaledOutArgs);
1202 invRowSums =
rcp(
new Epetra_Vector(ermW->OperatorRangeMap()));
1206 ermW->InvRowSums(*invRowSums);
1210 <<
"\nComputed inverse row sum scaling from W that"
1211 " will be used to scale f(...) and its derivatives:\n";
1212 double minVal = 0, maxVal = 0, avgVal = 0;
1213 invRowSums->MinValue(&minVal);
1214 invRowSums->MaxValue(&maxVal);
1215 invRowSums->MeanValue(&avgVal);
1218 <<
"min(invRowSums) = " << minVal <<
"\n"
1219 <<
"max(invRowSums) = " << maxVal <<
"\n"
1220 <<
"avg(invRowSums) = " << avgVal <<
"\n";
1243 const EpetraExt::ModelEvaluator::OutArgs &,
1244 RCP<LinearOpBase<double> > &W_op,
1247 const ModelEvaluatorBase::OutArgs<double> &
1251 using Teuchos::rcp_dynamic_cast;
1255 efwdW->setFullyInitialized(
true);
1260 if (W_op.shares_resource(efwdW)) {
1264 rcp_dynamic_cast<
EpetraLinearOp>(W_op,
true)->setFullyInitialized(
true);
1277 typedef EpetraExt::ModelEvaluator EME;
1283 EME::InArgs epetraOrigNominalValues;
1284 EpetraExt::gatherModelNominalValues(
1287 EME::InArgs epetraOrigLowerBounds;
1288 EME::InArgs epetraOrigUpperBounds;
1289 EpetraExt::gatherModelBounds(
1290 *
epetraModel_, &epetraOrigLowerBounds, &epetraOrigUpperBounds );
1299 if( epetraOrigNominalValues.supports(EME::IN_ARG_x_dot) ) {
1302 if( epetraOrigNominalValues.supports(EME::IN_ARG_x) ) {
1309 EME::InArgs epetraScaledNominalValues =
epetraModel_->createInArgs();
1310 EpetraExt::scaleModelVars(
1314 EME::InArgs epetraScaledLowerBounds =
epetraModel_->createInArgs();
1315 EME::InArgs epetraScaledUpperBounds =
epetraModel_->createInArgs();
1316 EpetraExt::scaleModelBounds(
1317 epetraOrigLowerBounds, epetraOrigUpperBounds,
epetraModel_->getInfBound(),
1319 &epetraScaledLowerBounds, &epetraScaledUpperBounds
1347 typedef EpetraExt::ModelEvaluator EME;
1349 const EpetraExt::ModelEvaluator &epetraModel = *
epetraModel_;
1350 EME::InArgs epetraInArgs = epetraModel.createInArgs();
1351 EME::OutArgs epetraOutArgs = epetraModel.createOutArgs();
1352 const int l_Np = epetraOutArgs.Np();
1353 const int l_Ng = epetraOutArgs.Ng();
1359 InArgsSetup<double> inArgs;
1360 inArgs.setModelEvalDescription(this->
description());
1361 inArgs.set_Np(epetraInArgs.Np());
1362 inArgs.setSupports(IN_ARG_x_dot, epetraInArgs.supports(EME::IN_ARG_x_dot));
1363 inArgs.setSupports(IN_ARG_x, epetraInArgs.supports(EME::IN_ARG_x));
1364 inArgs.setSupports(IN_ARG_x_dot_mp, epetraInArgs.supports(EME::IN_ARG_x_dot_mp));
1365 inArgs.setSupports(IN_ARG_x_mp, epetraInArgs.supports(EME::IN_ARG_x_mp));
1366 #ifdef HAVE_THYRA_ME_POLYNOMIAL
1367 inArgs.setSupports(IN_ARG_x_dot_poly,
1368 epetraInArgs.supports(EME::IN_ARG_x_dot_poly));
1369 inArgs.setSupports(IN_ARG_x_poly, epetraInArgs.supports(EME::IN_ARG_x_poly));
1370 #endif // HAVE_THYRA_ME_POLYNOMIAL
1371 inArgs.setSupports(IN_ARG_t, epetraInArgs.supports(EME::IN_ARG_t));
1372 inArgs.setSupports(IN_ARG_alpha, epetraInArgs.supports(EME::IN_ARG_alpha));
1373 inArgs.setSupports(IN_ARG_beta, epetraInArgs.supports(EME::IN_ARG_beta));
1374 inArgs.setSupports(IN_ARG_step_size, epetraInArgs.supports(EME::IN_ARG_step_size));
1375 inArgs.setSupports(IN_ARG_stage_number, epetraInArgs.supports(EME::IN_ARG_stage_number));
1376 for(
int l=0; l<l_Np; ++l) {
1377 inArgs.setSupports(IN_ARG_p_mp, l, epetraInArgs.supports(EME::IN_ARG_p_mp, l));
1385 OutArgsSetup<double> outArgs;
1386 outArgs.setModelEvalDescription(this->
description());
1387 outArgs.set_Np_Ng(l_Np, l_Ng);
1389 outArgs.setSupports(OUT_ARG_f, epetraOutArgs.supports(EME::OUT_ARG_f));
1390 outArgs.setSupports(OUT_ARG_f_mp, epetraOutArgs.supports(EME::OUT_ARG_f_mp));
1391 if (outArgs.supports(OUT_ARG_f)) {
1393 outArgs.setSupports(OUT_ARG_W_op, epetraOutArgs.supports(EME::OUT_ARG_W));
1394 outArgs.set_W_properties(
convert(epetraOutArgs.get_W_properties()));
1396 for(
int l=0; l<l_Np; ++l) {
1397 outArgs.setSupports(OUT_ARG_DfDp, l,
1398 convert(epetraOutArgs.supports(EME::OUT_ARG_DfDp, l)));
1399 if(!outArgs.supports(OUT_ARG_DfDp, l).none())
1400 outArgs.set_DfDp_properties(l,
1401 convert(epetraOutArgs.get_DfDp_properties(l)));
1402 if (outArgs.supports(OUT_ARG_f_mp))
1404 outArgs.setSupports(OUT_ARG_DfDp_mp, l,
1405 convert(epetraOutArgs.supports(EME::OUT_ARG_DfDp_mp, l)));
1406 if(!outArgs.supports(OUT_ARG_DfDp_mp, l).none())
1407 outArgs.set_DfDp_mp_properties(l,
1408 convert(epetraOutArgs.get_DfDp_mp_properties(l)));
1413 for(
int j=0; j<l_Ng; ++j) {
1414 if (inArgs.supports(IN_ARG_x_dot))
1415 outArgs.setSupports(OUT_ARG_DgDx_dot, j,
1416 convert(epetraOutArgs.supports(EME::OUT_ARG_DgDx_dot, j)));
1417 if(!outArgs.supports(OUT_ARG_DgDx_dot, j).none())
1418 outArgs.set_DgDx_dot_properties(j,
1419 convert(epetraOutArgs.get_DgDx_dot_properties(j)));
1420 if (inArgs.supports(IN_ARG_x))
1421 outArgs.setSupports(OUT_ARG_DgDx, j,
1422 convert(epetraOutArgs.supports(EME::OUT_ARG_DgDx, j)));
1423 if(!outArgs.supports(OUT_ARG_DgDx, j).none())
1424 outArgs.set_DgDx_properties(j,
1425 convert(epetraOutArgs.get_DgDx_properties(j)));
1426 if (inArgs.supports(IN_ARG_x_dot_mp))
1427 outArgs.setSupports(OUT_ARG_DgDx_dot_mp, j,
1428 convert(epetraOutArgs.supports(EME::OUT_ARG_DgDx_dot_mp, j)));
1429 if(!outArgs.supports(OUT_ARG_DgDx_dot_mp, j).none())
1430 outArgs.set_DgDx_dot_mp_properties(j,
1431 convert(epetraOutArgs.get_DgDx_dot_mp_properties(j)));
1432 if (inArgs.supports(IN_ARG_x_mp))
1433 outArgs.setSupports(OUT_ARG_DgDx_mp, j,
1434 convert(epetraOutArgs.supports(EME::OUT_ARG_DgDx_mp, j)));
1435 if(!outArgs.supports(OUT_ARG_DgDx_mp, j).none())
1436 outArgs.set_DgDx_mp_properties(j,
1437 convert(epetraOutArgs.get_DgDx_mp_properties(j)));
1440 for(
int j=0; j < l_Ng; ++j)
for(
int l=0; l < l_Np; ++l) {
1441 const EME::DerivativeSupport epetra_DgDp_j_l_support =
1442 epetraOutArgs.supports(EME::OUT_ARG_DgDp, j, l);
1443 outArgs.setSupports(OUT_ARG_DgDp, j, l,
1444 convert(epetra_DgDp_j_l_support));
1445 if(!outArgs.supports(OUT_ARG_DgDp, j, l).none())
1446 outArgs.set_DgDp_properties(j, l,
1447 convert(epetraOutArgs.get_DgDp_properties(j, l)));
1448 const EME::DerivativeSupport epetra_DgDp_mp_j_l_support =
1449 epetraOutArgs.supports(EME::OUT_ARG_DgDp_mp, j, l);
1450 outArgs.setSupports(OUT_ARG_DgDp_mp, j, l,
1451 convert(epetra_DgDp_mp_j_l_support));
1452 if(!outArgs.supports(OUT_ARG_DgDp_mp, j, l).none())
1453 outArgs.set_DgDp_mp_properties(j, l,
1454 convert(epetraOutArgs.get_DgDp_mp_properties(j, l)));
1456 for(
int j=0; j<l_Ng; ++j) {
1457 outArgs.setSupports(OUT_ARG_g_mp, j, epetraOutArgs.supports(EME::OUT_ARG_g_mp, j));
1459 #ifdef HAVE_THYRA_ME_POLYNOMIAL
1460 outArgs.setSupports(OUT_ARG_f_poly,
1461 epetraOutArgs.supports(EME::OUT_ARG_f_poly));
1462 #endif // HAVE_THYRA_ME_POLYNOMIAL
1474 return Thyra::partialNonconstEpetraLinearOp(
1490 Thyra::epetraModelEvaluator(
1491 const RCP<const EpetraExt::ModelEvaluator> &epetraModel,
1492 const RCP<LinearOpWithSolveFactoryBase<double> > &W_factory
1495 return Teuchos::rcp(
new EpetraModelEvaluator(epetraModel,W_factory));
1499 Thyra::ModelEvaluatorBase::EDerivativeMultiVectorOrientation
1501 const EpetraExt::ModelEvaluator::EDerivativeMultiVectorOrientation &mvOrientation
1504 switch(mvOrientation) {
1505 case EpetraExt::ModelEvaluator::DERIV_MV_BY_COL :
1506 return ModelEvaluatorBase::DERIV_MV_BY_COL;
1507 case EpetraExt::ModelEvaluator::DERIV_TRANS_MV_BY_ROW :
1508 return ModelEvaluatorBase::DERIV_TRANS_MV_BY_ROW;
1516 EpetraExt::ModelEvaluator::EDerivativeMultiVectorOrientation
1518 const ModelEvaluatorBase::EDerivativeMultiVectorOrientation &mvOrientation
1521 switch(mvOrientation) {
1522 case ModelEvaluatorBase::DERIV_MV_BY_COL :
1523 return EpetraExt::ModelEvaluator::DERIV_MV_BY_COL;
1524 case ModelEvaluatorBase::DERIV_TRANS_MV_BY_ROW :
1525 return EpetraExt::ModelEvaluator::DERIV_TRANS_MV_BY_ROW;
1533 Thyra::ModelEvaluatorBase::DerivativeProperties
1535 const EpetraExt::ModelEvaluator::DerivativeProperties &derivativeProperties
1538 ModelEvaluatorBase::EDerivativeLinearity linearity;
1539 switch(derivativeProperties.linearity) {
1540 case EpetraExt::ModelEvaluator::DERIV_LINEARITY_UNKNOWN:
1541 linearity = ModelEvaluatorBase::DERIV_LINEARITY_UNKNOWN;
1543 case EpetraExt::ModelEvaluator::DERIV_LINEARITY_CONST:
1544 linearity = ModelEvaluatorBase::DERIV_LINEARITY_CONST;
1546 case EpetraExt::ModelEvaluator::DERIV_LINEARITY_NONCONST:
1547 linearity = ModelEvaluatorBase::DERIV_LINEARITY_NONCONST;
1552 ModelEvaluatorBase::ERankStatus rank;
1553 switch(derivativeProperties.rank) {
1554 case EpetraExt::ModelEvaluator::DERIV_RANK_UNKNOWN:
1555 rank = ModelEvaluatorBase::DERIV_RANK_UNKNOWN;
1557 case EpetraExt::ModelEvaluator::DERIV_RANK_FULL:
1558 rank = ModelEvaluatorBase::DERIV_RANK_FULL;
1560 case EpetraExt::ModelEvaluator::DERIV_RANK_DEFICIENT:
1561 rank = ModelEvaluatorBase::DERIV_RANK_DEFICIENT;
1566 return ModelEvaluatorBase::DerivativeProperties(
1567 linearity,rank,derivativeProperties.supportsAdjoint);
1571 Thyra::ModelEvaluatorBase::DerivativeSupport
1573 const EpetraExt::ModelEvaluator::DerivativeSupport &derivativeSupport
1576 ModelEvaluatorBase::DerivativeSupport ds;
1577 if(derivativeSupport.supports(EpetraExt::ModelEvaluator::DERIV_LINEAR_OP))
1578 ds.plus(ModelEvaluatorBase::DERIV_LINEAR_OP);
1579 if(derivativeSupport.supports(EpetraExt::ModelEvaluator::DERIV_MV_BY_COL))
1580 ds.plus(ModelEvaluatorBase::DERIV_MV_BY_COL);
1581 if(derivativeSupport.supports(EpetraExt::ModelEvaluator::DERIV_TRANS_MV_BY_ROW))
1582 ds.plus(ModelEvaluatorBase::DERIV_TRANS_MV_BY_ROW);
1587 EpetraExt::ModelEvaluator::Derivative
1589 const ModelEvaluatorBase::Derivative<double> &derivative,
1590 const RCP<const Epetra_Map> &fnc_map,
1591 const RCP<const Epetra_Map> &var_map
1594 typedef ModelEvaluatorBase MEB;
1595 if(derivative.getLinearOp().get()) {
1596 return EpetraExt::ModelEvaluator::Derivative(
1597 Teuchos::rcp_const_cast<Epetra_Operator>(
1598 Teuchos::dyn_cast<const EpetraLinearOp>(*derivative.getLinearOp()).epetra_op()
1602 else if(derivative.getDerivativeMultiVector().getMultiVector().get()) {
1603 return EpetraExt::ModelEvaluator::Derivative(
1604 EpetraExt::ModelEvaluator::DerivativeMultiVector(
1606 ( derivative.getDerivativeMultiVector().getOrientation() == MEB::DERIV_MV_BY_COL
1610 ,derivative.getDerivativeMultiVector().getMultiVector()
1612 ,
convert(derivative.getDerivativeMultiVector().getOrientation())
1616 return EpetraExt::ModelEvaluator::Derivative();
1618 EpetraExt::ModelEvaluator::MPDerivative
1620 const ModelEvaluatorBase::MPDerivative &derivative,
1626 if(derivative.getLinearOp().get()) {
1627 return EpetraExt::ModelEvaluator::MPDerivative(
1628 derivative.getLinearOp()
1631 else if(derivative.getDerivativeMultiVector().getMultiVector().get()) {
1632 return EpetraExt::ModelEvaluator::MPDerivative(
1633 EpetraExt::ModelEvaluator::MPDerivativeMultiVector(
1634 derivative.getDerivativeMultiVector().getMultiVector()
1635 ,
convert(derivative.getDerivativeMultiVector().getOrientation())
1639 return EpetraExt::ModelEvaluator::MPDerivative();
RCP< Epetra_MultiVector > get_Epetra_MultiVector(const Epetra_Map &map, const RCP< MultiVectorBase< double > > &mv)
Get a non-const Epetra_MultiVector view from a non-const MultiVectorBase object if possible...
const ModelEvaluatorBase::InArgs< double > & getFinalPoint() const
ModelEvaluatorBase::InArgs< double > getLowerBounds() const
void setNominalValues(const ModelEvaluatorBase::InArgs< double > &nominalValues)
Set the nominal values.
void updateInArgsOutArgs() const
void convertInArgsFromEpetraToThyra(const EpetraExt::ModelEvaluator::InArgs &epetraInArgs, ModelEvaluatorBase::InArgs< double > *inArgs) const
Concrete LinearOpBase adapter subclass for Epetra_Operator object.
RCP< const EpetraExt::ModelEvaluator > getEpetraModel() const
RCP< const Epetra_Map > x_map_
std::string typeName(const T &t)
void reportFinalPoint(const ModelEvaluatorBase::InArgs< double > &finalPoint, const bool wasSolved)
bool is_null(const boost::shared_ptr< T > &p)
void preEvalScalingSetup(EpetraExt::ModelEvaluator::InArgs *epetraInArgs, EpetraExt::ModelEvaluator::OutArgs *epetraUnscaledOutArgs, const RCP< Teuchos::FancyOStream > &out, const Teuchos::EVerbosityLevel verbLevel) const
ModelEvaluatorBase::OutArgs< double > createOutArgsImpl() const
RCP< LinearOpBase< double > > create_DgDx_dot_op_impl(int j) const
RCP< const Teuchos::ParameterList > getValidParameters() const
RCP< LinearOpBase< double > > create_W_op() const
ModelEvaluatorBase::InArgs< double > getNominalValues() const
RCP< Teuchos::ParameterList > getNonconstParameterList()
p_map_is_local_t g_map_is_local_
RCP< LinearOpWithSolveFactoryBase< double > > W_factory_
RCP< const VectorSpaceBase< double > > create_VectorSpace(const RCP< const Epetra_Map > &epetra_map)
Create an VectorSpaceBase object given an Epetra_Map object.
EStateFunctionScaling stateFunctionScaling_
void convertOutArgsFromThyraToEpetra(const ModelEvaluatorBase::OutArgs< double > &outArgs, EpetraExt::ModelEvaluator::OutArgs *epetraUnscaledOutArgs, RCP< LinearOpBase< double > > *W_op, RCP< EpetraLinearOp > *efwdW, RCP< Epetra_Operator > *eW) const
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
RCP< PreconditionerBase< double > > create_W_prec() const
Returns null currently.
void setStateVariableScalingVec(const RCP< const Epetra_Vector > &stateVariableScalingVec)
Set the state variable scaling vector s_x (see above).
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
ModelEvaluatorBase::InArgs< double > createInArgs() const
basic_OSTab< char > OSTab
EpetraExt::ModelEvaluator::OutArgs epetraOutArgsScaling_
void convertInArgsFromThyraToEpetra(const ModelEvaluatorBase::InArgs< double > &inArgs, EpetraExt::ModelEvaluator::InArgs *epetraInArgs) const
RCP< const VectorSpaceBase< double > > x_space_
RCP< const Epetra_Vector > getStateVariableScalingVec() const
Get the inverse state variable scaling vector inv_s_x (see above).
ModelEvaluatorBase::OutArgs< double > prototypeOutArgs_
#define TEUCHOS_ASSERT_IN_RANGE_UPPER_EXCLUSIVE(index, lower_inclusive, upper_exclusive)
void start(bool reset=false)
RCP< const Teuchos::ParameterList > getParameterList() const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
RCP< const LinearOpWithSolveFactoryBase< double > > get_W_factory() const
EpetraExt::ModelEvaluator::MPDerivative convert(const ModelEvaluatorBase::MPDerivative &derivative, const RCP< const Epetra_Map > &fnc_map, const RCP< const Epetra_Map > &var_map)
RCP< const Epetra_Vector > getStateVariableInvScalingVec() const
Get the state variable scaling vector s_x (see above).
RCP< VectorBase< double > > create_Vector(const RCP< Epetra_Vector > &epetra_v, const RCP< const VectorSpaceBase< double > > &space=Teuchos::null)
Create a non-const VectorBase object from a non-const Epetra_Vector object.
RCP< const VectorSpaceBase< double > > f_space_
ModelEvaluatorBase::InArgs< double > getUpperBounds() const
ModelEvaluatorBase::InArgs< double > finalPoint_
void setStateFunctionScalingVec(const RCP< const Epetra_Vector > &stateFunctionScalingVec)
Set the state function scaling vector s_f (see above).
RCP< LinearOpBase< double > > create_DgDx_op_impl(int j) const
bool finalPointWasSolved() const
ModelEvaluatorBase::InArgs< double > lowerBounds_
RCP< const Epetra_Vector > getStateFunctionScalingVec() const
Get the state function scaling vector s_f (see above).
RCP< Epetra_Vector > get_Epetra_Vector(const Epetra_Map &map, const RCP< VectorBase< double > > &v)
Get a non-const Epetra_Vector view from a non-const VectorBase object if possible.
RCP< EpetraLinearOp > create_epetra_W_op() const
void initialize(const RCP< const EpetraExt::ModelEvaluator > &epetraModel, const RCP< LinearOpWithSolveFactoryBase< double > > &W_factory)
RCP< const VectorSpaceBase< double > > get_f_space() const
void updateNominalValuesAndBounds() const
RCP< LinearOpBase< double > > create_DfDp_op_impl(int l) const
RCP< const VectorSpaceBase< double > > get_x_space() const
RCP< const VectorSpaceBase< double > > get_g_space(int j) const
EpetraExt::ModelEvaluator::InArgs epetraInArgsScaling_
RCP< const VectorSpaceBase< double > > get_p_space(int l) const
bool currentInArgsOutArgs_
void resize(size_type new_size, const value_type &x=value_type())
bool finalPointWasSolved_
RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
void setParameterList(RCP< Teuchos::ParameterList > const ¶mList)
TEUCHOSCORE_LIB_DLL_EXPORT bool includesVerbLevel(const EVerbosityLevel verbLevel, const EVerbosityLevel requestedVerbLevel, const bool isDefaultLevel=false)
TypeTo implicit_cast(const TypeFrom &t)
const RCP< T > & assert_not_null() const
RCP< LinearOpBase< double > > create_DgDp_op_impl(int j, int l) const
std::string description() const
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
RCP< const Epetra_Vector > stateFunctionScalingVec_
ModelEvaluatorBase::EDerivativeMultiVectorOrientation convert(const EpetraExt::ModelEvaluator::EDerivativeMultiVectorOrientation &mvOrientation)
void postEvalScalingSetup(const EpetraExt::ModelEvaluator::OutArgs &epetraUnscaledOutArgs, const RCP< Teuchos::FancyOStream > &out, const Teuchos::EVerbosityLevel verbLevel) const
void uninitialize(RCP< const EpetraExt::ModelEvaluator > *epetraModel=NULL, RCP< LinearOpWithSolveFactoryBase< double > > *W_factory=NULL)
bool nonnull(const boost::shared_ptr< T > &p)
Teuchos::ArrayView< const std::string > get_g_names(int j) const
void evalModelImpl(const ModelEvaluatorBase::InArgs< double > &inArgs, const ModelEvaluatorBase::OutArgs< double > &outArgs) const
bool nominalValuesAndBoundsAreUpdated_
void finishConvertingOutArgsFromEpetraToThyra(const EpetraExt::ModelEvaluator::OutArgs &epetraOutArgs, RCP< LinearOpBase< double > > &W_op, RCP< EpetraLinearOp > &efwdW, RCP< Epetra_Operator > &eW, const ModelEvaluatorBase::OutArgs< double > &outArgs) const
ModelEvaluatorBase::InArgs< double > prototypeInArgs_
RCP< const Epetra_Map > f_map_
double totalElapsedTime(bool readCurrentTime=false) const
RCP< const Epetra_Vector > stateVariableScalingVec_
RCP< const EpetraExt::ModelEvaluator > epetraModel_
RCP< Teuchos::ParameterList > paramList_
ModelEvaluatorBase::InArgs< double > nominalValues_
#define TEUCHOS_ASSERT(assertion_test)
p_map_is_local_t p_map_is_local_
RCP< Teuchos::ParameterList > unsetParameterList()
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
RCP< const Epetra_Vector > invStateVariableScalingVec_
ModelEvaluatorBase::InArgs< double > upperBounds_