]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add function MappingCartesian::fill_jacobian_derivatives(..) 8773/head
authorSimon Sticko <simon@sticko.se>
Fri, 13 Sep 2019 08:54:02 +0000 (10:54 +0200)
committerSimon Sticko <simon@sticko.se>
Tue, 17 Sep 2019 14:31:02 +0000 (16:31 +0200)
Presently all of the fill_fe_*values functions in MappingCartesian
contain loops that fill the derivatives of the Jacobian with zeros.
In order to remove duplication, add a function which does this instead.

include/deal.II/fe/mapping_cartesian.h
source/fe/mapping_cartesian.cc

index 7b18a9e518cf06ed91e6d2193dcf76c7afb9c730..89f115096c49db54709028004778938d4008b490 100644 (file)
@@ -328,6 +328,18 @@ private:
     const unsigned int           face_no,
     const InternalData &         data,
     std::vector<Tensor<1, dim>> &normal_vectors) const;
+
+  /**
+   * Since the Jacobian is constant for this mapping all derivatives of the
+   * Jacobian are identically zero. Fill these quantities with zeros if the
+   * corresponding update flags say that they should be updated.
+   */
+  void
+  maybe_update_jacobian_derivatives(
+    const InternalData &             data,
+    const CellSimilarity::Similarity cell_similarity,
+    internal::FEValuesImplementation::MappingRelatedData<dim, spacedim>
+      &output_data) const;
 };
 
 /*@}*/
index 96aa3d35741239674b481a9c3fd544c33b89a6e4..f26f1b1aff6926fc242c085ec54d5e989764adcc 100644 (file)
@@ -315,6 +315,58 @@ MappingCartesian<dim, spacedim>::maybe_update_normal_vectors(
 
 
 
+template <int dim, int spacedim>
+void
+MappingCartesian<dim, spacedim>::maybe_update_jacobian_derivatives(
+  const InternalData &             data,
+  const CellSimilarity::Similarity cell_similarity,
+  internal::FEValuesImplementation::MappingRelatedData<dim, spacedim>
+    &output_data) const
+{
+  if (cell_similarity != CellSimilarity::translation)
+    {
+      if (data.update_each & update_jacobian_grads)
+        for (unsigned int i = 0; i < output_data.jacobian_grads.size(); ++i)
+          output_data.jacobian_grads[i] = DerivativeForm<2, dim, spacedim>();
+
+      if (data.update_each & update_jacobian_pushed_forward_grads)
+        for (unsigned int i = 0;
+             i < output_data.jacobian_pushed_forward_grads.size();
+             ++i)
+          output_data.jacobian_pushed_forward_grads[i] = Tensor<3, spacedim>();
+
+      if (data.update_each & update_jacobian_2nd_derivatives)
+        for (unsigned int i = 0;
+             i < output_data.jacobian_2nd_derivatives.size();
+             ++i)
+          output_data.jacobian_2nd_derivatives[i] =
+            DerivativeForm<3, dim, spacedim>();
+
+      if (data.update_each & update_jacobian_pushed_forward_2nd_derivatives)
+        for (unsigned int i = 0;
+             i < output_data.jacobian_pushed_forward_2nd_derivatives.size();
+             ++i)
+          output_data.jacobian_pushed_forward_2nd_derivatives[i] =
+            Tensor<4, spacedim>();
+
+      if (data.update_each & update_jacobian_3rd_derivatives)
+        for (unsigned int i = 0;
+             i < output_data.jacobian_3rd_derivatives.size();
+             ++i)
+          output_data.jacobian_3rd_derivatives[i] =
+            DerivativeForm<4, dim, spacedim>();
+
+      if (data.update_each & update_jacobian_pushed_forward_3rd_derivatives)
+        for (unsigned int i = 0;
+             i < output_data.jacobian_pushed_forward_3rd_derivatives.size();
+             ++i)
+          output_data.jacobian_pushed_forward_3rd_derivatives[i] =
+            Tensor<5, spacedim>();
+    }
+}
+
+
+
 template <int dim, int spacedim>
 CellSimilarity::Similarity
 MappingCartesian<dim, spacedim>::fill_fe_values(
@@ -361,51 +413,8 @@ MappingCartesian<dim, spacedim>::fill_fe_values(
           for (unsigned int j = 0; j < dim; ++j)
             output_data.jacobians[i][j][j] = data.cell_extents[j];
         }
-  // "compute" the derivative of the Jacobian at the quadrature
-  // points, which are all zero of course
-  if (data.update_each & update_jacobian_grads)
-    if (cell_similarity != CellSimilarity::translation)
-      for (unsigned int i = 0; i < output_data.jacobian_grads.size(); ++i)
-        output_data.jacobian_grads[i] = DerivativeForm<2, dim, spacedim>();
-
-  if (data.update_each & update_jacobian_pushed_forward_grads)
-    if (cell_similarity != CellSimilarity::translation)
-      for (unsigned int i = 0;
-           i < output_data.jacobian_pushed_forward_grads.size();
-           ++i)
-        output_data.jacobian_pushed_forward_grads[i] = Tensor<3, spacedim>();
-
-  // "compute" the hessian of the Jacobian at the quadrature points,
-  // which are all also zero of course
-  if (data.update_each & update_jacobian_2nd_derivatives)
-    if (cell_similarity != CellSimilarity::translation)
-      for (unsigned int i = 0; i < output_data.jacobian_2nd_derivatives.size();
-           ++i)
-        output_data.jacobian_2nd_derivatives[i] =
-          DerivativeForm<3, dim, spacedim>();
 
-  if (data.update_each & update_jacobian_pushed_forward_2nd_derivatives)
-    if (cell_similarity != CellSimilarity::translation)
-      for (unsigned int i = 0;
-           i < output_data.jacobian_pushed_forward_2nd_derivatives.size();
-           ++i)
-        output_data.jacobian_pushed_forward_2nd_derivatives[i] =
-          Tensor<4, spacedim>();
-
-  if (data.update_each & update_jacobian_3rd_derivatives)
-    if (cell_similarity != CellSimilarity::translation)
-      for (unsigned int i = 0; i < output_data.jacobian_3rd_derivatives.size();
-           ++i)
-        output_data.jacobian_3rd_derivatives[i] =
-          DerivativeForm<4, dim, spacedim>();
-
-  if (data.update_each & update_jacobian_pushed_forward_3rd_derivatives)
-    if (cell_similarity != CellSimilarity::translation)
-      for (unsigned int i = 0;
-           i < output_data.jacobian_pushed_forward_3rd_derivatives.size();
-           ++i)
-        output_data.jacobian_pushed_forward_3rd_derivatives[i] =
-          Tensor<5, spacedim>();
+  maybe_update_jacobian_derivatives(data, cell_similarity, output_data);
 
   // "compute" inverse Jacobian at the quadrature points, which are
   // all the same
@@ -481,41 +490,7 @@ MappingCartesian<dim, spacedim>::fill_fe_face_values(
           output_data.jacobians[i][d][d] = data.cell_extents[d];
       }
 
-  if (data.update_each & update_jacobian_grads)
-    for (unsigned int i = 0; i < output_data.jacobian_grads.size(); ++i)
-      output_data.jacobian_grads[i] = DerivativeForm<2, dim, spacedim>();
-
-  if (data.update_each & update_jacobian_pushed_forward_grads)
-    for (unsigned int i = 0;
-         i < output_data.jacobian_pushed_forward_grads.size();
-         ++i)
-      output_data.jacobian_pushed_forward_grads[i] = Tensor<3, spacedim>();
-
-  if (data.update_each & update_jacobian_2nd_derivatives)
-    for (unsigned int i = 0; i < output_data.jacobian_2nd_derivatives.size();
-         ++i)
-      output_data.jacobian_2nd_derivatives[i] =
-        DerivativeForm<3, dim, spacedim>();
-
-  if (data.update_each & update_jacobian_pushed_forward_2nd_derivatives)
-    for (unsigned int i = 0;
-         i < output_data.jacobian_pushed_forward_2nd_derivatives.size();
-         ++i)
-      output_data.jacobian_pushed_forward_2nd_derivatives[i] =
-        Tensor<4, spacedim>();
-
-  if (data.update_each & update_jacobian_3rd_derivatives)
-    for (unsigned int i = 0; i < output_data.jacobian_3rd_derivatives.size();
-         ++i)
-      output_data.jacobian_3rd_derivatives[i] =
-        DerivativeForm<4, dim, spacedim>();
-
-  if (data.update_each & update_jacobian_pushed_forward_3rd_derivatives)
-    for (unsigned int i = 0;
-         i < output_data.jacobian_pushed_forward_3rd_derivatives.size();
-         ++i)
-      output_data.jacobian_pushed_forward_3rd_derivatives[i] =
-        Tensor<5, spacedim>();
+  maybe_update_jacobian_derivatives(data, CellSimilarity::none, output_data);
 
   if (data.update_each & update_inverse_jacobians)
     for (unsigned int i = 0; i < output_data.inverse_jacobians.size(); ++i)
@@ -593,41 +568,7 @@ MappingCartesian<dim, spacedim>::fill_fe_subface_values(
           output_data.jacobians[i][d][d] = data.cell_extents[d];
       }
 
-  if (data.update_each & update_jacobian_grads)
-    for (unsigned int i = 0; i < output_data.jacobian_grads.size(); ++i)
-      output_data.jacobian_grads[i] = DerivativeForm<2, dim, spacedim>();
-
-  if (data.update_each & update_jacobian_pushed_forward_grads)
-    for (unsigned int i = 0;
-         i < output_data.jacobian_pushed_forward_grads.size();
-         ++i)
-      output_data.jacobian_pushed_forward_grads[i] = Tensor<3, spacedim>();
-
-  if (data.update_each & update_jacobian_2nd_derivatives)
-    for (unsigned int i = 0; i < output_data.jacobian_2nd_derivatives.size();
-         ++i)
-      output_data.jacobian_2nd_derivatives[i] =
-        DerivativeForm<3, dim, spacedim>();
-
-  if (data.update_each & update_jacobian_pushed_forward_2nd_derivatives)
-    for (unsigned int i = 0;
-         i < output_data.jacobian_pushed_forward_2nd_derivatives.size();
-         ++i)
-      output_data.jacobian_pushed_forward_2nd_derivatives[i] =
-        Tensor<4, spacedim>();
-
-  if (data.update_each & update_jacobian_3rd_derivatives)
-    for (unsigned int i = 0; i < output_data.jacobian_3rd_derivatives.size();
-         ++i)
-      output_data.jacobian_3rd_derivatives[i] =
-        DerivativeForm<4, dim, spacedim>();
-
-  if (data.update_each & update_jacobian_pushed_forward_3rd_derivatives)
-    for (unsigned int i = 0;
-         i < output_data.jacobian_pushed_forward_3rd_derivatives.size();
-         ++i)
-      output_data.jacobian_pushed_forward_3rd_derivatives[i] =
-        Tensor<5, spacedim>();
+  maybe_update_jacobian_derivatives(data, CellSimilarity::none, output_data);
 
   if (data.update_each & update_inverse_jacobians)
     for (unsigned int i = 0; i < output_data.inverse_jacobians.size(); ++i)

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.