From: Wolfgang Bangerth Date: Fri, 15 Dec 2023 05:04:07 +0000 (-0700) Subject: Use std::variant instead of a hand-rolled version. X-Git-Tag: relicensing~240^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F16355%2Fhead;p=dealii.git Use std::variant instead of a hand-rolled version. --- diff --git a/include/deal.II/fe/fe_values_base.h b/include/deal.II/fe/fe_values_base.h index 3d56f8b895..7fe4802204 100644 --- a/include/deal.II/fe/fe_values_base.h +++ b/include/deal.II/fe/fe_values_base.h @@ -48,6 +48,7 @@ #include #include #include +#include DEAL_II_NAMESPACE_OPEN @@ -1571,23 +1572,29 @@ protected: "FEValues::reinit() with an iterator type that allows to extract " "degrees of freedom, such as DoFHandler::cell_iterator."); + /** + * Constructor. Creates an unusable object that is not associated with + * any cell at all. + */ + CellIteratorContainer() = default; + /** * Constructor. */ - CellIteratorContainer(); + CellIteratorContainer( + const typename Triangulation::cell_iterator &cell); /** * Constructor. */ - template CellIteratorContainer( - const TriaIterator> &cell); + const typename DoFHandler::cell_iterator &cell); /** * Constructor. */ CellIteratorContainer( - const typename Triangulation::cell_iterator &cell); + const typename DoFHandler::level_cell_iterator &cell); /** * Indicate whether FEValues::reinit() was called. @@ -1628,9 +1635,16 @@ protected: Vector &out) const; private: - std::optional::cell_iterator> cell; - const DoFHandler *dof_handler; - bool level_dof_access; + /** + * The cell in question, if one has been assigned to this object. The + * concrete data type can either be a Triangulation cell iterator, a + * DoFHandler cell iterator, or a DoFHandler level cell iterator. + */ + std::optional< + std::variant::cell_iterator, + typename DoFHandler::cell_iterator, + typename DoFHandler::level_cell_iterator>> + cell; }; /** @@ -1784,18 +1798,6 @@ private: /*---------------------- Inline functions: FEValuesBase ---------------------*/ -template -template -inline FEValuesBase::CellIteratorContainer:: - CellIteratorContainer( - const TriaIterator> &cell) - : cell(cell) - , dof_handler(&cell->get_dof_handler()) - , level_dof_access(lda) -{} - - - template inline const FEValuesViews::Scalar & FEValuesBase::operator[]( diff --git a/source/fe/fe_values_base.cc b/source/fe/fe_values_base.cc index 6120a4f9b0..c05c23618b 100644 --- a/source/fe/fe_values_base.cc +++ b/source/fe/fe_values_base.cc @@ -107,21 +107,27 @@ namespace internal /* ------------ FEValuesBase::CellIteratorContainer ----------- */ + template -FEValuesBase::CellIteratorContainer::CellIteratorContainer() - : cell() - , dof_handler(nullptr) - , level_dof_access(false) +FEValuesBase::CellIteratorContainer::CellIteratorContainer( + const typename Triangulation::cell_iterator &cell) + : cell(cell) {} template FEValuesBase::CellIteratorContainer::CellIteratorContainer( - const typename Triangulation::cell_iterator &cell) + const typename DoFHandler::cell_iterator &cell) + : cell(cell) +{} + + + +template +FEValuesBase::CellIteratorContainer::CellIteratorContainer( + const typename DoFHandler::level_cell_iterator &cell) : cell(cell) - , dof_handler(nullptr) - , level_dof_access(false) {} @@ -141,7 +147,14 @@ operator typename Triangulation::cell_iterator() const { Assert(is_initialized(), ExcNotReinited()); - return cell.value(); + // We can always convert to a tria iterator, regardless of which of + // the three types of cell we store. + return std::visit( + [](auto &cell_iterator) -> + typename Triangulation::cell_iterator { + return cell_iterator; + }, + cell.value()); } @@ -152,9 +165,17 @@ FEValuesBase::CellIteratorContainer::n_dofs_for_dof_handler() const { Assert(is_initialized(), ExcNotReinited()); - Assert(dof_handler != nullptr, ExcNeedsDoFHandler()); - return dof_handler->n_dofs(); + switch (cell.value().index()) + { + case 1: + return std::get<1>(cell.value())->get_dof_handler().n_dofs(); + case 2: + return std::get<2>(cell.value())->get_dof_handler().n_dofs(); + default: + Assert(false, ExcNeedsDoFHandler()); + return numbers::invalid_dof_index; + } } @@ -167,20 +188,21 @@ FEValuesBase::CellIteratorContainer::get_interpolated_dof_values( Vector &out) const { Assert(is_initialized(), ExcNotReinited()); - Assert(dof_handler != nullptr, ExcNeedsDoFHandler()); - - if (level_dof_access) - DoFCellAccessor(&cell.value()->get_triangulation(), - cell.value()->level(), - cell.value()->index(), - dof_handler) - .get_interpolated_dof_values(in, out); - else - DoFCellAccessor(&cell.value()->get_triangulation(), - cell.value()->level(), - cell.value()->index(), - dof_handler) - .get_interpolated_dof_values(in, out); + + switch (cell.value().index()) + { + case 1: + std::get<1>(cell.value())->get_interpolated_dof_values(in, out); + break; + + case 2: + std::get<2>(cell.value())->get_interpolated_dof_values(in, out); + break; + + default: + Assert(false, ExcNeedsDoFHandler()); + break; + } } @@ -192,21 +214,32 @@ FEValuesBase::CellIteratorContainer::get_interpolated_dof_values( Vector &out) const { Assert(is_initialized(), ExcNotReinited()); - Assert(dof_handler != nullptr, ExcNeedsDoFHandler()); - Assert(level_dof_access == false, ExcNotImplemented()); - const DoFCellAccessor cell_dofs( - &cell.value()->get_triangulation(), - cell.value()->level(), - cell.value()->index(), - dof_handler); + switch (cell.value().index()) + { + case 1: + { + const typename DoFHandler::cell_iterator cell = + std::get<1>(this->cell.value()); + + std::vector dof_indices( + cell->get_fe().n_dofs_per_cell()); + + cell->get_dof_indices(dof_indices); - std::vector dof_indices( - cell_dofs.get_fe().n_dofs_per_cell()); - cell_dofs.get_dof_indices(dof_indices); + for (unsigned int i = 0; i < cell->get_fe().n_dofs_per_cell(); ++i) + out[i] = (in.is_element(dof_indices[i]) ? 1 : 0); - for (unsigned int i = 0; i < cell_dofs.get_fe().n_dofs_per_cell(); ++i) - out[i] = (in.is_element(dof_indices[i]) ? 1 : 0); + break; + } + case 2: + Assert(false, ExcNotImplemented()); + break; + + default: + Assert(false, ExcNeedsDoFHandler()); + break; + } }