]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Use identity<...>::type in the right place. 8987/head
authorWolfgang Bangerth <bangerth@colostate.edu>
Fri, 1 Nov 2019 15:50:54 +0000 (09:50 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Fri, 1 Nov 2019 22:45:27 +0000 (16:45 -0600)
include/deal.II/hp/refinement.h
source/hp/refinement.cc
source/hp/refinement.inst.in

index c87db67890f3701979dc823f36dc79b9a623f987..8cd5ae50ebed91b55a403e49d2fd1501719fb8f5 100644 (file)
@@ -138,8 +138,7 @@ namespace hp
      */
     template <typename Number>
     using ComparisonFunction =
-      std::function<bool(const typename identity<Number>::type &,
-                         const typename identity<Number>::type &)>;
+      std::function<bool(const Number &, const Number &)>;
 
     /**
      * @name Setting p-adaptivity flags
@@ -207,8 +206,10 @@ namespace hp
       const Vector<Number> &               criteria,
       const Number                         p_refine_threshold,
       const Number                         p_coarsen_threshold,
-      const ComparisonFunction<Number> &compare_refine = std::greater<Number>(),
-      const ComparisonFunction<Number> &compare_coarsen = std::less<Number>());
+      const ComparisonFunction<typename identity<Number>::type>
+        &compare_refine = std::greater<Number>(),
+      const ComparisonFunction<typename identity<Number>::type>
+        &compare_coarsen = std::less<Number>());
 
     /**
      * Adapt which finite element to use on cells whose criteria meet a certain
@@ -247,8 +248,10 @@ namespace hp
       const Vector<Number> &               criteria,
       const double                         p_refine_fraction  = 0.5,
       const double                         p_coarsen_fraction = 0.5,
-      const ComparisonFunction<Number> &compare_refine = std::greater<Number>(),
-      const ComparisonFunction<Number> &compare_coarsen = std::less<Number>());
+      const ComparisonFunction<typename identity<Number>::type>
+        &compare_refine = std::greater<Number>(),
+      const ComparisonFunction<typename identity<Number>::type>
+        &compare_coarsen = std::less<Number>());
 
     /**
      * Adapt which finite element to use on cells based on the regularity of the
@@ -317,11 +320,12 @@ namespace hp
     template <int dim, typename Number, int spacedim>
     void
     p_adaptivity_from_reference(
-      const hp::DoFHandler<dim, spacedim> &dof_handler,
-      const Vector<Number> &               criteria,
-      const Vector<Number> &               references,
-      const ComparisonFunction<Number> &   compare_refine,
-      const ComparisonFunction<Number> &   compare_coarsen);
+      const hp::DoFHandler<dim, spacedim> &                      dof_handler,
+      const Vector<Number> &                                     criteria,
+      const Vector<Number> &                                     references,
+      const ComparisonFunction<typename identity<Number>::type> &compare_refine,
+      const ComparisonFunction<typename identity<Number>::type>
+        &compare_coarsen);
     /**
      * @}
      */
index 2d60af474a5524006c52ed3baa77eaeaccaec575..c160347b1ce7167698fdbe336e8066fcc2089e3e 100644 (file)
@@ -91,8 +91,9 @@ namespace hp
       const Vector<Number> &               criteria,
       const Number                         p_refine_threshold,
       const Number                         p_coarsen_threshold,
-      const ComparisonFunction<Number> &   compare_refine,
-      const ComparisonFunction<Number> &   compare_coarsen)
+      const ComparisonFunction<typename identity<Number>::type> &compare_refine,
+      const ComparisonFunction<typename identity<Number>::type>
+        &compare_coarsen)
     {
       AssertDimension(dof_handler.get_triangulation().n_active_cells(),
                       criteria.size());
@@ -122,8 +123,9 @@ namespace hp
       const Vector<Number> &               criteria,
       const double                         p_refine_fraction,
       const double                         p_coarsen_fraction,
-      const ComparisonFunction<Number> &   compare_refine,
-      const ComparisonFunction<Number> &   compare_coarsen)
+      const ComparisonFunction<typename identity<Number>::type> &compare_refine,
+      const ComparisonFunction<typename identity<Number>::type>
+        &compare_coarsen)
     {
       AssertDimension(dof_handler.get_triangulation().n_active_cells(),
                       criteria.size());
@@ -259,11 +261,12 @@ namespace hp
     template <int dim, typename Number, int spacedim>
     void
     p_adaptivity_from_reference(
-      const hp::DoFHandler<dim, spacedim> &dof_handler,
-      const Vector<Number> &               criteria,
-      const Vector<Number> &               references,
-      const ComparisonFunction<Number> &   compare_refine,
-      const ComparisonFunction<Number> &   compare_coarsen)
+      const hp::DoFHandler<dim, spacedim> &                      dof_handler,
+      const Vector<Number> &                                     criteria,
+      const Vector<Number> &                                     references,
+      const ComparisonFunction<typename identity<Number>::type> &compare_refine,
+      const ComparisonFunction<typename identity<Number>::type>
+        &compare_coarsen)
     {
       AssertDimension(dof_handler.get_triangulation().n_active_cells(),
                       criteria.size());
index f0ffb981957dc92b1dc75fa441c352eebb1ab6b9..fdbf312fe831b42421fb8f8dd3ba47eba905584f 100644 (file)
@@ -59,8 +59,8 @@ for (deal_II_dimension : DIMENSIONS; S : REAL_SCALARS;
           const Vector<S> &,
           const S,
           const S,
-          const std::function<bool(const S &, const S &)> &,
-          const std::function<bool(const S &, const S &)> &);
+          const ComparisonFunction<typename identity<S>::type> &,
+          const ComparisonFunction<typename identity<S>::type> &);
 
         template void
         p_adaptivity_from_relative_threshold<deal_II_dimension,
@@ -70,8 +70,8 @@ for (deal_II_dimension : DIMENSIONS; S : REAL_SCALARS;
           const Vector<S> &,
           const double,
           const double,
-          const std::function<bool(const S &, const S &)> &,
-          const std::function<bool(const S &, const S &)> &);
+          const ComparisonFunction<typename identity<S>::type> &,
+          const ComparisonFunction<typename identity<S>::type> &);
 
         template void
         p_adaptivity_from_regularity<deal_II_dimension,
@@ -87,8 +87,8 @@ for (deal_II_dimension : DIMENSIONS; S : REAL_SCALARS;
           const hp::DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
           const Vector<S> &,
           const Vector<S> &,
-          const std::function<bool(const S &, const S &)> &,
-          const std::function<bool(const S &, const S &)> &);
+          const ComparisonFunction<typename identity<S>::type> &,
+          const ComparisonFunction<typename identity<S>::type> &);
 
         template void
         predict_error<deal_II_dimension, S, deal_II_space_dimension>(

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.