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