]> https://gitweb.dealii.org/ - dealii.git/commitdiff
fe/fe_values: convert internal::do_function_values to ArrayView
authorMatthias Maier <maier@tamu.edu>
Thu, 22 Jun 2023 05:49:43 +0000 (00:49 -0500)
committerMatthias Maier <tamiko@43-1.org>
Thu, 22 Jun 2023 05:59:29 +0000 (00:59 -0500)
source/fe/fe_values.cc

index 5f5a84e663c04193c7abae859d919e2ea26efc77..562c636d53c6230304e698c25c9522f4cb50173d 100644 (file)
@@ -2817,7 +2817,7 @@ namespace internal
 
   template <typename Number, typename Number2>
   void
-  do_function_values(const Number2 *                 dof_values_ptr,
+  do_function_values(const ArrayView<Number2> &      dof_values,
                      const dealii::Table<2, double> &shape_values,
                      std::vector<Number> &           values)
   {
@@ -2839,7 +2839,7 @@ namespace internal
     // the shape_values data stored contiguously
     for (unsigned int shape_func = 0; shape_func < dofs_per_cell; ++shape_func)
       {
-        const Number2 value = dof_values_ptr[shape_func];
+        const Number2 value = dof_values[shape_func];
         // For auto-differentiable numbers, the fact that a DoF value is zero
         // does not imply that its derivatives are zero as well. So we
         // can't filter by value for these number types.
@@ -2858,13 +2858,13 @@ namespace internal
   template <int dim, int spacedim, typename VectorType>
   void
   do_function_values(
-    const typename VectorType::value_type *dof_values_ptr,
-    const dealii::Table<2, double> &       shape_values,
-    const FiniteElement<dim, spacedim> &   fe,
-    const std::vector<unsigned int> &      shape_function_to_row_table,
-    ArrayView<VectorType>                  values,
-    const bool                             quadrature_points_fastest = false,
-    const unsigned int                     component_multiple        = 1)
+    const ArrayView<typename VectorType::value_type> &dof_values,
+    const dealii::Table<2, double> &                  shape_values,
+    const FiniteElement<dim, spacedim> &              fe,
+    const std::vector<unsigned int> &shape_function_to_row_table,
+    ArrayView<VectorType>            values,
+    const bool                       quadrature_points_fastest = false,
+    const unsigned int               component_multiple        = 1)
   {
     using Number = typename VectorType::value_type;
     // initialize with zero
@@ -2906,7 +2906,7 @@ namespace internal
       for (unsigned int shape_func = 0; shape_func < dofs_per_cell;
            ++shape_func)
         {
-          const Number &value = dof_values_ptr[shape_func + mc * dofs_per_cell];
+          const Number &value = dof_values[shape_func + mc * dofs_per_cell];
           // For auto-differentiable numbers, the fact that a DoF value is zero
           // does not imply that its derivatives are zero as well. So we
           // can't filter by value for these number types.
@@ -3281,7 +3281,8 @@ FEValuesBase<dim, spacedim>::get_function_values(
   // get function values of dofs on this cell
   Vector<Number> dof_values(dofs_per_cell);
   present_cell.get_interpolated_dof_values(fe_function, dof_values);
-  internal::do_function_values(dof_values.begin(),
+  internal::do_function_values(make_array_view(dof_values.begin(),
+                                               dof_values.end()),
                                this->finite_element_output.shape_values,
                                values);
 }
@@ -3305,7 +3306,8 @@ FEValuesBase<dim, spacedim>::get_function_values(
   boost::container::small_vector<Number, 200> dof_values(dofs_per_cell);
   for (unsigned int i = 0; i < dofs_per_cell; ++i)
     dof_values[i] = internal::get_vector_element(fe_function, indices[i]);
-  internal::do_function_values(dof_values.data(),
+  internal::do_function_values(make_array_view(dof_values.begin(),
+                                               dof_values.end()),
                                this->finite_element_output.shape_values,
                                values);
 }
@@ -3330,7 +3332,7 @@ FEValuesBase<dim, spacedim>::get_function_values(
   Vector<Number> dof_values(dofs_per_cell);
   present_cell.get_interpolated_dof_values(fe_function, dof_values);
   internal::do_function_values(
-    dof_values.begin(),
+    make_array_view(dof_values.begin(), dof_values.end()),
     this->finite_element_output.shape_values,
     *fe,
     this->finite_element_output.shape_function_to_row_table,
@@ -3359,7 +3361,7 @@ FEValuesBase<dim, spacedim>::get_function_values(
   for (unsigned int i = 0; i < dofs_per_cell; ++i)
     dof_values[i] = internal::get_vector_element(fe_function, indices[i]);
   internal::do_function_values(
-    dof_values.data(),
+    make_array_view(dof_values.begin(), dof_values.end()),
     this->finite_element_output.shape_values,
     *fe,
     this->finite_element_output.shape_function_to_row_table,
@@ -3392,7 +3394,7 @@ FEValuesBase<dim, spacedim>::get_function_values(
   for (unsigned int i = 0; i < indices.size(); ++i)
     dof_values[i] = internal::get_vector_element(fe_function, indices[i]);
   internal::do_function_values(
-    dof_values.data(),
+    make_array_view(dof_values.begin(), dof_values.end()),
     this->finite_element_output.shape_values,
     *fe,
     this->finite_element_output.shape_function_to_row_table,

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.