]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Make std_cxx20::ranges::iota_view compatible with std::ranges::iota_view
authorDaniel Arndt <arndtd@ornl.gov>
Mon, 4 May 2020 14:20:45 +0000 (10:20 -0400)
committerWolfgang Bangerth <bangerth@colostate.edu>
Thu, 7 May 2020 12:42:40 +0000 (06:42 -0600)
examples/step-69/step-69.cc
include/deal.II/base/geometry_info.h
include/deal.II/base/std_cxx20/iota_view.h
include/deal.II/fe/fe_values.h

index 5bf418310fbd31fd1e15c60a7e466af7484d06e4..b2ebd869d2d4bedd81051f735556338a8dbea3cc 100644 (file)
@@ -1254,14 +1254,17 @@ namespace Step69
       TimerOutput::Scope scope(computing_timer,
                                "offline_data - compute |c_ij|, and n_ij");
 
-      const std_cxx20::ranges::iota_view<unsigned int> indices(
+      const std_cxx20::ranges::iota_view<unsigned int, unsigned int> indices(
         0, n_locally_relevant);
 
       const auto on_subranges = //
-        [&](std_cxx20::ranges::iota_view<unsigned int>::iterator       i1,
-            const std_cxx20::ranges::iota_view<unsigned int>::iterator i2) {
+        [&](
+          std_cxx20::ranges::iota_view<unsigned int, unsigned int>::iterator i1,
+          const std_cxx20::ranges::iota_view<unsigned int,
+                                             unsigned int>::iterator i2) {
           for (const auto row_index :
-               std_cxx20::ranges::iota_view<unsigned int>(*i1, *i2))
+               std_cxx20::ranges::iota_view<unsigned int, unsigned int>(*i1,
+                                                                        *i2))
             {
               // First column-loop: we compute and store the entries of the
               // matrix norm_matrix and write normalized entries into the
@@ -1863,10 +1866,10 @@ namespace Step69
     const auto &n_locally_owned    = offline_data->n_locally_owned;
     const auto &n_locally_relevant = offline_data->n_locally_relevant;
 
-    const std_cxx20::ranges::iota_view<unsigned int> indices_owned(
-      0, n_locally_owned);
-    const std_cxx20::ranges::iota_view<unsigned int> indices_relevant(
-      0, n_locally_relevant);
+    const std_cxx20::ranges::iota_view<unsigned int, unsigned int>
+      indices_owned(0, n_locally_owned);
+    const std_cxx20::ranges::iota_view<unsigned int, unsigned int>
+      indices_relevant(0, n_locally_relevant);
 
     const auto &sparsity = offline_data->sparsity_pattern;
 
@@ -1920,10 +1923,13 @@ namespace Step69
                                "time_stepping - 1 compute d_ij");
 
       const auto on_subranges = //
-        [&](std_cxx20::ranges::iota_view<unsigned int>::iterator       i1,
-            const std_cxx20::ranges::iota_view<unsigned int>::iterator i2) {
+        [&](
+          std_cxx20::ranges::iota_view<unsigned int, unsigned int>::iterator i1,
+          const std_cxx20::ranges::iota_view<unsigned int,
+                                             unsigned int>::iterator i2) {
           for (const auto i :
-               std_cxx20::ranges::iota_view<unsigned int>(*i1, *i2))
+               std_cxx20::ranges::iota_view<unsigned int, unsigned int>(*i1,
+                                                                        *i2))
             {
               const auto U_i = gather(U, i);
 
@@ -2016,12 +2022,15 @@ namespace Step69
       // locally.
 
       const auto on_subranges = //
-        [&](std_cxx20::ranges::iota_view<unsigned int>::iterator       i1,
-            const std_cxx20::ranges::iota_view<unsigned int>::iterator i2) {
+        [&](
+          std_cxx20::ranges::iota_view<unsigned int, unsigned int>::iterator i1,
+          const std_cxx20::ranges::iota_view<unsigned int,
+                                             unsigned int>::iterator i2) {
           double tau_max_on_subrange = std::numeric_limits<double>::infinity();
 
           for (const auto i :
-               std_cxx20::ranges::iota_view<unsigned int>(*i1, *i2))
+               std_cxx20::ranges::iota_view<unsigned int, unsigned int>(*i1,
+                                                                        *i2))
             {
               double d_sum = 0.;
 
@@ -2102,8 +2111,10 @@ namespace Step69
                                "time_stepping - 3 perform update");
 
       const auto on_subranges =
-        [&](std_cxx20::ranges::iota_view<unsigned int>::iterator       i1,
-            const std_cxx20::ranges::iota_view<unsigned int>::iterator i2) {
+        [&](
+          std_cxx20::ranges::iota_view<unsigned int, unsigned int>::iterator i1,
+          const std_cxx20::ranges::iota_view<unsigned int,
+                                             unsigned int>::iterator i2) {
           for (const auto i : boost::make_iterator_range(i1, i2))
             {
               Assert(i < n_locally_owned, ExcInternalError());
@@ -2341,7 +2352,8 @@ namespace Step69
     const auto &n_locally_owned     = offline_data->n_locally_owned;
 
     const auto indices =
-      std_cxx20::ranges::iota_view<unsigned int>(0, n_locally_owned);
+      std_cxx20::ranges::iota_view<unsigned int, unsigned int>(0,
+                                                               n_locally_owned);
 
     // We define the r_i_max and r_i_min in the current MPI process as
     // atomic doubles in order to avoid race conditions between threads:
@@ -2352,8 +2364,10 @@ namespace Step69
     // global maxima and minima of the gradients.
     {
       const auto on_subranges = //
-        [&](std_cxx20::ranges::iota_view<unsigned int>::iterator       i1,
-            const std_cxx20::ranges::iota_view<unsigned int>::iterator i2) {
+        [&](
+          std_cxx20::ranges::iota_view<unsigned int, unsigned int>::iterator i1,
+          const std_cxx20::ranges::iota_view<unsigned int,
+                                             unsigned int>::iterator i2) {
           double r_i_max_on_subrange = 0.;
           double r_i_min_on_subrange = std::numeric_limits<double>::infinity();
 
@@ -2437,8 +2451,10 @@ namespace Step69
 
     {
       const auto on_subranges = //
-        [&](std_cxx20::ranges::iota_view<unsigned int>::iterator       i1,
-            const std_cxx20::ranges::iota_view<unsigned int>::iterator i2) {
+        [&](
+          std_cxx20::ranges::iota_view<unsigned int, unsigned int>::iterator i1,
+          const std_cxx20::ranges::iota_view<unsigned int,
+                                             unsigned int>::iterator i2) {
           for (const auto i : boost::make_iterator_range(i1, i2))
             {
               Assert(i < n_locally_owned, ExcInternalError());
index d8ea7c5e4e0f8030821882fbf408cde8d5a173f6..04f6646c5e3ef9585249d0f492d34b23ad5b5e57 100644 (file)
@@ -1946,7 +1946,7 @@ struct GeometryInfo
    * taking on all valid indices for faces (zero and one in 1d, zero
    * through three in 2d, and zero through 5 in 3d).
    */
-  static std_cxx20::ranges::iota_view<unsigned int>
+  static std_cxx20::ranges::iota_view<unsigned int, unsigned int>
   face_indices();
 
   /**
@@ -1977,7 +1977,7 @@ struct GeometryInfo
    * Here, we are looping over all vertices of all cells, with `vertex_index`
    * taking on all valid indices.
    */
-  static std_cxx20::ranges::iota_view<unsigned int>
+  static std_cxx20::ranges::iota_view<unsigned int, unsigned int>
   vertex_indices();
 
   /**
@@ -2832,19 +2832,19 @@ GeometryInfo<0>::vertex_indices()
 
 
 template <int dim>
-inline std_cxx20::ranges::iota_view<unsigned int>
+inline std_cxx20::ranges::iota_view<unsigned int, unsigned int>
 GeometryInfo<dim>::face_indices()
 {
-  return std_cxx20::ranges::iota_view<unsigned int>(0U, faces_per_cell);
+  return {0U, faces_per_cell};
 }
 
 
 
 template <int dim>
-inline std_cxx20::ranges::iota_view<unsigned int>
+inline std_cxx20::ranges::iota_view<unsigned int, unsigned int>
 GeometryInfo<dim>::vertex_indices()
 {
-  return std_cxx20::ranges::iota_view<unsigned int>(0U, vertices_per_cell);
+  return {0U, vertices_per_cell};
 }
 
 
index 8c2b34606b6224db4bdff8d5962c144166fa8f68..8c1e471943f7fa3024cc03b83730a155bf2cc212 100644 (file)
 
 #include <deal.II/base/config.h>
 
-#include <boost/range/irange.hpp>
-
+#if defined __cpp_lib_ranges && __cpp_lib_ranges >= 201911
+#  include <ranges>
+#else
+#  include <boost/range/irange.hpp>
+#endif
 
 DEAL_II_NAMESPACE_OPEN
 
@@ -26,6 +29,10 @@ namespace std_cxx20
 {
   namespace ranges
   {
+#if defined __cpp_lib_ranges && __cpp_lib_ranges >= 201911
+    template <typename IncrementableType, typename BoundType>
+    using iota_view = std::ranges::iota_view<IncrementableType, BoundType>;
+#else
     /**
      * A poor-man's implementation of std::ranges::iota_view using
      * boost's integer_range class. The two classes are not completely
@@ -38,8 +45,9 @@ namespace std_cxx20
      * class can be found at
      * https://en.cppreference.com/w/cpp/ranges/iota_view .
      */
-    template <typename IncrementableType>
+    template <typename IncrementableType, typename /*BoundType*/>
     using iota_view = boost::integer_range<IncrementableType>;
+#endif
   } // namespace ranges
 } // namespace std_cxx20
 
index 59aa92b46c8213e4e2cfe602b78d9876ae4abf4f..109a2c036a473d5cfa5bcc4667c39208ca672174 100644 (file)
@@ -3030,7 +3030,7 @@ public:
    * `q_point` taking on all valid indices for quadrature points, as defined
    * by the quadrature rule passed to `fe_values`.
    */
-  std_cxx20::ranges::iota_view<unsigned int>
+  std_cxx20::ranges::iota_view<unsigned int, unsigned int>
   quadrature_point_indices() const;
 
   /**
@@ -5561,10 +5561,10 @@ FEValuesBase<dim, spacedim>::dof_indices_ending_at(
 
 
 template <int dim, int spacedim>
-inline std_cxx20::ranges::iota_view<unsigned int>
+inline std_cxx20::ranges::iota_view<unsigned int, unsigned int>
 FEValuesBase<dim, spacedim>::quadrature_point_indices() const
 {
-  return std_cxx20::ranges::iota_view<unsigned int>(0U, n_quadrature_points);
+  return {0U, n_quadrature_points};
 }
 
 

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.