]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Deprecate old functions and use the new one
authorMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Tue, 7 Apr 2020 07:17:25 +0000 (09:17 +0200)
committerMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Wed, 8 Apr 2020 10:25:31 +0000 (12:25 +0200)
include/deal.II/fe/fe_tools.h
include/deal.II/fe/fe_tools.templates.h
include/deal.II/matrix_free/cuda_hanging_nodes_internal.h
source/fe/fe_bernstein.cc
source/fe/fe_trace.cc

index 860f58817cb22a6a8b772bef932b25698a41dcfb..1872d6f2ae2ab5588bd7d936b58b201f010c9a7b 100644 (file)
@@ -895,75 +895,72 @@ namespace FETools
    * in y-direction and finally in z-direction. Discontinuous elements of
    * class FE_DGQ() are numbered in this way, for example.
    *
-   * This function constructs a table which lexicographic index each degree of
-   * freedom in the hierarchic numbering would have. It operates on the
-   * continuous finite element given as first argument, and outputs the
-   * lexicographic indices in the second.
-   *
-   * Note that since this function uses specifics of the continuous finite
-   * elements, it can only operate on FiniteElementData<dim> objects inherent
-   * in FE_Q(). However, this function does not take a FE_Q object as it is
-   * also invoked by the FE_Q() constructor.
-   *
-   * It is assumed that the size of the output argument already matches the
-   * correct size, which is equal to the number of degrees of freedom in the
-   * finite element.
+   * This function returns a vector containing information about the
+   * lexicographic index each degree of freedom in the hierarchic numbering
+   * would have to a given degree of a continuous finite element.
    */
   template <int dim>
-  void
-  hierarchic_to_lexicographic_numbering(unsigned int               degree,
-                                        std::vector<unsigned int> &h2l);
+  std::vector<unsigned int>
+  hierarchic_to_lexicographic_numbering(unsigned int degree);
 
   /**
-   * Like the previous function but instead of returning its result through
-   * the last argument return it as a value.
+   * Like the previous function but instead of returning its result as a value
+   * return it through the last argument.
+   *
+   * @deprecated Use the function that returns the renumbering in a vector
+   * instead.
    */
   template <int dim>
-  std::vector<unsigned int>
-  hierarchic_to_lexicographic_numbering(unsigned int degree);
+  DEAL_II_DEPRECATED void
+  hierarchic_to_lexicographic_numbering(unsigned int               degree,
+                                        std::vector<unsigned int> &h2l);
 
   /**
    * Like the previous functions but using a FiniteElementData instead of the
    * polynomial degree.
+   *
+   * @deprecated Use the function that returns the renumbering in a vector and
+   * uses the degree of the basis as an argument instead.
    */
   template <int dim>
-  void
+  DEAL_II_DEPRECATED void
   hierarchic_to_lexicographic_numbering(const FiniteElementData<dim> &fe_data,
                                         std::vector<unsigned int> &   h2l);
 
   /**
-   * Like the previous function but instead of returning its result through
-   * the last argument return it as a value.
+   * @deprecated Use the function that uses the degree of the basis as an
+   * argument instead.
    */
   template <int dim>
-  std::vector<unsigned int>
-  hierarchic_to_lexicographic_numbering(const FiniteElementData<dim> &fe_data);
+  DEAL_II_DEPRECATED std::vector<unsigned int>
+                     hierarchic_to_lexicographic_numbering(const FiniteElementData<dim> &fe_data);
 
   /**
    * This is the reverse function to the above one, generating the map from
-   * the lexicographic to the hierarchical numbering. All the remarks made
-   * about the above function are also valid here.
+   * the lexicographic to the hierarchical numbering for a given polynomial
+   * degree of a continuous finite element. All the remarks made about the
+   * above function are also valid here.
    */
   template <int dim>
   std::vector<unsigned int>
   lexicographic_to_hierarchic_numbering(unsigned int degree);
 
   /**
-   * Like the previous functions but using a FiniteElementData instead of the
-   * polynomial degree.
+   * @deprecated Use the function that returns the renumbering in a vector and
+   * uses the degree of the basis as an argument instead.
    */
   template <int dim>
-  void
+  DEAL_II_DEPRECATED void
   lexicographic_to_hierarchic_numbering(const FiniteElementData<dim> &fe_data,
                                         std::vector<unsigned int> &   l2h);
 
   /**
-   * Like the previous function but instead of returning its result through
-   * the last argument return it as a value.
+   * @deprecated Use the function that uses the degree of the basis as an
+   * argument instead.
    */
   template <int dim>
-  std::vector<unsigned int>
-  lexicographic_to_hierarchic_numbering(const FiniteElementData<dim> &fe_data);
+  DEAL_II_DEPRECATED std::vector<unsigned int>
+                     lexicographic_to_hierarchic_numbering(const FiniteElementData<dim> &fe_data);
 
 
   /**
index b4a594fa3aefb0839d9e682a609eca82b191a342..a20409b117aaf7f2292858b1c18e8c921aa38350 100644 (file)
@@ -3123,17 +3123,15 @@ namespace FETools
 
 
   template <int dim>
-  void
-  hierarchic_to_lexicographic_numbering(const unsigned int         degree,
-                                        std::vector<unsigned int> &h2l)
+  std::vector<unsigned int>
+  hierarchic_to_lexicographic_numbering(const unsigned int degree)
   {
     // number of support points in each direction
     const unsigned int n = degree + 1;
 
     const unsigned int dofs_per_cell = Utilities::fixed_power<dim>(n);
 
-    // Assert size matches degree
-    AssertDimension(h2l.size(), dofs_per_cell);
+    std::vector<unsigned int> h2l(dofs_per_cell);
 
     // polynomial degree
     const unsigned int dofs_per_line = degree - 1;
@@ -3289,17 +3287,19 @@ namespace FETools
         default:
           Assert(false, ExcNotImplemented());
       }
+
+    return h2l;
   }
 
 
 
   template <int dim>
-  std::vector<unsigned int>
-  hierarchic_to_lexicographic_numbering(const unsigned int degree)
+  void
+  hierarchic_to_lexicographic_numbering(const unsigned int         degree,
+                                        std::vector<unsigned int> &h2l)
   {
-    std::vector<unsigned int> h2l(Utilities::pow(degree + 1, dim));
-    hierarchic_to_lexicographic_numbering<dim>(degree, h2l);
-    return h2l;
+    AssertDimension(h2l.size(), Utilities::fixed_power<dim>(degree + 1));
+    h2l = hierarchic_to_lexicographic_numbering<dim>(degree);
   }
 
 
@@ -3321,9 +3321,7 @@ namespace FETools
   hierarchic_to_lexicographic_numbering(const FiniteElementData<dim> &fe)
   {
     Assert(fe.n_components() == 1, ExcInvalidFE());
-    std::vector<unsigned int> h2l(fe.dofs_per_cell);
-    hierarchic_to_lexicographic_numbering<dim>(fe.dofs_per_line + 1, h2l);
-    return h2l;
+    return hierarchic_to_lexicographic_numbering<dim>(fe.dofs_per_line + 1);
   }
 
 
index a4c87d51d1938bcede6f279d875266438eb26908..605313d40464857418f1cf13c2b40d1ed9dca66a 100644 (file)
@@ -147,7 +147,7 @@ namespace CUDAWrappers
         fe_q.get_subface_interpolation_matrix(fe_q, 0, interpolation_matrix);
 
         std::vector<unsigned int> mapping =
-          FETools::lexicographic_to_hierarchic_numbering<1>(FE_Q<1>(fe_degree));
+          FETools::lexicographic_to_hierarchic_numbering<1>(fe_degree);
 
         FullMatrix<double> mapped_matrix(fe_q.dofs_per_face,
                                          fe_q.dofs_per_face);
@@ -282,8 +282,7 @@ namespace CUDAWrappers
       std::vector<types::global_dof_index> neighbor_dofs(dofs_per_face);
 
       const auto lex_face_mapping =
-        FETools::lexicographic_to_hierarchic_numbering<dim - 1>(
-          FE_Q<dim - 1>(fe_degree));
+        FETools::lexicographic_to_hierarchic_numbering<dim - 1>(fe_degree);
 
       for (const unsigned int face : GeometryInfo<dim>::face_indices())
         {
index 99294631a2eeac5d7290dffc59eefefffdb0f5c3..db2860bf1d046058172886b2c52458afdcaf526f 100644 (file)
@@ -367,10 +367,7 @@ FE_Bernstein<dim, spacedim>::renumber_bases(const unsigned int deg)
 {
   TensorProductPolynomials<dim> tpp(
     dealii::generate_complete_bernstein_basis<double>(deg));
-  std::vector<unsigned int>    renumber(Utilities::fixed_power<dim>(deg + 1));
-  const FiniteElementData<dim> fe(this->get_dpo_vector(deg), 1, deg);
-  FETools::hierarchic_to_lexicographic_numbering(fe, renumber);
-  tpp.set_numbering(renumber);
+  tpp.set_numbering(FETools::hierarchic_to_lexicographic_numbering<dim>(deg));
   return tpp;
 }
 
index f844990475593eee829796312fac3af45d099f67..d40fef3a92f326ab133d7c036e42b40a727dbc94 100644 (file)
@@ -50,9 +50,8 @@ FE_TraceQ<dim, spacedim>::FE_TraceQ(const unsigned int degree)
   Assert(degree > 0,
          ExcMessage("FE_Trace can only be used for polynomial degrees "
                     "greater than zero"));
-  std::vector<unsigned int> renumber(this->dofs_per_face);
-  FETools::hierarchic_to_lexicographic_numbering<dim - 1>(degree, renumber);
-  this->poly_space.set_numbering(renumber);
+  this->poly_space.set_numbering(
+    FETools::hierarchic_to_lexicographic_numbering<dim - 1>(degree));
 
   // Initialize face support points
   this->unit_face_support_points = fe_q.get_unit_face_support_points();

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.