24 #ifndef TIMEITERATOR_H_
25 #define TIMEITERATOR_H_
27 #include <deal.II/dofs/dof_handler.h>
28 #include <dofs/dof_accessor.h>
30 #include <grid/tria_iterator.h>
34 namespace IteratorState
82 present_index_ = present_index;
90 element_ = it.element_;
91 present_index_ = it.present_index_;
100 if (present_index_ >= 0 && present_index_
101 < static_cast<int>(element_->get_triangulation().n_active_cells()))
105 if (present_index_ == -1)
109 if (present_index_ == -2)
125 return present_index_;
131 element_ = element.element_;
132 present_index_ = element.present_index_;
140 present_index_ = present_index;
168 return element_->center()(0);
178 return element_->face(0)->center()(0);
188 return element_->face(1)->center()(0);
198 return element_->diameter();
211 element_->get_dof_indices(local_dof_indices);
222 if (present_index_ < static_cast<int>(element_->get_triangulation().n_active_cells()))
224 element_ = element_->neighbor(1);
244 if (present_index_ < static_cast<int>(element_->get_triangulation().n_active_cells()))
246 element_ = element_->neighbor(1);
263 if (present_index_ >= 0)
265 element_ = element_->neighbor(0);
283 if (present_index_ >= 0)
285 element_ = element_->neighbor(0);
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
const active_cell_it & getelement_() const
Definition: timeiterator.h:202
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
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