]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Deprecate MatrixFree::get_dof_handler. 11181/head
authorMarc Fehling <mafehling.git@gmail.com>
Tue, 10 Nov 2020 06:55:08 +0000 (23:55 -0700)
committerMarc Fehling <mafehling.git@gmail.com>
Sat, 14 Nov 2020 06:27:22 +0000 (23:27 -0700)
Removing the instantiation will throw an error if someone writes
matrix_free.template get_dof_handler<DoFHandler<dim>>();
which will probably never happen as the default template is way more convenient.

doc/news/changes/incompatibilities/20201113Fehling [new file with mode: 0644]
include/deal.II/matrix_free/matrix_free.h
include/deal.II/matrix_free/matrix_free.templates.h
source/matrix_free/matrix_free.inst.in

diff --git a/doc/news/changes/incompatibilities/20201113Fehling b/doc/news/changes/incompatibilities/20201113Fehling
new file mode 100644 (file)
index 0000000..75bff19
--- /dev/null
@@ -0,0 +1,5 @@
+Deprecated: The variant of the function MatrixFree::get_dof_handler
+expecting a DoFHandlerType template has been deprecated. Use the 
+template-less variant returning a DoFHandler instead.
+<br>
+(Marc Fehling, 2020/11/13)
index 25ad250d0f337e8595262ef459fd5b1bf4803005..b4f1b3a8219bb6fb21d38805f4d40ba714fe6b85 100644 (file)
@@ -1689,6 +1689,13 @@ public:
   get_faces_by_cells_boundary_id(const unsigned int macro_cell,
                                  const unsigned int face_number) const;
 
+  /**
+   * Return the DoFHandler with the index as given to the respective
+   * `std::vector` argument in the reinit() function.
+   */
+  const DoFHandler<dim> &
+  get_dof_handler(const unsigned int dof_handler_index = 0) const;
+
   /**
    * Return the DoFHandler with the index as given to the respective
    * `std::vector` argument in the reinit() function. Note that if you want to
@@ -1696,10 +1703,12 @@ public:
    * one, you will need to use the `template` before the function call, i.e.,
    * you will have something like `matrix_free.template
    * get_dof_handler<hp::DoFHandler<dim>>()`.
+   *
+   * @deprecated Use the non-templated equivalent of this function.
    */
-  template <typename DoFHandlerType = DoFHandler<dim>>
-  const DoFHandlerType &
-  get_dof_handler(const unsigned int dof_handler_index = 0) const;
+  template <typename DoFHandlerType>
+  DEAL_II_DEPRECATED const DoFHandlerType &
+                           get_dof_handler(const unsigned int dof_handler_index = 0) const;
 
   /**
    * Return the cell iterator in deal.II speak to a given cell in the
index fdd4febe56aeb804801b02e14a1d315634de9c3e..bbb211e0b95f62a3120e61901e62fe0b18a52a60 100644 (file)
@@ -134,6 +134,18 @@ MatrixFree<dim, Number, VectorizedArrayType>::renumber_dofs(
 
 
 
+template <int dim, typename Number, typename VectorizedArrayType>
+const DoFHandler<dim> &
+MatrixFree<dim, Number, VectorizedArrayType>::get_dof_handler(
+  const unsigned int dof_handler_index) const
+{
+  AssertIndexRange(dof_handler_index, n_components());
+
+  return *(dof_handlers[dof_handler_index]);
+}
+
+
+
 template <int dim, typename Number, typename VectorizedArrayType>
 template <typename DoFHandlerType>
 const DoFHandlerType &
index b256093dcb7f438abaea9dd1cc965ceaa3a65fc3..aaacd0001c7697d771d5e2c4e3bc30c1d7c57f99 100644 (file)
@@ -42,12 +42,6 @@ for (deal_II_dimension : DIMENSIONS;
         const std::vector<IndexSet> &,
         const std::vector<hp::QCollection<deal_II_dimension>> &,
         const AdditionalData &);
-
-    template const DoFHandler<deal_II_dimension> &
-    MatrixFree<deal_II_dimension,
-               deal_II_scalar_vectorized::value_type,
-               deal_II_scalar_vectorized>::
-      get_dof_handler<DoFHandler<deal_II_dimension>>(const unsigned int) const;
   }
 
 

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.