Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_TimeEventListIndex_impl.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_TimeEventListIndex_impl_hpp
11#define Tempus_TimeEventListIndex_impl_hpp
12
13namespace Tempus {
14
15template <class Scalar>
17{
18 this->setType("List Index");
19 this->setName("TimeEventListIndex");
20}
21
22template <class Scalar>
24 std::string name)
25{
26 this->setType("List Index");
27 if (name == "" && !indexList.empty()) {
28 std::ostringstream oss;
29 oss << "TimeEventListIndex (" << indexList_.front() << ", ... ,"
30 << indexList_.back() << ")";
31 this->setName(oss.str());
32 }
33 else {
34 this->setName(name);
35 }
36
37 this->setIndexList(indexList);
38}
39
40template <class Scalar>
41void TimeEventListIndex<Scalar>::setIndexList(std::vector<int> indexList,
42 bool sort)
43{
44 indexList_ = indexList;
45 if (sort) {
46 std::sort(indexList_.begin(), indexList_.end());
47 indexList_.erase(std::unique(indexList_.begin(), indexList_.end()),
48 indexList_.end());
49 }
50}
51
52template <class Scalar>
54{
55 if (indexList_.size() == 0) {
56 indexList_.push_back(index);
57 return;
58 }
59
60 std::vector<int>::iterator it;
61 it = std::find(indexList_.begin(), indexList_.end(), index);
62 // Check if index is already in list.
63 if (it != indexList_.end()) return;
64
65 it = std::upper_bound(indexList_.begin(), indexList_.end(), index);
66 indexList_.insert(it, index);
67}
68
69template <class Scalar>
71{
72 return (std::find(indexList_.begin(), indexList_.end(), index) !=
73 indexList_.end());
74}
75
76template <class Scalar>
78{
79 return indexOfNextEvent(index) - index; // Neg. indicating in the past.
80}
81
82template <class Scalar>
84{
85 if (indexList_.size() == 0) return this->getDefaultIndex();
86
87 // Check if before first event.
88 if (index < indexList_.front()) return indexList_.front();
89
90 // Check if after last event.
91 if (index >= indexList_.back()) return this->getDefaultIndex();
92
93 std::vector<int>::const_iterator it =
94 std::upper_bound(indexList_.begin(), indexList_.end(), index);
95
96 return int(*it);
97}
98
99template <class Scalar>
100bool TimeEventListIndex<Scalar>::eventInRangeIndex(int index1, int index2) const
101{
102 if (index1 > index2) {
103 int tmp = index1;
104 index1 = index2;
105 index2 = tmp;
106 }
107
108 if (indexList_.size() == 0) return false;
109
110 // Check if range is completely outside index events.
111 if (index2 < indexList_.front() || indexList_.back() < index1) return false;
112
113 Scalar indexEvent1 = indexOfNextEvent(index1);
114 Scalar indexEvent2 = indexOfNextEvent(index2);
115 // Check if the next index event is different for the two indices.
116 if (indexEvent1 != indexEvent2) return true;
117
118 // Check if indices bracket index event.
119 if (index1 < indexEvent1 && indexEvent1 <= index2) return true;
120
121 return false;
122}
123
124template <class Scalar>
126 Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
127{
128 auto l_out = Teuchos::fancyOStream(out.getOStream());
129 Teuchos::OSTab ostab(*l_out, 2, "TimeEventListIndex");
130 l_out->setOutputToRootOnly(0);
131
132 *l_out << "TimeEventListIndex:"
133 << "\n"
134 << " name = " << this->getName() << "\n"
135 << " Type = " << this->getType() << "\n"
136 << " IndexList_ = ";
137 if (!indexList_.empty()) {
138 for (auto it = indexList_.begin(); it != indexList_.end() - 1; ++it)
139 *l_out << *it << ", ";
140 *l_out << *(indexList_.end() - 1) << std::endl;
141 }
142 else {
143 *l_out << "<empty>" << std::endl;
144 }
145}
146
147template <class Scalar>
148Teuchos::RCP<const Teuchos::ParameterList>
150{
151 Teuchos::RCP<Teuchos::ParameterList> pl =
152 Teuchos::parameterList("Time Event List Index");
153
154 pl->setName(this->getName());
155 pl->set("Name", this->getName());
156 pl->set("Type", this->getType());
157
158 std::ostringstream list;
159 if (!indexList_.empty()) {
160 for (std::size_t i = 0; i < indexList_.size() - 1; ++i)
161 list << indexList_[i] << ", ";
162 list << indexList_[indexList_.size() - 1];
163 }
164 pl->set<std::string>("Index List", list.str(),
165 "Comma deliminated list of indices");
166
167 return pl;
168}
169
170// Nonmember constructors.
171// ------------------------------------------------------------------------
172
173template <class Scalar>
174Teuchos::RCP<TimeEventListIndex<Scalar> > createTimeEventListIndex(
175 Teuchos::RCP<Teuchos::ParameterList> pl)
176{
177 auto teli = Teuchos::rcp(new TimeEventListIndex<Scalar>());
178 if (pl == Teuchos::null) return teli; // Return default TimeEventListIndex.
179
180 TEUCHOS_TEST_FOR_EXCEPTION(
181 pl->get<std::string>("Type", "List Index") != "List Index",
182 std::logic_error,
183 "Error - Time Event Type != 'List Index'. (='" +
184 pl->get<std::string>("Type") + "')\n");
185
186 pl->validateParametersAndSetDefaults(*teli->getValidParameters());
187
188 teli->setName(pl->get("Name", "From createTimeEventListIndex"));
189
190 std::vector<int> indexList;
191 indexList.clear();
192 std::string str = pl->get<std::string>("Index List");
193 std::string delimiters(",");
194 std::string::size_type lastPos = str.find_first_not_of(delimiters, 0);
195 std::string::size_type pos = str.find_first_of(delimiters, lastPos);
196 while ((pos != std::string::npos) || (lastPos != std::string::npos)) {
197 std::string token = str.substr(lastPos, pos - lastPos);
198 indexList.push_back(int(std::stoi(token)));
199 if (pos == std::string::npos) break;
200
201 lastPos = str.find_first_not_of(delimiters, pos);
202 pos = str.find_first_of(delimiters, lastPos);
203 }
204 teli->setIndexList(indexList);
205
206 return teli;
207}
208
209} // namespace Tempus
210#endif // Tempus_TimeEventListIndex_impl_hpp
TimeEventListIndex specifies a list of index events.
virtual bool eventInRangeIndex(int index1, int index2) const
Test if an event occurs within the index range.
virtual void setIndexList(std::vector< int > indexList, bool sort=true)
Set the vector of event indices.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
Describe member data.
virtual void addIndex(int index)
Add the index to event vector.
virtual bool isIndex(int index) const
Test if index is a time event.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a valid ParameterList with current settings.
virtual int indexToNextEvent(int index) const
How many indices until the next event.
virtual int indexOfNextEvent(int index) const
Return the index of the next event following the input index.
Teuchos::RCP< TimeEventListIndex< Scalar > > createTimeEventListIndex(Teuchos::RCP< Teuchos::ParameterList > pList)
Nonmember Constructor via ParameterList.