]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove another ArrayView 15182/head
authorMartin Kronbichler <martin.kronbichler@uni-a.de>
Fri, 5 May 2023 08:45:50 +0000 (10:45 +0200)
committerMartin Kronbichler <martin.kronbichler@uni-a.de>
Mon, 8 May 2023 06:43:15 +0000 (08:43 +0200)
include/deal.II/matrix_free/fe_point_evaluation.h
include/deal.II/matrix_free/tensor_product_kernels.h

index b22152cd7c8d2b41678a0a2d0c61641485772792..9efda1086b6014fcb843586f80a42d3d836aa7c8 100644 (file)
@@ -1494,6 +1494,7 @@ FEPointEvaluation<n_components, dim, spacedim, Number>::do_reinit()
     mapping_info->get_n_q_points_unvectorized(current_cell_index,
                                               current_face_number);
 
+  // round up n_points_scalar / n_lanes_user_interface
   const_cast<unsigned int &>(n_q_points) =
     (n_q_points_scalar + n_lanes_user_interface - 1) / n_lanes_user_interface;
 
@@ -1533,6 +1534,7 @@ FEPointEvaluation<n_components, dim, spacedim, Number>::do_reinit()
 
   if (fast_path && !polynomials_are_hat_functions)
     {
+      // round up n_q_points_scalar / n_lanes_internal
       const std::size_t n_batches =
         (n_q_points_scalar + n_lanes_internal - 1) / n_lanes_internal;
       const std::size_t n_shapes = poly.size();
@@ -1592,7 +1594,7 @@ FEPointEvaluation<n_components, dim, spacedim, Number>::evaluate_fast(
       if (evaluation_flags & EvaluationFlags::values)
         {
           for (unsigned int v = 0;
-               v < stride && (stride > 1 ? q + v < n_q_points_scalar : true);
+               v < stride && (stride == 1 || q + v < n_q_points_scalar);
                ++v)
             ETT::set_value(val_and_grad.first, v, values[qb * stride + v]);
         }
@@ -1603,7 +1605,7 @@ FEPointEvaluation<n_components, dim, spacedim, Number>::evaluate_fast(
                  ExcNotInitialized());
 
           for (unsigned int v = 0;
-               v < stride && (stride > 1 ? q + v < n_q_points_scalar : true);
+               v < stride && (stride == 1 || q + v < n_q_points_scalar);
                ++v)
             {
               const unsigned int offset = qb * stride + v;
@@ -1764,7 +1766,7 @@ FEPointEvaluation<n_components, dim, spacedim, Number>::integrate_fast(
             }
 
           for (unsigned int v = 0;
-               v < stride && (stride > 1 ? q + v < n_q_points_scalar : true);
+               v < stride && (stride == 1 || q + v < n_q_points_scalar);
                ++v)
             ETT::get_value(value, v, values[qb * stride + v]);
         }
@@ -1782,7 +1784,7 @@ FEPointEvaluation<n_components, dim, spacedim, Number>::integrate_fast(
             }
 
           for (unsigned int v = 0;
-               v < stride && (stride > 1 ? q + v < n_q_points_scalar : true);
+               v < stride && (stride == 1 || q + v < n_q_points_scalar);
                ++v)
             {
               const unsigned int offset = qb * stride + v;
@@ -1806,9 +1808,7 @@ FEPointEvaluation<n_components, dim, spacedim, Number>::integrate_fast(
         internal::integrate_add_tensor_product_value_and_gradient_shapes<
           dim,
           VectorizedArrayType,
-          vectorized_value_type>(make_array_view(shapes,
-                                                 qb * n_shapes,
-                                                 n_shapes),
+          vectorized_value_type>(shapes.data() + qb * n_shapes,
                                  n_shapes,
                                  value,
                                  gradient,
index db7732e636b28783f6212cfbaf95b0ebe06e87ac..20ae6672ee1a21355f2b90a6e490a3005485147f 100644 (file)
@@ -3477,12 +3477,11 @@ namespace internal
     DEAL_II_ALWAYS_INLINE
 #endif
     void
-    do_apply_test_functions_xy(
-      AlignedVector<Number2> &                          values,
-      const ArrayView<dealii::ndarray<Number, 2, dim>> &shapes,
-      const std::array<Number2, 3> &                    test_grads_value,
-      const int                                         n_shapes_runtime,
-      int &                                             i)
+    do_apply_test_functions_xy(AlignedVector<Number2> &               values,
+                               const dealii::ndarray<Number, 2, dim> *shapes,
+                               const std::array<Number2, 3> &test_grads_value,
+                               const int                     n_shapes_runtime,
+                               int &                         i)
   {
     if (length > 0)
       {
@@ -3544,11 +3543,11 @@ namespace internal
   template <int dim, typename Number, typename Number2>
   inline void
   integrate_add_tensor_product_value_and_gradient_shapes(
-    const ArrayView<dealii::ndarray<Number, 2, dim>> &shapes,
-    const int                                         n_shapes,
-    const Number2 &                                   value,
-    const Tensor<1, dim, Number2> &                   gradient,
-    AlignedVector<Number2> &                          values)
+    const dealii::ndarray<Number, 2, dim> *shapes,
+    const int                              n_shapes,
+    const Number2 &                        value,
+    const Tensor<1, dim, Number2> &        gradient,
+    AlignedVector<Number2> &               values)
   {
     static_assert(dim >= 1 && dim <= 3, "Only dim=1,2,3 implemented");
 

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.