]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Improve loop nesting.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 6 Nov 2015 22:13:43 +0000 (16:13 -0600)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Sat, 7 Nov 2015 21:01:42 +0000 (15:01 -0600)
We had a number of places of the form
  for (...tight loop...)
    if (constant condition)
      data update;

The compiler almost certainly can hoist the condition out of the loop,
but why make it this complicated for the compiler. I also find the code
easier to read because the loop is really only over a counting index
and does not carry any particular meaning at a level higher than the
if-statement.

include/deal.II/fe/fe_poly.templates.h
include/deal.II/fe/fe_poly_face.templates.h
source/fe/fe_dgp_nonparametric.cc

index da43087d1cf38cc5cc744b108f98877c5273b4c8..cd3e28dbd0eb3de0e040b22b89362fbdc4779fd1 100644 (file)
@@ -249,41 +249,44 @@ fill_fe_values (const Mapping<dim,spacedim>                                  &ma
 
   const UpdateFlags flags(fe_data.update_each);
 
-  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+  if (flags & update_values)
+    for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+      for (unsigned int i=0; i<quadrature.size(); ++i)
+        output_data.shape_values(k,i) = fe_data.shape_values[k][i];
+
+  if (flags & update_gradients && cell_similarity != CellSimilarity::translation)
+    for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+      mapping.transform (fe_data.shape_gradients[k],
+                         mapping_covariant,
+                         mapping_internal,
+                         output_data.shape_gradients[k]);
+
+  if (flags & update_hessians && cell_similarity != CellSimilarity::translation)
     {
-      if (flags & update_values)
+      for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+        mapping.transform (fe_data.shape_hessians[k],
+                           mapping_covariant_gradient,
+                           mapping_internal,
+                           output_data.shape_hessians[k]);
+
+      for (unsigned int k=0; k<this->dofs_per_cell; ++k)
         for (unsigned int i=0; i<quadrature.size(); ++i)
-          output_data.shape_values(k,i) = fe_data.shape_values[k][i];
+          for (unsigned int j=0; j<spacedim; ++j)
+            output_data.shape_hessians[k][i] -=
+              mapping_data.jacobian_pushed_forward_grads[i][j]
+              * output_data.shape_gradients[k][i][j];
+    }
 
-      if (flags & update_gradients && cell_similarity != CellSimilarity::translation)
-        mapping.transform (fe_data.shape_gradients[k],
-                           mapping_covariant,
+  if (flags & update_3rd_derivatives && cell_similarity != CellSimilarity::translation)
+    {
+      for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+        mapping.transform (fe_data.shape_3rd_derivatives[k],
+                           mapping_covariant_hessian,
                            mapping_internal,
-                           output_data.shape_gradients[k]);
-
-      if (flags & update_hessians && cell_similarity != CellSimilarity::translation)
-        {
-          mapping.transform (fe_data.shape_hessians[k],
-                             mapping_covariant_gradient,
-                             mapping_internal,
-                             output_data.shape_hessians[k]);
-
-          for (unsigned int i=0; i<quadrature.size(); ++i)
-            for (unsigned int j=0; j<spacedim; ++j)
-              output_data.shape_hessians[k][i] -=
-                mapping_data.jacobian_pushed_forward_grads[i][j]
-                * output_data.shape_gradients[k][i][j];
-        }
-
-      if (flags & update_3rd_derivatives && cell_similarity != CellSimilarity::translation)
-        {
-          mapping.transform (fe_data.shape_3rd_derivatives[k],
-                             mapping_covariant_hessian,
-                             mapping_internal,
-                             output_data.shape_3rd_derivatives[k]);
-
-          correct_third_derivatives(output_data, mapping_data, quadrature.size(), k);
-        }
+                           output_data.shape_3rd_derivatives[k]);
+
+      for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+        correct_third_derivatives(output_data, mapping_data, quadrature.size(), k);
     }
 }
 
@@ -321,45 +324,48 @@ fill_fe_face_values (const Mapping<dim,spacedim>
 
   const UpdateFlags flags(fe_data.update_each);
 
-  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+  if (flags & update_values)
+    for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+      for (unsigned int i=0; i<quadrature.size(); ++i)
+        output_data.shape_values(k,i) = fe_data.shape_values[k][i+offset];
+
+  if (flags & update_gradients)
+    for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+      mapping.transform (make_slice(fe_data.shape_gradients[k], offset, quadrature.size()),
+                         mapping_covariant,
+                         mapping_internal,
+                         output_data.shape_gradients[k]);
+
+  if (flags & update_hessians)
     {
-      if (flags & update_values)
+      for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+        mapping.transform (make_slice(fe_data.shape_hessians[k],
+                                      offset,
+                                      quadrature.size()),
+                           mapping_covariant_gradient,
+                           mapping_internal,
+                           output_data.shape_hessians[k]);
+
+      for (unsigned int k=0; k<this->dofs_per_cell; ++k)
         for (unsigned int i=0; i<quadrature.size(); ++i)
-          output_data.shape_values(k,i) = fe_data.shape_values[k][i+offset];
+          for (unsigned int j=0; j<spacedim; ++j)
+            output_data.shape_hessians[k][i] -=
+              mapping_data.jacobian_pushed_forward_grads[i][j]
+              * output_data.shape_gradients[k][i][j];
+    }
 
-      if (flags & update_gradients)
-        mapping.transform (make_slice(fe_data.shape_gradients[k], offset, quadrature.size()),
-                           mapping_covariant,
+  if (flags & update_3rd_derivatives)
+    {
+      for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+        mapping.transform (make_slice(fe_data.shape_3rd_derivatives[k],
+                                      offset,
+                                      quadrature.size()),
+                           mapping_covariant_hessian,
                            mapping_internal,
-                           output_data.shape_gradients[k]);
-
-      if (flags & update_hessians)
-        {
-          mapping.transform (make_slice(fe_data.shape_hessians[k],
-                                        offset,
-                                        quadrature.size()),
-                             mapping_covariant_gradient,
-                             mapping_internal,
-                             output_data.shape_hessians[k]);
-
-          for (unsigned int i=0; i<quadrature.size(); ++i)
-            for (unsigned int j=0; j<spacedim; ++j)
-              output_data.shape_hessians[k][i] -=
-                mapping_data.jacobian_pushed_forward_grads[i][j]
-                * output_data.shape_gradients[k][i][j];
-        }
-
-      if (flags & update_3rd_derivatives)
-        {
-          mapping.transform (make_slice(fe_data.shape_3rd_derivatives[k],
-                                        offset,
-                                        quadrature.size()),
-                             mapping_covariant_hessian,
-                             mapping_internal,
-                             output_data.shape_3rd_derivatives[k]);
-
-          correct_third_derivatives(output_data, mapping_data, quadrature.size(), k);
-        }
+                           output_data.shape_3rd_derivatives[k]);
+
+      for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+        correct_third_derivatives(output_data, mapping_data, quadrature.size(), k);
     }
 }
 
@@ -397,48 +403,53 @@ fill_fe_subface_values (const Mapping<dim,spacedim>
 
   const UpdateFlags flags(fe_data.update_each);
 
-  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+  if (flags & update_values)
+    for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+      for (unsigned int i=0; i<quadrature.size(); ++i)
+        output_data.shape_values(k,i) = fe_data.shape_values[k][i+offset];
+
+  if (flags & update_gradients)
+    for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+      mapping.transform (make_slice(fe_data.shape_gradients[k], offset, quadrature.size()),
+                         mapping_covariant,
+                         mapping_internal,
+                         output_data.shape_gradients[k]);
+
+  if (flags & update_hessians)
     {
-      if (flags & update_values)
+      for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+        mapping.transform (make_slice(fe_data.shape_hessians[k],
+                                      offset,
+                                      quadrature.size()),
+                           mapping_covariant_gradient,
+                           mapping_internal,
+                           output_data.shape_hessians[k]);
+
+      for (unsigned int k=0; k<this->dofs_per_cell; ++k)
         for (unsigned int i=0; i<quadrature.size(); ++i)
-          output_data.shape_values(k,i) = fe_data.shape_values[k][i+offset];
+          for (unsigned int j=0; j<spacedim; ++j)
+            output_data.shape_hessians[k][i] -=
+              mapping_data.jacobian_pushed_forward_grads[i][j]
+              * output_data.shape_gradients[k][i][j];
+    }
 
-      if (flags & update_gradients)
-        mapping.transform (make_slice(fe_data.shape_gradients[k], offset, quadrature.size()),
-                           mapping_covariant,
+  if (flags & update_3rd_derivatives)
+    {
+      for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+        mapping.transform (make_slice(fe_data.shape_3rd_derivatives[k],
+                                      offset,
+                                      quadrature.size()),
+                           mapping_covariant_hessian,
                            mapping_internal,
-                           output_data.shape_gradients[k]);
-
-      if (flags & update_hessians)
-        {
-          mapping.transform (make_slice(fe_data.shape_hessians[k],
-                                        offset,
-                                        quadrature.size()),
-                             mapping_covariant_gradient,
-                             mapping_internal,
-                             output_data.shape_hessians[k]);
-
-          for (unsigned int i=0; i<quadrature.size(); ++i)
-            for (unsigned int j=0; j<spacedim; ++j)
-              output_data.shape_hessians[k][i] -=
-                mapping_data.jacobian_pushed_forward_grads[i][j]
-                * output_data.shape_gradients[k][i][j];
-        }
-
-      if (flags & update_3rd_derivatives)
-        {
-          mapping.transform (make_slice(fe_data.shape_3rd_derivatives[k],
-                                        offset,
-                                        quadrature.size()),
-                             mapping_covariant_hessian,
-                             mapping_internal,
-                             output_data.shape_3rd_derivatives[k]);
-
-          correct_third_derivatives(output_data, mapping_data, quadrature.size(), k);
-        }
+                           output_data.shape_3rd_derivatives[k]);
+
+      for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+        correct_third_derivatives(output_data, mapping_data, quadrature.size(), k);
     }
 }
 
+
+
 template <class POLY, int dim, int spacedim>
 inline void
 FE_Poly<POLY,dim,spacedim>::
index 972118fcd84115902cd5859d8b9c1dda47f0b72d..2e5a78a164ec61a87e91136bd6796831e680f4ec 100644 (file)
@@ -197,13 +197,14 @@ fill_fe_subface_values (const Mapping<dim,spacedim> &,
   const unsigned int offset = subface*quadrature.size();
 
   if (fe_data.update_each & update_values)
-    for (unsigned int i=0; i<quadrature.size(); ++i)
-      {
-        for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+    {
+      for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+        for (unsigned int i=0; i<quadrature.size(); ++i)
           output_data.shape_values(k,i) = 0.;
-        for (unsigned int k=0; k<fe_data.shape_values.size(); ++k)
+      for (unsigned int k=0; k<fe_data.shape_values.size(); ++k)
+        for (unsigned int i=0; i<quadrature.size(); ++i)
           output_data.shape_values(foffset+k,i) = fe_data.shape_values[k][i+offset];
-      }
+    }
 
   Assert (!(fe_data.update_each & update_gradients), ExcNotImplemented());
   Assert (!(fe_data.update_each & update_hessians), ExcNotImplemented());
index 8bc0ff473962bd23b81bc942f0bc51ba34f5a8e1..9cb3f39a4df13a5454b3c1e797b76fe827cbda58 100644 (file)
@@ -308,15 +308,18 @@ fill_fe_values (const Mapping<dim,spacedim> &,
                                  values, grads, grad_grads,
                                  empty_vector_of_3rd_order_tensors,
                                  empty_vector_of_4th_order_tensors);
-        for (unsigned int k=0; k<this->dofs_per_cell; ++k)
-          {
-            if (fe_data.update_each & update_values)
-              output_data.shape_values[k][i] = values[k];
-            if (fe_data.update_each & update_gradients)
-              output_data.shape_gradients[k][i] = grads[k];
-            if (fe_data.update_each & update_hessians)
-              output_data.shape_hessians[k][i] = grad_grads[k];
-          }
+
+        if (fe_data.update_each & update_values)
+          for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+            output_data.shape_values[k][i] = values[k];
+
+        if (fe_data.update_each & update_gradients)
+          for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+            output_data.shape_gradients[k][i] = grads[k];
+
+        if (fe_data.update_each & update_hessians)
+          for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+            output_data.shape_hessians[k][i] = grad_grads[k];
       }
 }
 
@@ -351,15 +354,18 @@ fill_fe_face_values (const Mapping<dim,spacedim> &,
                                  values, grads, grad_grads,
                                  empty_vector_of_3rd_order_tensors,
                                  empty_vector_of_4th_order_tensors);
-        for (unsigned int k=0; k<this->dofs_per_cell; ++k)
-          {
-            if (fe_data.update_each & update_values)
-              output_data.shape_values[k][i] = values[k];
-            if (fe_data.update_each & update_gradients)
-              output_data.shape_gradients[k][i] = grads[k];
-            if (fe_data.update_each & update_hessians)
-              output_data.shape_hessians[k][i] = grad_grads[k];
-          }
+
+        if (fe_data.update_each & update_values)
+          for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+            output_data.shape_values[k][i] = values[k];
+
+        if (fe_data.update_each & update_gradients)
+          for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+            output_data.shape_gradients[k][i] = grads[k];
+
+        if (fe_data.update_each & update_hessians)
+          for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+            output_data.shape_hessians[k][i] = grad_grads[k];
       }
 }
 
@@ -395,15 +401,18 @@ fill_fe_subface_values (const Mapping<dim,spacedim> &,
                                  values, grads, grad_grads,
                                  empty_vector_of_3rd_order_tensors,
                                  empty_vector_of_4th_order_tensors);
-        for (unsigned int k=0; k<this->dofs_per_cell; ++k)
-          {
-            if (fe_data.update_each & update_values)
-              output_data.shape_values[k][i] = values[k];
-            if (fe_data.update_each & update_gradients)
-              output_data.shape_gradients[k][i] = grads[k];
-            if (fe_data.update_each & update_hessians)
-              output_data.shape_hessians[k][i] = grad_grads[k];
-          }
+
+        if (fe_data.update_each & update_values)
+          for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+            output_data.shape_values[k][i] = values[k];
+
+        if (fe_data.update_each & update_gradients)
+          for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+            output_data.shape_gradients[k][i] = grads[k];
+
+        if (fe_data.update_each & update_hessians)
+          for (unsigned int k=0; k<this->dofs_per_cell; ++k)
+            output_data.shape_hessians[k][i] = grad_grads[k];
       }
 }
 

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.