Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_Interpolator.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ****************************************************************************
3 // Tempus: Copyright (2017) Sandia Corporation
4 //
5 // Distributed under BSD 3-clause license (See accompanying file Copyright.txt)
6 // ****************************************************************************
7 // @HEADER
8 
9 #ifndef Tempus_Interpolator_hpp
10 #define Tempus_Interpolator_hpp
11 
12 // Teuchos
13 #include "Teuchos_VerboseObject.hpp"
14 #include "Teuchos_Describable.hpp"
15 #include "Teuchos_ParameterListAcceptorDefaultBase.hpp"
16 #include "Teuchos_RCP.hpp"
17 #include "Teuchos_ScalarTraits.hpp"
18 
19 // Tempus
20 #include "Tempus_config.hpp"
21 #include "Tempus_SolutionState.hpp"
22 
23 #include <vector>
24 
25 namespace Tempus {
26 
29 template<class Scalar>
31  : virtual public Teuchos::Describable,
32  virtual public Teuchos::ParameterListAcceptor,
33  virtual public Teuchos::VerboseObject<Interpolator<Scalar> >
34 {
35 public:
36 
55  virtual void setNodes(
56  const Teuchos::RCP<const std::vector<Teuchos::RCP<SolutionState<Scalar> > > > & nodes) = 0;
57 
60  virtual void interpolate(const Scalar& t,
61  SolutionState<Scalar>* state_out) const = 0;
62 
65  virtual int order() const = 0;
66 
67 };
68 
70 template<class Scalar>
71 void interpolate(const Interpolator<Scalar>& interpolator,
72  const Scalar& t,
73  SolutionState<Scalar>* state_out)
74 {
75  interpolator.interpolate(t, state_out);
76 }
77 
79 template<class Scalar>
80 void interpolate(Interpolator<Scalar>& interpolator,
81  const Teuchos::RCP<const std::vector<Teuchos::RCP<SolutionState<Scalar> > > >& nodes,
82  const Scalar& t,
83  SolutionState<Scalar>* state_out)
84 {
85  interpolator.setNodes(nodes);
86  interpolator.interpolate(t, state_out);
87 }
88 
90 template <typename Scalar>
91 bool floating_compare_equals(const Scalar& a, const Scalar& b,
92  const Scalar& scale)
93 {
95  typedef typename ST::magnitudeType mag_type;
96 
97  const mag_type tol = ST::magnitude(scale)*ST::eps()*mag_type(1000.0);
98  const mag_type diff = ST::magnitude(a-b);
99 
100  return diff <= tol || diff <= ST::sfmin();
101 }
102 
103 } // namespace Tempus
104 
105 #endif //Tempus_Interpolator_hpp
virtual void setNodes(const Teuchos::RCP< const std::vector< Teuchos::RCP< SolutionState< Scalar > > > > &nodes)=0
Store pointer to interpolation nodes.
Base strategy class for interpolation functionality.
void interpolate(const Interpolator< Scalar > &interpolator, const Scalar &t, SolutionState< Scalar > *state_out)
Nonmember functions.
virtual void interpolate(const Scalar &t, SolutionState< Scalar > *state_out) const =0
Perform an interpolation.
bool floating_compare_equals(const Scalar &a, const Scalar &b, const Scalar &scale)
Helper function for comparing times.
virtual int order() const =0
Return the order of the interpolation.
Solution state for integrators and steppers. SolutionState contains the metadata for solutions and th...