template <int dim, int spacedim>
DoFHandler<dim, spacedim>::~DoFHandler()
{
- // unsubscribe all attachments to refinement signals of the underlying
- // triangulation
+ // unsubscribe all attachments to signals of the underlying triangulation
for (auto &connection : this->tria_listeners)
connection.disconnect();
this->tria_listeners.clear();
+ for (auto &connection : this->tria_listeners_for_transfer)
+ connection.disconnect();
+ this->tria_listeners_for_transfer.clear();
+
// release allocated memory
// virtual functions called in constructors and destructors never use the
// override in a derived class
connection.disconnect();
this->tria_listeners.clear();
+ for (auto &connection : this->tria_listeners_for_transfer)
+ connection.disconnect();
+ this->tria_listeners_for_transfer.clear();
+
// release allocated memory and policy
DoFHandler<dim, spacedim>::clear();
this->policy.reset();
{
hp_capability_enabled = false;
- // unsubscribe connections to signals
- for (auto &connection : this->tria_listeners)
+ // unsubscribe connections to signals that are only relevant for
+ // hp-mode, since we only have a single element here
+ for (auto &connection : this->tria_listeners_for_transfer)
connection.disconnect();
- this->tria_listeners.clear();
+ this->tria_listeners_for_transfer.clear();
// release active and future finite element tables
this->hp_cell_active_fe_indices.clear();
{
hp_capability_enabled = false;
- // unsubscribe connections to signals
- for (auto &connection : this->tria_listeners)
+ // unsubscribe connections to signals that are only relevant for
+ // hp-mode, since we only have a single element here
+ for (auto &connection : this->tria_listeners_for_transfer)
connection.disconnect();
- this->tria_listeners.clear();
+ this->tria_listeners_for_transfer.clear();
// release active and future finite element tables
this->hp_cell_active_fe_indices.clear();
void
DoFHandler<dim, spacedim>::connect_to_triangulation_signals()
{
+ // make sure this is called during initialization in hp-mode
+ Assert(hp_capability_enabled, ExcOnlyAvailableWithHP());
+
// connect functions to signals of the underlying triangulation
this->tria_listeners.push_back(this->tria->signals.create.connect(
[this]() { this->reinit(*(this->tria)); }));
&this->get_triangulation()))
{
// repartitioning signals
- this->tria_listeners.push_back(
+ this->tria_listeners_for_transfer.push_back(
this->tria->signals.pre_distributed_repartition.connect([this]() {
internal::hp::DoFHandlerImplementation::Implementation::
ensure_absence_of_future_fe_indices<dim, spacedim>(*this);
}));
- this->tria_listeners.push_back(
+ this->tria_listeners_for_transfer.push_back(
this->tria->signals.pre_distributed_repartition.connect(
[this]() { this->pre_distributed_transfer_action(); }));
- this->tria_listeners.push_back(
+ this->tria_listeners_for_transfer.push_back(
this->tria->signals.post_distributed_repartition.connect(
[this] { this->post_distributed_transfer_action(); }));
// refinement signals
- this->tria_listeners.push_back(
+ this->tria_listeners_for_transfer.push_back(
this->tria->signals.pre_distributed_refinement.connect(
[this]() { this->pre_distributed_transfer_action(); }));
- this->tria_listeners.push_back(
+ this->tria_listeners_for_transfer.push_back(
this->tria->signals.post_distributed_refinement.connect(
[this]() { this->post_distributed_transfer_action(); }));
// serialization signals
- this->tria_listeners.push_back(
+ this->tria_listeners_for_transfer.push_back(
this->tria->signals.post_distributed_save.connect(
[this]() { this->active_fe_index_transfer.reset(); }));
- this->tria_listeners.push_back(
+ this->tria_listeners_for_transfer.push_back(
this->tria->signals.post_distributed_load.connect(
[this]() { this->update_active_fe_table(); }));
}
&this->get_triangulation()) != nullptr)
{
// partitioning signals
- this->tria_listeners.push_back(
+ this->tria_listeners_for_transfer.push_back(
this->tria->signals.pre_partition.connect([this]() {
internal::hp::DoFHandlerImplementation::Implementation::
ensure_absence_of_future_fe_indices(*this);
}));
// refinement signals
- this->tria_listeners.push_back(this->tria->signals.pre_refinement.connect(
- [this] { this->pre_transfer_action(); }));
- this->tria_listeners.push_back(
+ this->tria_listeners_for_transfer.push_back(
+ this->tria->signals.pre_refinement.connect(
+ [this] { this->pre_transfer_action(); }));
+ this->tria_listeners_for_transfer.push_back(
this->tria->signals.post_refinement.connect(
[this] { this->post_transfer_action(); }));
}
else
{
// refinement signals
- this->tria_listeners.push_back(this->tria->signals.pre_refinement.connect(
- [this] { this->pre_transfer_action(); }));
- this->tria_listeners.push_back(
+ this->tria_listeners_for_transfer.push_back(
+ this->tria->signals.pre_refinement.connect(
+ [this] { this->pre_transfer_action(); }));
+ this->tria_listeners_for_transfer.push_back(
this->tria->signals.post_refinement.connect(
[this] { this->post_transfer_action(); }));
}
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+
+// check that DoFHandler::clear() will be called
+// whenever Triangulation::Signals::clear is triggered
+
+
+#include <deal.II/distributed/fully_distributed_tria.h>
+#include <deal.II/distributed/shared_tria.h>
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/dofs/dof_handler.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include "../tests.h"
+
+
+void
+test_serial()
+{
+ constexpr const unsigned dim = 2;
+
+ Triangulation<dim> tria;
+ GridGenerator::hyper_cube(tria);
+
+ DoFHandler<dim> dh(tria);
+ dh.distribute_dofs(FE_Q<dim>(1));
+
+ tria.clear();
+ deallog << "ndofs:" << dh.n_dofs() << std::endl;
+}
+
+
+void
+test_parallel_shared()
+{
+ constexpr const unsigned dim = 2;
+
+ parallel::shared::Triangulation<dim> tria(MPI_COMM_WORLD);
+ GridGenerator::hyper_cube(tria);
+
+ DoFHandler<dim> dh(tria);
+ dh.distribute_dofs(FE_Q<dim>(1));
+
+ tria.clear();
+ deallog << "ndofs:" << dh.n_dofs() << std::endl;
+}
+
+
+void
+test_parallel_distributed()
+{
+ constexpr const unsigned dim = 2;
+
+ parallel::distributed::Triangulation<dim> tria(MPI_COMM_WORLD);
+ GridGenerator::hyper_cube(tria);
+
+ DoFHandler<dim> dh(tria);
+ dh.distribute_dofs(FE_Q<dim>(1));
+
+ tria.clear();
+ deallog << "ndofs:" << dh.n_dofs() << std::endl;
+}
+
+
+void
+test_parallel_fullydistributed()
+{
+ constexpr const unsigned dim = 2;
+
+ parallel::distributed::Triangulation<dim> tria_base(MPI_COMM_WORLD);
+ GridGenerator::hyper_cube(tria_base);
+
+ const auto description =
+ TriangulationDescription::Utilities::create_description_from_triangulation(
+ tria_base, MPI_COMM_WORLD);
+
+ parallel::fullydistributed::Triangulation<dim> tria(MPI_COMM_WORLD);
+ tria.create_triangulation(description);
+
+ DoFHandler<dim> dh(tria);
+ dh.distribute_dofs(FE_Q<dim>(1));
+
+ tria.clear();
+ deallog << "ndofs:" << dh.n_dofs() << std::endl;
+}
+
+
+int
+main(int argc, char *argv[])
+{
+ initlog();
+
+ Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+
+ test_serial();
+ test_parallel_shared();
+ test_parallel_distributed();
+ test_parallel_fullydistributed();
+}