]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Reduce template contortion level by one in hopes that this helps MS VC++.
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Tue, 4 Feb 2014 21:45:30 +0000 (21:45 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Tue, 4 Feb 2014 21:45:30 +0000 (21:45 +0000)
git-svn-id: https://svn.dealii.org/trunk@32407 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/matrix_free/mapping_info.h
deal.II/include/deal.II/matrix_free/mapping_info.templates.h

index 692cf4de579286a6c819c600e0630b06accb32fd..ee1f8c1961ce6a980d9dfd3a4a12c6224967201d 100644 (file)
@@ -57,6 +57,11 @@ namespace internal
        */
       static const unsigned int n_cell_types = 1U<<n_cell_type_bits;
 
+      /**
+       * An abbreviation for the length of vector lines of the current data type.
+       */
+      static const unsigned int n_vector_elements = VectorizedArray<Number>::n_array_elements;
+
       /**
        * Empty constructor.
        */
@@ -279,7 +284,6 @@ namespace internal
 
       /**
        * Contains all the stuff that depends on the quadrature formula
-
        */
       std::vector<MappingInfoDependent> mapping_data_gen;
 
@@ -320,8 +324,8 @@ namespace internal
                              const std::pair<unsigned int,unsigned int> *cells,
                              const unsigned int  cell,
                              const unsigned int  my_q,
-                             CellType (&cell_t_prev)[VectorizedArray<Number>::n_array_elements],
-                             CellType (&cell_t)[VectorizedArray<Number>::n_array_elements],
+                             CellType (&cell_t_prev)[n_vector_elements],
+                             CellType (&cell_t)[n_vector_elements],
                              FEValues<dim,dim> &fe_values,
                              CellData          &cell_data) const;
     };
index f19e7be6c54a8247217059f1b34ff4f79fc58248..5d093bc971a398ea8a6acbcece26c34effeede5b 100644 (file)
@@ -606,13 +606,11 @@ end_set:
                                                const std::pair<unsigned int,unsigned int> *cells,
                                                const unsigned int  cell,
                                                const unsigned int  my_q,
-                                               CellType (&cell_t_prev)[VectorizedArray<Number>::n_array_elements],
-                                               CellType (&cell_t)[VectorizedArray<Number>::n_array_elements],
+                                               CellType (&cell_t_prev)[n_vector_elements],
+                                               CellType (&cell_t)[n_vector_elements],
                                                FEValues<dim,dim> &fe_val,
                                                CellData          &data) const
     {
-      const unsigned int vectorization_length =
-        VectorizedArray<Number>::n_array_elements;
       const unsigned int n_q_points = fe_val.n_quadrature_points;
       const UpdateFlags update_flags = fe_val.get_update_flags();
 
@@ -620,7 +618,7 @@ end_set:
       // not have that field here)
       const double zero_tolerance_double = data.jac_size *
                                            std::numeric_limits<double>::epsilon() * 1024.;
-      for (unsigned int j=0; j<vectorization_length; ++j)
+      for (unsigned int j=0; j<n_vector_elements; ++j)
         {
           typename dealii::Triangulation<dim>::cell_iterator
           cell_it (&tria, cells[j].first, cells[j].second);
@@ -653,7 +651,7 @@ end_set:
               if (j==0)
                 {
                   Assert (cell>0, ExcInternalError());
-                  cell_t[j] = cell_t_prev[vectorization_length-1];
+                  cell_t[j] = cell_t_prev[n_vector_elements-1];
                 }
               else
                 cell_t[j] = cell_t[j-1];
@@ -791,11 +789,11 @@ end_set:
                         data.general_jac_grad[q][d][e][f][j] = jacobian_grad[d][e][f];
                 }
             }
-        } // end loop over all entries in vectorization (vectorization_length
+        } // end loop over all entries in vectorization (n_vector_elements
       // cells)
 
       // set information for next cell
-      for (unsigned int j=0; j<vectorization_length; ++j)
+      for (unsigned int j=0; j<n_vector_elements; ++j)
         cell_t_prev[j] = cell_t[j];
     }
 

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.