10 #ifndef Tempus_InterpolatorLagrange_impl_hpp
11 #define Tempus_InterpolatorLagrange_impl_hpp
17 #include "Thyra_VectorStdOps.hpp"
21 template <
class Scalar>
25 int n = (*nodes_).size();
28 lagrange(order_, t, state_out);
30 lagrange(n - 1, t, state_out);
33 template <
class Scalar>
37 pl_ = Teuchos::parameterList();
38 if (paramList == Teuchos::null)
39 *pl_ = *(this->getValidParameters());
43 order_ = pl_->get<
int>(
"Order");
46 template <
class Scalar>
53 template <
class Scalar>
62 template <
class Scalar>
67 tmp->
set<std::string>(
"Interpolator Type",
"Lagrange");
68 tmp->
set<
int>(
"Order", 0);
72 template <
class Scalar>
78 using Thyra::V_StVpStV;
82 int n = nodes_->size();
85 const Scalar t_begin = (*nodes_)[0]->getTime();
86 const Scalar t_final = (*nodes_)[n - 1]->getTime();
93 else if (t >= t_final)
96 auto it = std::find_if(nodes_->begin(), nodes_->end(),
97 [=](
const RCP<SolutionState<Scalar> >& s) {
98 return t <= s->getTime();
100 i = std::distance(nodes_->begin(), it) - 1;
106 state_out->
copy((*nodes_)[i]);
111 state_out->
copy((*nodes_)[i + 1]);
118 if (k + p + 1 > n) k = n - p - 1;
121 RCP<VectorBase<Scalar> > x = state_out->
getX();
122 RCP<VectorBase<Scalar> > xdot = state_out->
getXDot();
123 RCP<VectorBase<Scalar> > xdotdot = state_out->
getXDotDot();
126 Thyra::assign(x.ptr(), zero);
127 if (!
is_null(xdot)) Thyra::assign(xdot.ptr(), zero);
128 if (!
is_null(xdotdot)) Thyra::assign(xdotdot.ptr(), zero);
130 for (
int j = 0; j <= p; ++j) {
131 const Scalar tj = (*nodes_)[k + j]->getTime();
132 RCP<const VectorBase<Scalar> > xj = (*nodes_)[k + j]->getX();
133 RCP<const VectorBase<Scalar> > xdotj = (*nodes_)[k + j]->getXDot();
134 RCP<const VectorBase<Scalar> > xdotdotj = (*nodes_)[k + j]->getXDotDot();
138 for (
int l = 0; l <= p; ++l) {
140 const Scalar tl = (*nodes_)[k + l]->getTime();
145 const Scalar a = num / den;
146 Thyra::Vp_StV(x.ptr(), a, *xj);
147 if (!
is_null(xdot)) Thyra::Vp_StV(xdot.ptr(), a, *xdotj);
148 if (!
is_null(xdotdot)) Thyra::Vp_StV(xdotdot.ptr(), a, *xdotdotj);
152 state_out->
getMetaData()->copy((*nodes_)[i]->getMetaData());
154 state_out->
getMetaData()->setDt(t - (*nodes_)[i]->getTime());
161 #endif // Tempus_InterpolatorLagrange_impl_hpp
bool is_null(const boost::shared_ptr< T > &p)
virtual Teuchos::RCP< Thyra::VectorBase< Scalar > > getXDot()
bool is_null(const std::shared_ptr< T > &p)
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
virtual Teuchos::RCP< Thyra::VectorBase< Scalar > > getXDotDot()
void validateParametersAndSetDefaults(ParameterList const &validParamList, int const depth=1000)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
void lagrange(const int p, const Scalar &t, SolutionState< Scalar > *state_out) const
virtual Teuchos::RCP< const SolutionStateMetaData< Scalar > > getMetaData() const
bool floating_compare_equals(const Scalar &a, const Scalar &b, const Scalar &scale)
Helper function for comparing times.
void interpolate(const Scalar &t, SolutionState< Scalar > *state_out) const
Perform an interpolation.
#define TEUCHOS_ASSERT(assertion_test)
virtual Teuchos::RCP< Thyra::VectorBase< Scalar > > getX()
Solution state for integrators and steppers.
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > ¶mList)
virtual void copy(const Teuchos::RCP< const SolutionState< Scalar > > &ss)
This is a deep copy.