DOpE
timeiterator.h
Go to the documentation of this file.
1 
24 #ifndef TIMEITERATOR_H_
25 #define TIMEITERATOR_H_
26 
27 #include <deal.II/dofs/dof_handler.h>
28 #include <deal.II/dofs/dof_accessor.h>
29 #include <deal.II/grid/tria_iterator.h>
30 
31 namespace DOpE
32 {
33  namespace IteratorState
34  {
35 
41  {
50  };
51  }
52 
53  typedef dealii::DoFHandler<1>::active_cell_iterator active_cell_it;
54 
55 
63  {
64  public:
65  //constructors
66 
71  {
72  present_index_ = -3;//i.e. invalid
73  }
74 
78  TimeIterator(const active_cell_it& element, int present_index) :
79  element_(element)
80  {
81  present_index_ = present_index;
82  }
83 
88  {
89  element_ = it.element_;
90  present_index_ = it.present_index_;
91  }
92 
97  GetState() const
98  {
99  if (present_index_ >= 0 && present_index_
100  < static_cast<int>(element_->get_triangulation().n_active_cells()))
101  return IteratorState::valid;
102  else
103  {
104  if (present_index_ == -1)
106  else
107  {
108  if (present_index_ == -2)
110  else
111  return IteratorState::invalid;
112  }
113  }
114  }
115 
121  int
122  GetIndex() const
123  {
124  return present_index_;
125  }
126 
127  TimeIterator&
128  operator=(const TimeIterator& element)
129  {
130  element_ = element.element_;
131  present_index_ = element.present_index_;
132  return *this;
133  }
134 
135  void
136  Initialize(const active_cell_it& element, int present_index)
137  {
138  element_ = element;
139  present_index_ = present_index;
140  }
141 
142  bool
143  operator==(const TimeIterator& element)
144  {
145  if (GetState() == element.GetState())
146  return (element_ == element.getelement_());
147  else
148  return false;
149  }
150 
151  bool
152  operator!=(const TimeIterator& element)
153  {
154  if (GetState() == element.GetState())
155  return (element_ != element.getelement_());
156  else
157  return true;
158  }
159 
163  double
164  get_center() const
165  {
166  assert(GetState()==IteratorState::valid);
167  return element_->center()(0);
168  }
169 
173  double
174  get_left() const
175  {
176  assert(GetState()==IteratorState::valid);
177  return element_->face(0)->center()(0);
178  }
179 
183  double
184  get_right() const
185  {
186  assert(GetState()==IteratorState::valid);
187  return element_->face(1)->center()(0);
188  }
189 
193  double
194  get_k() const
195  {
196  assert(GetState()==IteratorState::valid);
197  return element_->diameter();
198  }
199 
200  const active_cell_it&
201  getelement_() const
202  {
203  return element_;
204  }
205 
206  void
207  get_time_dof_indices(std::vector<unsigned int>& local_dof_indices) const
208  {
209  assert(GetState()==IteratorState::valid);
210  element_->get_dof_indices(local_dof_indices);
211  }
212 
216  TimeIterator &
218  {
219  assert(GetState()==IteratorState::valid);
220  ++present_index_;
221  if (present_index_ < static_cast<int>(element_->get_triangulation().n_active_cells()))
222  {
223  element_ = element_->neighbor(1);
224  }
225  else
226  {
227  present_index_ = -1;
228  }
229 
230  return *this;
231  }
232 
238  {
239  assert(GetState()==IteratorState::valid);
240  TimeIterator tmp(element_, present_index_);
241 
242  ++present_index_;
243  if (present_index_ < static_cast<int>(element_->get_triangulation().n_active_cells()))
244  {
245  element_ = element_->neighbor(1);
246  }
247  else
248  {
249  present_index_ = -1;
250  }
251  return tmp;
252  }
253 
257  TimeIterator&
259  {
260  assert(GetState()==IteratorState::valid);
261  --present_index_;
262  if (present_index_ >= 0)
263  {
264  element_ = element_->neighbor(0);
265  }
266  else
267  {
268  present_index_ = -2;
269  }
270 
271  return *this;
272  }
273 
279  {
280  assert(GetState()==IteratorState::valid);
281  TimeIterator tmp(element_, present_index_);
282  if (present_index_ >= 0)
283  {
284  element_ = element_->neighbor(0);
285  }
286  else
287  {
288  present_index_ = -2;
289  }
290  return tmp;
291  }
292  private:
293  active_cell_it element_;
294  int present_index_;
295  };
296 }
297 
298 #endif /* TIMEITERATOR_H_ */
IteratorState::IteratorStates GetState() const
Definition: timeiterator.h:97
TimeIterator(const active_cell_it &element, int present_index)
Definition: timeiterator.h:78
Iterator reached end of container.
Definition: timeiterator.h:45
TimeIterator & operator=(const TimeIterator &element)
Definition: timeiterator.h:128
bool operator==(const TimeIterator &element)
Definition: timeiterator.h:143
TimeIterator operator--(int)
Definition: timeiterator.h:278
dealii::DoFHandler< 1 >::active_cell_iterator active_cell_it
Definition: timeiterator.h:53
const active_cell_it & getelement_() const
Definition: timeiterator.h:201
TimeIterator & operator++()
Definition: timeiterator.h:217
Iterator is invalid, probably due to an error.
Definition: timeiterator.h:49
double get_center() const
Definition: timeiterator.h:164
TimeIterator operator++(int)
Definition: timeiterator.h:237
Definition: timeiterator.h:47
Definition: timeiterator.h:62
double get_k() const
Definition: timeiterator.h:194
bool operator!=(const TimeIterator &element)
Definition: timeiterator.h:152
TimeIterator()
Definition: timeiterator.h:70
void get_time_dof_indices(std::vector< unsigned int > &local_dof_indices) const
Definition: timeiterator.h:207
void Initialize(const active_cell_it &element, int present_index)
Definition: timeiterator.h:136
IteratorStates
Definition: timeiterator.h:40
double get_right() const
Definition: timeiterator.h:184
TimeIterator(const TimeIterator &it)
Definition: timeiterator.h:87
double get_left() const
Definition: timeiterator.h:174
Iterator points to a valid object.
Definition: timeiterator.h:43
int GetIndex() const
Definition: timeiterator.h:122
TimeIterator & operator--()
Definition: timeiterator.h:258