]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Make MappingFEField use only update_each. 1339/head
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Wed, 12 Aug 2015 14:39:16 +0000 (09:39 -0500)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Wed, 12 Aug 2015 15:45:00 +0000 (10:45 -0500)
This is the exact equivalent to #1337 (commit 04430b5) but
for MappingFEField. The logic there seems to simply have been
copied from MappingQ1, so this was easy (and passes the testsuite).

source/fe/mapping_fe_field.cc

index e942f7a839ed6e46bab6728cf87a8a12e038dcf5..2a51d304cbd146b9b17104147f41298e4cab7deb 100644 (file)
@@ -125,27 +125,10 @@ MappingFEField<dim,spacedim,VECTOR,DH>::compute_shapes_virtual (
 
 template<int dim, int spacedim, class VECTOR, class DH>
 UpdateFlags
-MappingFEField<dim,spacedim,VECTOR,DH>::update_once (const UpdateFlags in) const
+MappingFEField<dim,spacedim,VECTOR,DH>::update_once (const UpdateFlags) const
 {
-  UpdateFlags out = UpdateFlags(in & (update_transformation_values
-                                      | update_transformation_gradients));
-
-  // Shape function values
-  if (in & update_quadrature_points)
-    out |= update_transformation_values;
-
-  // Shape function gradients
-  if (in & (update_covariant_transformation
-            | update_contravariant_transformation
-            | update_JxW_values
-            | update_boundary_forms
-            | update_normal_vectors
-            | update_jacobians
-            | update_jacobian_grads
-            | update_inverse_jacobians))
-    out |= update_transformation_gradients;
-
-  return out;
+  // nothing here
+  return update_default;
 }
 
 
@@ -239,10 +222,19 @@ MappingFEField<dim,spacedim,VECTOR,DH>::compute_data (const UpdateFlags      upd
 
   const UpdateFlags flags(data.update_flags);
 
-  if (flags & update_transformation_values)
+  // see if we need the (transformation) shape function values
+  // and/or gradients and resize the necessary arrays
+  if (flags & update_quadrature_points)
     data.shape_values.resize(data.n_shape_functions * n_q_points);
 
-  if (flags & update_transformation_gradients)
+  if (flags & (update_covariant_transformation
+               | update_contravariant_transformation
+               | update_JxW_values
+               | update_boundary_forms
+               | update_normal_vectors
+               | update_jacobians
+               | update_jacobian_grads
+               | update_inverse_jacobians))
     data.shape_derivatives.resize(data.n_shape_functions * n_q_points);
 
   if (flags & update_covariant_transformation)
@@ -390,7 +382,7 @@ namespace internal
                             const std::vector<unsigned int>      &fe_to_real,
                             std::vector<Point<spacedim> >        &quadrature_points)
     {
-      const UpdateFlags update_flags(data.current_update_flags());
+      const UpdateFlags update_flags = data.update_each;
 
       if (update_flags & update_quadrature_points)
         {
@@ -427,7 +419,7 @@ namespace internal
                             const ComponentMask                &fe_mask,
                             const std::vector<unsigned int>    &fe_to_real)
     {
-      const UpdateFlags update_flags(data.current_update_flags());
+      const UpdateFlags update_flags = data.update_each;
 
       // then Jacobians
       if (update_flags & update_contravariant_transformation)
@@ -497,7 +489,7 @@ namespace internal
                                  const std::vector<unsigned int>              &fe_to_real,
                                  std::vector<DerivativeForm<2,dim,spacedim> > &jacobian_grads)
     {
-      const UpdateFlags update_flags(data.current_update_flags());
+      const UpdateFlags update_flags = data.update_each;
       if (update_flags & update_jacobian_grads)
         {
           const unsigned int n_q_points = jacobian_grads.size();
@@ -551,7 +543,7 @@ namespace internal
                              const typename dealii::MappingFEField<dim,spacedim>::InternalData &data,
                              internal::FEValues::MappingRelatedData<dim,spacedim>         &output_data)
     {
-      const UpdateFlags update_flags(data.current_update_flags());
+      const UpdateFlags update_flags = data.update_each;
 
       if (update_flags & update_boundary_forms)
         {
@@ -749,7 +741,7 @@ fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
                                    QProjector<dim>::DataSetDescriptor::cell (),
                                    data, *fe, fe_mask, fe_to_real);
 
-  const UpdateFlags update_flags(data.current_update_flags());
+  const UpdateFlags update_flags = data.update_each;
   const std::vector<double> &weights=quadrature.get_weights();
 
   // Multiply quadrature weights by absolute value of Jacobian determinants or
@@ -1080,8 +1072,8 @@ transform_unit_to_real_cell (const typename Triangulation<dim,spacedim>::cell_it
 //  the right size and transformation shape values already computed at point
 //  p.
   const Quadrature<dim> point_quadrature(p);
-  std_cxx11::unique_ptr<InternalData>
-  mdata (get_data(update_transformation_values, point_quadrature));
+  std_cxx11::unique_ptr<InternalData> mdata (get_data(update_quadrature_points | update_jacobians,
+                                                      point_quadrature));
 
   update_internal_dofs(cell, *mdata);
 
@@ -1138,7 +1130,7 @@ transform_real_to_unit_cell (const typename Triangulation<dim,spacedim>::cell_it
 
   const Quadrature<dim> point_quadrature(initial_p_unit);
 
-  UpdateFlags update_flags = update_transformation_values|update_transformation_gradients;
+  UpdateFlags update_flags = update_quadrature_points | update_jacobians;
   if (spacedim>dim)
     update_flags |= update_jacobian_grads;
   std_cxx11::unique_ptr<InternalData>
@@ -1259,7 +1251,6 @@ failure:
 }
 
 
-
 template<int dim, int spacedim, class VECTOR, class DH>
 unsigned int
 MappingFEField<dim,spacedim,VECTOR,DH>::get_degree() const

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.