]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Adjust return parameter. Add inverse mass matrix 10922/head
authorMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Thu, 17 Sep 2020 10:59:57 +0000 (12:59 +0200)
committerMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Thu, 17 Sep 2020 10:59:57 +0000 (12:59 +0200)
include/deal.II/matrix_free/evaluation_kernels.h
include/deal.II/matrix_free/operators.h

index da6a399cc88e13b7f7001dface2a816f9176bcd7..8f40cb85d37a6310d639052a90bdbd180d6f04ef 100644 (file)
@@ -1396,7 +1396,7 @@ namespace internal
   struct FEEvaluationImplEvaluateSelector
   {
     template <int fe_degree, int n_q_points_1d>
-    static void
+    static bool
     run(const unsigned int                                      n_components,
         const EvaluationFlags::EvaluationFlags                  evaluation_flag,
         const internal::MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
@@ -1519,6 +1519,8 @@ namespace internal
         }
       else
         AssertThrow(false, ExcNotImplemented());
+
+      return false;
     }
   };
 
@@ -1543,7 +1545,7 @@ namespace internal
   struct FEEvaluationImplIntegrateSelector
   {
     template <int fe_degree, int n_q_points_1d>
-    static void
+    static bool
     run(const unsigned int                     n_components,
         const EvaluationFlags::EvaluationFlags integration_flag,
         const internal::MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
@@ -1666,6 +1668,8 @@ namespace internal
         }
       else
         AssertThrow(false, ExcNotImplemented());
+
+      return false;
     }
   };
 
@@ -2313,7 +2317,7 @@ namespace internal
   struct FEFaceEvaluationImplEvaluateSelector
   {
     template <int fe_degree, int n_q_points_1d>
-    static void
+    static bool
     run(const unsigned int                                         n_components,
         const MatrixFreeFunctions::ShapeInfo<VectorizedArrayType> &data,
         const VectorizedArrayType *                                values_array,
@@ -2392,6 +2396,8 @@ namespace internal
                                     scratch_data,
                                     values_quad,
                                     gradients_quad);
+
+      return false;
     }
   };
 
@@ -2401,7 +2407,7 @@ namespace internal
   struct FEFaceEvaluationImplIntegrateSelector
   {
     template <int fe_degree, int n_q_points_1d>
-    static void
+    static bool
     run(const unsigned int                                         n_components,
         const MatrixFreeFunctions::ShapeInfo<VectorizedArrayType> &data,
         VectorizedArrayType *                                      values_array,
@@ -2484,6 +2490,7 @@ namespace internal
                                            values_array,
                                            integrate_gradients,
                                            face_no);
+      return false;
     }
   };
 
@@ -3515,18 +3522,20 @@ namespace internal
 
   /**
    * This struct implements the action of the inverse mass matrix operation
+   * using an FEEvaluationBaseData argument.
    */
-  template <int dim, int fe_degree, typename Number>
-  struct CellwiseInverseMassMatrixImpl
+  template <int dim, typename Number>
+  struct CellwiseInverseMassMatrixImplBasic
   {
-    static void
-    apply(const unsigned int                  n_components,
-          const FEEvaluationBaseData<dim,
-                                     typename Number::value_type,
-                                     false,
-                                     Number> &fe_eval,
-          const Number *                      in_array,
-          Number *                            out_array)
+    template <int fe_degree, int = 0>
+    static bool
+    run(const unsigned int                  n_components,
+        const FEEvaluationBaseData<dim,
+                                   typename Number::value_type,
+                                   false,
+                                   Number> &fe_eval,
+        const Number *                      in_array,
+        Number *                            out_array)
     {
       constexpr unsigned int dofs_per_component =
         Utilities::pow(fe_degree + 1, dim);
@@ -3573,14 +3582,26 @@ namespace internal
             evaluator.template hessians<1, false, false>(out, out);
           evaluator.template hessians<0, false, false>(out, out);
         }
+      return false;
     }
+  };
 
-    static void
-    apply(const unsigned int           n_desired_components,
-          const AlignedVector<Number> &inverse_shape,
-          const AlignedVector<Number> &inverse_coefficients,
-          const Number *               in_array,
-          Number *                     out_array)
+
+
+  /**
+   * This struct implements the action of the inverse mass matrix operation
+   * using an FEEvaluationBaseData argument.
+   */
+  template <int dim, typename Number>
+  struct CellwiseInverseMassMatrixImplFlexible
+  {
+    template <int fe_degree, int = 0>
+    static bool
+    run(const unsigned int           n_desired_components,
+        const AlignedVector<Number> &inverse_shape,
+        const AlignedVector<Number> &inverse_coefficients,
+        const Number *               in_array,
+        Number *                     out_array)
     {
       constexpr unsigned int dofs_per_component =
         Utilities::pow(fe_degree + 1, dim);
@@ -3630,13 +3651,25 @@ namespace internal
 
           inv_coefficient += shift_coefficient;
         }
+      return false;
     }
+  };
 
-    static void
-    transform_from_q_points_to_basis(const unsigned int n_desired_components,
-                                     const AlignedVector<Number> &inverse_shape,
-                                     const Number *               in_array,
-                                     Number *                     out_array)
+
+
+  /**
+   * This struct implements the action of the inverse mass matrix operation
+   * using an FEEvaluationBaseData argument.
+   */
+  template <int dim, typename Number>
+  struct CellwiseInverseMassMatrixImplTransformFromQPoints
+  {
+    template <int fe_degree, int = 0>
+    static bool
+    run(const unsigned int           n_desired_components,
+        const AlignedVector<Number> &inverse_shape,
+        const Number *               in_array,
+        Number *                     out_array)
     {
       constexpr unsigned int dofs_per_cell = Utilities::pow(fe_degree + 1, dim);
       internal::EvaluatorTensorProduct<internal::evaluate_evenodd,
@@ -3667,6 +3700,7 @@ namespace internal
           if (dim == 1)
             evaluator.template hessians<0, false, false>(in, out);
         }
+      return false;
     }
   };
 
index 9e68e74c9809fec2d0b339ab445a047a4db95086..cd5cb3e7cae5fb107f76405447bc07dc1a5b6ded 100644 (file)
@@ -1016,9 +1016,8 @@ namespace MatrixFreeOperators
     VectorizedArrayType>::apply(const VectorizedArrayType *in_array,
                                 VectorizedArrayType *      out_array) const
   {
-    internal::
-      CellwiseInverseMassMatrixImpl<dim, fe_degree, VectorizedArrayType>::apply(
-        n_components, fe_eval, in_array, out_array);
+    internal::CellwiseInverseMassMatrixImplBasic<dim, VectorizedArrayType>::
+      template run<fe_degree>(n_components, fe_eval, in_array, out_array);
   }
 
 
@@ -1039,8 +1038,8 @@ namespace MatrixFreeOperators
           const VectorizedArrayType *               in_array,
           VectorizedArrayType *                     out_array) const
   {
-    internal::
-      CellwiseInverseMassMatrixImpl<dim, fe_degree, VectorizedArrayType>::apply(
+    internal::CellwiseInverseMassMatrixImplFlexible<dim, VectorizedArrayType>::
+      template run<fe_degree>(
         n_actual_components,
         fe_eval.get_shape_info().data.front().inverse_shape_values_eo,
         inverse_coefficients,
@@ -1065,13 +1064,14 @@ namespace MatrixFreeOperators
                                      const VectorizedArrayType *in_array,
                                      VectorizedArrayType *      out_array) const
   {
-    internal::
-      CellwiseInverseMassMatrixImpl<dim, fe_degree, VectorizedArrayType>::
-        transform_from_q_points_to_basis(
-          n_actual_components,
-          fe_eval.get_shape_info().data.front().inverse_shape_values_eo,
-          in_array,
-          out_array);
+    internal::CellwiseInverseMassMatrixImplTransformFromQPoints<
+      dim,
+      VectorizedArrayType>::template run<fe_degree>(n_actual_components,
+                                                    fe_eval.get_shape_info()
+                                                      .data.front()
+                                                      .inverse_shape_values_eo,
+                                                    in_array,
+                                                    out_array);
   }
 
 

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.