Tempus
Version of the Day
Time Integration
Toggle main menu visibility
Loading...
Searching...
No Matches
src
Tempus_TimeEventListIndex_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_TimeEventListIndex_decl_hpp
10
#define Tempus_TimeEventListIndex_decl_hpp
11
12
#include <vector>
13
14
#include "Teuchos_Time.hpp"
15
#include "Teuchos_ParameterList.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
TimeEventListIndex
:
virtual
public
TimeEventBase
<Scalar>
30
{
31
public
:
32
34
TimeEventListIndex
();
35
37
TimeEventListIndex
(std::vector<int> indexList,
38
std::string name =
"TimeEventListIndex"
);
39
41
virtual
~TimeEventListIndex
() {}
42
44
45
52
virtual
bool
isIndex
(
int
index)
const
;
53
59
virtual
int
indexToNextEvent
(
int
index)
const
;
60
72
virtual
int
indexOfNextEvent
(
int
index)
const
;
73
83
virtual
bool
eventInRangeIndex
(
int
index1,
int
index2)
const
;
84
86
virtual
void
describe
(Teuchos::FancyOStream &out,
87
const
Teuchos::EVerbosityLevel verbLevel)
const
;
89
91
92
93
virtual
std::vector<int>
getIndexList
()
const
{
return
indexList_
; }
94
102
virtual
void
setIndexList
(std::vector<int> indexList,
bool
sort =
true
);
103
112
virtual
void
addIndex
(
int
index);
113
115
virtual
void
clearIndexList
() {
indexList_
.clear(); }
117
126
Teuchos::RCP<const Teuchos::ParameterList>
getValidParameters
()
const
;
127
128
129
protected
:
130
131
std::vector<int>
indexList_
;
// Sorted and unique list of index events.
132
133
};
134
135
136
// Nonmember Contructors
137
// ------------------------------------------------------------------------
138
148
template
<
class
Scalar>
149
Teuchos::RCP<TimeEventListIndex<Scalar> >
150
createTimeEventListIndex
(Teuchos::RCP<Teuchos::ParameterList> pList);
151
152
153
}
// namespace Tempus
154
155
#endif
// Tempus_TimeEventListIndex_decl_hpp
Tempus_TimeEventBase.hpp
Tempus::TimeEventBase::TimeEventBase
TimeEventBase()
Constructor.
Definition
Tempus_TimeEventBase.hpp:40
Tempus::TimeEventListIndex::TimeEventListIndex
TimeEventListIndex()
Default constructor.
Definition
Tempus_TimeEventListIndex_impl.hpp:16
Tempus::TimeEventListIndex::eventInRangeIndex
virtual bool eventInRangeIndex(int index1, int index2) const
Test if an event occurs within the index range.
Definition
Tempus_TimeEventListIndex_impl.hpp:104
Tempus::TimeEventListIndex::indexList_
std::vector< int > indexList_
Definition
Tempus_TimeEventListIndex_decl.hpp:131
Tempus::TimeEventListIndex::~TimeEventListIndex
virtual ~TimeEventListIndex()
Destructor.
Definition
Tempus_TimeEventListIndex_decl.hpp:41
Tempus::TimeEventListIndex::setIndexList
virtual void setIndexList(std::vector< int > indexList, bool sort=true)
Set the vector of event indices.
Definition
Tempus_TimeEventListIndex_impl.hpp:42
Tempus::TimeEventListIndex::describe
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
Describe member data.
Definition
Tempus_TimeEventListIndex_impl.hpp:130
Tempus::TimeEventListIndex::addIndex
virtual void addIndex(int index)
Add the index to event vector.
Definition
Tempus_TimeEventListIndex_impl.hpp:54
Tempus::TimeEventListIndex::isIndex
virtual bool isIndex(int index) const
Test if index is a time event.
Definition
Tempus_TimeEventListIndex_impl.hpp:72
Tempus::TimeEventListIndex::clearIndexList
virtual void clearIndexList()
Clear the vector of all events.
Definition
Tempus_TimeEventListIndex_decl.hpp:115
Tempus::TimeEventListIndex::getValidParameters
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a valid ParameterList with current settings.
Definition
Tempus_TimeEventListIndex_impl.hpp:153
Tempus::TimeEventListIndex::getIndexList
virtual std::vector< int > getIndexList() const
Return a vector of event indices.
Definition
Tempus_TimeEventListIndex_decl.hpp:93
Tempus::TimeEventListIndex::indexToNextEvent
virtual int indexToNextEvent(int index) const
How many indices until the next event.
Definition
Tempus_TimeEventListIndex_impl.hpp:79
Tempus::TimeEventListIndex::indexOfNextEvent
virtual int indexOfNextEvent(int index) const
Return the index of the next event following the input index.
Definition
Tempus_TimeEventListIndex_impl.hpp:86
Tempus
Definition
Tempus_AdjointAuxSensitivityModelEvaluator_decl.hpp:21
Tempus::createTimeEventListIndex
Teuchos::RCP< TimeEventListIndex< Scalar > > createTimeEventListIndex(Teuchos::RCP< Teuchos::ParameterList > pList)
Nonmember Constructor via ParameterList.
Definition
Tempus_TimeEventListIndex_impl.hpp:179
Generated by
1.17.0