From: Marc Fehling Date: Fri, 30 Aug 2019 21:08:10 +0000 (+0200) Subject: Introduced DoFCellAccessor:get_future_fe(). X-Git-Tag: v9.2.0-rc1~1158^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=ed552d14dd0345c360a3d45f9d083cf63758c856;p=dealii.git Introduced DoFCellAccessor:get_future_fe(). --- diff --git a/include/deal.II/dofs/dof_accessor.h b/include/deal.II/dofs/dof_accessor.h index d0458fdfaa..e1bd0a72c3 100644 --- a/include/deal.II/dofs/dof_accessor.h +++ b/include/deal.II/dofs/dof_accessor.h @@ -1946,6 +1946,28 @@ public: * @{ */ + /** + * Return the finite element that will be assigned to this cell next time the + * triangulation gets refined and coarsened. If no future finite element has + * been specified for this cell via the set_future_fe_index() function, the + * active one will remain unchanged, in which case the active finite element + * will be returned. + * + * For non-hp DoF handlers, this is of course always the same element, + * independent of the cell we are presently on, but for hp DoF handlers, this + * may change from cell to cell. + * + * @note Since degrees of freedom only exist on active cells for + * hp::DoFHandler (i.e., there is currently no implementation of multilevel + * hp::DoFHandler objects), it does not make sense to query the finite + * element on non-active cells since they do not have finite element spaces + * associated with them without having any degrees of freedom. Consequently, + * this function will produce an exception when called on non-active cells. + */ + const FiniteElement & + get_future_fe() const; + /** * Return the fe_index of the finite element that will be assigned to this * cell next time the triangulation gets refined and coarsened. If no future diff --git a/include/deal.II/dofs/dof_accessor.templates.h b/include/deal.II/dofs/dof_accessor.templates.h index bb90ee4087..ce84e10cce 100644 --- a/include/deal.II/dofs/dof_accessor.templates.h +++ b/include/deal.II/dofs/dof_accessor.templates.h @@ -3974,6 +3974,7 @@ inline const FiniteElement & DoFCellAccessor::get_fe() const { + Assert(this->dof_handler != nullptr, typename BaseClass::ExcInvalidObject()); Assert( (dynamic_cast *>( @@ -3982,7 +3983,6 @@ DoFCellAccessor::get_fe() const ExcMessage("In hp::DoFHandler objects, finite elements are only associated " "with active cells. Consequently, you can not ask for the " "active finite element on cells with children.")); - Assert(this->dof_handler != nullptr, typename BaseClass::ExcInvalidObject()); return this->dof_handler->get_fe(active_fe_index()); } @@ -4048,6 +4048,26 @@ DoFCellAccessor::set_active_fe_index( +template +inline const FiniteElement & +DoFCellAccessor::get_future_fe() const +{ + Assert(this->dof_handler != nullptr, typename BaseClass::ExcInvalidObject()); + Assert( + (dynamic_cast *>( + this->dof_handler) != nullptr) || + (this->has_children() == false), + ExcMessage("In hp::DoFHandler objects, finite elements are only associated " + "with active cells. Consequently, you can not ask for the " + "future finite element on cells with children.")); + + return this->dof_handler->get_fe(future_fe_index()); +} + + + template inline unsigned int DoFCellAccessor::future_fe_index() const diff --git a/tests/hp/future_fe_indices.cc b/tests/hp/future_fe_indices.cc index 4408536725..6e4ba26084 100644 --- a/tests/hp/future_fe_indices.cc +++ b/tests/hp/future_fe_indices.cc @@ -52,20 +52,15 @@ test() auto cell = dh.begin_active(); cell->set_future_fe_index(1); - // try to set future index to an invalid one - try - { - (++cell)->set_future_fe_index(2); - } - catch (const ExcIndexRange &) - { - deallog << "Set to 2 failed" << std::endl; - } - // verify flags for (const auto &cell : dh.active_cell_iterators()) - deallog << "cell:" << cell->id().to_string() - << ", future_fe:" << cell->future_fe_index() << std::endl; + { + deallog << "cell:" << cell->id().to_string() + << ", future_fe:" << cell->future_fe_index() << std::endl; + Assert(&(dh.get_fe(cell->future_fe_index())) == &(cell->get_future_fe()), + ExcMessage( + "DoFCellAccessor::get_future_fe() returns the wrong object.")); + } // clear all flags and check if all were cleared for (const auto &cell : dh.active_cell_iterators()) diff --git a/tests/hp/future_fe_indices.output b/tests/hp/future_fe_indices.output index d5f4992a09..a77d9096ba 100644 --- a/tests/hp/future_fe_indices.output +++ b/tests/hp/future_fe_indices.output @@ -1,12 +1,12 @@ DEAL:1D::cell:0_1:0, future_fe:1 -DEAL:1D::cell:0_1:1, future_fe:2 +DEAL:1D::cell:0_1:1, future_fe:0 DEAL:2D::cell:0_1:0, future_fe:1 -DEAL:2D::cell:0_1:1, future_fe:2 +DEAL:2D::cell:0_1:1, future_fe:0 DEAL:2D::cell:0_1:2, future_fe:0 DEAL:2D::cell:0_1:3, future_fe:0 DEAL:3D::cell:0_1:0, future_fe:1 -DEAL:3D::cell:0_1:1, future_fe:2 +DEAL:3D::cell:0_1:1, future_fe:0 DEAL:3D::cell:0_1:2, future_fe:0 DEAL:3D::cell:0_1:3, future_fe:0 DEAL:3D::cell:0_1:4, future_fe:0