]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove DoFHandler's initialize and set_fe 14047/head
authorDaniel Arndt <arndtd@ornl.gov>
Fri, 24 Jun 2022 17:46:24 +0000 (13:46 -0400)
committerDaniel Arndt <arndtd@ornl.gov>
Fri, 24 Jun 2022 18:31:36 +0000 (14:31 -0400)
doc/news/changes/incompatibilities/20220624Arndt-2 [new file with mode: 0644]
include/deal.II/dofs/dof_handler.h
source/dofs/dof_handler.cc

diff --git a/doc/news/changes/incompatibilities/20220624Arndt-2 b/doc/news/changes/incompatibilities/20220624Arndt-2
new file mode 100644 (file)
index 0000000..0276fa2
--- /dev/null
@@ -0,0 +1,4 @@
+Removed: The deprecated member functions DoFHandler::initialize() and
+DoFHandler::set_fe() have been removed.
+<br>
+(Daniel Arndt, 2022/06/24)
index f850febf210e12192479e147b52367a789a5cd8b..87926405a4c696cafe7a78bf3c01e85558c2698b 100644 (file)
@@ -573,63 +573,6 @@ public:
   DoFHandler &
   operator=(const DoFHandler &) = delete;
 
-  /**
-   * Assign a Triangulation and a FiniteElement to the DoFHandler and compute
-   * the distribution of degrees of freedom over the mesh.
-   *
-   * @deprecated Use reinit() and distribute_dofs() instead.
-   */
-  DEAL_II_DEPRECATED
-  void
-  initialize(const Triangulation<dim, spacedim> &tria,
-             const FiniteElement<dim, spacedim> &fe);
-
-  /**
-   * Same as above but taking an hp::FECollection object.
-   *
-   * @deprecated Use reinit() and distribute_dofs() instead.
-   */
-  DEAL_II_DEPRECATED
-  void
-  initialize(const Triangulation<dim, spacedim> &   tria,
-             const hp::FECollection<dim, spacedim> &fe);
-
-  /**
-   * Assign a FiniteElement @p fe to this object.
-   *
-   * @note This function makes a copy of the finite element given as
-   * argument, and stores it as a member variable. Consequently, it is
-   * possible to write code such as
-   * @code
-   *   dof_handler.set_fe(FE_Q<dim>(2));
-   * @endcode
-   * You can then access the finite element later on by calling
-   * DoFHandler::get_fe(). However, it is often more convenient to
-   * keep a named finite element object as a member variable in your
-   * main class and refer to it directly whenever you need to access
-   * properties of the finite element (such as
-   * FiniteElementData::dofs_per_cell). This is what all tutorial programs do.
-   *
-   * @warning This function only sets a FiniteElement. Degrees of freedom have
-   * either not been distributed yet, or are distributed using a previously set
-   * element. In both cases, accessing degrees of freedom will lead to invalid
-   * results. To restore consistency, call distribute_dofs().
-   *
-   * @deprecated Use distribute_dofs() instead.
-   */
-  DEAL_II_DEPRECATED
-  void
-  set_fe(const FiniteElement<dim, spacedim> &fe);
-
-  /**
-   * Same as above but taking an hp::FECollection object.
-   *
-   * @deprecated Use distribute_dofs() instead.
-   */
-  DEAL_II_DEPRECATED
-  void
-  set_fe(const hp::FECollection<dim, spacedim> &fe);
-
   /**
    * Go through the triangulation and set the active FE indices of all
    * active cells to the values given in @p active_fe_indices.
index 5726e3f1d16e91406277845f61a38c824111437e..6a3093f98c3ef8cb5f9a5638088eee5816eed6e2 100644 (file)
@@ -1807,27 +1807,6 @@ DoFHandler<dim, spacedim>::~DoFHandler()
 
 
 
-template <int dim, int spacedim>
-void
-DoFHandler<dim, spacedim>::initialize(const Triangulation<dim, spacedim> &tria,
-                                      const FiniteElement<dim, spacedim> &fe)
-{
-  this->initialize(tria, hp::FECollection<dim, spacedim>(fe));
-}
-
-
-
-template <int dim, int spacedim>
-void
-DoFHandler<dim, spacedim>::initialize(const Triangulation<dim, spacedim> &tria,
-                                      const hp::FECollection<dim, spacedim> &fe)
-{
-  this->reinit(tria);
-  this->distribute_dofs(fe);
-}
-
-
-
 template <int dim, int spacedim>
 void
 DoFHandler<dim, spacedim>::reinit(const Triangulation<dim, spacedim> &tria)
@@ -2189,81 +2168,6 @@ DoFHandler<dim, spacedim>::memory_consumption() const
 
 
 
-template <int dim, int spacedim>
-void
-DoFHandler<dim, spacedim>::set_fe(const FiniteElement<dim, spacedim> &fe)
-{
-  this->set_fe(hp::FECollection<dim, spacedim>(fe));
-}
-
-
-
-template <int dim, int spacedim>
-void
-DoFHandler<dim, spacedim>::set_fe(const hp::FECollection<dim, spacedim> &ff)
-{
-  Assert(
-    this->tria != nullptr,
-    ExcMessage(
-      "You need to set the Triangulation in the DoFHandler using reinit() or "
-      "in the constructor before you can distribute DoFs."));
-  Assert(this->tria->n_levels() > 0,
-         ExcMessage("The Triangulation you are using is empty!"));
-  Assert(ff.size() > 0, ExcMessage("The hp::FECollection given is empty!"));
-
-  // don't create a new object if the one we have is already appropriate
-  if (this->fe_collection != ff)
-    {
-      this->fe_collection = hp::FECollection<dim, spacedim>(ff);
-
-      const bool contains_multiple_fes = (this->fe_collection.size() > 1);
-
-      // disable hp-mode if only a single finite element has been registered
-      if (hp_capability_enabled && !contains_multiple_fes)
-        {
-          hp_capability_enabled = false;
-
-          // 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_for_transfer.clear();
-
-          // release active and future finite element tables
-          this->hp_cell_active_fe_indices.clear();
-          this->hp_cell_active_fe_indices.shrink_to_fit();
-          this->hp_cell_future_fe_indices.clear();
-          this->hp_cell_future_fe_indices.shrink_to_fit();
-        }
-
-      // re-enabling hp-mode is not permitted since the active and future FE
-      // tables are no longer available
-      AssertThrow(
-        hp_capability_enabled || !contains_multiple_fes,
-        ExcMessage(
-          "You cannot re-enable hp-capabilities after you registered a single "
-          "finite element. Please create a new DoFHandler object instead."));
-    }
-
-  if (hp_capability_enabled)
-    {
-      // make sure every processor knows the active FE indices
-      // on both its own cells and all ghost cells
-      dealii::internal::hp::DoFHandlerImplementation::Implementation::
-        communicate_active_fe_indices(*this);
-
-      // make sure that the FE collection is large enough to
-      // cover all FE indices presently in use on the mesh
-      for (const auto &cell : this->active_cell_iterators())
-        if (!cell->is_artificial())
-          Assert(cell->active_fe_index() < this->fe_collection.size(),
-                 ExcInvalidFEIndex(cell->active_fe_index(),
-                                   this->fe_collection.size()));
-    }
-}
-
-
-
 template <int dim, int spacedim>
 void
 DoFHandler<dim, spacedim>::distribute_dofs(

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.