42#ifndef THYRA_MODEL_EVALUATOR_HELPERS_HPP
43#define THYRA_MODEL_EVALUATOR_HELPERS_HPP
46#include "Thyra_ModelEvaluator.hpp"
103template<
class Scalar>
114template<
class Scalar>
118 ,
const std::string &derivName
123template<
class Scalar>
127 ,
const std::string &derivName
137template<
class Scalar>
139 const std::string &modelEvalDescription,
141 const std::string &deriv_name,
143 const std::string &fnc_space_name,
145 const std::string &var_space_name
153template<
class Scalar>
155 const std::string &modelEvalDescription,
165template<
class Scalar>
176template<
class Scalar>
185template<
class Scalar>
194template<
class Scalar>
204template<
class Scalar>
214template<
class Scalar>
225template<
class Scalar>
237template<
class Scalar>
249template<
class Scalar>
261template<
class Scalar>
274#ifdef HAVE_THYRA_ME_POLYNOMIAL
278template<
class Scalar>
288template<
class Scalar>
309#include "Thyra_AssertOp.hpp"
313template<
class Scalar>
324template<
class Scalar>
326Thyra::derivativeGradient(
327 const RCP<MultiVectorBase<Scalar> > &grad
330 return ModelEvaluatorBase::Derivative<Scalar>(
332 ModelEvaluatorBase::DERIV_MV_GRADIENT_FORM
337template<
class Scalar>
339Thyra::create_DfDp_mv(
340 const ModelEvaluator<Scalar>& model,
342 ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation
346 return createMembers( model.get_f_space(), model.get_p_space(l)->dim() );
350template<
class Scalar>
352Thyra::create_DgDx_dot_mv(
353 const ModelEvaluator<Scalar>& model,
355 ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation
358 typedef ModelEvaluatorBase MEB;
359 switch(orientation) {
360 case MEB::DERIV_MV_BY_COL:
362 MEB::DerivativeMultiVector<Scalar>(
363 createMembers( model.get_g_space(j), model.get_x_space()->dim() )
364 ,MEB::DERIV_MV_BY_COL
366 case MEB::DERIV_TRANS_MV_BY_ROW:
368 MEB::DerivativeMultiVector<Scalar>(
369 createMembers( model.get_x_space(), model.get_g_space(j)->dim() )
370 ,MEB::DERIV_TRANS_MV_BY_ROW
379template<
class Scalar>
381Thyra::create_DgDx_mv(
382 const ModelEvaluator<Scalar>& model,
384 ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation
387 return create_DgDx_dot_mv(model,j,orientation);
391template<
class Scalar>
393Thyra::create_DgDp_mv(
394 const ModelEvaluator<Scalar>& model,
397 ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation
400 typedef ModelEvaluatorBase MEB;
401 switch(orientation) {
402 case MEB::DERIV_MV_BY_COL:
404 MEB::DerivativeMultiVector<Scalar>(
405 createMembers( model.get_g_space(j), model.get_p_space(l)->dim() )
406 ,MEB::DERIV_MV_BY_COL
408 case MEB::DERIV_TRANS_MV_BY_ROW:
410 MEB::DerivativeMultiVector<Scalar>(
411 createMembers( model.get_p_space(l), model.get_g_space(j)->dim() )
412 ,MEB::DERIV_TRANS_MV_BY_ROW
421template<
class Scalar>
424 const ModelEvaluatorBase::Derivative<Scalar> &deriv
425 ,
const std::string &derivName
429 deriv.getLinearOp().get()!=NULL, std::logic_error
430 ,
"Error, LinearOpBase type not expected for " << derivName <<
"!"
432 return deriv.getDerivativeMultiVector();
436template<
class Scalar>
439 const ModelEvaluatorBase::Derivative<Scalar> &deriv
440 ,
const std::string &derivName
441 ,ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation
444 typedef ModelEvaluatorBase MEB;
446 deriv.getLinearOp().get()!=NULL, std::logic_error
447 ,
"Error, LinearOpBase type not expected for " << derivName <<
"!"
449 MEB::DerivativeMultiVector<Scalar>
450 dmv = deriv.getDerivativeMultiVector();
451 RCP<MultiVectorBase<Scalar> >
452 mv = dmv.getMultiVector();
455 dmv.getOrientation() != orientation, std::logic_error
456 ,
"Error, the orientation " <<
toString(dmv.getOrientation()) <<
" is not the"
457 " expected orientation of " <<
toString(orientation)
458 <<
" for " << derivName <<
"!"
465template<
class Scalar>
467 const std::string &modelEvalDescription,
469 const std::string &deriv_name,
471 const std::string &fnc_space_name,
473 const std::string &var_space_name
479 if (!is_null(lo->range())) {
481 modelEvalDescription,
482 *lo->range(), deriv_name +
".range()",
483 fnc_space, fnc_space_name
486 modelEvalDescription,
487 *lo->domain(), deriv_name +
".domain()",
488 var_space, var_space_name
495 case MEB::DERIV_MV_BY_COL: {
497 modelEvalDescription,
498 *mv->range(), deriv_name +
".range()",
499 fnc_space, fnc_space_name
502 modelEvalDescription,
503 *mv->domain(), deriv_name +
".domain()",
504 var_space, var_space_name
508 case MEB::DERIV_TRANS_MV_BY_ROW: {
510 modelEvalDescription,
511 *mv->range(), deriv_name +
"^T.range()",
512 var_space, var_space_name
515 modelEvalDescription,
516 *mv->domain(), deriv_name +
"^T.domain()",
517 fnc_space, fnc_space_name
530template<
class Scalar>
532 const std::string &modelEvalDescription,
540 const int Ng = outArgs.
Ng();
541 const int Np = outArgs.
Np();
549 inArgs.
Np() != outArgs.
Np(), std::logic_error,
550 "Error: The underlying model " << modelEvalDescription <<
" incorrectly\n"
551 "set inArgs.Np() = "<<inArgs.
Np()<<
" != outArgs.Np() = "
559 "Error: The underlying model " << modelEvalDescription <<
" supports\n"
560 "x_dot but does not support x!"
567 "Error: The underlying model " << modelEvalDescription <<
" supports\n"
568 "x_dot but does not support t!"
579 "Error: The underlying model " << modelEvalDescription <<
" says that\n"
580 "it supports W and/or W_op but it does not support x!"
588 !( inArgs.
supports(MEB::IN_ARG_alpha) && inArgs.
supports(MEB::IN_ARG_beta) )
591 "Error: The underlying model " << modelEvalDescription <<
" supports W and/or W_op\n"
592 "and x_dot but it does not support alpha and beta as InArgs! \n"
593 "If the model supports W and x_dot, then it can be interpreted as an implicit \n"
594 "ODE/DAE and therefore D(f)/D(x_dot) must be non-zero and therefore alpha must \n"
595 "be supported. If, however, the model can be interpreted as an explicit ODE \n"
596 "then x_dot should not be supported at all."
599 for (
int l = 0; l <
Np; ++l ) {
603 for (
int j = 0; j <
Ng; ++j ) {
607 ( !outArgs.
supports(MEB::OUT_ARG_DgDx_dot,j).none()
608 && !inArgs.
supports(MEB::IN_ARG_x_dot) ),
610 "Error: The underlying model " << modelEvalDescription <<
" says that\n"
611 "it supports DgDx_dot("<<j<<
") but it does not support x_dot!"
616 ( !outArgs.
supports(MEB::OUT_ARG_DgDx,j).none()
617 && !inArgs.
supports(MEB::IN_ARG_x) ),
619 "Error: The underlying model " << modelEvalDescription <<
" says that\n"
620 "it supports DgDx("<<j<<
") but it does not support x!"
632template<
class Scalar>
642 const int Np = inArgs.
Np();
653 if ( inArgs.
supports(MEB::IN_ARG_x) && !is_null(inArgs.
get_x()) ) {
659 for (
int l = 0; l <
Np; ++l ) {
660 if (!is_null(inArgs.
get_p(l))) {
669template<
class Scalar>
682 const int Ng = outArgs.
Ng();
683 const int Np = outArgs.
Np();
692 if ( outArgs.
supports(MEB::OUT_ARG_f) && !is_null(outArgs.
get_f()) ) {
698 if ( outArgs.
supports(MEB::OUT_ARG_W) && !is_null(outArgs.
get_W()) ) {
699 if (!is_null(outArgs.
get_W()->range())) {
708 if ( outArgs.
supports(MEB::OUT_ARG_W_op) && !is_null(outArgs.
get_W_op()) ) {
709 if (!is_null(outArgs.
get_W_op()->range())) {
723 ( outArgs.
supports(MEB::OUT_ARG_W) && !is_null(outArgs.
get_W()) )
729 if ( inArgs->
supports(MEB::IN_ARG_alpha) && inArgs->
supports(MEB::IN_ARG_beta) ) {
737 else if ( inArgs->
supports(MEB::IN_ARG_beta) ) {
743 if (outArgs.
supports(MEB::OUT_ARG_f)) {
744 for (
int l = 0; l <
Np; ++l ) {
745 if (!outArgs.
supports(MEB::OUT_ARG_DfDp,l).none()) {
748 outArgs.
get_DfDp(l),
"DfDp("+TU::toString(l)+
")",
750 *model.
get_p_space(l),
"p_space("+TU::toString(l)+
")"
757 for (
int j = 0; j <
Ng; ++j ) {
758 if (!is_null(outArgs.
get_g(j))) {
765 for (
int j = 0; j <
Ng; ++j ) {
766 if (!outArgs.
supports(MEB::OUT_ARG_DgDx_dot,j).none()) {
769 outArgs.
get_DgDx_dot(j),
"DgDx_dot("+TU::toString(j)+
")",
770 *model.
get_g_space(j),
"g_space("+TU::toString(j)+
")",
777 for (
int j = 0; j <
Ng; ++j ) {
778 if (!outArgs.
supports(MEB::OUT_ARG_DgDx,j).none()) {
781 outArgs.
get_DgDx(j),
"DgDx("+TU::toString(j)+
")",
782 *model.
get_g_space(j),
"g_space("+TU::toString(j)+
")",
789 for (
int j = 0; j <
Ng; ++j ) {
790 for (
int l = 0; l <
Np; ++l ) {
791 if (!outArgs.
supports(MEB::OUT_ARG_DgDp,j,l).none()) {
792 const std::string j_str = TU::toString(j);
793 const std::string l_str = TU::toString(l);
796 outArgs.
get_DgDp(j,l),
"DgDp("+j_str+
","+l_str+
")",
807template<
class Scalar>
823template<
class Scalar>
825 const ModelEvaluator<Scalar> &model
826 ,
const VectorBase<Scalar> &x
827 ,VectorBase<Scalar> *f
828 ,LinearOpWithSolveBase<Scalar> *W
832 typedef ModelEvaluatorBase MEB;
834 MEB::InArgs<Scalar> inArgs = model.createInArgs();
835 MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
842 model.evalModel(inArgs,outArgs);
847template<
class Scalar>
849 const ModelEvaluator<Scalar> &model
850 ,
const VectorBase<Scalar> &x
852 ,VectorBase<Scalar> *f
855 typedef ModelEvaluatorBase MEB;
856 MEB::InArgs<Scalar> inArgs = model.createInArgs();
857 MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
859 if(inArgs.supports(MEB::IN_ARG_t)) inArgs.set_t(t);
861 model.evalModel(inArgs,outArgs);
865template<
class Scalar>
867 const ModelEvaluator<Scalar> &model,
869 const VectorBase<Scalar> &p_l,
871 const Ptr<VectorBase<Scalar> > &g_j
874 typedef ModelEvaluatorBase MEB;
875 MEB::InArgs<Scalar> inArgs = model.createInArgs();
876 MEB::OutArgs<Scalar> outArgs= model.createOutArgs();
877 inArgs.set_p(l, Teuchos::rcpFromRef(p_l));
878 outArgs.set_g(j, Teuchos::rcpFromRef(*g_j));
879 model.evalModel(inArgs,outArgs);
883template<
class Scalar>
885 const ModelEvaluator<Scalar> &model,
887 const VectorBase<Scalar> &p_l,
890 VectorBase<Scalar> *g_j
893 typedef ModelEvaluatorBase MEB;
894 MEB::InArgs<Scalar> inArgs = model.createInArgs();
895 MEB::OutArgs<Scalar> outArgs= model.createOutArgs();
899 model.evalModel(inArgs,outArgs);
903template<
class Scalar>
904void Thyra::eval_g_DgDp(
905 const ModelEvaluator<Scalar> &model,
907 const VectorBase<Scalar> &p_l,
909 const Ptr<VectorBase<Scalar> > &g_j,
910 const ModelEvaluatorBase::Derivative<Scalar> &DgDp_j_l
913 typedef ModelEvaluatorBase MEB;
914 MEB::InArgs<Scalar> inArgs = model.createInArgs();
915 MEB::OutArgs<Scalar> outArgs= model.createOutArgs();
916 inArgs.set_p(l, Teuchos::rcpFromRef(p_l));
918 outArgs.set_g(j, Teuchos::rcpFromPtr(g_j));
920 if (!DgDp_j_l.isEmpty()) {
921 outArgs.set_DgDp(j, l, DgDp_j_l);
923 model.evalModel(inArgs,outArgs);
927template<
class Scalar>
929 const ModelEvaluator<Scalar> &model
930 ,
const VectorBase<Scalar> &x_dot
931 ,
const VectorBase<Scalar> &x
932 ,
const typename ModelEvaluatorBase::InArgs<Scalar>::ScalarMag &t
933 ,VectorBase<Scalar> *f
937 typedef ModelEvaluatorBase MEB;
939 MEB::InArgs<Scalar> inArgs = model.createInArgs();
940 MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
944 if(inArgs.supports(MEB::IN_ARG_t))
949 model.evalModel(inArgs,outArgs);
954template<
class Scalar>
956 const ModelEvaluator<Scalar> &model
957 ,
const VectorBase<Scalar> &x_dot
958 ,
const VectorBase<Scalar> &x
959 ,
const typename ModelEvaluatorBase::InArgs<Scalar>::ScalarMag &t
962 ,VectorBase<Scalar> *f
963 ,LinearOpWithSolveBase<Scalar> *W
967 typedef ModelEvaluatorBase MEB;
969 MEB::InArgs<Scalar> inArgs = model.createInArgs();
970 MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
974 if(inArgs.supports(MEB::IN_ARG_t))
976 inArgs.set_alpha(alpha);
977 inArgs.set_beta(beta);
982 model.evalModel(inArgs,outArgs);
987#ifdef HAVE_THYRA_ME_POLYNOMIAL
990template<
class Scalar>
991void Thyra::eval_f_poly(
992 const ModelEvaluator<Scalar> &model
994 ,
const typename ModelEvaluatorBase::InArgs<Scalar>::ScalarMag &t
999 typedef ModelEvaluatorBase MEB;
1001 MEB::InArgs<Scalar> inArgs = model.createInArgs();
1002 MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
1005 if(inArgs.supports(MEB::IN_ARG_t))
1010 model.evalModel(inArgs,outArgs);
1015template<
class Scalar>
1016void Thyra::eval_f_poly(
1017 const ModelEvaluator<Scalar> &model
1019 ,
const VectorBase<Scalar> &x_poly
1020 ,
const typename ModelEvaluatorBase::InArgs<Scalar>::ScalarMag &t
1025 typedef ModelEvaluatorBase MEB;
1027 MEB::InArgs<Scalar> inArgs = model.createInArgs();
1028 MEB::OutArgs<Scalar> outArgs = model.createOutArgs();
1030 inArgs.set_x_dot_poly(
Teuchos::rcp(&x_dot_poly,
false));
1032 if(inArgs.supports(MEB::IN_ARG_t))
1037 model.evalModel(inArgs,outArgs);
virtual std::string description() const
Base class for all linear operators that can support a high-level solve operation.
Simple aggregate class for a derivative object represented as a column-wise multi-vector or its trans...
Simple aggregate class that stores a derivative object as a general linear operator or as a multi-vec...
RCP< MultiVectorBase< Scalar > > getMultiVector() const
RCP< LinearOpBase< Scalar > > getLinearOp() const
EDerivativeMultiVectorOrientation getMultiVectorOrientation() const
Concrete aggregate class for all input arguments computable by a ModelEvaluator subclass object.
std::string modelEvalDescription() const
RCP< const VectorBase< Scalar > > get_p(int l) const
Get p(l) where 0 <= l && l < this->Np().
Scalar get_beta() const
Precondition: supports(IN_ARG_beta)==true.
RCP< const VectorBase< Scalar > > get_x() const
Precondition: supports(IN_ARG_x)==true.
Teuchos::ScalarTraits< Scalar >::magnitudeType ScalarMag
RCP< const VectorBase< Scalar > > get_x_dot() const
Precondition: supports(IN_ARG_x_dot)==true.
bool supports(EInArgsMembers arg) const
Determines if an input argument is supported or not.
int Np() const
Return the number of parameter subvectors p(l) supported (Np >= 0).
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< LinearOpWithSolveBase< Scalar > > get_W() const
Precondition: supports(OUT_ARG_W)==true.
int Ng() const
Return the number of axillary response functions g(j)(...) supported (Ng >= 0).
std::string modelEvalDescription() const
Derivative< Scalar > get_DgDp(int j, int l) const
Precondition: supports(OUT_ARG_DgDp,j,l)==true.
Evaluation< VectorBase< Scalar > > get_f() const
Precondition: supports(OUT_ARG_f)==true.
RCP< LinearOpBase< Scalar > > get_W_op() const
Precondition: supports(OUT_ARG_W_op)==true.
Derivative< Scalar > get_DgDx_dot(int j) const
Precondition: supports(OUT_ARG_DgDx_dot,j)==true.
Derivative< Scalar > get_DgDx(int j) const
Precondition: supports(OUT_ARG_DgDx,j)==true.
bool supports(EOutArgsMembers arg) const
Determine if an input argument is supported or not.
int Np() const
Return the number of parameter subvectors p(l) supported (Np >= 0).
Evaluation< VectorBase< Scalar > > get_g(int j) const
Precondition: supports(OUT_ARG_g)==true..
Base subclass for ModelEvaluator that defines some basic types.
EDerivativeMultiVectorOrientation
void assertInArgsEvalObjects(const ModelEvaluator< Scalar > &model, const ModelEvaluatorBase::InArgs< Scalar > &inArgs)
Assert that the objects in an InArgs object match a given model.
void assertInArgsOutArgsSetup(const std::string &modelEvalDescription, const ModelEvaluatorBase::InArgs< Scalar > &inArgs, const ModelEvaluatorBase::OutArgs< Scalar > &outArgs)
Assert that an InArgs and OutArgs object are setup consistently.
void assertOutArgsEvalObjects(const ModelEvaluator< Scalar > &model, const ModelEvaluatorBase::OutArgs< Scalar > &outArgs, const ModelEvaluatorBase::InArgs< Scalar > *inArgs=0)
Assert that the objects in an OutArgs object match a given model.
void assertDerivSpaces(const std::string &modelEvalDescription, const ModelEvaluatorBase::Derivative< Scalar > &deriv, const std::string &deriv_name, const VectorSpaceBase< Scalar > &fnc_space, const std::string &fnc_space_name, const VectorSpaceBase< Scalar > &var_space, const std::string &var_space_name)
Assert that that Thyra objects imbedded in a Derivative object matches its function and variable spac...
RCP< ModelEvaluatorBase::InArgs< Scalar > > clone(const ModelEvaluatorBase::InArgs< Scalar > &inArgs)
Create a clone of an InArgs object.
Pure abstract base interface for evaluating a stateless "model" that can be mapped into a number of d...
virtual RCP< const VectorSpaceBase< Scalar > > get_f_space() const =0
Return the vector space for the state function f(...) <: RE^n_x.
virtual ModelEvaluatorBase::OutArgs< Scalar > createOutArgs() const =0
Create an empty output functions/derivatives object that can be set up and passed to evalModel().
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.
virtual ModelEvaluatorBase::InArgs< Scalar > createInArgs() const =0
Create an empty input arguments object that can be set up and passed to evalModel().
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.
virtual RCP< const VectorSpaceBase< Scalar > > get_x_space() const =0
Return the vector space for the state variables x <: RE^n_x.
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.
Interface for a collection of column vectors called a multi-vector.
Abstract interface for finite-dimensional dense vectors.
Abstract interface for objects that represent a space for vectors.
RCP< MultiVectorBase< Scalar > > createMembers(const RCP< const VectorSpaceBase< Scalar > > &vs, int numMembers, const std::string &label="")
Create a set of vector members (a MultiVectorBase) from the vector space.
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
#define TEUCHOS_ASSERT_EQUALITY(val1, val2)
bool is_null(const std::shared_ptr< T > &p)
#define THYRA_ASSERT_VEC_SPACES(FUNC_NAME, VS1, VS2)
This is a very useful macro that should be used to validate that two vector spaces are compatible.
#define THYRA_ASSERT_VEC_SPACES_NAMES(FUNC_NAME, VS1, VS1_NAME, VS2, VS2_NAME)
Helper assertion macro.
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
std::string toString(const T &t)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)