42 #ifndef THYRA_DIRECTIONAL_FINITE_DIFF_CALCULATOR_DEF_HPP
43 #define THYRA_DIRECTIONAL_FINITE_DIFF_CALCULATOR_DEF_HPP
46 #include "Thyra_DirectionalFiniteDiffCalculator_decl.hpp"
47 #include "Thyra_ModelEvaluatorHelpers.hpp"
48 #include "Thyra_DetachedVectorView.hpp"
49 #include "Thyra_DetachedMultiVectorView.hpp"
50 #include "Thyra_StateFuncModelEvaluatorBase.hpp"
51 #include "Thyra_MultiVectorStdOps.hpp"
52 #include "Thyra_VectorStdOps.hpp"
53 #include "Teuchos_TimeMonitor.hpp"
54 #include "Teuchos_VerboseObjectParameterListHelpers.hpp"
60 namespace DirectionalFiniteDiffCalculatorTypes {
73 template<
class Scalar>
74 class OutArgsCreator :
public StateFuncModelEvaluatorBase<Scalar>
78 RCP<const VectorSpaceBase<Scalar> > get_x_space()
const
80 RCP<const VectorSpaceBase<Scalar> > get_f_space()
const
82 ModelEvaluatorBase::InArgs<Scalar> createInArgs()
const
84 ModelEvaluatorBase::OutArgs<Scalar> createOutArgs()
const
87 const ModelEvaluatorBase::InArgs<Scalar> &inArgs,
88 const ModelEvaluatorBase::OutArgs<Scalar> &outArgs
92 static ModelEvaluatorBase::OutArgs<Scalar> createOutArgs(
93 const ModelEvaluator<Scalar> &model,
94 const SelectedDerivatives &fdDerivatives
100 const MEB::OutArgs<Scalar> wrappedOutArgs = model.createOutArgs();
101 const int Np = wrappedOutArgs.Np(),
Ng = wrappedOutArgs.Ng();
102 MEB::OutArgsSetup<Scalar> outArgs;
104 outArgs.setModelEvalDescription(
105 "DirectionalFiniteDiffCalculator: " + model.description()
111 outArgs.set_Np_Ng(Np,Ng);
112 outArgs.setSupports(wrappedOutArgs);
116 const SelectedDerivatives::supports_DfDp_t
117 &supports_DfDp = fdDerivatives.supports_DfDp_;
119 SelectedDerivatives::supports_DfDp_t::const_iterator
120 itr = supports_DfDp.begin();
121 itr != supports_DfDp.end();
126 assert_p_space(model,l);
127 outArgs.setSupports(MEB::OUT_ARG_DfDp,l,MEB::DERIV_MV_BY_COL);
132 const SelectedDerivatives::supports_DgDp_t
133 &supports_DgDp = fdDerivatives.supports_DgDp_;
135 SelectedDerivatives::supports_DgDp_t::const_iterator
136 itr = supports_DgDp.begin();
137 itr != supports_DgDp.end();
141 const int j = itr->first;
142 const int l = itr->second;
143 assert_p_space(model,l);
144 outArgs.setSupports(MEB::OUT_ARG_DgDp,j,l,MEB::DERIV_MV_BY_COL);
153 static void assert_p_space(
const ModelEvaluator<Scalar> &model,
const int l )
156 const bool p_space_l_is_in_core = model.get_p_space(l)->hasInCoreView();
158 !p_space_l_is_in_core, std::logic_error,
159 "Error, for the model " << model.description()
160 <<
", the space p_space("<<l<<
") must be in-core so that they can"
161 " act as the domain space for the multi-vector derivative!"
178 template<
class Scalar>
179 const std::string& DirectionalFiniteDiffCalculator<Scalar>::FDMethod_name()
181 static std::string loc_FDMethod_name =
"FD Method";
182 return loc_FDMethod_name;
186 template<
class Scalar>
189 Thyra::DirectionalFiniteDiffCalculatorTypes::EFDMethodType
192 DirectionalFiniteDiffCalculator<Scalar>::fdMethodValidator()
194 static RCP<Teuchos::StringToIntegralParameterEntryValidator<EFDMethodType> >
195 loc_fdMethodValidator
198 Teuchos::tuple<std::string>(
204 ,
"order-four-central"
207 ,Teuchos::tuple<std::string>(
208 "Use O(eps) one sided finite differences (cramped bounds)"
209 ,
"Use O(eps^2) one sided finite differences (cramped bounds)"
210 ,
"Use O(eps^2) two sided central finite differences"
211 ,
"Use \"order-two-central\" when not cramped by bounds, otherwise use \"order-two\""
212 ,
"Use O(eps^4) one sided finite differences (cramped bounds)"
213 ,
"Use O(eps^4) two sided central finite differences"
214 ,
"Use \"order-four-central\" when not cramped by bounds, otherwise use \"order-four\""
216 ,Teuchos::tuple<Thyra::DirectionalFiniteDiffCalculatorTypes::EFDMethodType>(
217 Thyra::DirectionalFiniteDiffCalculatorTypes::FD_ORDER_ONE
218 ,Thyra::DirectionalFiniteDiffCalculatorTypes::FD_ORDER_TWO
219 ,Thyra::DirectionalFiniteDiffCalculatorTypes::FD_ORDER_TWO_CENTRAL
220 ,Thyra::DirectionalFiniteDiffCalculatorTypes::FD_ORDER_TWO_AUTO
221 ,Thyra::DirectionalFiniteDiffCalculatorTypes::FD_ORDER_FOUR
222 ,Thyra::DirectionalFiniteDiffCalculatorTypes::FD_ORDER_FOUR_CENTRAL
223 ,Thyra::DirectionalFiniteDiffCalculatorTypes::FD_ORDER_FOUR_AUTO
228 return loc_fdMethodValidator;
232 template<
class Scalar>
234 DirectionalFiniteDiffCalculator<Scalar>::FDMethod_default()
236 static std::string loc_FDMethod_default =
"order-one";
237 return loc_FDMethod_default;
241 template<
class Scalar>
243 DirectionalFiniteDiffCalculator<Scalar>::FDStepSelectType_name()
245 static std::string loc_FDStepSelectType_name =
"FD Step Select Type";
246 return loc_FDStepSelectType_name;
250 template<
class Scalar>
253 Thyra::DirectionalFiniteDiffCalculatorTypes::EFDStepSelectType
256 DirectionalFiniteDiffCalculator<Scalar>::fdStepSelectTypeValidator()
258 static const RCP<Teuchos::StringToIntegralParameterEntryValidator<EFDStepSelectType> >
259 loc_fdStepSelectTypeValidator
262 Teuchos::tuple<std::string>(
266 ,Teuchos::tuple<std::string>(
267 "Use absolute step size \""+FDStepLength_name()+
"\""
268 ,
"Use relative step size \""+FDStepLength_name()+
"\"*||xo||inf"
270 ,Teuchos::tuple<Thyra::DirectionalFiniteDiffCalculatorTypes::EFDStepSelectType>(
271 Thyra::DirectionalFiniteDiffCalculatorTypes::FD_STEP_ABSOLUTE
272 ,Thyra::DirectionalFiniteDiffCalculatorTypes::FD_STEP_RELATIVE
277 return loc_fdStepSelectTypeValidator;
281 template<
class Scalar>
283 DirectionalFiniteDiffCalculator<Scalar>::FDStepSelectType_default()
285 static std::string loc_FDStepSelectType_default =
"Absolute";
286 return loc_FDStepSelectType_default;
290 template<
class Scalar>
292 DirectionalFiniteDiffCalculator<Scalar>::FDStepLength_name()
294 static std::string loc_FDStepLength_name =
"FD Step Length";
295 return loc_FDStepLength_name;
299 template<
class Scalar>
301 DirectionalFiniteDiffCalculator<Scalar>::FDStepLength_default()
303 static double loc_FDStepLength_default = -1.0;
304 return loc_FDStepLength_default;
311 template<
class Scalar>
318 :fd_method_type_(fd_method_type_in)
319 ,fd_step_select_type_(fd_step_select_type_in)
320 ,fd_step_size_(fd_step_size_in)
321 ,fd_step_size_min_(fd_step_size_min_in)
328 template<
class Scalar>
335 paramList_ = paramList;
337 *paramList_, FDMethod_name(), FDMethod_default());
338 fd_step_select_type_ = fdStepSelectTypeValidator()->getIntegralValue(
339 *paramList_, FDStepSelectType_name(), FDStepSelectType_default());
340 fd_step_size_ = paramList_->get(
341 FDStepLength_name(),FDStepLength_default());
342 Teuchos::readVerboseObjectSublist(&*paramList_,
this);
346 template<
class Scalar>
354 template<
class Scalar>
359 paramList_ = Teuchos::null;
364 template<
class Scalar>
372 template<
class Scalar>
376 using Teuchos::rcp_implicit_cast;
380 pl = Teuchos::parameterList();
382 FDMethod_name(), FDMethod_default(),
383 "The method used to compute the finite differences.",
384 rcp_implicit_cast<const PEV>(fdMethodValidator())
387 FDStepSelectType_name(), FDStepSelectType_default(),
388 "Method used to select the finite difference step length.",
389 rcp_implicit_cast<const PEV>(fdStepSelectTypeValidator())
392 FDStepLength_name(), FDStepLength_default()
393 ,
"The length of the finite difference step to take.\n"
394 "A value of < 0.0 means that the step length will be determined automatically."
396 Teuchos::setupVerboseObjectSublist(&*pl);
402 template<
class Scalar>
409 return DirectionalFiniteDiffCalculatorTypes::OutArgsCreator<Scalar>::createOutArgs(
410 model, fdDerivatives );
414 template<
class Scalar>
426 THYRA_FUNC_TIME_MONITOR(
427 string(
"Thyra::DirectionalFiniteDiffCalculator<")+ST::name()+
">::calcVariations(...)"
435 namespace DFDCT = DirectionalFiniteDiffCalculatorTypes;
443 if(out.
get() && trace)
444 *out <<
"\nEntering DirectionalFiniteDiffCalculator<Scalar>::calcVariations(...)\n";
446 if(out.
get() && trace)
448 <<
"\nbasePoint=\n" << describe(bp,verbLevel)
449 <<
"\ndirections=\n" << describe(dir,verbLevel)
450 <<
"\nbaseFunctionValues=\n" << describe(bfunc,verbLevel)
455 var.
isEmpty(), std::logic_error,
456 "Error, all of the variations can not be null!"
488 switch(this->fd_method_type()) {
489 case DFDCT::FD_ORDER_ONE:
490 if(out.
get()&&trace) *out<<
"\nUsing one-sided, first-order finite differences ...\n";
492 case DFDCT::FD_ORDER_TWO:
493 if(out.
get()&&trace) *out<<
"\nUsing one-sided, second-order finite differences ...\n";
495 case DFDCT::FD_ORDER_TWO_CENTRAL:
496 if(out.
get()&&trace) *out<<
"\nUsing second-order central finite differences ...\n";
498 case DFDCT::FD_ORDER_TWO_AUTO:
499 if(out.
get()&&trace) *out<<
"\nUsing auto selection of some second-order finite difference method ...\n";
501 case DFDCT::FD_ORDER_FOUR:
502 if(out.
get()&&trace) *out<<
"\nUsing one-sided, fourth-order finite differences ...\n";
504 case DFDCT::FD_ORDER_FOUR_CENTRAL:
505 if(out.
get()&&trace) *out<<
"\nUsing fourth-order central finite differences ...\n";
507 case DFDCT::FD_ORDER_FOUR_AUTO:
508 if(out.
get()&&trace) *out<<
"\nUsing auto selection of some fourth-order finite difference method ...\n";
522 sqrt_epsilon = SMT::squareroot(SMT::eps()),
523 u_optimal_1 = sqrt_epsilon,
524 u_optimal_2 = SMT::squareroot(sqrt_epsilon),
525 u_optimal_4 = SMT::squareroot(u_optimal_2);
528 bp_norm = SMT::zero();
533 switch(this->fd_method_type()) {
534 case DFDCT::FD_ORDER_ONE:
535 uh_opt = u_optimal_1 * ( fd_step_select_type() == DFDCT::FD_STEP_ABSOLUTE ? 1.0 : bp_norm + 1.0 );
537 case DFDCT::FD_ORDER_TWO:
538 case DFDCT::FD_ORDER_TWO_CENTRAL:
539 case DFDCT::FD_ORDER_TWO_AUTO:
540 uh_opt = u_optimal_2 * ( fd_step_select_type() == DFDCT::FD_STEP_ABSOLUTE ? 1.0 : bp_norm + 1.0 );
542 case DFDCT::FD_ORDER_FOUR:
543 case DFDCT::FD_ORDER_FOUR_CENTRAL:
544 case DFDCT::FD_ORDER_FOUR_AUTO:
545 uh_opt = u_optimal_4 * ( fd_step_select_type() == DFDCT::FD_STEP_ABSOLUTE ? 1.0 : bp_norm + 1.0 );
551 if(out.
get()&&trace) *out<<
"\nDefault optimal step length uh_opt = " << uh_opt <<
" ...\n";
558 uh = this->fd_step_size();
562 else if( fd_step_select_type() == DFDCT::FD_STEP_RELATIVE )
563 uh *= (bp_norm + 1.0);
565 if(out.
get()&&trace) *out<<
"\nStep size to be used uh="<<uh<<
"\n";
579 DFDCT::EFDMethodType l_fd_method_type = this->fd_method_type();
580 switch(l_fd_method_type) {
581 case DFDCT::FD_ORDER_TWO_AUTO:
582 l_fd_method_type = DFDCT::FD_ORDER_TWO_CENTRAL;
584 case DFDCT::FD_ORDER_FOUR_AUTO:
585 l_fd_method_type = DFDCT::FD_ORDER_FOUR_CENTRAL;
595 p_saved = out->precision();
600 const int Np = var.
Np(), Ng = var.
Ng();
603 VectorPtr per_x, per_x_dot, per_x_dot_dot;
604 std::vector<VectorPtr> per_p(Np);
608 if( dir.
get_x().get() )
611 pp.set_x(bp.
get_x());
613 if( bp.
supports(MEB::IN_ARG_x_dot) ) {
615 pp.set_x_dot(per_x_dot=createMember(model.
get_x_space()));
619 if( bp.
supports(MEB::IN_ARG_x_dot_dot) ) {
621 pp.set_x_dot_dot(per_x_dot_dot=createMember(model.
get_x_space()));
625 for(
int l = 0; l < Np; ++l ) {
626 if( dir.
get_p(l).get() )
627 pp.set_p(l,per_p[l]=createMember(model.
get_p_space(l)));
629 pp.set_p(l,bp.
get_p(l));
631 if(out.
get() && trace)
633 <<
"\nperturbedPoint after initial setup (with some uninitialized vectors) =\n"
634 << describe(pp,verbLevel);
637 bool all_funcs_at_base_computed =
true;
641 if( var.
supports(MEB::OUT_ARG_f) && (f=var.
get_f()).
get() ) {
643 assign(f.ptr(),ST::zero());
644 if(!bfunc.
get_f().get()) all_funcs_at_base_computed =
false;
646 for(
int j = 0; j < Ng; ++j ) {
648 if( (g_j=var.
get_g(j)).
get() ) {
650 assign(g_j.ptr(),ST::zero());
651 if(!bfunc.
get_g(j).get()) all_funcs_at_base_computed =
false;
655 if(out.
get() && trace)
657 <<
"\nperturbedFunctions after initial setup (with some uninitialized vectors) =\n"
658 << describe(pfunc,verbLevel);
660 const int dbl_p = 15;
662 *out << std::setprecision(dbl_p);
670 switch(l_fd_method_type) {
671 case DFDCT::FD_ORDER_ONE:
675 case DFDCT::FD_ORDER_TWO:
679 case DFDCT::FD_ORDER_TWO_CENTRAL:
683 case DFDCT::FD_ORDER_FOUR:
687 case DFDCT::FD_ORDER_FOUR_CENTRAL:
694 for(
int eval_i = 1; eval_i <= num_evals; ++eval_i ) {
699 switch(l_fd_method_type) {
700 case DFDCT::FD_ORDER_ONE: {
713 case DFDCT::FD_ORDER_TWO: {
730 case DFDCT::FD_ORDER_TWO_CENTRAL: {
743 case DFDCT::FD_ORDER_FOUR: {
768 case DFDCT::FD_ORDER_FOUR_CENTRAL: {
789 case DFDCT::FD_ORDER_TWO_AUTO:
790 case DFDCT::FD_ORDER_FOUR_AUTO:
796 if(out.
get() && trace)
797 *out <<
"\neval_i="<<eval_i<<
", uh_i="<<uh_i<<
", wgt_i="<<wgt_i<<
"\n";
801 if(uh_i == ST::zero()) {
802 MEB::OutArgs<Scalar> bfuncall;
803 if(!all_funcs_at_base_computed) {
806 if( pfunc.supports(MEB::OUT_ARG_f) && pfunc.get_f().get() && !bfunc.
get_f().get() ) {
809 for(
int j = 0; j < Ng; ++j ) {
810 if( pfunc.get_g(j).get() && !bfunc.
get_g(j).get() ) {
811 bfuncall.set_g(j,createMember(model.
get_g_space(j)));
815 bfuncall.setArgs(bfunc);
821 if(out.
get() && trace)
822 *out <<
"\nSetting variations = wgt_i * basePoint ...\n";
824 if( pfunc.supports(MEB::OUT_ARG_f) && (f=var.
get_f()).
get() ) {
825 V_StV<Scalar>(f.ptr(), wgt_i, *bfuncall.get_f());
827 for(
int j = 0; j < Ng; ++j ) {
829 if( (g_j=var.
get_g(j)).
get() ) {
830 V_StV<Scalar>(g_j.ptr(), wgt_i, *bfuncall.get_g(j));
835 if(out.
get() && trace)
836 *out <<
"\nSetting perturbedPoint = basePoint + uh_i*uh*direction ...\n";
840 V_StVpV(per_x.ptr(),as<Scalar>(uh_i*uh),*dir.
get_x(),*bp.
get_x());
845 for (
int l = 0; l < Np; ++l ) {
846 if( dir.
get_p(l).get() )
847 V_StVpV(per_p[l].ptr(), as<Scalar>(uh_i*uh), *dir.
get_p(l), *bp.
get_p(l));
850 if(out.
get() && trace)
851 *out <<
"\nperturbedPoint=\n" << describe(pp,verbLevel);
853 if(out.
get() && trace)
854 *out <<
"\nCompute perturbedFunctions at perturbedPoint...\n";
856 if(out.
get() && trace)
857 *out <<
"\nperturbedFunctions=\n" << describe(pfunc,verbLevel);
861 if(out.
get() && trace)
862 *out <<
"\nComputing variations += wgt_i*perturbedfunctions ...\n";
864 if( pfunc.supports(MEB::OUT_ARG_f) && (f=pfunc.get_f()).
get() )
865 Vp_StV<Scalar>(var.
get_f().ptr(), wgt_i, *f);
866 for(
int j = 0; j < Ng; ++j ) {
868 if( (g_j=pfunc.get_g(j)).
get() )
869 Vp_StV<Scalar>(var.
get_g(j).ptr(), wgt_i, *g_j);
873 if(out.
get() && trace)
874 *out <<
"\nvariations=\n" << describe(var,verbLevel);
883 const Scalar alpha = ST::one()/(dwgt*uh);
884 if(out.
get() && trace)
886 <<
"\nComputing variations *= (1.0)/(dwgt*uh),"
887 <<
" where (1.0)/(dwgt*uh) = (1.0)/("<<dwgt<<
"*"<<uh<<
") = "<<alpha<<
" ...\n";
891 for(
int j = 0; j < Ng; ++j ) {
893 if( (g_j=var.
get_g(j)).
get() )
894 Vt_S(g_j.ptr(),alpha);
896 if(out.
get() && trace)
897 *out <<
"\nFinal variations=\n" << describe(var,verbLevel);
901 *out << std::setprecision(p_saved);
903 if(out.
get() && trace)
904 *out <<
"\nLeaving DirectionalFiniteDiffCalculator<Scalar>::calcVariations(...)\n";
909 template<
class Scalar>
921 THYRA_FUNC_TIME_MONITOR(
922 string(
"Thyra::DirectionalFiniteDiffCalculator<")+ST::name()+
">::calcDerivatives(...)"
934 if(out.
get() && trace)
935 *out <<
"\nEntering DirectionalFiniteDiffCalculator<Scalar>::calcDerivatives(...)\n";
937 if(out.
get() && trace)
939 <<
"\nbasePoint=\n" << describe(bp,verbLevel)
940 <<
"\nbaseFunctionValues=\n" << describe(bfunc,verbLevel)
946 const int Np = bp.
Np(), Ng = bfunc.
Ng();
949 MultiVectorPtr DfDp_l;
950 std::vector<MEB::DerivativeMultiVector<Scalar> > DgDp_l(Ng);
951 std::vector<VectorPtr> var_g(Ng);
952 for(
int l = 0; l < Np; ++l ) {
953 if(out.
get() && trace)
954 *out <<
"\nComputing derivatives for parameter subvector p("<<l<<
") ...\n";
959 bool hasDerivObject =
false;
962 !deriv.
supports(MEB::OUT_ARG_DfDp,l).none()
966 hasDerivObject =
true;
967 std::ostringstream name; name <<
"DfDp("<<l<<
")";
968 DfDp_l = get_mv(deriv.
get_DfDp(l),name.str(),MEB::DERIV_MV_BY_COL);
971 DfDp_l = Teuchos::null;
974 for (
int j = 0; j < Ng; ++j ) {
976 !deriv.
supports(MEB::OUT_ARG_DgDp,j,l).none()
981 hasDerivObject =
true;
982 std::ostringstream name; name <<
"DgDp("<<j<<
","<<l<<
")";
983 DgDp_l[j] = get_dmv(deriv.
get_DgDp(j,l),name.str());
984 if( DgDp_l[j].getMultiVector().get() && !var_g[j].get() )
991 DgDp_l[j] = MEB::DerivativeMultiVector<Scalar>();
992 var_g[j] = Teuchos::null;
998 if (hasDerivObject) {
999 VectorPtr e_i = createMember(model.
get_p_space(l));
1001 assign(e_i.ptr(),ST::zero());
1002 const int np_l = e_i->space()->dim();
1003 for(
int i = 0 ; i < np_l; ++ i ) {
1004 if(out.
get() && trace)
1005 *out <<
"\nComputing derivatives for single variable p("<<l<<
")("<<i<<
") ...\n";
1007 if(DfDp_l.get()) var.set_f(DfDp_l->col(i));
1008 for(
int j = 0; j < Ng; ++j) {
1009 MultiVectorPtr DgDp_j_l;
1010 if( (DgDp_j_l=DgDp_l[j].getMultiVector()).get() ) {
1011 var.set_g(j,var_g[j]);
1014 set_ele(i,ST::one(),e_i.ptr());
1015 this->calcVariations(
1016 model,bp,dir,bfunc,var
1018 set_ele(i,ST::zero(),e_i.ptr());
1019 if (DfDp_l.get()) var.set_f(Teuchos::null);
1020 for (
int j = 0; j < Ng; ++j) {
1021 MultiVectorPtr DgDp_j_l;
1022 if ( !
is_null(DgDp_j_l=DgDp_l[j].getMultiVector()) ) {
1023 assign( DgDp_j_l->col(i).ptr(), *var_g[j] );
1027 dir.set_p(l,Teuchos::null);
1031 if(out.
get() && trace)
1033 <<
"\nderivatives=\n" << describe(deriv,verbLevel);
1035 if(out.
get() && trace)
1036 *out <<
"\nLeaving DirectionalFiniteDiffCalculator<Scalar>::calcDerivatives(...)\n";
1044 #endif // THYRA_DIRECTIONAL_FINITE_DIFF_CALCULATOR_DEF_HPP
virtual ModelEvaluatorBase::OutArgs< Scalar > createOutArgs() const =0
Create an empty output functions/derivatives object that can be set up and passed to evalModel()...
int Ng() const
Return the number of axillary response functions g(j)(...) supported (Ng >= 0).
virtual RCP< const VectorSpaceBase< Scalar > > get_x_space() const =0
Return the vector space for the state variables x <: RE^n_x.
ModelEvaluatorBase::OutArgs< Scalar > createOutArgs(const ModelEvaluator< Scalar > &model, const SelectedDerivatives &fdDerivatives)
Create an augmented out args object for holding finite difference objects.
Pure abstract base interface for evaluating a stateless "model" that can be mapped into a number of d...
Simple utility class used to select finite difference derivatives for OutArgs object.
void setParameterList(RCP< ParameterList > const ¶mList)
bool is_null(const boost::shared_ptr< T > &p)
DirectionalFiniteDiffCalculatorTypes::EFDStepSelectType EFDStepSelectType
Concrete aggregate class for all output arguments computable by a ModelEvaluator subclass object...
Derivative< Scalar > get_DfDp(int l) const
Precondition: supports(OUT_ARG_DfDp,l)==true.
RCP< const VectorBase< Scalar > > get_x_dot() const
Precondition: supports(IN_ARG_x_dot)==true.
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual void evalModel(const ModelEvaluatorBase::InArgs< Scalar > &inArgs, const ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const =0
Compute all of the requested functions/derivatives at the given point.
Evaluation< VectorBase< Scalar > > get_g(int j) const
Precondition: supports(OUT_ARG_g)==true..
Evaluation< VectorBase< Scalar > > get_f() const
Precondition: supports(OUT_ARG_f)==true.
virtual ModelEvaluatorBase::InArgs< Scalar > createInArgs() const =0
Create an empty input arguments object that can be set up and passed to evalModel().
bool supports(EOutArgsMembers arg) const
Determine if an input argument is supported or not.
Derivative< Scalar > get_DgDp(int j, int l) const
Precondition: supports(OUT_ARG_DgDp,j,l)==true.
IntegralType getIntegralValue(ParameterList const ¶mList, std::string const ¶mName)
DirectionalFiniteDiffCalculatorTypes::EFDMethodType EFDMethodType
void calcVariations(const ModelEvaluator< Scalar > &model, const ModelEvaluatorBase::InArgs< Scalar > &basePoint, const ModelEvaluatorBase::InArgs< Scalar > &directions, const ModelEvaluatorBase::OutArgs< Scalar > &baseFunctionValues, const ModelEvaluatorBase::OutArgs< Scalar > &variations) const
Compute variations using directional finite differences..
int Np() const
Return the number of parameter subvectors p(l) supported (Np >= 0).
RCP< const ParameterList > getValidParameters() const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual RCP< const VectorSpaceBase< Scalar > > get_p_space(int l) const =0
Return the vector space for the auxiliary parameters p(l) <: RE^n_p_l.
void calcDerivatives(const ModelEvaluator< Scalar > &model, const ModelEvaluatorBase::InArgs< Scalar > &basePoint, const ModelEvaluatorBase::OutArgs< Scalar > &baseFunctionValues, const ModelEvaluatorBase::OutArgs< Scalar > &derivatives) const
Compute entire derivative objects using finite differences.
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
Base subclass for ModelEvaluator that defines some basic types.
virtual RCP< const VectorSpaceBase< Scalar > > get_f_space() const =0
Return the vector space for the state function f(...) <: RE^n_x.
DirectionalFiniteDiffCalculator(EFDMethodType fd_method_type=DirectionalFiniteDiffCalculatorTypes::FD_ORDER_FOUR_AUTO, EFDStepSelectType fd_step_select_type=DirectionalFiniteDiffCalculatorTypes::FD_STEP_ABSOLUTE, ScalarMag fd_step_size=-1.0, ScalarMag fd_step_size_min=-1.0)
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
TypeTo as(const TypeFrom &t)
RCP< const ParameterList > getParameterList() const
virtual RCP< const VectorSpaceBase< Scalar > > get_g_space(int j) const =0
Return the vector space for the auxiliary response functions g(j) <: RE^n_g_j.
ST::magnitudeType ScalarMag
RCP< ParameterList > unsetParameterList()
ModelEvaluatorBase()
constructor
RCP< ParameterList > getNonconstParameterList()
int Np() const
Return the number of parameter subvectors p(l) supported (Np >= 0).
RCP< const VectorBase< Scalar > > get_x_dot_dot() const
Precondition: supports(IN_ARG_x_dot_dot)==true.
bool supports(EInArgsMembers arg) const
Determines if an input argument is supported or not.
RCP< const VectorBase< Scalar > > get_x() const
Precondition: supports(IN_ARG_x)==true.
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
RCP< const VectorBase< Scalar > > get_p(int l) const
Get p(l) where 0 <= l && l < this->Np().
Concrete aggregate class for all input arguments computable by a ModelEvaluator subclass object...