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 // Teuchos
15 #include "Teuchos_Time.hpp"
16 
17 // Tempus
18 #include "Tempus_config.hpp"
19 #include "Tempus_TimeEventBase.hpp"
20 
21 
22 namespace Tempus {
23 
24 
29 template<class Scalar>
30 class TimeEventList : virtual public TimeEventBase<Scalar>
31 {
32 public:
33 
35  TimeEventList();
36 
38  TimeEventList(std::string name, std::vector<Scalar> timeList,
39  Scalar relTol, bool landOnExactly);
40 
42  virtual ~TimeEventList() {}
43 
45 
46  virtual bool isTime(Scalar time) const;
48 
50  virtual Scalar timeToNextEvent(Scalar time) const;
51 
53  virtual Scalar timeOfNextEvent(Scalar time) const;
54 
56  virtual bool eventInRange(Scalar time1, Scalar time2) const;
57 
59  virtual void describe() const;
61 
63 
64  virtual std::vector<Scalar> getTimeList() const { return timeList_; }
65  virtual void setTimeList(std::vector<Scalar> timeList);
66  virtual void addTime(Scalar time);
67  virtual void clearTimeList() { timeList_.clear(); }
68 
69  virtual Scalar getRelTol() const { return relTol_; }
70  virtual void setRelTol(Scalar relTol);
71 
72  virtual Scalar getAbsTol() const { return absTol_; }
73 
74  virtual bool getLandOnExactly() const { return landOnExactly_; }
75  virtual void setLandOnExactly(bool LOE) { landOnExactly_ = LOE; }
77 
78 
79 protected:
80 
81  virtual void setTimeScale();
82 
83  std::vector<Scalar> timeList_;
84 
85  Scalar timeScale_;
86  Scalar relTol_;
87  Scalar absTol_;
89 };
90 
91 
92 } // namespace Tempus
93 
94 #endif // Tempus_TimeEventList_decl_hpp
virtual void setTimeList(std::vector< Scalar > timeList)
virtual bool isTime(Scalar time) const
Test if time is near a TimeEvent (within tolerance).
Scalar relTol_
Relative time tolerance for matching time events.
virtual ~TimeEventList()
Destructor.
Scalar absTol_
Absolute time tolerance, relTol_*timeScale_.
virtual void addTime(Scalar time)
virtual Scalar timeOfNextEvent(Scalar time) const
Time of the next event. Negative indicating the last event is in the past.
virtual bool eventInRange(Scalar time1, Scalar time2) const
Test if an event occurs within the time range.
virtual void describe() const
Describe member data.
std::vector< Scalar > timeList_
Sorted and unique list of time events.
virtual void setRelTol(Scalar relTol)
This class defines time events which can be used to &quot;trigger&quot; an action. Time events are points in ti...
virtual bool getLandOnExactly() const
Scalar timeScale_
A reference time scale, max(abs(start_,stop_)).
virtual Scalar getAbsTol() const
TimeEventList()
Default constructor.
virtual Scalar timeToNextEvent(Scalar time) const
How much time until the next event. Negative indicating the last event is in the past.
virtual Scalar getRelTol() const
virtual void setLandOnExactly(bool LOE)
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...
virtual std::vector< Scalar > getTimeList() const