9 #ifndef TEMPUS_TEST_DAHLQUIST_TEST_MODEL_IMPL_HPP
10 #define TEMPUS_TEST_DAHLQUIST_TEST_MODEL_IMPL_HPP
12 #include "Thyra_DefaultSpmdVectorSpace.hpp"
13 #include "Thyra_DetachedVectorView.hpp"
14 #include "Thyra_DetachedMultiVectorView.hpp"
15 #include "Thyra_DefaultSerialDenseLinearOpWithSolveFactory.hpp"
16 #include "Thyra_DefaultMultiVectorLinearOpWithSolve.hpp"
17 #include "Thyra_DefaultLinearOpSource.hpp"
18 #include "Thyra_VectorStdOps.hpp"
19 #include "Thyra_MultiVectorStdOps.hpp"
20 #include "Thyra_DefaultMultiVectorProductVector.hpp"
25 namespace Tempus_Test {
27 template<
class Scalar>
29 { constructDahlquistTestModel(-1.0,
false); }
32 template<
class Scalar>
35 { constructDahlquistTestModel(lambda, includeXDot); }
38 template<
class Scalar>
44 includeXDot_ = includeXDot;
45 isInitialized_ =
false;
47 xDotIC_ = Scalar(lambda);
54 acceptModelParams_ =
false;
57 x_space_ = Thyra::defaultSpmdVectorSpace<Scalar>(dim_);
58 f_space_ = Thyra::defaultSpmdVectorSpace<Scalar>(dim_);
64 MEB::InArgsSetup<Scalar> inArgs;
65 inArgs.setModelEvalDescription(this->description());
66 inArgs.setSupports( MEB::IN_ARG_t );
67 inArgs.setSupports( MEB::IN_ARG_x );
68 inArgs.setSupports( MEB::IN_ARG_x_dot );
70 inArgs.setSupports( MEB::IN_ARG_beta );
71 inArgs.setSupports( MEB::IN_ARG_alpha );
72 if (acceptModelParams_) inArgs.set_Np( Np_ );
79 MEB::OutArgsSetup<Scalar> outArgs;
80 outArgs.setModelEvalDescription(this->description());
81 outArgs.setSupports( MEB::OUT_ARG_f );
83 outArgs.setSupports( MEB::OUT_ARG_W_op );
84 if (acceptModelParams_) {
85 outArgs.set_Np_Ng(Np_,Ng_);
86 outArgs.setSupports( MEB::OUT_ARG_DfDp,0,
87 MEB::DERIV_MV_JACOBIAN_FORM );
88 outArgs.setSupports( MEB::OUT_ARG_DgDp,0,0,
89 MEB::DERIV_MV_JACOBIAN_FORM );
90 outArgs.setSupports( MEB::OUT_ARG_DgDx,0,
91 MEB::DERIV_MV_GRADIENT_FORM );
97 nominalValues_ = inArgs_;
99 nominalValues_.set_t(Scalar(0.0));
105 nominalValues_.set_x(x_ic);
112 x_dot_ic_view[0] = xDotIC_;
114 nominalValues_.set_x_dot(x_dot_ic);
117 isInitialized_ =
true;
121 template<
class Scalar>
128 inArgs.
set_t(exact_t);
134 exact_x_view[0] = exp(lambda_*exact_t);
136 inArgs.
set_x(exact_x);
143 exact_x_dot_view[0] = lambda_ * exp(lambda_*exact_t);
152 template<
class Scalar>
161 template<
class Scalar>
170 template<
class Scalar>
175 return nominalValues_;
179 template<
class Scalar>
188 template<
class Scalar>
197 template<
class Scalar>
208 using Teuchos::rcp_dynamic_cast;
210 "Error, setupInOutArgs_ must be called first!\n");
222 f_out_view[0] = lambda_*x_in_view[0];
225 "Error -- Dahlquist Test Model requires f_out!\n");
232 matrix_view(0,0) = lambda_;
239 x_dot_in = inArgs.
get_x_dot().assert_not_null();
246 f_out_view[0] = x_dot_in_view[0] - lambda_*x_in_view[0];
253 matrix_view(0,0) = alpha - beta*lambda_;
261 template<
class Scalar>
275 vec_view[0] = lambda_;
277 V_V(multivec->col(0).
ptr(),*vec);
281 Thyra::linearOpWithSolve<Scalar>(
288 template<
class Scalar>
299 template<
class Scalar>
305 Thyra::defaultSerialDenseLinearOpWithSolveFactory<Scalar>();
309 template<
class Scalar>
314 if (!acceptModelParams_) {
315 return Teuchos::null;
320 else if (l == 1 || l == 2)
322 return Teuchos::null;
326 template<
class Scalar>
331 if (!acceptModelParams_) {
332 return Teuchos::null;
338 p_strings->push_back(
"Model Coefficient: a");
343 p_strings->push_back(
"DxDp");
345 p_strings->push_back(
"Dx_dotDp");
349 template<
class Scalar>
359 #endif // TEMPUS_TEST_DAHLQUIST_TEST_MODEL_IMPL_HPP
Thyra::ModelEvaluatorBase::InArgs< Scalar > getNominalValues() const
void evalModelImpl(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs_bar, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs_bar) const
bool is_null(const boost::shared_ptr< T > &p)
RCP< const VectorBase< Scalar > > get_x_dot() const
Teuchos::RCP< Thyra::LinearOpBase< Scalar > > create_W_op() const
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Evaluation< VectorBase< Scalar > > get_f() const
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_x_space() const
#define TEUCHOS_ASSERT_IN_RANGE_UPPER_EXCLUSIVE(index, lower_inclusive, upper_exclusive)
Thyra::ModelEvaluatorBase::OutArgs< Scalar > createOutArgsImpl() const
Thyra::ModelEvaluatorBase::InArgs< Scalar > createInArgs() const
void set_x(const RCP< const VectorBase< Scalar > > &x)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Thyra::ModelEvaluatorBase::InArgs< Scalar > getExactSolution(double t) const
Exact solution.
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_f_space() const
Teuchos::RCP< const Thyra::LinearOpWithSolveFactoryBase< Scalar > > get_W_factory() const
void set_x_dot(const RCP< const VectorBase< Scalar > > &x_dot)
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_g_space(int j) const
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_p_space(int l) const
Teuchos::RCP< Thyra::LinearOpWithSolveBase< Scalar > > create_W() const
void constructDahlquistTestModel(Scalar lambda, bool includeXDot)
RCP< LinearOpBase< Scalar > > get_W_op() const
RCP< const VectorBase< Scalar > > get_x() const