10#ifndef Tempus_Interpolator_hpp
11#define Tempus_Interpolator_hpp
14#include "Teuchos_VerboseObject.hpp"
15#include "Teuchos_Describable.hpp"
16#include "Teuchos_ParameterListAcceptorDefaultBase.hpp"
17#include "Teuchos_RCP.hpp"
18#include "Teuchos_ScalarTraits.hpp"
21#include "Tempus_config.hpp"
22#include "Tempus_SolutionState.hpp"
30template <
class Scalar>
32 :
virtual public Teuchos::Describable,
33 virtual public Teuchos::ParameterListAcceptor,
34 virtual public Teuchos::VerboseObject<Interpolator<Scalar> > {
71template <
class Scalar>
79template <
class Scalar>
91template <
typename Scalar>
95 typedef Teuchos::ScalarTraits<Scalar> ST;
96 typedef typename ST::magnitudeType mag_type;
98 const mag_type tol = ST::magnitude(scale) * ST::eps() * mag_type(1000.0);
99 const mag_type diff = ST::magnitude(a - b);
101 return diff <= tol || diff <= ST::sfmin();
Base strategy class for interpolation functionality.
virtual void setNodes(const Teuchos::RCP< const std::vector< Teuchos::RCP< SolutionState< Scalar > > > > &nodes)=0
Store pointer to interpolation nodes.
virtual void interpolate(const Scalar &t, SolutionState< Scalar > *state_out) const =0
Perform an interpolation.
virtual int order() const =0
Return the order of the interpolation.
Solution state for integrators and steppers.
void interpolate(const Interpolator< Scalar > &interpolator, const Scalar &t, SolutionState< Scalar > *state_out)
Nonmember functions.
bool floating_compare_equals(const Scalar &a, const Scalar &b, const Scalar &scale)
Helper function for comparing times.