10 #ifndef Tempus_TimeEventList_decl_hpp
11 #define Tempus_TimeEventList_decl_hpp
18 #include "Tempus_config.hpp"
27 template <
class Scalar>
35 std::string name =
"TimeEventList",
bool landOnExactly =
true,
36 Scalar relTol = std::numeric_limits<Scalar>::epsilon() *
54 virtual bool isTime(Scalar time)
const;
93 virtual bool eventInRange(Scalar time1, Scalar time2)
const;
114 virtual void setTimeList(std::vector<Scalar> timeList,
bool sort =
true);
125 virtual void addTime(Scalar time);
219 template <
class Scalar>
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 "trigger" 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.