From bd498e0df6343cc8581293fbbd2fba85f7fac3f5 Mon Sep 17 00:00:00 2001
From: Markus Blatt <markus@dr-blatt.de>
Date: Wed, 16 Apr 2014 16:14:24 +0200
Subject: [PATCH] [g++-4.4] Fixes missing std::initializer_list::const_iterator
 for g++-4.4.
MIME-Version: 1.0
Content-Type: text/plain; charset=UTF-8
Content-Transfer-Encoding: 8bit

With g++-4.4 I got the following errors:

fmatrix.hh:106: error: no type named ‘const_iterator’ in ‘class std::initializer_list<std::initializer_list<int> >’

Looking into the header there is indeed the const_iterator (demanded by the standard)
missing. Therefore I refactored the code to use std::copy_n and thus am omitting the iterator
type.

At the same time I made the actual copying safer for the case that the initializer_list
is of different size and people compile with -DNDEBUG. The assertion for the size is still
left in there. Remove the assert or the size computation at your discretion.
---
 dune/common/fmatrix.hh | 19 +++++++++----------
 dune/common/fvector.hh | 15 +++++++--------
 2 files changed, 16 insertions(+), 18 deletions(-)

diff --git a/dune/common/fmatrix.hh b/dune/common/fmatrix.hh
index 43d715ed2..d57af2dba 100644
--- a/dune/common/fmatrix.hh
+++ b/dune/common/fmatrix.hh
@@ -7,6 +7,7 @@
 #include <cmath>
 #include <cstddef>
 #include <iostream>
+#include <algorithm>
 #include <initializer_list>
 
 #include <dune/common/exceptions.hh>
@@ -100,21 +101,19 @@ namespace Dune
      */
     FieldMatrix (std::initializer_list<std::initializer_list<K> > const &ll)
     {
-      assert(ll.size() == rows);
-      size_t i = 0;
-      for (typename std::initializer_list<std::initializer_list<K> >::
-             const_iterator lit = ll.begin(); lit != ll.end(); ++lit)
-        _data[i++] = *lit;
+      assert(ll.size() == rows); // Actually, this is not needed any more!
+      std::copy_n(ll.begin(), std::min(static_cast<std::size_t>(ROWS),
+                                       ll.size()),
+                 _data.begin());
     }
 
     /** \brief Constructor initializing the matrix from a list of vector
      */
     FieldMatrix(std::initializer_list<Dune::FieldVector<K, cols> > const &l) {
-      assert(l.size() == rows);
-      size_t i = 0;
-      for (typename std::initializer_list<Dune::FieldVector<K, cols> >::
-             const_iterator lit = l.begin(); lit != l.end(); ++lit)
-        _data[i++] = *lit;
+      assert(l.size() == rows); // Actually, this is not needed any more!
+      std::copy_n(l.begin(), std::min(static_cast<std::size_t>(ROWS),
+                                      l.size()),
+                 _data.begin());
     }
 
     //===== assignment
diff --git a/dune/common/fvector.hh b/dune/common/fvector.hh
index bc212c1e5..606e6be5c 100644
--- a/dune/common/fvector.hh
+++ b/dune/common/fvector.hh
@@ -11,6 +11,7 @@
 #include <cstring>
 #include <utility>
 #include <initializer_list>
+#include <algorithm>
 
 #include <dune/common/std/constexpr.hh>
 
@@ -119,11 +120,10 @@ namespace Dune {
 
     FieldVector (std::initializer_list<K> const &l)
     {
-      assert(l.size() == dimension);
-      size_t i = 0;
-      for (typename std::initializer_list<K>::const_iterator it = l.begin();
-           it != l.end(); ++it)
-        _data[i++] = *it;
+      assert(l.size() == dimension);// Actually, this is not needed any more!
+      std::copy_n(l.begin(), std::min(static_cast<std::size_t>(dimension),
+                                      l.size()),
+                 _data.begin());
     }
 
     /**
@@ -142,9 +142,8 @@ namespace Dune {
     {
       DUNE_UNUSED_PARAMETER(dummy);
       // do a run-time size check, for the case that x is not a FieldVector
-      assert(x.size() == SIZE);
-      for (size_type i = 0; i<SIZE; i++)
-        _data[i] = x[i];
+      assert(x.size() == SIZE); // Actually this is not needed any more!
+      std::copy_n(x.begin(), std::min(static_cast<std::size_t>(SIZE),x.size()), _data.begin());
     }
 
     //! Constructor making vector with identical coordinates
-- 
GitLab