Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_InterpolatorLagrange_impl.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_impl_hpp
11#define Tempus_InterpolatorLagrange_impl_hpp
12
13#include <algorithm>
14#include <iterator>
15
16#include "Teuchos_ScalarTraits.hpp"
17#include "Thyra_VectorStdOps.hpp"
18
19namespace Tempus {
20
21template <class Scalar>
23 const Scalar& t, SolutionState<Scalar>* state_out) const
24{
25 int n = (*nodes_).size();
26 TEUCHOS_ASSERT(n > 0);
27 if (n >= order_ + 1)
28 lagrange(order_, t, state_out);
29 else
30 lagrange(n - 1, t, state_out);
31}
32
33template <class Scalar>
35 const Teuchos::RCP<Teuchos::ParameterList>& paramList)
36{
37 pl_ = Teuchos::parameterList();
38 if (paramList == Teuchos::null)
39 *pl_ = *(this->getValidParameters());
40 else
41 *pl_ = *paramList;
42 pl_->validateParametersAndSetDefaults(*this->getValidParameters());
43 order_ = pl_->get<int>("Order");
44}
45
46template <class Scalar>
47Teuchos::RCP<Teuchos::ParameterList>
52
53template <class Scalar>
54Teuchos::RCP<Teuchos::ParameterList>
56{
57 Teuchos::RCP<Teuchos::ParameterList> tmp = pl_;
58 pl_ = Teuchos::null;
59 return tmp;
60}
61
62template <class Scalar>
63Teuchos::RCP<const Teuchos::ParameterList>
65{
66 Teuchos::RCP<Teuchos::ParameterList> tmp = Teuchos::parameterList();
67 tmp->set<std::string>("Interpolator Type", "Lagrange");
68 tmp->set<int>("Order", 0);
69 return tmp;
70}
71
72template <class Scalar>
74 const int p, const Scalar& t, SolutionState<Scalar>* state_out) const
75{
76 using Teuchos::is_null;
77 using Teuchos::RCP;
78 using Thyra::V_StVpStV;
80
81 // Here we assume we have at least p nodes
82 int n = nodes_->size();
83 TEUCHOS_ASSERT(n >= p);
84
85 const Scalar t_begin = (*nodes_)[0]->getTime();
86 const Scalar t_final = (*nodes_)[n - 1]->getTime();
87
88 // Find largest index i such that
89 // (*nodes)[i]->getTime() <= t
90 int i;
91 if (t <= t_begin)
92 i = 0;
93 else if (t >= t_final)
94 i = n - 1;
95 else {
96 auto it = std::find_if(nodes_->begin(), nodes_->end(),
97 [=](const RCP<SolutionState<Scalar> >& s) {
98 return t <= s->getTime();
99 });
100 i = std::distance(nodes_->begin(), it) - 1;
101 }
102 TEUCHOS_ASSERT(i >= 0 && i <= n - 1);
103
104 // First we check for exact node matches:
105 if (floating_compare_equals((*nodes_)[i]->getTime(), t, t_final)) {
106 state_out->copy((*nodes_)[i]);
107 return;
108 }
109 if (i < n - 1 &&
110 floating_compare_equals((*nodes_)[i + 1]->getTime(), t, t_final)) {
111 state_out->copy((*nodes_)[i + 1]);
112 return;
113 }
114
115 // Put t roughly in the middle of the interpolation window
116 int k = i - p / 2;
117 if (k < 0) k = 0;
118 if (k + p + 1 > n) k = n - p - 1;
119 TEUCHOS_ASSERT(k >= 0 && k + p + 1 <= n);
120
121 RCP<VectorBase<Scalar> > x = state_out->getX();
122 RCP<VectorBase<Scalar> > xdot = state_out->getXDot();
123 RCP<VectorBase<Scalar> > xdotdot = state_out->getXDotDot();
124
125 Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
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);
129
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();
135
136 Scalar num = 1.0;
137 Scalar den = 1.0;
138 for (int l = 0; l <= p; ++l) {
139 if (l != j) {
140 const Scalar tl = (*nodes_)[k + l]->getTime();
141 num *= t - tl;
142 den *= tj - tl;
143 }
144 }
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);
149 }
150
151 // Set meta-data for interpolated state
152 state_out->getMetaData()->copy((*nodes_)[i]->getMetaData());
153 state_out->getMetaData()->setTime(t);
154 state_out->getMetaData()->setDt(t - (*nodes_)[i]->getTime());
155 state_out->getMetaData()->setIsInterpolated(true);
156 // What else should we set??
157}
158
159} // namespace Tempus
160
161#endif // Tempus_InterpolatorLagrange_impl_hpp
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
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.
Solution state for integrators and steppers.
virtual void copy(const Teuchos::RCP< const SolutionState< Scalar > > &ss)
This is a deep copy.
virtual Teuchos::RCP< Thyra::VectorBase< Scalar > > getX()
virtual Teuchos::RCP< Thyra::VectorBase< Scalar > > getXDotDot()
virtual Teuchos::RCP< Thyra::VectorBase< Scalar > > getXDot()
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.