Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_TimeEventBase.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_TimeEventBase_decl_hpp
11 #define Tempus_TimeEventBase_decl_hpp
12 
13 // Teuchos
14 #include "Teuchos_Time.hpp"
15 #include "Teuchos_VerboseObject.hpp"
16 
17 #include "Tempus_config.hpp"
18 
19 namespace Tempus {
20 
35 template <class Scalar>
37  public:
40  : timeEventType_("Base"),
41  name_("TimeEventBase"),
42  defaultTime_(std::numeric_limits<Scalar>::max()),
43  defaultTol_(std::numeric_limits<Scalar>::epsilon() * Scalar(100.0)),
44  defaultIndex_(std::numeric_limits<int>::max())
45  {
46  }
47 
49  virtual ~TimeEventBase() {}
50 
52 
53 
63  virtual bool isTime(Scalar time) const { return false; }
64 
75  virtual Scalar timeToNextEvent(Scalar time) const { return defaultTime_; }
76 
93  virtual Scalar timeOfNextEvent(Scalar time) const { return defaultTime_; }
94 
112  virtual bool eventInRange(Scalar time1, Scalar time2) const { return false; }
113 
124  virtual bool isIndex(int index) const { return false; }
125 
135  virtual int indexToNextEvent(int index) const { return defaultIndex_; }
136 
153  virtual int indexOfNextEvent(int index) const { return defaultIndex_; }
154 
170  virtual bool eventInRangeIndex(int index1, int index2) const { return false; }
171 
173  virtual void describe(Teuchos::FancyOStream &out,
174  const Teuchos::EVerbosityLevel verbLevel) const
175  {
176  auto l_out = Teuchos::fancyOStream(out.getOStream());
177  Teuchos::OSTab ostab(*l_out, 2, "TimeEventBase");
178  l_out->setOutputToRootOnly(0);
179 
180  *l_out << "TimeEventBase name = " << getName() << std::endl;
181  }
182 
195  virtual Scalar getAbsTol() const { return defaultTol_; }
196 
214  virtual bool getLandOnExactly() const { return false; }
215 
217 
219 
220 
229  virtual std::string getName() const { return name_; }
231  virtual void setName(std::string name) { name_ = name; }
237  virtual std::string getType() const { return timeEventType_; }
239  virtual Scalar getDefaultTime() const { return defaultTime_; }
241  virtual Scalar getDefaultTol() const { return defaultTol_; }
243  virtual int getDefaultIndex() const { return defaultIndex_; }
245 
248  {
250  Teuchos::parameterList("Time Event Base");
251 
252  pl->setName(this->getName());
253  pl->set("Name", this->getName());
254  pl->set("Type", this->getType());
255 
256  return pl;
257  }
258 
259  protected:
260  virtual void setType(std::string s) { timeEventType_ = s; }
261 
262  private:
263  std::string timeEventType_;
264  std::string name_;
265  const Scalar defaultTime_;
266  const Scalar defaultTol_;
267  const int defaultIndex_;
268 };
269 
270 } // namespace Tempus
271 
272 #endif // Tempus_TimeEventBase_decl_hpp
virtual bool getLandOnExactly() const
Return if the time events need to be landed on exactly.
virtual std::string getType() const
Return the type of TimeEvent.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
Describe member data.
virtual Scalar getDefaultTime() const
Return the default time used for TimeEvents.
const Scalar defaultTol_
Default tolerance.
virtual void setType(std::string s)
const int defaultIndex_
Default index.
virtual void setName(std::string name)
Set the name of the TimeEvent.
virtual int indexOfNextEvent(int index) const
Return the index of the next event following the input index.
virtual bool eventInRange(Scalar time1, Scalar time2) const
Test if an event occurs within the time range.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
virtual Scalar timeToNextEvent(Scalar time) const
How much time until the next event.
std::string name_
Name to identify the TimeEvent.
virtual bool isTime(Scalar time) const
Test if time is near an event (within tolerance).
virtual int getDefaultIndex() const
Return the default index used by TimeEvents.
const Scalar defaultTime_
Default time.
virtual Scalar getDefaultTol() const
Return the default tolerance used by TimeEvents.
virtual bool isIndex(int index) const
Test if index is an event.
std::string timeEventType_
Time Event type.
virtual Scalar timeOfNextEvent(Scalar time) const
Return the time of the next event following the input time.
virtual Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return ParameterList with current values.
This class defines time events which can be used to &quot;trigger&quot; an action.
virtual bool eventInRangeIndex(int index1, int index2) const
Test if an event occurs within the index range.
virtual std::string getName() const
Return the name of the TimeEvent.
RCP< std::basic_ostream< char_type, traits_type > > getOStream()
virtual Scalar getAbsTol() const
Return the absolute tolerance.
virtual ~TimeEventBase()
Destructor.
virtual int indexToNextEvent(int index) const
How many indices until the next event.
ParameterList & setName(const std::string &name)