]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove VectorSlice from some internal functions. 5099/head
authorDavid Wells <wellsd2@rpi.edu>
Sat, 16 Sep 2017 20:36:18 +0000 (16:36 -0400)
committerDavid Wells <wellsd2@rpi.edu>
Sun, 17 Sep 2017 21:34:00 +0000 (17:34 -0400)
source/fe/fe_values.cc

index 4892d82e31973ea956969503582a6410d175bb40..87d9337dcedb35e205d4c5e17b72d7fe7d2eff18 100644 (file)
@@ -2700,16 +2700,17 @@ namespace internal
       }
   }
 
-  template <int dim, int spacedim, typename VectorType, typename Number>
+  template <int dim, int spacedim, typename VectorType>
   void
-  do_function_values (const Number                      *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,
-                      VectorSlice<std::vector<VectorType> > &values,
-                      const bool quadrature_points_fastest  = false,
-                      const unsigned int component_multiple = 1)
+  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)
   {
+    typedef typename VectorType::value_type Number;
     // initialize with zero
     for (unsigned int i=0; i<values.size(); ++i)
       std::fill_n (values[i].begin(), values[i].size(),
@@ -2835,13 +2836,13 @@ namespace internal
 
   template <int order, int dim, int spacedim, typename Number>
   void
-  do_function_derivatives (const Number                      *dof_values_ptr,
-                           const dealii::Table<2,Tensor<order,spacedim> > &shape_derivatives,
-                           const FiniteElement<dim,spacedim> &fe,
-                           const std::vector<unsigned int> &shape_function_to_row_table,
-                           VectorSlice<std::vector<std::vector<Tensor<order,spacedim,Number> > > > &derivatives,
-                           const bool quadrature_points_fastest  = false,
-                           const unsigned int component_multiple = 1)
+  do_function_derivatives (const Number                                            *dof_values_ptr,
+                           const dealii::Table<2,Tensor<order,spacedim> >          &shape_derivatives,
+                           const FiniteElement<dim,spacedim>                       &fe,
+                           const std::vector<unsigned int>                         &shape_function_to_row_table,
+                           ArrayView<std::vector<Tensor<order,spacedim,Number> > >  derivatives,
+                           const bool quadrature_points_fastest                   = false,
+                           const unsigned int component_multiple                  = 1)
   {
     // initialize with zero
     for (unsigned int i=0; i<derivatives.size(); ++i)
@@ -3123,9 +3124,11 @@ void 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);
-  VectorSlice<std::vector<Vector<Number> > > val(values);
-  internal::do_function_values(dof_values.begin(), this->finite_element_output.shape_values, *fe,
-                               this->finite_element_output.shape_function_to_row_table, val);
+  internal::do_function_values(dof_values.begin(),
+                               this->finite_element_output.shape_values,
+                               *fe,
+                               this->finite_element_output.shape_function_to_row_table,
+                               make_array_view(values.begin(), values.end()));
 }
 
 
@@ -3145,13 +3148,16 @@ void FEValuesBase<dim,spacedim>::get_function_values (
   Assert (this->update_flags & update_values,
           ExcAccessToUninitializedField("update_values"));
 
-  VectorSlice<std::vector<Vector<Number> > > val(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] = get_vector_element (fe_function, indices[i]);
-  internal::do_function_values(dof_values.data(), this->finite_element_output.shape_values, *fe,
-                               this->finite_element_output.shape_function_to_row_table, val,
-                               false, indices.size()/dofs_per_cell);
+  internal::do_function_values(dof_values.data(),
+                               this->finite_element_output.shape_values,
+                               *fe,
+                               this->finite_element_output.shape_function_to_row_table,
+                               make_array_view(values.begin(), values.end()),
+                               false,
+                               indices.size()/dofs_per_cell);
 }
 
 
@@ -3176,8 +3182,11 @@ void FEValuesBase<dim,spacedim>::get_function_values (
   boost::container::small_vector<Number, 200> dof_values(indices.size());
   for (unsigned int i=0; i<indices.size(); ++i)
     dof_values[i] = get_vector_element (fe_function, indices[i]);
-  internal::do_function_values(dof_values.data(), this->finite_element_output.shape_values, *fe,
-                               this->finite_element_output.shape_function_to_row_table, values,
+  internal::do_function_values(dof_values.data(),
+                               this->finite_element_output.shape_values,
+                               *fe,
+                               this->finite_element_output.shape_function_to_row_table,
+                               make_array_view(values.begin(), values.end()),
                                quadrature_points_fastest,
                                indices.size()/dofs_per_cell);
 }
@@ -3248,10 +3257,11 @@ FEValuesBase<dim,spacedim>::get_function_gradients (
   // 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);
-  VectorSlice<std::vector<std::vector<Tensor<1,spacedim,Number> > > > grads(gradients);
-  internal::do_function_derivatives(dof_values.begin(), this->finite_element_output.shape_gradients,
-                                    *fe, this->finite_element_output.shape_function_to_row_table,
-                                    grads);
+  internal::do_function_derivatives(dof_values.begin(),
+                                    this->finite_element_output.shape_gradients,
+                                    *fe,
+                                    this->finite_element_output.shape_function_to_row_table,
+                                    make_array_view(gradients.begin(), gradients.end()));
 }
 
 
@@ -3275,9 +3285,12 @@ void FEValuesBase<dim,spacedim>::get_function_gradients (
   boost::container::small_vector<Number, 200> dof_values(indices.size());
   for (unsigned int i=0; i<indices.size(); ++i)
     dof_values[i] = get_vector_element (fe_function, indices[i]);
-  internal::do_function_derivatives(dof_values.data(), this->finite_element_output.shape_gradients,
-                                    *fe, this->finite_element_output.shape_function_to_row_table,
-                                    gradients, quadrature_points_fastest,
+  internal::do_function_derivatives(dof_values.data(),
+                                    this->finite_element_output.shape_gradients,
+                                    *fe,
+                                    this->finite_element_output.shape_function_to_row_table,
+                                    make_array_view(gradients.begin(), gradients.end()),
+                                    quadrature_points_fastest,
                                     indices.size()/dofs_per_cell);
 }
 
@@ -3348,10 +3361,12 @@ get_function_hessians (const InputVector                         &fe_function,
   // 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);
-  VectorSlice<std::vector<std::vector<Tensor<2,spacedim,Number> > > > hes(hessians);
-  internal::do_function_derivatives(dof_values.begin(), this->finite_element_output.shape_hessians,
-                                    *fe, this->finite_element_output.shape_function_to_row_table,
-                                    hes, quadrature_points_fastest);
+  internal::do_function_derivatives(dof_values.begin(),
+                                    this->finite_element_output.shape_hessians,
+                                    *fe,
+                                    this->finite_element_output.shape_function_to_row_table,
+                                    make_array_view(hessians.begin(), hessians.end()),
+                                    quadrature_points_fastest);
 }
 
 
@@ -3373,9 +3388,12 @@ void FEValuesBase<dim, spacedim>::get_function_hessians (
   boost::container::small_vector<Number, 200> dof_values(indices.size());
   for (unsigned int i=0; i<indices.size(); ++i)
     dof_values[i] = get_vector_element (fe_function, indices[i]);
-  internal::do_function_derivatives(dof_values.data(), this->finite_element_output.shape_hessians,
-                                    *fe, this->finite_element_output.shape_function_to_row_table,
-                                    hessians, quadrature_points_fastest,
+  internal::do_function_derivatives(dof_values.data(),
+                                    this->finite_element_output.shape_hessians,
+                                    *fe,
+                                    this->finite_element_output.shape_function_to_row_table,
+                                    make_array_view(hessians.begin(), hessians.end()),
+                                    quadrature_points_fastest,
                                     indices.size()/dofs_per_cell);
 }
 
@@ -3565,10 +3583,12 @@ get_function_third_derivatives (const InputVector                         &fe_fu
   // 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);
-  VectorSlice<std::vector<std::vector<Tensor<3,spacedim,Number> > > > third(third_derivatives);
-  internal::do_function_derivatives(dof_values.begin(), this->finite_element_output.shape_3rd_derivatives,
-                                    *fe, this->finite_element_output.shape_function_to_row_table,
-                                    third, quadrature_points_fastest);
+  internal::do_function_derivatives(dof_values.begin(),
+                                    this->finite_element_output.shape_3rd_derivatives,
+                                    *fe,
+                                    this->finite_element_output.shape_function_to_row_table,
+                                    make_array_view(third_derivatives.begin(), third_derivatives.end()),
+                                    quadrature_points_fastest);
 }
 
 
@@ -3590,9 +3610,12 @@ void FEValuesBase<dim, spacedim>::get_function_third_derivatives (
   boost::container::small_vector<Number, 200> dof_values(indices.size());
   for (unsigned int i=0; i<indices.size(); ++i)
     dof_values[i] = get_vector_element (fe_function, indices[i]);
-  internal::do_function_derivatives(dof_values.data(), this->finite_element_output.shape_3rd_derivatives,
-                                    *fe, this->finite_element_output.shape_function_to_row_table,
-                                    third_derivatives, quadrature_points_fastest,
+  internal::do_function_derivatives(dof_values.data(),
+                                    this->finite_element_output.shape_3rd_derivatives,
+                                    *fe,
+                                    this->finite_element_output.shape_function_to_row_table,
+                                    make_array_view(third_derivatives.begin(), third_derivatives.end()),
+                                    quadrature_points_fastest,
                                     indices.size()/dofs_per_cell);
 }
 

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.