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: Copyright (2017) Sandia Corporation
4 //
5 // Distributed under BSD 3-clause license (See accompanying file Copyright.txt)
6 // ****************************************************************************
7 // @HEADER
8 
9 #ifndef Tempus_TimeEventList_decl_hpp
10 #define Tempus_TimeEventList_decl_hpp
11 
12 #include <vector>
13 
14 #include "Teuchos_Time.hpp"
16 
17 #include "Tempus_config.hpp"
18 #include "Tempus_TimeEventBase.hpp"
19 
20 
21 namespace Tempus {
22 
23 
28 template<class Scalar>
29 class TimeEventList : virtual public TimeEventBase<Scalar>
30 {
31 public:
32 
34  TimeEventList();
35 
38  std::vector<Scalar> timeList,
39  std::string name = "TimeEventList",
40  bool landOnExactly = true,
41  Scalar relTol = std::numeric_limits<Scalar>::epsilon()*Scalar(100.0));
42 
44  virtual ~TimeEventList() {}
45 
47 
48 
58  virtual bool isTime(Scalar time) const;
59 
70  virtual Scalar timeToNextEvent(Scalar time) const;
71 
83  virtual Scalar timeOfNextEvent(Scalar time) const;
84 
97  virtual bool eventInRange(Scalar time1, Scalar time2) const;
98 
100  virtual void describe(Teuchos::FancyOStream &out,
101  const Teuchos::EVerbosityLevel verbLevel) const;
103 
105 
106  virtual std::vector<Scalar> getTimeList() const { return timeList_; }
108 
118  virtual void setTimeList(std::vector<Scalar> timeList, bool sort = true);
119 
129  virtual void addTime(Scalar time);
130 
132  virtual void clearTimeList() { timeList_.clear(); }
133 
135  virtual Scalar getRelTol() const { return relTol_; }
136 
146  virtual void setRelTol(Scalar relTol);
147 
158  virtual Scalar getAbsTol() const { return absTol_; }
159 
171  virtual bool getLandOnExactly() const { return landOnExactly_; }
172 
174  virtual void setLandOnExactly(bool LOE) { landOnExactly_ = LOE; }
176 
188 
189 
190 protected:
191 
200  virtual void setTimeScale();
201 
202  std::vector<Scalar> timeList_;
203 
204  Scalar timeScale_;
205  Scalar relTol_;
206  Scalar absTol_;
208 };
209 
210 
211 // Nonmember Contructors
212 // ------------------------------------------------------------------------
213 
223 template<class Scalar>
226 
227 
228 } // namespace Tempus
229 
230 #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.
bool landOnExactly_
Should these time events be landed on exactly, i.e, adjust the timestep to hit time event...
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.