]> https://gitweb.dealii.org/ - dealii.git/commitdiff
simplify
authorMatthias Maier <tamiko@43-1.org>
Fri, 15 Sep 2017 15:32:41 +0000 (10:32 -0500)
committerMatthias Maier <tamiko@43-1.org>
Sun, 17 Sep 2017 20:26:40 +0000 (15:26 -0500)
include/deal.II/fe/fe_tools.templates.h

index cd064e4b859a7f265bff03f71c56cfbd3fb2a6fb..ac7b14fa98fcc3330ad3b51538348338a2a0f229 100644 (file)
@@ -2912,55 +2912,50 @@ namespace FETools
 
   namespace
   {
-    // internal helper class to partially specialize the number type
+    // Helper functions for
+    // FETools::convert_generalized_support_point_values_to_dof_values
 
     template <int dim, int spacedim, typename number>
-    struct ConvertHelper
+    static void convert_helper(
+      const FiniteElement<dim, spacedim> &finite_element,
+      const std::vector<Vector<number> > &support_point_values,
+      std::vector<number>                &dof_values)
     {
-      static void convert_generalized_support_point_values_to_dof_values(
-        const FiniteElement<dim, spacedim> &finite_element,
-        const std::vector<Vector<number> > &support_point_values,
-        std::vector<number>                &dof_values)
-      {
-        static Threads::ThreadLocalStorage<std::vector<Vector<double> > >
-        double_support_point_values;
-        static Threads::ThreadLocalStorage<std::vector<double> >
-        double_dof_values;
+      static Threads::ThreadLocalStorage<std::vector<Vector<double> > >
+      double_support_point_values;
+      static Threads::ThreadLocalStorage<std::vector<double> >
+      double_dof_values;
 
-        double_support_point_values.get().resize(support_point_values.size());
-        double_dof_values.get().resize(dof_values.size());
+      double_support_point_values.get().resize(support_point_values.size());
+      double_dof_values.get().resize(dof_values.size());
 
-        for (unsigned int i = 0; i < support_point_values.size(); ++i)
-          {
-            double_support_point_values.get()[i].reinit(
-              finite_element.n_components(), false);
-            std::copy(std::begin(support_point_values[i]),
-                      std::end(support_point_values[i]),
-                      std::begin(double_support_point_values.get()[i]));
-          }
+      for (unsigned int i = 0; i < support_point_values.size(); ++i)
+        {
+          double_support_point_values.get()[i].reinit(
+            finite_element.n_components(), false);
+          std::copy(std::begin(support_point_values[i]),
+                    std::end(support_point_values[i]),
+                    std::begin(double_support_point_values.get()[i]));
+        }
 
-        finite_element.convert_generalized_support_point_values_to_dof_values(
-          double_support_point_values.get(), double_dof_values.get());
+      finite_element.convert_generalized_support_point_values_to_dof_values(
+        double_support_point_values.get(), double_dof_values.get());
+
+      std::copy(std::begin(double_dof_values.get()),
+                std::end(double_dof_values.get()),
+                std::begin(dof_values));
+    }
 
-        std::copy(std::begin(double_dof_values.get()),
-                  std::end(double_dof_values.get()),
-                  std::begin(dof_values));
-      }
-    };
 
-    // In case of double simply call the finite element directly
     template <int dim, int spacedim>
-    struct ConvertHelper<dim, spacedim, double>
+    static void convert_helper(
+      const FiniteElement<dim, spacedim> &finite_element,
+      const std::vector<Vector<double> > &support_point_values,
+      std::vector<double>                &dof_values)
     {
-      static void convert_generalized_support_point_values_to_dof_values(
-        const FiniteElement<dim, spacedim> &finite_element,
-        const std::vector<Vector<double> > &support_point_values,
-        std::vector<double>                &dof_values)
-      {
-        finite_element.convert_generalized_support_point_values_to_dof_values(
-          support_point_values, dof_values);
-      }
-    };
+      finite_element.convert_generalized_support_point_values_to_dof_values(
+        support_point_values, dof_values);
+    }
 
   } /* anonymous namespace */
 
@@ -2976,8 +2971,7 @@ namespace FETools
                     finite_element.get_generalized_support_points().size());
     AssertDimension(dof_values.size(), finite_element.dofs_per_cell);
 
-    ConvertHelper<dim, spacedim, number>::
-    convert_generalized_support_point_values_to_dof_values(
+    convert_helper<dim, spacedim>(
       finite_element, support_point_values, dof_values);
   }
 

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.