10#include "Thyra_EpetraModelEvaluator.hpp"
11#include "Thyra_LinearOpWithSolveFactoryHelpers.hpp"
13#include "Thyra_EpetraLinearOp.hpp"
14#include "Thyra_DetachedMultiVectorView.hpp"
15#include "Thyra_ModelEvaluatorDelegatorBase.hpp"
16#include "EpetraExt_ModelEvaluatorScalingTools.h"
17#include "Epetra_RowMatrix.h"
18#include "Teuchos_Time.hpp"
19#include "Teuchos_implicit_cast.hpp"
20#include "Teuchos_Assert.hpp"
21#include "Teuchos_StandardParameterEntryValidators.hpp"
22#include "Teuchos_VerboseObjectParameterListHelpers.hpp"
28const std::string StateFunctionScaling_name =
"State Function Scaling";
31 Thyra::EpetraModelEvaluator::EStateFunctionScaling
34stateFunctionScalingValidator;
35const 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>
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)
101 :nominalValuesAndBoundsAreUpdated_(false), stateFunctionScaling_(STATE_FUNC_SCALING_NONE),
102 currentInArgsOutArgs_(false), finalPointWasSolved_(false)
116 epetraModel_ = epetraModel;
118 W_factory_ = W_factory;
120 x_map_ = epetraModel_->get_x_map();
121 f_map_ = epetraModel_->get_f_map();
122 if (!is_null(x_map_)) {
127 EpetraExt::ModelEvaluator::InArgs inArgs = epetraModel_->createInArgs();
128 p_map_.
resize(inArgs.Np()); p_space_.
resize(inArgs.Np());
129 p_map_is_local_.resize(inArgs.Np(),
false);
130 for(
int l = 0; l < implicit_cast<int>(p_space_.
size()); ++l ) {
132 p_map_l = ( p_map_[l] = epetraModel_->get_p_map(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!");
140 p_map_is_local_[l] = !p_map_l->DistributedGlobal();
144 EpetraExt::ModelEvaluator::OutArgs outArgs = epetraModel_->createOutArgs();
145 g_map_.
resize(outArgs.Ng()); g_space_.
resize(outArgs.Ng());
146 g_map_is_local_.resize(outArgs.Ng(),
false);
147 for(
int j = 0; j < implicit_cast<int>(g_space_.
size()); ++j ) {
149 g_map_j = ( g_map_[j] = epetraModel_->get_g_map(j) );
150 g_map_is_local_[j] = !g_map_j->DistributedGlobal();
154 epetraInArgsScaling_ = epetraModel_->createInArgs();
155 epetraOutArgsScaling_ = epetraModel_->createOutArgs();
156 nominalValuesAndBoundsAreUpdated_ =
false;
157 finalPointWasSolved_ =
false;
160 currentInArgsOutArgs_ =
false;
175 nominalValues_.
setArgs(nominalValues);
190 nominalValuesAndBoundsAreUpdated_ =
false;
197 return stateVariableScalingVec_;
204 updateNominalValuesAndBounds();
205 return invStateVariableScalingVec_;
213 stateFunctionScalingVec_ = stateFunctionScalingVec;
220 return stateFunctionScalingVec_;
229 if(epetraModel) *epetraModel = epetraModel_;
230 if(W_factory) *W_factory = W_factory_;
236 currentInArgsOutArgs_ =
false;
249 return finalPointWasSolved_;
258 std::ostringstream oss;
259 oss <<
"Thyra::EpetraModelEvaluator{";
260 oss <<
"epetraModel=";
261 if(epetraModel_.
get())
262 oss <<
"\'"<<epetraModel_->description()<<
"\'";
265 oss <<
",W_factory=";
267 oss <<
"\'"<<W_factory_->description()<<
"\'";
284 paramList_ = paramList;
285 const EStateFunctionScaling stateFunctionScaling_old = stateFunctionScaling_;
286 stateFunctionScaling_ = stateFunctionScalingValidator->getIntegralValue(
287 *paramList_, StateFunctionScaling_name, StateFunctionScaling_default
289 if( stateFunctionScaling_ != stateFunctionScaling_old )
291 Teuchos::readVerboseObjectSublist(&*paramList_,
this);
326 using Teuchos::tuple;
327 using Teuchos::rcp_implicit_cast;
330 if(is_null(validPL)) {
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>(
348 STATE_FUNC_SCALING_NONE,
349 STATE_FUNC_SCALING_ROW_SUM
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);
373 return p_space_.
size();
379 return g_space_.
size();
413 return epetraModel_->get_p_names(l);
431 return epetraModel_->get_g_names(j);
438 updateNominalValuesAndBounds();
439 return nominalValues_;
446 updateNominalValuesAndBounds();
454 updateNominalValuesAndBounds();
462 return this->create_epetra_W_op();
482 if (!currentInArgsOutArgs_)
483 updateInArgsOutArgs();
484 return prototypeInArgs_;
494 finalPoint_.
setArgs(finalPoint);
495 finalPointWasSolved_ = wasSolved;
503EpetraModelEvaluator::create_DfDp_op_impl(
int )
const
510RCP<LinearOpBase<double> >
511EpetraModelEvaluator::create_DgDx_dot_op_impl(
int )
const
518RCP<LinearOpBase<double> >
519EpetraModelEvaluator::create_DgDx_op_impl(
int )
const
526RCP<LinearOpBase<double> >
527EpetraModelEvaluator::create_DgDp_op_impl(
int ,
int )
const
534ModelEvaluatorBase::OutArgs<double>
535EpetraModelEvaluator::createOutArgsImpl()
const
537 if (!currentInArgsOutArgs_) {
538 updateInArgsOutArgs();
540 return prototypeOutArgs_;
544void EpetraModelEvaluator::evalModelImpl(
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;
562 this->updateNominalValuesAndBounds();
573 inArgs.setArgs(inArgs_in);
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
593 stateFunctionScaling_ != STATE_FUNC_SCALING_NONE
594 &&
is_null(stateFunctionScalingVec_)
598 VOTSLOWSF W_factory_outputTempState(W_factory_,out,verbLevel);
612 *out <<
"\nSetting-up/creating input arguments ...\n";
616 EME::InArgs epetraScaledInArgs = epetraModel_->createInArgs();
617 convertInArgsFromThyraToEpetra( inArgs, &epetraScaledInArgs );
621 EME::InArgs epetraInArgs = epetraModel_->createInArgs();
622 EpetraExt::unscaleModelVars(
623 epetraScaledInArgs, epetraInArgsScaling_, &epetraInArgs,
629 OSTab(out).o() <<
"\nTime to setup InArgs = "<<timer.totalElapsedTime()<<
" sec\n";
636 *out <<
"\nSetting-up/creating output arguments ...\n";
641 EME::OutArgs epetraUnscaledOutArgs = epetraModel_->createOutArgs();
645 RCP<LinearOpBase<double> > W_op;
646 RCP<EpetraLinearOp> efwdW;
647 RCP<Epetra_Operator> eW;
651 convertOutArgsFromThyraToEpetra(
653 &epetraUnscaledOutArgs,
661 if (firstTimeStateFuncScaling) {
662 preEvalScalingSetup(&epetraInArgs,&epetraUnscaledOutArgs,out,verbLevel);
668 <<
"\nTime to setup OutArgs = "
669 << timer.totalElapsedTime() <<
" sec\n";
676 *out <<
"\nEvaluating the Epetra output functions ...\n";
679 epetraModel_->evalModel(epetraInArgs,epetraUnscaledOutArgs);
684 <<
"\nTime to evaluate Epetra output functions = "
685 << timer.totalElapsedTime() <<
" sec\n";
696 *out <<
"\nCompute scale factors if needed ...\n";
699 if (firstTimeStateFuncScaling) {
700 postEvalScalingSetup(epetraUnscaledOutArgs,out,verbLevel);
706 <<
"\nTime to compute scale factors = "
707 << timer.totalElapsedTime() <<
" sec\n";
714 *out <<
"\nScale the output objects ...\n";
717 EME::OutArgs epetraOutArgs = epetraModel_->createOutArgs();
718 bool allFuncsWhereScaled =
false;
719 EpetraExt::scaleModelFuncs(
720 epetraUnscaledOutArgs, epetraInArgsScaling_, epetraOutArgsScaling_,
721 &epetraOutArgs, &allFuncsWhereScaled,
725 if (stateFunctionScaling_ != STATE_FUNC_SCALING_NONE)
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 = "
738 << timer.totalElapsedTime() <<
" sec\n";
746 *out <<
"\nFinish processing and wrapping the output objects ...\n";
749 finishConvertingOutArgsFromEpetraToThyra(
750 epetraOutArgs, W_op, efwdW, eW,
757 <<
"\nTime to finish processing and wrapping the output objects = "
758 << timer.totalElapsedTime() <<
" sec\n";
764 THYRA_MODEL_EVALUATOR_DECORATOR_EVAL_MODEL_END();
772void EpetraModelEvaluator::convertInArgsFromEpetraToThyra(
773 const EpetraExt::ModelEvaluator::InArgs &epetraInArgs,
774 ModelEvaluatorBase::InArgs<double> *inArgs
783 if(inArgs->supports(MEB::IN_ARG_x)) {
784 inArgs->set_x(
create_Vector( epetraInArgs.get_x(), x_space_ ) );
787 if(inArgs->supports(MEB::IN_ARG_x_dot)) {
788 inArgs->set_x_dot(
create_Vector( epetraInArgs.get_x_dot(), x_space_ ) );
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 ) {
801 inArgs->set_p( l,
create_Vector( epetraInArgs.get_p(l), p_space_[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());
815void EpetraModelEvaluator::convertInArgsFromThyraToEpetra(
816 const ModelEvaluatorBase::InArgs<double> &inArgs,
817 EpetraExt::ModelEvaluator::InArgs *epetraInArgs
822 using Teuchos::rcp_const_cast;
823#ifdef HAVE_THYRA_ME_POLYNOMIAL
830 RCP<const VectorBase<double> > x_dot;
831 if( inArgs.supports(
IN_ARG_x_dot) && (x_dot = inArgs.get_x_dot()).get() ) {
833 epetraInArgs->set_x_dot(e_x_dot);
835 RCP<const Stokhos::ProductEpetraVector > x_dot_mp;
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);
840 RCP<const VectorBase<double> > x;
841 if( inArgs.supports(
IN_ARG_x) && (x = inArgs.get_x()).get() ) {
843 epetraInArgs->set_x(e_x);
845 RCP<const Stokhos::ProductEpetraVector > x_mp;
846 if( inArgs.supports(
IN_ARG_x_mp) && (x_mp = inArgs.get_x_mp()).get() ) {
847 epetraInArgs->set_x_mp(x_mp);
850 RCP<const VectorBase<double> > p_l;
851 for(
int l = 0; l < inArgs.Np(); ++l ) {
852 p_l = inArgs.get_p(l);
855 RCP<const Stokhos::ProductEpetraVector > p_mp_l;
856 for(
int l = 0; l < inArgs.Np(); ++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
865 RCP<const Polynomial< VectorBase<double> > > x_dot_poly;
866 RCP<Epetra_Vector> epetra_ptr;
869 && (x_dot_poly = inArgs.get_x_dot_poly()).get()
872 RCP<Polynomial<Epetra_Vector> > epetra_x_dot_poly =
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);
882 RCP<const Polynomial< VectorBase<double> > > x_poly;
885 && (x_poly = inArgs.get_x_poly()).get()
888 RCP<Polynomial<Epetra_Vector> > epetra_x_poly =
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);
901 epetraInArgs->set_t(inArgs.get_t());
904 epetraInArgs->set_alpha(inArgs.get_alpha());
907 epetraInArgs->set_beta(inArgs.get_beta());
910 epetraInArgs->set_step_size(inArgs.get_step_size());
913 epetraInArgs->set_stage_number(inArgs.get_stage_number());
917void EpetraModelEvaluator::convertOutArgsFromThyraToEpetra(
918 const ModelEvaluatorBase::OutArgs<double> &outArgs,
919 EpetraExt::ModelEvaluator::OutArgs *epetraUnscaledOutArgs_inout,
920 RCP<LinearOpBase<double> > *W_op_out,
921 RCP<EpetraLinearOp> *efwdW_out,
922 RCP<Epetra_Operator> *eW_out
927 using Teuchos::rcp_const_cast;
928 using Teuchos::rcp_dynamic_cast;
943 EpetraExt::ModelEvaluator::OutArgs &epetraUnscaledOutArgs = *epetraUnscaledOutArgs_inout;
944 RCP<LinearOpBase<double> > &W_op = *W_op_out;
945 RCP<EpetraLinearOp> &efwdW = *efwdW_out;
946 RCP<Epetra_Operator> &eW = *eW_out;
950 RCP<VectorBase<double> > f;
951 if( outArgs.supports(
OUT_ARG_f) && (f = outArgs.get_f()).get() )
953 RCP<Stokhos::ProductEpetraVector > f_mp;
954 if( outArgs.supports(
OUT_ARG_f_mp) && (f_mp = outArgs.get_f_mp()).get() )
955 epetraUnscaledOutArgs.set_f_mp(f_mp);
960 RCP<VectorBase<double> > g_j;
961 for(
int j = 0; j < outArgs.Ng(); ++j ) {
962 g_j = outArgs.get_g(j);
965 RCP<Stokhos::ProductEpetraVector > g_mp_j;
966 for(
int j = 0; j < outArgs.Ng(); ++j ) {
968 g_mp_j = outArgs.get_g_mp(j);
969 if(g_mp_j.get()) epetraUnscaledOutArgs.set_g_mp(j,g_mp_j);
976 RCP<Stokhos::ProductEpetraOperator > W_mp;
977 if( outArgs.supports(
OUT_ARG_W_mp) && (W_mp = outArgs.get_W_mp()).get() )
978 epetraUnscaledOutArgs.set_W_mp(W_mp);
986 efwdW = rcp_const_cast<EpetraLinearOp>(
987 rcp_dynamic_cast<const EpetraLinearOp>(W_op,
true));
996 eW = efwdW->epetra_op();
997 epetraUnscaledOutArgs.set_W(eW);
1006 Derivative<double> DfDp_l;
1007 for(
int l = 0; l < outArgs.Np(); ++l ) {
1009 && !(DfDp_l = outArgs.get_DfDp(l)).isEmpty() )
1011 epetraUnscaledOutArgs.set_DfDp(l,
convert(DfDp_l,f_map_,p_map_[l]));
1014 MPDerivative DfDp_mp_l;
1015 for(
int l = 0; l < outArgs.Np(); ++l ) {
1017 && !(DfDp_mp_l = outArgs.get_DfDp_mp(l)).isEmpty() )
1019 epetraUnscaledOutArgs.set_DfDp_mp(l,
convert(DfDp_mp_l,f_map_,p_map_[l]));
1026 Derivative<double> DgDx_dot_j;
1027 for(
int j = 0; j < outArgs.Ng(); ++j ) {
1029 && !(DgDx_dot_j = outArgs.get_DgDx_dot(j)).isEmpty() )
1031 epetraUnscaledOutArgs.set_DgDx_dot(j,
convert(DgDx_dot_j,g_map_[j],x_map_));
1034 MPDerivative DgDx_dot_mp_j;
1035 for(
int j = 0; j < outArgs.Ng(); ++j ) {
1037 && !(DgDx_dot_mp_j = outArgs.get_DgDx_dot_mp(j)).isEmpty() )
1039 epetraUnscaledOutArgs.set_DgDx_dot_mp(j,
convert(DgDx_dot_mp_j,g_map_[j],x_map_));
1046 Derivative<double> DgDx_j;
1047 for(
int j = 0; j < outArgs.Ng(); ++j ) {
1049 && !(DgDx_j = outArgs.get_DgDx(j)).isEmpty() )
1051 epetraUnscaledOutArgs.set_DgDx(j,
convert(DgDx_j,g_map_[j],x_map_));
1054 MPDerivative DgDx_mp_j;
1055 for(
int j = 0; j < outArgs.Ng(); ++j ) {
1057 && !(DgDx_mp_j = outArgs.get_DgDx_mp(j)).isEmpty() )
1059 epetraUnscaledOutArgs.set_DgDx_mp(j,
convert(DgDx_mp_j,g_map_[j],x_map_));
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() )
1073 epetraUnscaledOutArgs.set_DgDp(j,l,
convert(DgDp_j_l,g_map_[j],p_map_[l]));
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() )
1084 epetraUnscaledOutArgs.set_DgDp_mp(j,l,
convert(DgDp_mp_j_l,g_map_[j],p_map_[l]));
1090#ifdef HAVE_THYRA_ME_POLYNOMIAL
1093 RCP<const Teuchos::Polynomial< VectorBase<double> > > f_poly;
1094 if( outArgs.supports(
OUT_ARG_f_poly) && (f_poly = outArgs.get_f_poly()).get() )
1096 RCP<Teuchos::Polynomial<Epetra_Vector> > epetra_f_poly =
1098 for (
unsigned int i=0; i<=f_poly->degree(); i++) {
1099 RCP<Epetra_Vector> epetra_ptr
1101 f_poly->getCoefficient(i)));
1102 epetra_f_poly->setCoefficientPtr(i,epetra_ptr);
1104 epetraUnscaledOutArgs.set_f_poly(epetra_f_poly);
1112void EpetraModelEvaluator::preEvalScalingSetup(
1113 EpetraExt::ModelEvaluator::InArgs *epetraInArgs_inout,
1114 EpetraExt::ModelEvaluator::OutArgs *epetraUnscaledOutArgs_inout,
1115 const RCP<Teuchos::FancyOStream> &out,
1120 typedef EpetraExt::ModelEvaluator EME;
1127 EpetraExt::ModelEvaluator::InArgs
1128 &epetraInArgs = *epetraInArgs_inout;
1129 EpetraExt::ModelEvaluator::OutArgs
1130 &epetraUnscaledOutArgs = *epetraUnscaledOutArgs_inout;
1133 ( stateFunctionScaling_ == STATE_FUNC_SCALING_ROW_SUM )
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);
1175void EpetraModelEvaluator::postEvalScalingSetup(
1176 const EpetraExt::ModelEvaluator::OutArgs &epetraUnscaledOutArgs,
1177 const RCP<Teuchos::FancyOStream> &out,
1184 using Teuchos::rcp_const_cast;
1188 switch(stateFunctionScaling_) {
1190 case STATE_FUNC_SCALING_ROW_SUM: {
1194 const RCP<Epetra_RowMatrix>
1195 ermW = get_Epetra_RowMatrix(epetraUnscaledOutArgs);
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";
1223 stateFunctionScalingVec_ = invRowSums;
1234 epetraOutArgsScaling_ = epetraModel_->createOutArgs();
1236 epetraOutArgsScaling_.set_f(
1237 rcp_const_cast<Epetra_Vector>(stateFunctionScalingVec_) );
1242void EpetraModelEvaluator::finishConvertingOutArgsFromEpetraToThyra(
1243 const EpetraExt::ModelEvaluator::OutArgs &,
1244 RCP<LinearOpBase<double> > &W_op,
1245 RCP<EpetraLinearOp> &efwdW,
1246 RCP<Epetra_Operator> &,
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);
1271void EpetraModelEvaluator::updateNominalValuesAndBounds()
const
1277 typedef EpetraExt::ModelEvaluator EME;
1279 if( !nominalValuesAndBoundsAreUpdated_ ) {
1283 EME::InArgs epetraOrigNominalValues;
1284 EpetraExt::gatherModelNominalValues(
1285 *epetraModel_, &epetraOrigNominalValues );
1287 EME::InArgs epetraOrigLowerBounds;
1288 EME::InArgs epetraOrigUpperBounds;
1289 EpetraExt::gatherModelBounds(
1290 *epetraModel_, &epetraOrigLowerBounds, &epetraOrigUpperBounds );
1294 epetraInArgsScaling_ = epetraModel_->createInArgs();
1296 if( !
is_null(stateVariableScalingVec_) ) {
1297 invStateVariableScalingVec_
1298 = EpetraExt::createInverseModelScalingVector(stateVariableScalingVec_);
1299 if( epetraOrigNominalValues.supports(EME::IN_ARG_x_dot) ) {
1300 epetraInArgsScaling_.set_x_dot(invStateVariableScalingVec_);
1302 if( epetraOrigNominalValues.supports(EME::IN_ARG_x) ) {
1303 epetraInArgsScaling_.set_x(invStateVariableScalingVec_);
1309 EME::InArgs epetraScaledNominalValues = epetraModel_->createInArgs();
1310 EpetraExt::scaleModelVars(
1311 epetraOrigNominalValues, epetraInArgsScaling_, &epetraScaledNominalValues
1314 EME::InArgs epetraScaledLowerBounds = epetraModel_->createInArgs();
1315 EME::InArgs epetraScaledUpperBounds = epetraModel_->createInArgs();
1316 EpetraExt::scaleModelBounds(
1317 epetraOrigLowerBounds, epetraOrigUpperBounds, epetraModel_->getInfBound(),
1318 epetraInArgsScaling_,
1319 &epetraScaledLowerBounds, &epetraScaledUpperBounds
1327 convertInArgsFromEpetraToThyra(epetraScaledNominalValues, &nominalValues_);
1328 convertInArgsFromEpetraToThyra(epetraScaledLowerBounds, &lowerBounds_);
1329 convertInArgsFromEpetraToThyra(epetraScaledUpperBounds, &upperBounds_);
1331 nominalValuesAndBoundsAreUpdated_ =
true;
1344void EpetraModelEvaluator::updateInArgsOutArgs()
const
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
1368 epetraInArgs.supports(EME::IN_ARG_x_dot_poly));
1369 inArgs.setSupports(
IN_ARG_x_poly, epetraInArgs.supports(EME::IN_ARG_x_poly));
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));
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));
1379 prototypeInArgs_ = inArgs;
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));
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) {
1398 convert(epetraOutArgs.supports(EME::OUT_ARG_DfDp, l)));
1400 outArgs.set_DfDp_properties(l,
1401 convert(epetraOutArgs.get_DfDp_properties(l)));
1405 convert(epetraOutArgs.supports(EME::OUT_ARG_DfDp_mp, l)));
1407 outArgs.set_DfDp_mp_properties(l,
1408 convert(epetraOutArgs.get_DfDp_mp_properties(l)));
1413 for(
int j=0; j<l_Ng; ++j) {
1416 convert(epetraOutArgs.supports(EME::OUT_ARG_DgDx_dot, j)));
1418 outArgs.set_DgDx_dot_properties(j,
1419 convert(epetraOutArgs.get_DgDx_dot_properties(j)));
1422 convert(epetraOutArgs.supports(EME::OUT_ARG_DgDx, j)));
1424 outArgs.set_DgDx_properties(j,
1425 convert(epetraOutArgs.get_DgDx_properties(j)));
1428 convert(epetraOutArgs.supports(EME::OUT_ARG_DgDx_dot_mp, j)));
1430 outArgs.set_DgDx_dot_mp_properties(j,
1431 convert(epetraOutArgs.get_DgDx_dot_mp_properties(j)));
1434 convert(epetraOutArgs.supports(EME::OUT_ARG_DgDx_mp, j)));
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);
1444 convert(epetra_DgDp_j_l_support));
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);
1451 convert(epetra_DgDp_mp_j_l_support));
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
1461 epetraOutArgs.supports(EME::OUT_ARG_f_poly));
1463 prototypeOutArgs_ = outArgs;
1466 currentInArgsOutArgs_ =
true;
1472EpetraModelEvaluator::create_epetra_W_op()
const
1474 return Thyra::partialNonconstEpetraLinearOp(
1476 create_and_assert_W(*epetraModel_)
1490Thyra::epetraModelEvaluator(
1491 const RCP<const EpetraExt::ModelEvaluator> &epetraModel,
1492 const RCP<LinearOpWithSolveFactoryBase<double> > &W_factory
1495 return Teuchos::rcp(
new EpetraModelEvaluator(epetraModel,W_factory));
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;
1516EpetraExt::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;
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);
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);
1587EpetraExt::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(
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();
1618EpetraExt::ModelEvaluator::MPDerivative
1620 const ModelEvaluatorBase::MPDerivative &derivative,
1621 const RCP<const Epetra_Map> &,
1622 const RCP<const Epetra_Map> &
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();
void resize(size_type new_size, const value_type &x=value_type())
const RCP< T > & assert_not_null() const
RCP< Teuchos::ParameterList > unsetParameterList()
RCP< Teuchos::ParameterList > getNonconstParameterList()
bool finalPointWasSolved() const
void setStateFunctionScalingVec(const RCP< const Epetra_Vector > &stateFunctionScalingVec)
Set the state function scaling vector s_f (see above).
RCP< const Epetra_Vector > getStateFunctionScalingVec() const
Get the state function scaling vector s_f (see above).
std::string description() const
void initialize(const RCP< const EpetraExt::ModelEvaluator > &epetraModel, const RCP< LinearOpWithSolveFactoryBase< double > > &W_factory)
RCP< const VectorSpaceBase< double > > get_f_space() const
void setStateVariableScalingVec(const RCP< const Epetra_Vector > &stateVariableScalingVec)
Set the state variable scaling vector s_x (see above).
RCP< LinearOpBase< double > > create_W_op() const
ModelEvaluatorBase::InArgs< double > getUpperBounds() const
void setNominalValues(const ModelEvaluatorBase::InArgs< double > &nominalValues)
Set the nominal values.
RCP< PreconditionerBase< double > > create_W_prec() const
Returns null currently.
RCP< const LinearOpWithSolveFactoryBase< double > > get_W_factory() const
ModelEvaluatorBase::InArgs< double > getNominalValues() const
RCP< const Teuchos::ParameterList > getParameterList() const
void reportFinalPoint(const ModelEvaluatorBase::InArgs< double > &finalPoint, const bool wasSolved)
ModelEvaluatorBase::InArgs< double > createInArgs() const
void setParameterList(RCP< Teuchos::ParameterList > const ¶mList)
RCP< const VectorSpaceBase< double > > get_g_space(int j) const
RCP< const EpetraExt::ModelEvaluator > getEpetraModel() const
void uninitialize(RCP< const EpetraExt::ModelEvaluator > *epetraModel=NULL, RCP< LinearOpWithSolveFactoryBase< double > > *W_factory=NULL)
RCP< const VectorSpaceBase< double > > get_x_space() const
RCP< const Epetra_Vector > getStateVariableInvScalingVec() const
Get the state variable scaling vector s_x (see above).
ModelEvaluatorBase::InArgs< double > getLowerBounds() const
RCP< const Epetra_Vector > getStateVariableScalingVec() const
Get the inverse state variable scaling vector inv_s_x (see above).
RCP< const Teuchos::ParameterList > getValidParameters() const
Teuchos::ArrayView< const std::string > get_g_names(int j) const
RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
RCP< const VectorSpaceBase< double > > get_p_space(int l) const
const ModelEvaluatorBase::InArgs< double > & getFinalPoint() const
ModelEvaluatorBase::EDerivativeMultiVectorOrientation convert(const EpetraExt::ModelEvaluator::EDerivativeMultiVectorOrientation &mvOrientation)
Factory interface for creating LinearOpWithSolveBase objects from compatible LinearOpBase objects.
Determines the forms of a general derivative that are supported.
Concrete aggregate class for all input arguments computable by a ModelEvaluator subclass object.
void setArgs(const InArgs< Scalar > &inArgs, bool ignoreUnsupported=false, bool cloneObjects=false)
Set non-null arguments (does not overwrite non-NULLs with NULLs) .
bool supports(EInArgsMembers arg) const
Determines if an input argument is supported or not.
Base subclass for ModelEvaluator that defines some basic types.
ModelEvaluatorBase()
constructor
EDerivativeMultiVectorOrientation
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 > > create_VectorSpace(const RCP< const Epetra_Map > &epetra_map)
Create an VectorSpaceBase object given an Epetra_Map object.
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< 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.
#define TEUCHOS_ASSERT(assertion_test)
#define TEUCHOS_ASSERT_IN_RANGE_UPPER_EXCLUSIVE(index, lower_inclusive, upper_exclusive)
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
bool is_null(const std::shared_ptr< T > &p)
bool nonnull(const std::shared_ptr< T > &p)
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
TypeTo implicit_cast(const TypeFrom &t)
std::string typeName(const T &t)
T_To & dyn_cast(T_From &from)
basic_OSTab< char > OSTab
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
TEUCHOSCORE_LIB_DLL_EXPORT bool includesVerbLevel(const EVerbosityLevel verbLevel, const EVerbosityLevel requestedVerbLevel, const bool isDefaultLevel=false)
Simple public strict containing properties of a derivative object.