template <typename CellIterator>
struct PeriodicFacePair;
}
+
+namespace internal
+{
+ namespace parallel
+ {
+ namespace distributed
+ {
+ template <int, int>
+ class TemporarilyMatchRefineFlags;
+ }
+ } // namespace parallel
+} // namespace internal
# endif
namespace parallel
template <int, int, class>
friend class dealii::FETools::internal::ExtrapolateImplementation;
+
+ template <int, int>
+ friend class dealii::internal::parallel::distributed::
+ TemporarilyMatchRefineFlags;
};
virtual types::coarse_cell_id
coarse_cell_index_to_coarse_cell_id(
const unsigned int coarse_cell_index) const override;
+
+ template <int, int>
+ friend class dealii::internal::parallel::distributed::
+ TemporarilyMatchRefineFlags;
};
} // namespace distributed
} // namespace parallel
*/
template <int dim, int spacedim = dim>
class Triangulation
- : public dealii::parallel::TriangulationBase<dim, spacedim>
+ : public dealii::parallel::DistributedTriangulationBase<dim, spacedim>
{
public:
/**
#endif
+namespace internal
+{
+ namespace parallel
+ {
+ namespace distributed
+ {
+ /**
+ * This class temporarily modifies the refine and coarsen flags of all
+ * active cells to match the p4est oracle.
+ *
+ * The modification only happens on parallel::distributed::Triangulation
+ * objects, and persists for the lifetime of an instantiation of this
+ * class.
+ *
+ * The TemporarilyMatchRefineFlags class should only be used in
+ * combination with the Triangulation::Signals::post_p4est_refinement
+ * signal. At this stage, the p4est orcale already has been refined, but
+ * the triangulation is still unchanged. After the modification, all
+ * refine and coarsen flags describe how the traingulation will acutally
+ * be refined.
+ */
+ template <int dim, int spacedim = dim>
+ class TemporarilyMatchRefineFlags : public Subscriptor
+ {
+ public:
+ /**
+ * Constructor.
+ *
+ * Stores the refine and coarsen flags of all active cells if the
+ * provided Triangulation is of type
+ * parallel::distributed::Triangulation.
+ *
+ * Adjusts them to be consistent with the p4est oracle.
+ */
+ TemporarilyMatchRefineFlags(Triangulation<dim, spacedim> &tria);
+
+ /**
+ * Destructor.
+ *
+ * Returns the refine and coarsen flags of all active cells on the
+ * parallel::distributed::Triangulation into their previous state.
+ */
+ ~TemporarilyMatchRefineFlags();
+
+ private:
+ /**
+ * The modified parallel::distributed::Triangulation.
+ */
+ const SmartPointer<
+ dealii::parallel::distributed::Triangulation<dim, spacedim>>
+ distributed_tria;
+
+ /**
+ * A vector that temporarily stores the refine flags before they have
+ * been modified on the parallel::distributed::Triangulation.
+ */
+ std::vector<bool> saved_refine_flags;
+
+ /**
+ * A vector that temporarily stores the coarsen flags before they have
+ * been modified on the parallel::distributed::Triangulation.
+ */
+ std::vector<bool> saved_coarsen_flags;
+ };
+ } // namespace distributed
+ } // namespace parallel
+} // namespace internal
+
+
DEAL_II_NAMESPACE_CLOSE
#endif
+namespace internal
+{
+ namespace parallel
+ {
+ namespace distributed
+ {
+ template <int dim, int spacedim>
+ TemporarilyMatchRefineFlags<dim, spacedim>::TemporarilyMatchRefineFlags(
+ Triangulation<dim, spacedim> &tria)
+ : distributed_tria(
+ dynamic_cast<
+ dealii::parallel::distributed::Triangulation<dim, spacedim> *>(
+ &tria))
+ {
+#ifdef DEAL_II_WITH_P4EST
+ if (distributed_tria != nullptr)
+ {
+ // Save the current set of refinement flags, and adjust the
+ // refinement flags to be consistent with the p4est oracle.
+ distributed_tria->save_coarsen_flags(saved_coarsen_flags);
+ distributed_tria->save_refine_flags(saved_refine_flags);
+
+ for (const auto &pair : distributed_tria->local_cell_relations)
+ {
+ const auto &cell = pair.first;
+ const auto &status = pair.second;
+
+ switch (status)
+ {
+ case dealii::Triangulation<dim, spacedim>::CELL_PERSIST:
+ // cell remains unchanged
+ cell->clear_refine_flag();
+ cell->clear_coarsen_flag();
+ break;
+
+ case dealii::Triangulation<dim, spacedim>::CELL_REFINE:
+ // cell will be refined
+ cell->clear_coarsen_flag();
+ cell->set_refine_flag();
+ break;
+
+ case dealii::Triangulation<dim, spacedim>::CELL_COARSEN:
+ // children of this cell will be coarsened
+ for (const auto &child : cell->child_iterators())
+ {
+ child->clear_refine_flag();
+ child->set_coarsen_flag();
+ }
+ break;
+
+ case dealii::Triangulation<dim, spacedim>::CELL_INVALID:
+ // do nothing as cell does not exist yet
+ break;
+
+ default:
+ Assert(false, ExcInternalError());
+ break;
+ }
+ }
+ }
+#endif
+ }
+
+
+
+ template <int dim, int spacedim>
+ TemporarilyMatchRefineFlags<dim, spacedim>::~TemporarilyMatchRefineFlags()
+ {
+#ifdef DEAL_II_WITH_P4EST
+ if (distributed_tria)
+ {
+ // Undo the refinement flags modification.
+ distributed_tria->load_coarsen_flags(saved_coarsen_flags);
+ distributed_tria->load_refine_flags(saved_refine_flags);
+ }
+#else
+ // pretend that this destructor does something to silence clang-tidy
+ (void)distributed_tria;
+#endif
+ }
+ } // namespace distributed
+ } // namespace parallel
+} // namespace internal
+
+
+
/*-------------- Explicit Instantiations -------------------------------*/
#include "tria.inst"
-for (deal_II_dimension : DIMENSIONS)
+for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
{
namespace parallel
\{
namespace distributed
\{
- template class Triangulation<deal_II_dimension>;
-#if deal_II_dimension < 3
- template class Triangulation<deal_II_dimension, deal_II_dimension + 1>;
+#if deal_II_dimension <= deal_II_space_dimension
+ template class Triangulation<deal_II_dimension,
+ deal_II_space_dimension>;
#endif
-#if deal_II_dimension < 2
- template class Triangulation<deal_II_dimension, deal_II_dimension + 2>;
+ \}
+ \}
+
+ namespace internal
+ \{
+ namespace parallel
+ \{
+ namespace distributed
+ \{
+#if deal_II_dimension <= deal_II_space_dimension
+ template class TemporarilyMatchRefineFlags<deal_II_dimension,
+ deal_II_space_dimension>;
#endif
+ \}
\}
\}
}
if (cell->is_locally_owned() && cell->center() == Point<dim>())
cell->set_active_fe_index(sequence.back());
- const bool fe_indices_changed =
- hp::Refinement::limit_p_level_difference(dofh,
- max_difference,
- contains_fe_index);
+ bool fe_indices_changed = false;
+ tria.signals.post_p4est_refinement.connect(
+ [&]() {
+ const internal::parallel::distributed::TemporarilyMatchRefineFlags<dim>
+ refine_modifier(tria);
+ fe_indices_changed =
+ hp::Refinement::limit_p_level_difference(dofh,
+ max_difference,
+ contains_fe_index);
+ },
+ boost::signals2::at_front);
+
tria.execute_coarsening_and_refinement();
- (void)fe_indices_changed;
Assert(fe_indices_changed, ExcInternalError());
// display number of cells for each FE index
DoFHandler<dim> dofh(tria);
dofh.distribute_dofs(fes);
+ bool fe_indices_changed;
+ tria.signals.post_p4est_refinement.connect(
+ [&]() {
+ const internal::parallel::distributed::TemporarilyMatchRefineFlags<dim>
+ refine_modifier(tria);
+ fe_indices_changed =
+ hp::Refinement::limit_p_level_difference(dofh,
+ max_difference,
+ contains_fe_index);
+ },
+ boost::signals2::at_front);
+
// increase p-level in center subsequently
for (unsigned int i = 0; i < sequence.size() - 1; ++i)
{
cell->set_future_fe_index(
fes.next_in_hierarchy(cell->active_fe_index()));
- const bool fe_indices_changed =
- hp::Refinement::limit_p_level_difference(dofh,
- max_difference,
- contains_fe_index);
+ fe_indices_changed = false;
tria.execute_coarsening_and_refinement();
- (void)fe_indices_changed;
if (i >= max_difference)
Assert(fe_indices_changed, ExcInternalError());
cell->set_active_fe_index(sequence.back());
dofh.distribute_dofs(fes);
- const bool fe_indices_changed =
- hp::Refinement::limit_p_level_difference(dofh,
- max_difference,
- contains_fe_index);
+ bool fe_indices_changed = false;
+ tria.signals.post_p4est_refinement.connect(
+ [&]() {
+ const internal::parallel::distributed::TemporarilyMatchRefineFlags<dim>
+ refine_modifier(tria);
+ fe_indices_changed =
+ hp::Refinement::limit_p_level_difference(dofh,
+ max_difference,
+ contains_fe_index);
+ },
+ boost::signals2::at_front);
+
tria.execute_coarsening_and_refinement();
- (void)fe_indices_changed;
Assert(fe_indices_changed, ExcInternalError());
#ifdef DEBUG
cell->set_active_fe_index(sequence.back());
dofh.distribute_dofs(fes);
- const bool fe_indices_changed =
- hp::Refinement::limit_p_level_difference(dofh,
- max_difference,
- contains_fe_index);
- (void)fe_indices_changed;
- Assert(fe_indices_changed, ExcInternalError());
+ bool fe_indices_changed = false;
+ tria.signals.post_p4est_refinement.connect(
+ [&]() {
+ const internal::parallel::distributed::TemporarilyMatchRefineFlags<dim>
+ refine_modifier(tria);
+ fe_indices_changed =
+ hp::Refinement::limit_p_level_difference(dofh,
+ max_difference,
+ contains_fe_index);
+ },
+ boost::signals2::at_front);
tria.execute_coarsening_and_refinement();
+ Assert(fe_indices_changed, ExcInternalError());
+
deallog << "active FE indices after adaptation:" << std::endl;
for (const auto &cell : dofh.active_cell_iterators())
if (cell->is_locally_owned())