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: Time Integration and Sensitivity Analysis Package
4 //
5 // Copyright 2017 NTESS and the Tempus contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 //@HEADER
9 
10 #ifndef Tempus_Interpolator_hpp
11 #define Tempus_Interpolator_hpp
12 
13 // Teuchos
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"
19 
20 // Tempus
21 #include "Tempus_config.hpp"
22 #include "Tempus_SolutionState.hpp"
23 
24 #include <vector>
25 
26 namespace Tempus {
27 
30 template <class Scalar>
32  : virtual public Teuchos::Describable,
33  virtual public Teuchos::ParameterListAcceptor,
34  virtual public Teuchos::VerboseObject<Interpolator<Scalar> > {
35  public:
55  virtual void setNodes(
56  const Teuchos::RCP<
57  const std::vector<Teuchos::RCP<SolutionState<Scalar> > > >&
58  nodes) = 0;
59 
62  virtual void interpolate(const Scalar& t,
63  SolutionState<Scalar>* state_out) const = 0;
64 
67  virtual int order() const = 0;
68 };
69 
71 template <class Scalar>
72 void interpolate(const Interpolator<Scalar>& interpolator, const Scalar& t,
73  SolutionState<Scalar>* state_out)
74 {
75  interpolator.interpolate(t, state_out);
76 }
77 
79 template <class Scalar>
81  Interpolator<Scalar>& interpolator,
82  const Teuchos::RCP<
83  const std::vector<Teuchos::RCP<SolutionState<Scalar> > > >& nodes,
84  const Scalar& t, SolutionState<Scalar>* state_out)
85 {
86  interpolator.setNodes(nodes);
87  interpolator.interpolate(t, state_out);
88 }
89 
91 template <typename Scalar>
92 bool floating_compare_equals(const Scalar& a, const Scalar& b,
93  const Scalar& scale)
94 {
96  typedef typename ST::magnitudeType mag_type;
97 
98  const mag_type tol = ST::magnitude(scale) * ST::eps() * mag_type(1000.0);
99  const mag_type diff = ST::magnitude(a - b);
100 
101  return diff <= tol || diff <= ST::sfmin();
102 }
103 
104 } // namespace Tempus
105 
106 #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.