]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Avoid workaround for zero-dimensional C-style arrays
authorDaniel Arndt <arndtd@ornl.gov>
Tue, 19 Nov 2019 23:49:12 +0000 (18:49 -0500)
committerDaniel Arndt <arndtd@ornl.gov>
Wed, 20 Nov 2019 05:38:18 +0000 (00:38 -0500)
include/deal.II/base/tensor_product_polynomials.h
source/base/tensor_product_polynomials.cc
source/base/tensor_product_polynomials_bubbles.cc
source/base/tensor_product_polynomials_const.cc

index 546904b53fa51f34849c56965551db186c530b4d..9423032102faff10bb748647928516de03001ddd 100644 (file)
@@ -231,10 +231,9 @@ protected:
    * dimensional polynomials for each space direction, given the index
    * <i>i</i>.
    */
-  // fix to avoid compiler warnings about zero length arrays
   void
-  compute_index(const unsigned int i,
-                unsigned int (&indices)[(dim > 0 ? dim : 1)]) const;
+  compute_index(const unsigned int             i,
+                std::array<unsigned int, dim> &indices) const;
 
   /**
    * TensorProductPolynomialsBubbles has a TensorProductPolynomials class
@@ -410,8 +409,8 @@ private:
    * <tt>i</tt>.
    */
   void
-  compute_index(const unsigned int i,
-                unsigned int (&indices)[(dim > 0 ? dim : 1)]) const;
+  compute_index(const unsigned int             i,
+                std::array<unsigned int, dim> &indices) const;
 
   /**
    * Given the input to the constructor, compute <tt>n_pols</tt>.
@@ -479,10 +478,10 @@ TensorProductPolynomials<dim, PolynomialType>::compute_derivative(
   const unsigned int i,
   const Point<dim> & p) const
 {
-  unsigned int indices[dim];
+  std::array<unsigned int, dim> indices;
   compute_index(i, indices);
 
-  double v[dim][5];
+  std::array<std::array<double, 5>, dim> v;
   {
     std::vector<double> tmp(5);
     for (unsigned int d = 0; d < dim; ++d)
@@ -607,7 +606,7 @@ Tensor<order, dim>
 AnisotropicPolynomials<dim>::compute_derivative(const unsigned int i,
                                                 const Point<dim> & p) const
 {
-  unsigned int indices[dim];
+  std::array<unsigned int, dim> indices;
   compute_index(i, indices);
 
   std::vector<std::vector<double>> v(dim, std::vector<double>(order + 1));
index d10d94b28c1a677a0e7ae50cf57a3b28b58f4150..521fd1732355a4f3b14a9deb4b1745db1ae1f067 100644 (file)
@@ -38,35 +38,38 @@ namespace internal
     compute_tensor_index(const unsigned int,
                          const unsigned int,
                          const unsigned int,
-                         unsigned int (&)[dim])
+                         std::array<unsigned int, dim> &)
     {
       Assert(false, ExcNotImplemented());
     }
 
+    template <>
     inline void
-    compute_tensor_index(const unsigned int n,
-                         const unsigned int,
-                         const unsigned int,
-                         unsigned int (&indices)[1])
+    compute_tensor_index<1>(const unsigned int n,
+                            const unsigned int,
+                            const unsigned int,
+                            std::array<unsigned int, 1> &indices)
     {
       indices[0] = n;
     }
 
+    template <>
     inline void
-    compute_tensor_index(const unsigned int n,
-                         const unsigned int n_pols_0,
-                         const unsigned int,
-                         unsigned int (&indices)[2])
+    compute_tensor_index<2>(const unsigned int n,
+                            const unsigned int n_pols_0,
+                            const unsigned int,
+                            std::array<unsigned int, 2> &indices)
     {
       indices[0] = n % n_pols_0;
       indices[1] = n / n_pols_0;
     }
 
+    template <>
     inline void
-    compute_tensor_index(const unsigned int n,
-                         const unsigned int n_pols_0,
-                         const unsigned int n_pols_1,
-                         unsigned int (&indices)[3])
+    compute_tensor_index<3>(const unsigned int           n,
+                            const unsigned int           n_pols_0,
+                            const unsigned int           n_pols_1,
+                            std::array<unsigned int, 3> &indices)
     {
       indices[0] = n % n_pols_0;
       indices[1] = (n / n_pols_0) % n_pols_1;
@@ -80,15 +83,15 @@ namespace internal
 template <int dim, typename PolynomialType>
 inline void
 TensorProductPolynomials<dim, PolynomialType>::compute_index(
-  const unsigned int i,
-  unsigned int (&indices)[(dim > 0 ? dim : 1)]) const
+  const unsigned int             i,
+  std::array<unsigned int, dim> &indices) const
 {
   Assert(i < Utilities::fixed_power<dim>(polynomials.size()),
          ExcInternalError());
-  internal::compute_tensor_index(index_map[i],
-                                 polynomials.size(),
-                                 polynomials.size(),
-                                 indices);
+  internal::compute_tensor_index<dim>(index_map[i],
+                                      polynomials.size(),
+                                      polynomials.size(),
+                                      indices);
 }
 
 
@@ -98,7 +101,7 @@ void
 TensorProductPolynomials<dim, PolynomialType>::output_indices(
   std::ostream &out) const
 {
-  unsigned int ix[(dim > 0) ? dim : 1];
+  std::array<unsigned int, dim> ix;
   for (unsigned int i = 0; i < this->n(); ++i)
     {
       compute_index(i, ix);
@@ -146,7 +149,7 @@ TensorProductPolynomials<dim, PolynomialType>::compute_value(
 {
   Assert(dim > 0, ExcNotImplemented());
 
-  unsigned int indices[dim];
+  std::array<unsigned int, dim> indices;
   compute_index(i, indices);
 
   double value = 1.;
@@ -164,14 +167,14 @@ TensorProductPolynomials<dim, PolynomialType>::compute_grad(
   const unsigned int i,
   const Point<dim> & p) const
 {
-  unsigned int indices[dim];
+  std::array<unsigned int, dim> indices;
   compute_index(i, indices);
 
   // compute values and
   // uni-directional derivatives at
   // the given point in each
   // co-ordinate direction
-  double v[dim][2];
+  std::array<std::array<double, 2>, dim> v;
   {
     std::vector<double> tmp(2);
     for (unsigned int d = 0; d < dim; ++d)
@@ -212,10 +215,10 @@ TensorProductPolynomials<dim, PolynomialType>::compute_grad_grad(
   const unsigned int i,
   const Point<dim> & p) const
 {
-  unsigned int indices[(dim > 0) ? dim : 1];
+  std::array<unsigned int, dim> indices;
   compute_index(i, indices);
 
-  double v[dim][3];
+  std::array<std::array<double, 3>, dim> v;
   {
     std::vector<double> tmp(3);
     for (unsigned int d = 0; d < dim; ++d)
@@ -456,8 +459,8 @@ AnisotropicPolynomials<dim>::AnisotropicPolynomials(
 template <int dim>
 void
 AnisotropicPolynomials<dim>::compute_index(
-  const unsigned int i,
-  unsigned int (&indices)[(dim > 0 ? dim : 1)]) const
+  const unsigned int             i,
+  std::array<unsigned int, dim> &indices) const
 {
 #ifdef DEBUG
   unsigned int n_poly = 1;
@@ -470,15 +473,15 @@ AnisotropicPolynomials<dim>::compute_index(
     {
     }
   else if (dim == 1)
-    internal::compute_tensor_index(i,
-                                   polynomials[0].size(),
-                                   0 /*not used*/,
-                                   indices);
+    internal::compute_tensor_index<dim>(i,
+                                        polynomials[0].size(),
+                                        0 /*not used*/,
+                                        indices);
   else
-    internal::compute_tensor_index(i,
-                                   polynomials[0].size(),
-                                   polynomials[1].size(),
-                                   indices);
+    internal::compute_tensor_index<dim>(i,
+                                        polynomials[0].size(),
+                                        polynomials[1].size(),
+                                        indices);
 }
 
 
@@ -488,7 +491,7 @@ double
 AnisotropicPolynomials<dim>::compute_value(const unsigned int i,
                                            const Point<dim> & p) const
 {
-  unsigned int indices[(dim > 0) ? dim : 1];
+  std::array<unsigned int, dim> indices;
   compute_index(i, indices);
 
   double value = 1.;
@@ -504,7 +507,7 @@ Tensor<1, dim>
 AnisotropicPolynomials<dim>::compute_grad(const unsigned int i,
                                           const Point<dim> & p) const
 {
-  unsigned int indices[(dim > 0) ? dim : 1];
+  std::array<unsigned int, dim> indices;
   compute_index(i, indices);
 
   // compute values and
@@ -532,7 +535,7 @@ Tensor<2, dim>
 AnisotropicPolynomials<dim>::compute_grad_grad(const unsigned int i,
                                                const Point<dim> & p) const
 {
-  unsigned int indices[(dim > 0) ? dim : 1];
+  std::array<unsigned int, dim> indices;
   compute_index(i, indices);
 
   std::vector<std::vector<double>> v(dim, std::vector<double>(3));
@@ -627,7 +630,7 @@ AnisotropicPolynomials<dim>::evaluate(
       // one-dimensional indices of
       // this particular tensor
       // product polynomial
-      unsigned int indices[(dim > 0) ? dim : 1];
+      std::array<unsigned int, dim> indices;
       compute_index(i, indices);
 
       if (update_values)
index 3912dd03c5856ccb42757976fd5e603850b88897..59551672c7a33b2892bc812de252db82e09d8562 100644 (file)
@@ -30,7 +30,7 @@ template <int dim>
 void
 TensorProductPolynomialsBubbles<dim>::output_indices(std::ostream &out) const
 {
-  unsigned int ix[dim];
+  std::array<unsigned int, dim> ix;
   for (unsigned int i = 0; i < tensor_polys.n(); ++i)
     {
       tensor_polys.compute_index(i, ix);
index 80efe0e01380853259c4fe271f963a493718cba0..f452a77d2f2b8d8cb7c0a456ae8cf4428e2d01d5 100644 (file)
@@ -30,7 +30,7 @@ template <int dim>
 void
 TensorProductPolynomialsConst<dim>::output_indices(std::ostream &out) const
 {
-  unsigned int ix[dim];
+  std::array<unsigned int, dim> ix;
   for (unsigned int i = 0; i < tensor_polys.n(); ++i)
     {
       tensor_polys.compute_index(i, ix);

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.