Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_InterpolatorLagrange_decl.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_InterpolatorLagrange_decl_hpp
11#define Tempus_InterpolatorLagrange_decl_hpp
12
13#include "Tempus_config.hpp"
15
16namespace Tempus {
17
21template <class Scalar>
22class InterpolatorLagrange : virtual public Interpolator<Scalar> {
23 public:
26
29
31
32
35 const Teuchos::RCP<
36 const std::vector<Teuchos::RCP<SolutionState<Scalar> > > >& nodes)
37 {
38 nodes_ = nodes;
39 }
40
42 void interpolate(const Scalar& t, SolutionState<Scalar>* state_out) const;
43
45 int order() const { return order_; }
46
48
50
51 std::string description() const { return "Tempus::InterpolatorLagrange"; }
52 void describe(Teuchos::FancyOStream& out,
53 const Teuchos::EVerbosityLevel /* verbLevel */) const
54 {
55 out.setOutputToRootOnly(0);
56 out << description() << "::describe" << std::endl;
57 }
59
61
62 void setParameterList(const Teuchos::RCP<Teuchos::ParameterList>& paramList);
63 Teuchos::RCP<Teuchos::ParameterList> getNonconstParameterList();
64 Teuchos::RCP<Teuchos::ParameterList> unsetParameterList();
65 Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
67
68 private:
69 void lagrange(const int p, const Scalar& t,
70 SolutionState<Scalar>* state_out) const;
71
72 Teuchos::RCP<const std::vector<Teuchos::RCP<SolutionState<Scalar> > > >
74 Teuchos::RCP<Teuchos::ParameterList> pl_;
75 int order_;
76};
77
78// Nonmember constructor
79template <class Scalar>
80Teuchos::RCP<InterpolatorLagrange<Scalar> > lagrangeInterpolator()
81{
82 return Teuchos::rcp(new InterpolatorLagrange<Scalar>);
83}
84
85} // namespace Tempus
86
87#endif // Tempus_InterpolatorLagrange_decl_hpp
Concrete implemenation of Interpolator that does simple lagrange interpolation.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel) const
int order() const
Return the order of the interpolation.
Teuchos::RCP< Teuchos::ParameterList > pl_
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
Teuchos::RCP< const std::vector< Teuchos::RCP< SolutionState< Scalar > > > > nodes_
void setNodes(const Teuchos::RCP< const std::vector< Teuchos::RCP< SolutionState< Scalar > > > > &nodes)
Store pointer to interpolation nodes.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &paramList)
void lagrange(const int p, const Scalar &t, SolutionState< Scalar > *state_out) const
void interpolate(const Scalar &t, SolutionState< Scalar > *state_out) const
Perform an interpolation.
Base strategy class for interpolation functionality.
Solution state for integrators and steppers.
Teuchos::RCP< InterpolatorLagrange< Scalar > > lagrangeInterpolator()