#include <deal.II/base/std_cxx14/memory.h>
#include <deal.II/base/subscriptor.h>
#include <deal.II/base/utilities.h>
-
+#include <deal.II/base/std_cxx14/algorithm.h>
#include <boost/archive/basic_archive.hpp>
#include <boost/core/demangle.hpp>
struct RankInfo<T, typename std::enable_if<is_map_compatible<T>::value>::type>
{
static constexpr int list_rank =
- std::max(internal::RankInfo<typename T::key_type>::list_rank,
- RankInfo<typename T::mapped_type>::list_rank) +
+ std_cxx14::max(internal::RankInfo<typename T::key_type>::list_rank,
+ RankInfo<typename T::mapped_type>::list_rank) +
1;
static constexpr int map_rank =
- std::max(internal::RankInfo<typename T::key_type>::map_rank,
- RankInfo<typename T::mapped_type>::map_rank) +
+ std_cxx14::max(internal::RankInfo<typename T::key_type>::map_rank,
+ RankInfo<typename T::mapped_type>::map_rank) +
1;
};
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+#ifndef dealii__cxx14_algorithm_h
+#define dealii__cxx14_algorithm_h
+
+#include <deal.II/base/config.h>
+
+#ifdef DEAL_II_WITH_CXX14
+#include <algorithm>
+#endif
+
+DEAL_II_NAMESPACE_OPEN
+namespace std_cxx14
+{
+#ifdef DEAL_II_WITH_CXX14
+ using std::max;
+ using std::min;
+#else
+
+ template<class T>
+ constexpr const T &max(const T &a, const T &b)
+ {
+ return (a < b) ? b : a;
+ }
+
+ template<class T, class Compare>
+ constexpr const T &max(const T &a, const T &b, Compare comp)
+ {
+ return (comp(a, b)) ? b : a;
+ }
+
+ template<class T>
+ constexpr const T &min(const T &a, const T &b)
+ {
+ return (b < a) ? b : a;
+ }
+
+ template<class T, class Compare>
+ constexpr const T &min(const T &a, const T &b, Compare comp)
+ {
+ return (comp(b, a)) ? b : a;
+ }
+
+#endif
+}
+DEAL_II_NAMESPACE_CLOSE
+
+#endif // dealii__cxx14_algorithm_h
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// test the implementation of theh constexpr min/max functions
+// in case deal.II is compiled with C++14 support, we only test
+// the STL.
+
+
+#include "../tests.h"
+#include <deal.II/base/std_cxx14/algorithm.h>
+#include <deal.II/base/logstream.h>
+#include <fstream>
+#include <iomanip>
+
+constexpr bool comp (const int &a, const int &b)
+{
+ return b<a;
+}
+
+
+int main ()
+{
+ initlog();
+ deallog.threshold_double(1.e-10);
+
+ constexpr int max_1 = std_cxx14::max(0,1);
+ deallog << max_1 << std::endl;
+ constexpr int max_2 = std_cxx14::max(3,2,comp);
+ deallog << max_2 << std::endl;
+
+ constexpr int min_1 = std_cxx14::min(1,2);
+ deallog << min_1 << std::endl;
+ constexpr int min_2 = std_cxx14::min(1,2,comp);
+ deallog << min_2 << std::endl;
+
+ deallog << "OK" << std::endl;
+}
+
--- /dev/null
+
+DEAL::1
+DEAL::2
+DEAL::1
+DEAL::2
+DEAL::OK