Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_TimeEventRange_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_TimeEventRange_impl_hpp
11#define Tempus_TimeEventRange_impl_hpp
12
14
15namespace Tempus {
16
17template <class Scalar>
19 : start_(0.0),
20 stop_(0.0),
21 stride_(0.0),
22 numEvents_(1),
23 timeScale_(1.0),
24 relTol_(this->getDefaultTol()),
25 absTol_(this->getDefaultTol()),
26 landOnExactly_(true)
27{
28 this->setType("Range");
29 setRelTol(this->getDefaultTol()), setTimeRange(0.0, 0.0, 0.0);
30 setLandOnExactly(true);
31 std::ostringstream oss;
32 oss << "TimeEventRange (" << start_ << "; " << stop_ << "; " << stride_
33 << ")";
34 this->setName(oss.str());
35}
36
37template <class Scalar>
38TimeEventRange<Scalar>::TimeEventRange(Scalar start, Scalar stop, Scalar stride,
39 std::string name, bool landOnExactly,
40 Scalar relTol)
41 : start_(start),
42 stop_(stop),
43 stride_(stride),
44 numEvents_(0),
45 timeScale_(std::max(std::abs(start_), std::abs(stop_))),
46 relTol_(relTol),
47 absTol_(relTol_ * timeScale_),
48 landOnExactly_(landOnExactly)
49{
50 this->setType("Range");
51 if (name == "") {
52 std::ostringstream oss;
53 oss << "TimeEventRange (" << start << "; " << stop << "; " << stride << ")";
54 this->setName(oss.str());
55 }
56 else {
57 this->setName(name);
58 }
59 setRelTol(relTol);
60 setTimeRange(start, stop, stride);
61 setLandOnExactly(landOnExactly);
62}
63
64template <class Scalar>
65TimeEventRange<Scalar>::TimeEventRange(Scalar start, Scalar stop, int numEvents,
66 std::string name, bool landOnExactly,
67 Scalar relTol)
68 : start_(start),
69 stop_(stop),
70 stride_(0.0),
71 numEvents_(numEvents),
72 timeScale_(std::max(std::abs(start_), std::abs(stop_))),
73 relTol_(relTol),
74 absTol_(relTol_ * timeScale_),
75 landOnExactly_(landOnExactly)
76{
77 if (name == "") {
78 std::ostringstream oss;
79 oss << "TimeEventRange (" << start << "; " << stop << "; " << numEvents
80 << ")";
81 this->setName(oss.str());
82 }
83 else {
84 this->setName(name);
85 }
86 setRelTol(relTol);
87 setTimeRange(start, stop, numEvents);
88 setLandOnExactly(landOnExactly);
89}
90
91template <class Scalar>
92void TimeEventRange<Scalar>::setTimeRange(Scalar start, Scalar stop,
93 Scalar stride)
94{
95 start_ = start;
96 stop_ = stop;
97 if (stop_ < start_) {
98 Scalar tmp = start_;
99 start_ = stop_;
100 stop_ = tmp;
101 }
102 setTimeScale();
103 setTimeStride(stride);
104}
105
106template <class Scalar>
107void TimeEventRange<Scalar>::setTimeRange(Scalar start, Scalar stop,
108 int numEvents)
109{
110 start_ = start;
111 stop_ = stop;
112 if (stop_ < start_) {
113 Scalar tmp = start_;
114 start_ = stop_;
115 stop_ = tmp;
116 }
117 setTimeScale();
118 setNumEvents(numEvents);
119}
120
121template <class Scalar>
123{
124 start_ = start;
125 if (stop_ < start_) stop_ = start_;
126 setTimeScale();
127 setTimeStride(stride_); // Reset numEvents with the current stride.
128}
129
130template <class Scalar>
132{
133 stop_ = stop;
134 if (start_ > stop_) start_ = stop_;
135 setTimeScale();
136 setTimeStride(stride_); // Reset numEvents with the current stride.
137}
138
139template <class Scalar>
141{
142 timeScale_ = std::max(std::abs(start_), std::abs(stop_));
143 absTol_ = relTol_ * timeScale_;
144
145 // Check if timeScale is near zero.
146 if (approxZero(timeScale_, absTol_)) {
147 timeScale_ = 1.0;
148 absTol_ = relTol_ * timeScale_;
149 }
150}
151
152template <class Scalar>
154{
155 stride_ = Teuchos::ScalarTraits<Scalar>::magnitude(stride);
156 if (approxEqualAbsTol(start_, stop_, absTol_)) {
157 stride_ = 0.0;
158 numEvents_ = 1;
159 return;
160 }
161
162 if ((stride_ > stop_ - start_) || (stride_ < 2 * absTol_)) {
163 stride_ = stop_ - start_;
164 }
165
166 numEvents_ = int((stop_ + absTol_ - start_) / stride_) + 1;
167}
168
169template <class Scalar>
171{
172 if (numEvents < 0) { // Do not use numEvents to set stride!
173 // Reset numEvents_ from stride_.
174 if (stride_ < 2 * absTol_) stride_ = 2 * absTol_;
175 numEvents_ = int((stop_ + absTol_ - start_) / stride_) + 1;
176 return;
177 }
178 else if (approxEqualAbsTol(start_, stop_, absTol_)) {
179 numEvents_ = 1;
180 stride_ = stop_ - start_;
181 }
182 else {
183 numEvents_ = numEvents;
184 if (numEvents_ < 2) numEvents_ = 2;
185 stride_ = (stop_ - start_) / Scalar(numEvents_ - 1);
186 }
187
188 // If stride_ is smaller than twice the absolute tolerance,
189 // the time steps cannot be differentiated!
190 if (stride_ <= 2 * absTol_) setTimeStride(2 * absTol_);
191}
192
193template <class Scalar>
195{
196 relTol_ = std::abs(relTol);
197 setTimeScale();
198}
199
200template <class Scalar>
201bool TimeEventRange<Scalar>::isTime(Scalar time) const
202{
203 // Check if before first event.
204 if (time < start_ - absTol_) return false;
205
206 // Check if after last event.
207 const Scalar timeOfLast = start_ + (numEvents_ - 1) * stride_;
208 if (time > timeOfLast + absTol_) return false;
209
210 int numStrides = 0;
211 if (!approxZero(stride_, 2 * absTol_)) numStrides = (time - start_) / stride_;
212
213 numStrides = std::min(std::max(0, numStrides), int(numEvents_ - 1));
214 const Scalar leftBracket = start_ + numStrides * stride_;
215
216 // Check if close to left bracket.
217 if (approxEqualAbsTol(time, leftBracket, absTol_)) return true;
218
219 // Check if close to right bracket.
220 const Scalar rightBracket = leftBracket + stride_;
221 if (approxEqualAbsTol(time, rightBracket, absTol_)) return true;
222
223 return false;
224}
225
226template <class Scalar>
228{
229 return timeOfNextEvent(time) - time; // Neg. time indicates in the past.
230}
231
232template <class Scalar>
234{
235 // Check if before first event.
236 if (time < start_ - absTol_) return start_;
237
238 const Scalar timeOfLast = start_ + (numEvents_ - 1) * stride_;
239 // Check if after or close to last event.
240 if (time > timeOfLast - absTol_) return std::numeric_limits<Scalar>::max();
241
242 int numStrides = 0;
243 if (!approxZero(stride_, 2 * absTol_))
244 numStrides = int((time - start_) / stride_) + 1;
245 const Scalar timeEvent = start_ + numStrides * stride_;
246
247 // Check timeEvent is near time. If so, return the next event.
248 if (approxEqualAbsTol(time, timeEvent, absTol_)) return timeEvent + stride_;
249
250 return timeEvent;
251}
252
253template <class Scalar>
254bool TimeEventRange<Scalar>::eventInRange(Scalar time1, Scalar time2) const
255{
256 if (time1 > time2) {
257 Scalar tmp = time1;
258 time1 = time2;
259 time2 = tmp;
260 }
261
262 // Check if range is completely outside time events.
263 const Scalar timeOfLast = start_ + (numEvents_ - 1) * stride_;
264 if (time2 < start_ - absTol_ || timeOfLast + absTol_ < time1) return false;
265
266 if (approxZero(stride_))
267 return (time1 < start_ - absTol_ && start_ - absTol_ <= time2);
268
269 const int strideJustBeforeTime1 = std::min(
270 int(numEvents_ - 1),
271 std::max(int(0), int((time1 - start_ + absTol_) / stride_ - 0.5)));
272
273 const int strideJustAfterTime2 = std::max(
274 int(0), std::min(int(numEvents_ - 1),
275 int((time2 - start_ + absTol_) / stride_ + 0.5)));
276
277 for (int i = strideJustBeforeTime1; i <= strideJustAfterTime2; i++) {
278 const Scalar timeEvent = start_ + i * stride_;
279 if (time1 < timeEvent - absTol_ && timeEvent - absTol_ <= time2)
280 return true;
281 }
282
283 return false;
284}
285
286template <class Scalar>
288 Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
289{
290 auto l_out = Teuchos::fancyOStream(out.getOStream());
291 Teuchos::OSTab ostab(*l_out, 2, "TimeEventRange");
292 l_out->setOutputToRootOnly(0);
293
294 *l_out << "TimeEventRange:"
295 << "\n"
296 << " name = " << this->getName() << "\n"
297 << " Type = " << this->getType() << "\n"
298 << " start_ = " << start_ << "\n"
299 << " stop_ = " << stop_ << "\n"
300 << " stride_ = " << stride_ << "\n"
301 << " numEvents_ = " << numEvents_ << "\n"
302 << " timeScale_ = " << timeScale_ << "\n"
303 << " relTol_ = " << relTol_ << "\n"
304 << " absTol_ = " << absTol_ << "\n"
305 << " landOnExactly_ = " << landOnExactly_ << std::endl;
306}
307
308template <class Scalar>
309Teuchos::RCP<const Teuchos::ParameterList>
311{
312 Teuchos::RCP<Teuchos::ParameterList> pl =
313 Teuchos::parameterList("Time Event Range");
314
315 pl->setName(this->getName());
316 pl->set("Name", this->getName());
317 pl->set("Type", this->getType());
318
319 pl->set("Start Time", getTimeStart(), "Start of time range");
320 pl->set("Stop Time", getTimeStop(), "Stop of time range");
321 pl->set("Stride Time", getTimeStride(), "Stride of time range");
322
323 if (getTimeStride() * Scalar(getNumEvents() - 1) -
324 (getTimeStop() - getTimeStart()) <
325 getAbsTol())
326 pl->set("Number of Events", getNumEvents(),
327 "Number of events in time range. If specified, 'Stride Time' is "
328 "reset.");
329
330 pl->set("Relative Tolerance", getRelTol(),
331 "Relative time tolerance for matching time events.");
332
333 pl->set("Land On Exactly", getLandOnExactly(),
334 "Should these time events be landed on exactly, i.e, adjust the "
335 "timestep to hit time event, versus stepping over and keeping the "
336 "time step unchanged.");
337
338 return pl;
339}
340
341// Nonmember constructors.
342// ------------------------------------------------------------------------
343
344template <class Scalar>
345Teuchos::RCP<TimeEventRange<Scalar> > createTimeEventRange(
346 Teuchos::RCP<Teuchos::ParameterList> pl)
347{
348 auto ter = Teuchos::rcp(new TimeEventRange<Scalar>());
349 if (pl == Teuchos::null) return ter; // Return default TimeEventRange.
350
351 TEUCHOS_TEST_FOR_EXCEPTION(pl->get<std::string>("Type", "Range") != "Range",
352 std::logic_error,
353 "Error - Time Event Type != 'Range'. (='" +
354 pl->get<std::string>("Type") + "')\n");
355
356 auto validPL = *ter->getValidParameters();
357 bool numEventsFound = pl->isParameter("Number of Events");
358 if (!numEventsFound) validPL.remove("Number of Events");
359
360 pl->validateParametersAndSetDefaults(validPL);
361
362 ter->setName(pl->get("Name", "From createTimeEventRange"));
363 ter->setTimeStart(pl->get("Start Time", ter->getTimeStart()));
364 ter->setTimeStop(pl->get("Stop Time", ter->getTimeStop()));
365 ter->setTimeStride(pl->get("Stride Time", ter->getTimeStride()));
366 if (numEventsFound)
367 ter->setNumEvents(pl->get("Number of Events", ter->getNumEvents()));
368 ter->setRelTol(pl->get("Relative Tolerance", ter->getRelTol()));
369 ter->setLandOnExactly(pl->get("Land On Exactly", ter->getLandOnExactly()));
370
371 return ter;
372}
373
374} // namespace Tempus
375#endif // Tempus_TimeEventRange_impl_hpp
virtual void setType(std::string s)
virtual void setName(std::string name)
Set the name of the TimeEvent.
virtual Scalar getDefaultTol() const
Return the default tolerance used by TimeEvents.
TimeEventRange specifies a start, stop and stride time.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a valid ParameterList with current settings.
virtual void setRelTol(Scalar relTol)
Set the relative tolerance.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
Describe member data.
virtual void setNumEvents(int numEvents)
Set the number of time events.
virtual void setTimeScale()
Set the time scale for the time events.
virtual void setLandOnExactly(bool LOE)
Set if the time event should be landed on exactly.
virtual void setTimeStop(Scalar stop)
Set the stop of the time range.
virtual void setTimeStride(Scalar stride)
Set the stride of the time range.
Scalar stride_
Stride of time range.
virtual Scalar timeOfNextEvent(Scalar time) const
Return the time of the next event following the input time.
virtual bool eventInRange(Scalar time1, Scalar time2) const
Test if an event occurs within the time range.
Scalar start_
Start of time range.
virtual void setTimeStart(Scalar start)
Set the start of the time range.
virtual void setTimeRange(Scalar start, Scalar stop, Scalar stride)
Set the range of time events from start, stop and stride.
virtual bool isTime(Scalar time) const
Test if time is near an event (within tolerance).
Scalar stop_
Stop of time range.
virtual Scalar timeToNextEvent(Scalar time) const
How much time until the next event.
bool approxZero(Scalar value, Scalar tol=Teuchos::ScalarTraits< Scalar >::sfmin())
Test if value is approximately zero within tolerance.
Teuchos::RCP< TimeEventRange< Scalar > > createTimeEventRange(Teuchos::RCP< Teuchos::ParameterList > pList)
Nonmember Constructor via ParameterList.
bool approxEqualAbsTol(Scalar a, Scalar b, Scalar absTol)
Test if values are approximately equal within the absolute tolerance.