]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Make FEEvaluationImpl work for ElementType::tensor_none 11188/head
authorPeter Munch <peterrmuench@gmail.com>
Sun, 15 Nov 2020 16:57:46 +0000 (17:57 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Tue, 17 Nov 2020 15:18:41 +0000 (16:18 +0100)
include/deal.II/matrix_free/evaluation_kernels.h
include/deal.II/matrix_free/shape_info.h

index eb9c18af22b76568cfab35fc73039e5b655ea268..036b628f597b4ab9e332041a0c9b42d2ee9d31e6 100644 (file)
@@ -137,6 +137,40 @@ namespace internal
 
 
 
+  /**
+   * Specialization for MatrixFreeFunctions::tensor_none, which cannot use the
+   * sum-factorization kernels.
+   */
+  template <int dim, int fe_degree, int n_q_points_1d, typename Number>
+  struct FEEvaluationImpl<MatrixFreeFunctions::tensor_none,
+                          dim,
+                          fe_degree,
+                          n_q_points_1d,
+                          Number>
+  {
+    static void
+    evaluate(const unsigned int                            n_components,
+             const EvaluationFlags::EvaluationFlags        evaluation_flag,
+             const MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
+             const Number *                                values_dofs_actual,
+             Number *                                      values_quad,
+             Number *                                      gradients_quad,
+             Number *                                      hessians_quad,
+             Number *                                      scratch_data);
+
+    static void
+    integrate(const unsigned int                            n_components,
+              const EvaluationFlags::EvaluationFlags        integration_flag,
+              const MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
+              Number *                                      values_dofs_actual,
+              Number *                                      values_quad,
+              Number *                                      gradients_quad,
+              Number *                                      scratch_data,
+              const bool add_into_values_array);
+  };
+
+
+
   template <MatrixFreeFunctions::ElementType type,
             int                              dim,
             int                              fe_degree,
@@ -641,6 +675,120 @@ namespace internal
 
 
 
+  template <int dim, int fe_degree, int n_q_points_1d, typename Number>
+  inline void
+  FEEvaluationImpl<
+    MatrixFreeFunctions::tensor_none,
+    dim,
+    fe_degree,
+    n_q_points_1d,
+    Number>::evaluate(const unsigned int                     n_components,
+                      const EvaluationFlags::EvaluationFlags evaluation_flag,
+                      const MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
+                      const Number *values_dofs_actual,
+                      Number *      values_quad,
+                      Number *      gradients_quad,
+                      Number *      hessians_quad,
+                      Number *      scratch_data)
+  {
+    (void)scratch_data;
+
+    const unsigned int n_dofs     = shape_info.dofs_per_component_on_cell;
+    const unsigned int n_q_points = shape_info.n_q_points;
+
+    if (evaluation_flag & EvaluationFlags::values)
+      {
+        const auto shape_values = shape_info.data.front().shape_values.data();
+        for (unsigned int c = 0; c < n_components; ++c)
+          for (unsigned int q = 0; q < n_q_points; ++q)
+            {
+              values_quad[n_q_points * c + q] = 0.0;
+              for (unsigned int i = 0; i < n_dofs; ++i)
+                values_quad[n_q_points * c + q] +=
+                  shape_values[i * n_q_points + q] *
+                  values_dofs_actual[n_dofs * c + i];
+            }
+      }
+
+    if (evaluation_flag & EvaluationFlags::gradients)
+      {
+        const auto shape_gradients =
+          shape_info.data.front().shape_gradients.data();
+        for (unsigned int c = 0; c < n_components; ++c)
+          for (unsigned int d = 0; d < dim; ++d)
+            for (unsigned int q = 0; q < n_q_points; ++q)
+              {
+                gradients_quad[n_q_points * dim * c + n_q_points * d + q] = 0.0;
+                for (unsigned int i = 0; i < n_dofs; ++i)
+                  gradients_quad[n_q_points * dim * c + n_q_points * d + q] +=
+                    shape_gradients[d * n_dofs * n_q_points + i * n_q_points +
+                                    q] *
+                    values_dofs_actual[n_dofs * dim * c + i];
+              }
+      }
+
+    if (evaluation_flag & EvaluationFlags::hessians)
+      {
+        Assert(false, ExcNotImplemented());
+        (void)hessians_quad;
+      }
+  }
+
+
+
+  template <int dim, int fe_degree, int n_q_points_1d, typename Number>
+  inline void
+  FEEvaluationImpl<
+    MatrixFreeFunctions::tensor_none,
+    dim,
+    fe_degree,
+    n_q_points_1d,
+    Number>::integrate(const unsigned int                     n_components,
+                       const EvaluationFlags::EvaluationFlags integration_flag,
+                       const MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
+                       Number *   values_dofs_actual,
+                       Number *   values_quad,
+                       Number *   gradients_quad,
+                       Number *   scratch_data,
+                       const bool add_into_values_array)
+  {
+    (void)scratch_data;
+
+    const unsigned int n_dofs     = shape_info.dofs_per_component_on_cell;
+    const unsigned int n_q_points = shape_info.n_q_points;
+
+    if (add_into_values_array == false)
+      for (unsigned int i = 0; i < n_components * n_dofs; ++i)
+        values_dofs_actual[i] = 0.0;
+
+    if (integration_flag & EvaluationFlags::values)
+      {
+        const auto shape_values = shape_info.data.front().shape_values.data();
+        for (unsigned int c = 0; c < n_components; ++c)
+          for (unsigned int q = 0; q < n_q_points; ++q)
+            for (unsigned int i = 0; i < n_dofs; ++i)
+              values_dofs_actual[n_dofs * c + i] +=
+                shape_values[i * n_q_points + q] *
+                values_quad[n_q_points * c + q];
+      }
+
+    if (integration_flag & EvaluationFlags::gradients)
+      {
+        const auto shape_gradients =
+          shape_info.data.front().shape_gradients.data();
+        for (unsigned int c = 0; c < n_components; ++c)
+          for (unsigned int d = 0; d < dim; ++d)
+            for (unsigned int q = 0; q < n_q_points; ++q)
+              for (unsigned int i = 0; i < n_dofs; ++i)
+                values_dofs_actual[n_dofs * dim * c + i] +=
+                  shape_gradients[d * n_dofs * n_q_points + i * n_q_points +
+                                  q] *
+                  gradients_quad[n_q_points * dim * c + n_q_points * d + q];
+      }
+  }
+
+
+
   /**
    * This struct implements the change between two different bases. This is an
    * ingredient in the FEEvaluationImplTransformToCollocation class where we
@@ -1499,6 +1647,22 @@ namespace internal
                               hessians_quad,
                               scratch_data);
         }
+      else if (shape_info.element_type ==
+               internal::MatrixFreeFunctions::tensor_none)
+        {
+          internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_none,
+                                     dim,
+                                     fe_degree,
+                                     n_q_points_1d,
+                                     Number>::evaluate(n_components,
+                                                       evaluation_flag,
+                                                       shape_info,
+                                                       values_dofs_actual,
+                                                       values_quad,
+                                                       gradients_quad,
+                                                       hessians_quad,
+                                                       scratch_data);
+        }
       else
         {
           internal::FEEvaluationImpl<
@@ -1645,6 +1809,22 @@ namespace internal
                                scratch_data,
                                sum_into_values_array);
         }
+      else if (shape_info.element_type ==
+               internal::MatrixFreeFunctions::tensor_none)
+        {
+          internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_none,
+                                     dim,
+                                     fe_degree,
+                                     n_q_points_1d,
+                                     Number>::integrate(n_components,
+                                                        integration_flag,
+                                                        shape_info,
+                                                        values_dofs_actual,
+                                                        values_quad,
+                                                        gradients_quad,
+                                                        scratch_data,
+                                                        sum_into_values_array);
+        }
       else
         {
           internal::FEEvaluationImpl<
index bb896a9c54b24bbd094c644e5de27274021759b7..12f3b6eea1a13f3c9de6ac70fe77882f6ed950b8 100644 (file)
@@ -91,7 +91,12 @@ namespace internal
        * of the unit interval 0.5 that additionally add a constant shape
        * function according to FE_Q_DG0.
        */
-      tensor_symmetric_plus_dg0 = 5
+      tensor_symmetric_plus_dg0 = 5,
+
+      /**
+       * Shape functions without an tensor product properties.
+       */
+      tensor_none = 6
     };
 
 

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.