Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_TimeEventList_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_TimeEventList_decl_hpp
11 #define Tempus_TimeEventList_decl_hpp
12 
13 #include <vector>
14 
15 #include "Teuchos_Time.hpp"
17 
18 #include "Tempus_config.hpp"
19 #include "Tempus_TimeEventBase.hpp"
20 
21 namespace Tempus {
22 
27 template <class Scalar>
28 class TimeEventList : virtual public TimeEventBase<Scalar> {
29  public:
31  TimeEventList();
32 
34  TimeEventList(std::vector<Scalar> timeList,
35  std::string name = "TimeEventList", bool landOnExactly = true,
36  Scalar relTol = std::numeric_limits<Scalar>::epsilon() *
37  Scalar(100.0));
38 
40  virtual ~TimeEventList() {}
41 
43 
44 
54  virtual bool isTime(Scalar time) const;
55 
66  virtual Scalar timeToNextEvent(Scalar time) const;
67 
79  virtual Scalar timeOfNextEvent(Scalar time) const;
80 
93  virtual bool eventInRange(Scalar time1, Scalar time2) const;
94 
96  virtual void describe(Teuchos::FancyOStream &out,
97  const Teuchos::EVerbosityLevel verbLevel) const;
99 
101 
102  virtual std::vector<Scalar> getTimeList() const { return timeList_; }
104 
114  virtual void setTimeList(std::vector<Scalar> timeList, bool sort = true);
115 
125  virtual void addTime(Scalar time);
126 
128  virtual void clearTimeList() { timeList_.clear(); }
129 
131  virtual Scalar getRelTol() const { return relTol_; }
132 
142  virtual void setRelTol(Scalar relTol);
143 
154  virtual Scalar getAbsTol() const { return absTol_; }
155 
168  virtual bool getLandOnExactly() const { return landOnExactly_; }
169 
171  virtual void setLandOnExactly(bool LOE) { landOnExactly_ = LOE; }
173 
185 
186  protected:
195  virtual void setTimeScale();
196 
197  std::vector<Scalar> timeList_;
198 
199  Scalar timeScale_;
200  Scalar relTol_;
201  Scalar absTol_;
203 };
206 
207 // Nonmember Contructors
208 // ------------------------------------------------------------------------
209 
219 template <class Scalar>
222 
223 } // namespace Tempus
224 
225 #endif // Tempus_TimeEventList_decl_hpp
virtual bool isTime(Scalar time) const
Test if time is near an event (within tolerance).
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a valid ParameterList with current settings.
Scalar relTol_
Relative time tolerance for matching time events.
virtual ~TimeEventList()
Destructor.
Scalar absTol_
Absolute time tolerance, relTol_*timeScale_.
virtual void addTime(Scalar time)
Add the time to event vector.
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.
virtual void setTimeScale()
Set the time scale for the time events.
std::vector< Scalar > timeList_
Sorted and unique list of time events.
virtual void clearTimeList()
Clear the vector of all events.
virtual void setRelTol(Scalar relTol)
Set the relative tolerance.
virtual void setTimeList(std::vector< Scalar > timeList, bool sort=true)
Set the list of time events.
This class defines time events which can be used to &quot;trigger&quot; an action.
virtual bool getLandOnExactly() const
Return if the time events need to be landed on exactly.
Scalar timeScale_
A reference time scale.
virtual Scalar getAbsTol() const
Return the absolute tolerance.
TimeEventList()
Default constructor.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
Describe member data.
virtual Scalar timeToNextEvent(Scalar time) const
How much time until the next event.
virtual Scalar getRelTol() const
Return the relative tolerance.
virtual void setLandOnExactly(bool LOE)
Set if the time event should be landed on exactly.
TimeEventList specifies a list of time events.
Teuchos::RCP< TimeEventList< Scalar > > createTimeEventList(Teuchos::RCP< Teuchos::ParameterList > pList)
Nonmember Constructor via ParameterList.
virtual std::vector< Scalar > getTimeList() const
Return the list of time events.