]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fixed support for arbitrary map types.
authorLuca Heltai <luca.heltai@sissa.it>
Thu, 20 Jul 2017 08:58:28 +0000 (10:58 +0200)
committerLuca Heltai <luca.heltai@sissa.it>
Sun, 23 Jul 2017 14:57:28 +0000 (16:57 +0200)
include/deal.II/base/patterns_tools.h
tests/base/patterns_tools_06.cc [new file with mode: 0644]
tests/base/patterns_tools_06.output [new file with mode: 0644]
tests/base/patterns_tools_07.cc [new file with mode: 0644]
tests/base/patterns_tools_07.output [new file with mode: 0644]

index 41a7af2c031bb0b2ecfa18a913a38f7dece11499..48ae17dfa549f80208ca8cb10659ce4f2fd06308 100644 (file)
 #ifndef dealii__patterns_tools_h
 #define dealii__patterns_tools_h
 
-
 #include <deal.II/base/config.h>
-#include <deal.II/base/utilities.h>
 #include <deal.II/base/exceptions.h>
-#include <deal.II/base/subscriptor.h>
-#include <deal.II/base/point.h>
 #include <deal.II/base/parameter_handler.h>
+#include <deal.II/base/point.h>
 #include <deal.II/base/std_cxx14/memory.h>
+#include <deal.II/base/subscriptor.h>
+#include <deal.II/base/utilities.h>
 
 #include <boost/core/demangle.hpp>
 
-#include <map>
-#include <vector>
 #include <array>
-#include <string>
-#include <memory>
-#include <sstream>
 #include <deque>
-#include <forward_list>
 #include <list>
 #include <map>
-#include <queue>
 #include <set>
-#include <stack>
+#include <sstream>
 #include <string>
-#include <tuple>
+#include <type_traits>
 #include <unordered_map>
 #include <unordered_set>
 #include <utility>
 #include <vector>
-#include <type_traits>
-
 
 DEAL_II_NAMESPACE_OPEN
 
 /**
- * Namespace for a few classes and functions that act on values and patterns, and
- * allow to convert from non elementary types to strings and vice versa.
+ * Namespace for a few classes and functions that act on values and patterns,
+ * and allow to convert from non elementary types to strings and vice versa.
  *
  * A typical usage of these tools is in the following example:
  *
@@ -85,8 +75,8 @@ DEAL_II_NAMESPACE_OPEN
  * namespace. Internally it uses the internal::RankInfo<T> class to decide how
  * many different separators are required to convert the given type to a string.
  *
- * For example, to write vectors of vectors, the default is to use "," for the first
- * (inner) separator, and ";" for the second (outer) separator, i.e.
+ * For example, to write vectors of vectors, the default is to use "," for the
+ * first (inner) separator, and ";" for the second (outer) separator, i.e.
  *
  * @code
  * std::vector<std::vector<unsigned int>> vec;
@@ -109,21 +99,21 @@ DEAL_II_NAMESPACE_OPEN
  * default_map_separator {{":"  ,  "="  ,  "@"  ,   "#"}};
  * @endcode
  *
- * When one needs a mixture of Patterns::List and Patterns::Map types, their RankInfo is
- * computed by taking
- * the maximum of the vector_rank of the Key and of the Value type, so that, for example,
- * it is possible to have the following
+ * When one needs a mixture of Patterns::List and Patterns::Map types, their
+ * RankInfo is computed by taking the maximum of the vector_rank of the Key and
+ * of the Value type, so that, for example, it is possible to have the following
  * @code
  * ... // Build compare class
  * std::map<std::vector<unsigned int>, std::vector<double>, compare> map;
  *
- * map = convert<decltype(map)>::to_value("1,2,3 : 5.0,6.0,7.0  ; 8,9,10 : 11.0,12.0,13.0");
+ * map = convert<decltype(map)>::to_value("1,2,3 : 5.0,6.0,7.0  ; 8,9,10 :
+ * 11.0,12.0,13.0");
  *
  * @endcode
  *
- * Some non elementary types are supported, like Point(), or std::complex<double>.
- * If you wish to support more types, you have to specialize the Convert struct as well
- * as the RankInfo struct.
+ * Some non elementary types are supported, like Point(), or
+ * std::complex<double>. If you wish to support more types, you have to
+ * specialize the Convert struct as well as the RankInfo struct.
  *
  * @ingroup input
  * @author Luca Heltai, 2017
@@ -131,14 +121,15 @@ DEAL_II_NAMESPACE_OPEN
 namespace PatternTools
 {
   /**
-   * Converter class. This class is used to generate strings and Patterns associated to
-   * the given type, and to convert from a string to the given type and viceversa.
+   * Converter class. This class is used to generate strings and Patterns
+   * associated to the given type, and to convert from a string to the given
+   * type and viceversa.
    *
-   * The second template parameter is used internally to allow for advanced SFINAE
-   * (substitution failure is not an error) tricks used to specialise this class for
-   * arbitrary STL containers and maps.
+   * The second template parameter is used internally to allow for advanced
+   * SFINAE (substitution failure is not an error) tricks used to specialise
+   * this class for arbitrary STL containers and maps.
    */
-  template <class T, class Enable=void>
+  template <class T, class Enable = void>
   struct Convert
   {
 
@@ -148,23 +139,24 @@ namespace PatternTools
      */
     static std::unique_ptr<Patterns::PatternBase> to_pattern() = delete;
 
-
     /**
-     * Return a string containing a textual version of the variable s. Use the pattern
-     * passed to perform the conversion, or create and use a default one.
+     * Return a string containing a textual version of the variable s. Use the
+     * pattern passed to perform the conversion, or create and use a default
+     * one.
      *
      * The created string is
      */
     static std::string to_string(const T &s,
-                                 const std::unique_ptr<Patterns::PatternBase> &p = Convert<T>::to_pattern()) = delete;
-
+                                 const std::unique_ptr<Patterns::PatternBase>
+                                 &p = Convert<T>::to_pattern()) = delete;
 
     /**
      * Convert a string to a value, using the given pattern. Use the pattern
      * passed to perform the conversion, or create and use a default one.
      */
     static T to_value(const std::string &s,
-                      const std::unique_ptr<Patterns::PatternBase> &p = Convert<T>::to_pattern()) = delete;
+                      const std::unique_ptr<Patterns::PatternBase> &p =
+                        Convert<T>::to_pattern()) = delete;
   };
 
   /**
@@ -177,11 +169,12 @@ namespace PatternTools
    * PatternTools::Convert<T>::to_pattern(), but a custom one can be used.
    */
   template <class ParameterType>
-  void add_parameter(const std::string           &entry,
-                     ParameterType               &parameter,
-                     ParameterHandler            &prm,
-                     const std::string           &documentation = std::string(),
-                     const Patterns::PatternBase &pattern = *Convert<ParameterType>::to_pattern());
+  void add_parameter(const std::string &entry,
+                     ParameterType &parameter,
+                     ParameterHandler &prm,
+                     const std::string &documentation = std::string(),
+                     const Patterns::PatternBase &pattern =
+                       *Convert<ParameterType>::to_pattern());
 
   /**
    * @addtogroup Exceptions
@@ -191,14 +184,15 @@ namespace PatternTools
   /**
    * Exception.
    */
-  DeclException2 (ExcNoMatch,
-                  std::string, Patterns::PatternBase &,
-                  << "The string " << arg1 << " does not match the pattern \""
-                  << arg2.description() << "\"");
+  DeclException2(ExcNoMatch,
+                 std::string,
+                 Patterns::PatternBase &,
+                 << "The string " << arg1 << " does not match the pattern \""
+                 << arg2.description()
+                 << "\"");
   //@}
 }
 
-
 // ---------------------- inline and template functions --------------------
 namespace PatternTools
 {
@@ -214,9 +208,10 @@ namespace PatternTools
      * with a Patterns::List pattern or with a Patterns::Map pattern.
      *
      * Objects like Point() or std::complex<double> are vector-likes, and
-     * have vector_rank 1. Elementary types, like `int`, `unsigned int`, `double`, etc.
-     * have vector_rank 0. `std::vector`, `std::list` and in general containers have rank
-     * equal to 1 + vector_rank of the contained type. Similarly for map types.
+     * have vector_rank 1. Elementary types, like `int`, `unsigned int`,
+     * `double`, etc. have vector_rank 0. `std::vector`, `std::list` and in
+     * general containers have rank equal to 1 + vector_rank of the contained
+     * type. Similarly for map types.
      *
      * A class with list_rank::value = 0 is either elementary or a
      * map. A class with map_rank::value = 0 is either a List compatible
@@ -229,7 +224,7 @@ namespace PatternTools
      *
      * @author Luca Heltai, 2017
      */
-    template <class T, class Enable=void>
+    template <class T, class Enable = void>
     struct RankInfo
     {
       static constexpr int list_rank = 0;
@@ -245,12 +240,16 @@ namespace PatternTools
     static std::unique_ptr<Patterns::PatternBase> to_pattern()
     {
       if (std::is_integral<T>::value)
-        return std_cxx14::make_unique<Patterns::Integer>(std::numeric_limits<T>::min(), std::numeric_limits<T>::max());
+        return std_cxx14::make_unique<Patterns::Integer>(
+                 std::numeric_limits<T>::min(), std::numeric_limits<T>::max());
       else if (std::is_floating_point<T>::value)
-        return std_cxx14::make_unique<Patterns::Double>(-std::numeric_limits<T>::max(), std::numeric_limits<T>::max());
+        return std_cxx14::make_unique<Patterns::Double>(
+                 -std::numeric_limits<T>::max(), std::numeric_limits<T>::max());
     }
 
-    static std::string to_string(const T &value, const std::unique_ptr<Patterns::PatternBase> &p = Convert<T>::to_pattern())
+    static std::string to_string(const T &value,
+                                 const std::unique_ptr<Patterns::PatternBase>
+                                 &p = Convert<T>::to_pattern())
     {
       std::string str;
       if (std::is_same<T, unsigned char>() || std::is_same<T, char>())
@@ -262,26 +261,28 @@ namespace PatternTools
     }
 
     static T to_value(const std::string &s,
-                      const std::unique_ptr<Patterns::PatternBase> &p = Convert<T>::to_pattern())
+                      const std::unique_ptr<Patterns::PatternBase> &p =
+                        Convert<T>::to_pattern())
     {
       AssertThrow(p->match(s), ExcNoMatch(s, *p));
       std::istringstream is(s);
       T value;
-      if (std::is_same<T,unsigned char>::value || std::is_same<T,char>::value)
+      if (std::is_same<T, unsigned char>::value || std::is_same<T, char>::value)
         {
           int i;
           is >> i;
-          value=i;
+          value = i;
         }
       else
         is >> value;
-      AssertThrow(!is.fail(), ExcMessage("Failed to convert from \"" + s + "\" to the type \""
-                                         +boost::core::demangle(typeid(T).name()) + "\""));
+      AssertThrow(!is.fail(),
+                  ExcMessage("Failed to convert from \"" + s +
+                             "\" to the type \"" +
+                             boost::core::demangle(typeid(T).name()) + "\""));
       return value;
     }
   };
 
-
   namespace internal
   {
     const std::array<std::string, 4> default_list_separator {{","  ,  ";"  ,  "|"  ,   "%"}};
@@ -293,14 +294,10 @@ namespace PatternTools
     template <typename... Args> struct is_stl_container<std::vector            <Args...>> : std::true_type {};
     template <typename... Args> struct is_stl_container<std::deque             <Args...>> : std::true_type {};
     template <typename... Args> struct is_stl_container<std::list              <Args...>> : std::true_type {};
-    template <typename... Args> struct is_stl_container<std::forward_list      <Args...>> : std::true_type {};
     template <typename... Args> struct is_stl_container<std::set               <Args...>> : std::true_type {};
     template <typename... Args> struct is_stl_container<std::multiset          <Args...>> : std::true_type {};
     template <typename... Args> struct is_stl_container<std::unordered_set     <Args...>> : std::true_type {};
     template <typename... Args> struct is_stl_container<std::unordered_multiset<Args...>> : std::true_type {};
-    template <typename... Args> struct is_stl_container<std::stack             <Args...>> : std::true_type {};
-    template <typename... Args> struct is_stl_container<std::queue             <Args...>> : std::true_type {};
-    template <typename... Args> struct is_stl_container<std::priority_queue    <Args...>> : std::true_type {};
 
     template <typename T>       struct is_stl_map : std::false_type {};
     template <typename... Args> struct is_stl_map<std::map               <Args...>> : std::true_type {};
@@ -309,38 +306,46 @@ namespace PatternTools
     template <typename... Args> struct is_stl_map<std::unordered_multimap<Args...>> : std::true_type {};
   }
 
-  //type trait to use the implementation type traits as well as decay the type
-  template <typename T> struct is_stl_container
+  // type trait to use the implementation type traits as well as decay the type
+  template <typename T>
+  struct is_stl_container
   {
-    static constexpr bool const value = internal::is_stl_container<typename std::decay<T>::type>::value;
+    static constexpr bool const value =
+      internal::is_stl_container<typename std::decay<T>::type>::value;
   };
 
-
-  template <typename T> struct is_stl_map
+  template <typename T>
+  struct is_stl_map
   {
-    static constexpr bool const value = internal::is_stl_map<typename std::decay<T>::type>::value;
+    static constexpr bool const value =
+      internal::is_stl_map<typename std::decay<T>::type>::value;
   };
 
-
   namespace internal
   {
     // Rank of vector types
     template <class T>
-    struct RankInfo<T, typename std::enable_if<is_stl_container<T>::value>::type >
+    struct RankInfo<T,
+      typename std::enable_if<is_stl_container<T>::value>::type>
     {
-      static constexpr int list_rank = RankInfo<typename T::value_type>::list_rank + 1;
-      static constexpr int map_rank = RankInfo<typename T::value_type>::map_rank;
+      static constexpr int list_rank =
+        RankInfo<typename T::value_type>::list_rank + 1;
+      static constexpr int map_rank =
+        RankInfo<typename T::value_type>::map_rank;
     };
 
-
     // Rank of map types
     template <class T>
-    struct RankInfo<T, typename std::enable_if<is_stl_map<T>::value>::type >
+    struct RankInfo<T, typename std::enable_if<is_stl_map<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) + 1;
-      static constexpr int map_rank = std::max(internal::RankInfo<typename T::key_type>::map_rank,
-                                               RankInfo<typename T::mapped_type>::map_rank) + 1;
+      static constexpr int list_rank =
+        std::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) +
+        1;
     };
 
     // Rank of Tensor types
@@ -352,47 +357,55 @@ namespace PatternTools
     };
 
     template <int dim, class Number>
-    struct RankInfo<Point<dim,Number>>:RankInfo<Tensor<1,dim,Number>> {};
+    struct RankInfo<Point<dim, Number>> : RankInfo<Tensor<1, dim, Number>>
+    {
+    };
 
     // Rank of complex types
     template <class Number>
     struct RankInfo<std::complex<Number>>
     {
-      static constexpr int list_rank = RankInfo<Number>::list_rank+1;
+      static constexpr int list_rank = RankInfo<Number>::list_rank + 1;
       static constexpr int map_rank = RankInfo<Number>::map_rank;
     };
   }
 
-
   // stl containers
   template <class T>
   struct Convert<T, typename std::enable_if<is_stl_container<T>::value>::type>
   {
     static std::unique_ptr<Patterns::PatternBase> to_pattern()
     {
-      static_assert(internal::RankInfo<T>::list_rank>0, "Cannot use this class for non List-compatible types.");
-      return std_cxx14::make_unique<Patterns::List>(*Convert<typename T::value_type>::to_pattern(),
-                                                    0, std::numeric_limits<unsigned int>::max(),
-                                                    internal::default_list_separator[internal::RankInfo<T>::list_rank-1]);
+      static_assert(internal::RankInfo<T>::list_rank > 0,
+                    "Cannot use this class for non List-compatible types.");
+      return std_cxx14::make_unique<Patterns::List>(
+               *Convert<typename T::value_type>::to_pattern(),
+               0,
+               std::numeric_limits<unsigned int>::max(),
+               internal::default_list_separator[internal::RankInfo<T>::list_rank - 1]);
     }
 
     static std::string to_string(const T &t,
-                                 const std::unique_ptr<Patterns::PatternBase> &pattern = Convert<T>::to_pattern())
+                                 const std::unique_ptr<Patterns::PatternBase>
+                                 &pattern = Convert<T>::to_pattern())
     {
 
       auto p = dynamic_cast<const Patterns::List *>(pattern.get());
-      AssertThrow(p, ExcMessage("I need a List pattern to convert a string to a List type."));
+      AssertThrow(
+        p,
+        ExcMessage(
+          "I need a List pattern to convert a string to a List type."));
       auto base_p = p->get_base_pattern().clone();
       std::vector<std::string> vec(t.size());
 
-      unsigned int i=0;
+      unsigned int i = 0;
       for (const auto &ti : t)
         vec[i++] = Convert<typename T::value_type>::to_string(ti, base_p);
 
       std::string s;
       if (vec.size() > 0)
         s = vec[0];
-      for (unsigned int i=1; i<vec.size(); ++i)
+      for (unsigned int i = 1; i < vec.size(); ++i)
         s += p->get_separator() + " " + vec[i];
 
       AssertThrow(pattern->match(s), ExcNoMatch(s, *p));
@@ -403,125 +416,146 @@ namespace PatternTools
      * Convert a string to a value, using the given pattern, or a default one.
      */
     static T to_value(const std::string &s,
-                      const std::unique_ptr<Patterns::PatternBase> &pattern = Convert<T>::to_pattern())
+                      const std::unique_ptr<Patterns::PatternBase> &pattern =
+                        Convert<T>::to_pattern())
     {
 
       AssertThrow(pattern->match(s), ExcNoMatch(s, *pattern));
 
       auto p = dynamic_cast<const Patterns::List *>(pattern.get());
-      AssertThrow(p, ExcMessage("I need a List pattern to convert a string to a List type."));
+      AssertThrow(
+        p,
+        ExcMessage(
+          "I need a List pattern to convert a string to a List type."));
 
       auto base_p = p->get_base_pattern().clone();
       T t;
 
-      auto v = Utilities::split_string_list(s,p->get_separator());
+      auto v = Utilities::split_string_list(s, p->get_separator());
       for (const auto &str : v)
-        t.insert(t.end(), Convert<typename T::value_type>::to_value(str, base_p));
+        t.insert(t.end(),
+                 Convert<typename T::value_type>::to_value(str, base_p));
 
       return t;
     }
   };
 
-
   // stl maps
   template <class T>
   struct Convert<T, typename std::enable_if<is_stl_map<T>::value>::type>
   {
     static std::unique_ptr<Patterns::PatternBase> to_pattern()
     {
-      static_assert(internal::RankInfo<T>::list_rank>0, "Cannot use this class for non List-compatible types.");
-      static_assert(internal::RankInfo<T>::map_rank>0, "Cannot use this class for non Map-compatible types.");
-      return std_cxx14::make_unique<Patterns::Map>(*Convert<typename T::key_type>::to_pattern(),
-                                                   *Convert<typename T::mapped_type>::to_pattern(),
-                                                   0, std::numeric_limits<unsigned int>::max(),
-                                                   internal::default_list_separator[internal::RankInfo<T>::list_rank-1],
-                                                   internal::default_map_separator[internal::RankInfo<T>::map_rank-1]
-                                                  );
+      static_assert(internal::RankInfo<T>::list_rank > 0,
+                    "Cannot use this class for non List-compatible types.");
+      static_assert(internal::RankInfo<T>::map_rank > 0,
+                    "Cannot use this class for non Map-compatible types.");
+      return std_cxx14::make_unique<Patterns::Map>(
+               *Convert<typename T::key_type>::to_pattern(),
+               *Convert<typename T::mapped_type>::to_pattern(),
+               0,
+               std::numeric_limits<unsigned int>::max(),
+               internal::default_list_separator[internal::RankInfo<T>::list_rank - 1],
+               internal::default_map_separator[internal::RankInfo<T>::map_rank - 1]);
     }
 
     static std::string to_string(const T &t,
-                                 const std::unique_ptr<Patterns::PatternBase> &pattern = Convert<T>::to_pattern())
+                                 const std::unique_ptr<Patterns::PatternBase>
+                                 &pattern = Convert<T>::to_pattern())
     {
 
       auto p = dynamic_cast<const Patterns::Map *>(pattern.get());
-      AssertThrow(p, ExcMessage("I need a Map pattern to convert a string to a List type."));
+      AssertThrow(
+        p,
+        ExcMessage("I need a Map pattern to convert a string to a List type."));
       auto key_p = p->get_key_pattern().clone();
       auto val_p = p->get_value_pattern().clone();
       std::vector<std::string> vec(t.size());
 
-      unsigned int i=0;
+      unsigned int i = 0;
       for (const auto &ti : t)
         vec[i++] =
           Convert<typename T::key_type>::to_string(ti.first, key_p) +
-          p->get_key_value_separator()+
+          p->get_key_value_separator() +
           Convert<typename T::mapped_type>::to_string(ti.second, val_p);
 
       std::string s;
       if (vec.size() > 0)
         s = vec[0];
-      for (unsigned int i=1; i<vec.size(); ++i)
+      for (unsigned int i = 1; i < vec.size(); ++i)
         s += p->get_separator() + " " + vec[i];
 
       AssertThrow(p->match(s), ExcNoMatch(s, *p));
       return s;
     }
 
-    static T  to_value(const std::string &s,
-                       const std::unique_ptr<Patterns::PatternBase> &pattern = Convert<T>::to_pattern())
+    static T to_value(const std::string &s,
+                      const std::unique_ptr<Patterns::PatternBase> &pattern =
+                        Convert<T>::to_pattern())
     {
 
       AssertThrow(pattern->match(s), ExcNoMatch(s, *pattern));
 
       auto p = dynamic_cast<const Patterns::Map *>(pattern.get());
-      AssertThrow(p, ExcMessage("I need a Map pattern to convert a string to a List type."));
+      AssertThrow(
+        p,
+        ExcMessage("I need a Map pattern to convert a string to a List type."));
 
       auto key_p = p->get_key_pattern().clone();
       auto val_p = p->get_value_pattern().clone();
       T t;
 
-      auto v = Utilities::split_string_list(s,p->get_separator());
+      auto v = Utilities::split_string_list(s, p->get_separator());
       for (const auto &str : v)
         {
-          auto key_val = Utilities::split_string_list(str, p->get_key_value_separator());
+          auto key_val =
+            Utilities::split_string_list(str, p->get_key_value_separator());
           AssertDimension(key_val.size(), 2);
-          t[Convert<typename T::key_type>::to_value(key_val[0], key_p)] =
-            Convert<typename T::mapped_type>::to_value(key_val[1]);
+          t.insert(std::make_pair(
+                     Convert<typename T::key_type>::to_value(key_val[0], key_p),
+                     Convert<typename T::mapped_type>::to_value(key_val[1])));
         }
 
       return t;
     }
   };
 
-
   // Tensors
   template <int rank, int dim, class Number>
-  struct Convert<Tensor<rank,dim,Number>>
+  struct Convert<Tensor<rank, dim, Number>>
   {
-    typedef Tensor<rank,dim,Number> T;
+    typedef Tensor<rank, dim, Number> T;
     static std::unique_ptr<Patterns::PatternBase> to_pattern()
     {
-      static_assert(internal::RankInfo<T>::list_rank>0, "Cannot use this class for non List-compatible types.");
-      return std_cxx14::make_unique<Patterns::List>(*Convert<typename T::value_type>::to_pattern(),
-                                                    dim, dim,
-                                                    internal::default_list_separator[internal::RankInfo<T>::list_rank-1]);
+      static_assert(internal::RankInfo<T>::list_rank > 0,
+                    "Cannot use this class for non List-compatible types.");
+      return std_cxx14::make_unique<Patterns::List>(
+               *Convert<typename T::value_type>::to_pattern(),
+               dim,
+               dim,
+               internal::default_list_separator[internal::RankInfo<T>::list_rank - 1]);
     }
 
     static std::string to_string(const T &t,
-                                 const std::unique_ptr<Patterns::PatternBase> &pattern = Convert<T>::to_pattern())
+                                 const std::unique_ptr<Patterns::PatternBase>
+                                 &pattern = Convert<T>::to_pattern())
     {
 
       auto p = dynamic_cast<const Patterns::List *>(pattern.get());
-      AssertThrow(p, ExcMessage("I need a List pattern to convert a string to a List type."));
+      AssertThrow(
+        p,
+        ExcMessage(
+          "I need a List pattern to convert a string to a List type."));
       auto base_p = p->get_base_pattern().clone();
       std::vector<std::string> vec(dim);
 
-      for (unsigned int i=0; i<dim; ++i)
+      for (unsigned int i = 0; i < dim; ++i)
         vec[i] = Convert<typename T::value_type>::to_string(t[i], base_p);
 
       std::string s;
       if (vec.size() > 0)
         s = vec[0];
-      for (unsigned int i=1; i<vec.size(); ++i)
+      for (unsigned int i = 1; i < vec.size(); ++i)
         s += p->get_separator() + " " + vec[i];
 
       AssertThrow(p->match(s), ExcNoMatch(s, *p));
@@ -529,53 +563,58 @@ namespace PatternTools
     }
 
     static T to_value(const std::string &s,
-                      const std::unique_ptr<Patterns::PatternBase> &pattern = Convert<T>::to_pattern())
+                      const std::unique_ptr<Patterns::PatternBase> &pattern =
+                        Convert<T>::to_pattern())
     {
 
       AssertThrow(pattern->match(s), ExcNoMatch(s, *pattern));
 
       auto p = dynamic_cast<const Patterns::List *>(pattern.get());
-      AssertThrow(p, ExcMessage("I need a List pattern to convert a string to a List type."));
+      AssertThrow(
+        p,
+        ExcMessage(
+          "I need a List pattern to convert a string to a List type."));
 
       auto base_p = p->get_base_pattern().clone();
       T t;
 
-      auto v = Utilities::split_string_list(s,p->get_separator());
-      unsigned int i=0;
+      auto v = Utilities::split_string_list(s, p->get_separator());
+      unsigned int i = 0;
       for (const auto &str : v)
-        t[i++] =  Convert<typename T::value_type>::to_value(str, base_p);
+        t[i++] = Convert<typename T::value_type>::to_value(str, base_p);
 
       return t;
     }
   };
 
-
   // Points
   template <int dim, class Number>
-  struct Convert<Point<dim,Number>>
+  struct Convert<Point<dim, Number>>
   {
 
-    typedef Point<dim,Number> T;
+    typedef Point<dim, Number> T;
 
     static std::unique_ptr<Patterns::PatternBase> to_pattern()
     {
-      return Convert<Tensor<1,dim,Number>>::to_pattern();
+      return Convert<Tensor<1, dim, Number>>::to_pattern();
     }
 
     static std::string to_string(const T &t,
-                                 const std::unique_ptr<Patterns::PatternBase> &pattern = Convert<T>::to_pattern())
+                                 const std::unique_ptr<Patterns::PatternBase>
+                                 &pattern = Convert<T>::to_pattern())
     {
-      return Convert<Tensor<1,dim,Number>>::to_string(Tensor<1,dim,Number>(t),pattern);
+      return Convert<Tensor<1, dim, Number>>::to_string(
+               Tensor<1, dim, Number>(t), pattern);
     }
 
     static T to_value(const std::string &s,
-                      const std::unique_ptr<Patterns::PatternBase> &pattern = Convert<T>::to_pattern())
+                      const std::unique_ptr<Patterns::PatternBase> &pattern =
+                        Convert<T>::to_pattern())
     {
-      return T(Convert<Tensor<1,dim,Number>>::to_value(s,pattern));
+      return T(Convert<Tensor<1, dim, Number>>::to_value(s, pattern));
     }
   };
 
-
   // Complex numbers
   template <class Number>
   struct Convert<std::complex<Number>>
@@ -584,18 +623,25 @@ namespace PatternTools
 
     static std::unique_ptr<Patterns::PatternBase> to_pattern()
     {
-      static_assert(internal::RankInfo<T>::list_rank>0, "Cannot use this class for non List-compatible types.");
-      return std_cxx14::make_unique<Patterns::List>(*Convert<typename T::value_type>::to_pattern(),
-                                                    2, 2,
-                                                    internal::default_list_separator[internal::RankInfo<T>::list_rank-1]);
+      static_assert(internal::RankInfo<T>::list_rank > 0,
+                    "Cannot use this class for non List-compatible types.");
+      return std_cxx14::make_unique<Patterns::List>(
+               *Convert<typename T::value_type>::to_pattern(),
+               2,
+               2,
+               internal::default_list_separator[internal::RankInfo<T>::list_rank - 1]);
     }
 
     static std::string to_string(const T &t,
-                                 const std::unique_ptr<Patterns::PatternBase> &pattern = Convert<T>::to_pattern())
+                                 const std::unique_ptr<Patterns::PatternBase>
+                                 &pattern = Convert<T>::to_pattern())
     {
 
       auto p = dynamic_cast<const Patterns::List *>(pattern.get());
-      AssertThrow(p, ExcMessage("I need a List pattern to convert a string to a List type."));
+      AssertThrow(
+        p,
+        ExcMessage(
+          "I need a List pattern to convert a string to a List type."));
       auto base_p = p->get_base_pattern().clone();
       std::string s =
         Convert<typename T::value_type>::to_string(t.real(), base_p) +
@@ -610,17 +656,21 @@ namespace PatternTools
      * Convert a string to a value, using the given pattern, or a default one.
      */
     static T to_value(const std::string &s,
-                      const std::unique_ptr<Patterns::PatternBase> &pattern = Convert<T>::to_pattern())
+                      const std::unique_ptr<Patterns::PatternBase> &pattern =
+                        Convert<T>::to_pattern())
     {
 
       AssertThrow(pattern->match(s), ExcNoMatch(s, *pattern));
 
       auto p = dynamic_cast<const Patterns::List *>(pattern.get());
-      AssertThrow(p, ExcMessage("I need a List pattern to convert a string to a List type."));
+      AssertThrow(
+        p,
+        ExcMessage(
+          "I need a List pattern to convert a string to a List type."));
 
       auto base_p = p->get_base_pattern().clone();
 
-      auto v = Utilities::split_string_list(s,p->get_separator());
+      auto v = Utilities::split_string_list(s, p->get_separator());
       AssertDimension(v.size(), 2);
       T t(Convert<typename T::value_type>::to_value(v[0], base_p),
           Convert<typename T::value_type>::to_value(v[1], base_p));
@@ -628,7 +678,6 @@ namespace PatternTools
     }
   };
 
-
   // Strings
   template <>
   struct Convert<std::string>
@@ -641,26 +690,27 @@ namespace PatternTools
     }
 
     static std::string to_string(const T &t,
-                                 const std::unique_ptr<Patterns::PatternBase> &pattern = Convert<T>::to_pattern())
+                                 const std::unique_ptr<Patterns::PatternBase>
+                                 &pattern = Convert<T>::to_pattern())
     {
-      AssertThrow(pattern->match(t), ExcNoMatch(t,*pattern));
+      AssertThrow(pattern->match(t), ExcNoMatch(t, *pattern));
       return t;
     }
 
     static T to_value(const std::string &s,
-                      const std::unique_ptr<Patterns::PatternBase> &pattern = Convert<T>::to_pattern())
+                      const std::unique_ptr<Patterns::PatternBase> &pattern =
+                        Convert<T>::to_pattern())
     {
-      AssertThrow(pattern->match(s), ExcNoMatch(s,*pattern));
+      AssertThrow(pattern->match(s), ExcNoMatch(s, *pattern));
       return s;
     }
   };
 
-
   template <class ParameterType>
-  void add_parameter(const std::string           &entry,
-                     ParameterType                &parameter,
-                     ParameterHandler            &prm,
-                     const std::string           &documentation,
+  void add_parameter(const std::string &entry,
+                     ParameterType &parameter,
+                     ParameterHandler &prm,
+                     const std::string &documentation,
                      const Patterns::PatternBase &pattern)
   {
 
@@ -668,18 +718,21 @@ namespace PatternTools
                   "You tried to add a parameter using a type "
                   "that is const. Use a non-const type.");
 
-    prm.declare_entry(entry, PatternTools::Convert<ParameterType>::to_string(parameter, pattern.clone()),
-                      pattern, documentation);
+    prm.declare_entry(entry,
+                      PatternTools::Convert<ParameterType>::to_string(
+                        parameter, pattern.clone()),
+                      pattern,
+                      documentation);
 
-    auto action = [&] (const std::string &val)
+    auto action = [&](const std::string &val)
     {
-      parameter = PatternTools::Convert<ParameterType>::to_value(val, pattern.clone());
+      parameter =
+        PatternTools::Convert<ParameterType>::to_value(val, pattern.clone());
     };
     prm.add_action(entry, action);
   }
 }
 
-
 DEAL_II_NAMESPACE_CLOSE
 
 #endif
diff --git a/tests/base/patterns_tools_06.cc b/tests/base/patterns_tools_06.cc
new file mode 100644 (file)
index 0000000..70deeba
--- /dev/null
@@ -0,0 +1,74 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2005 - 2015 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.
+//
+// ---------------------------------------------------------------------
+
+
+#include "../tests.h"
+#include <deal.II/base/patterns_tools.h>
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/point.h>
+#include <deal.II/base/numbers.h>
+#include <boost/core/demangle.hpp>
+
+#include <memory>
+
+using namespace dealii;
+using namespace PatternTools;
+
+// Try conversion on complex map types
+
+template<class T>
+void test(T t)
+{
+  auto p = Convert<T>::to_pattern();
+  deallog << "Pattern  : " << p->description() << std::endl;
+  auto s = Convert<T>::to_string(t);
+  deallog << "To String: " << s << std::endl;
+  deallog << "To value : " << Convert<T>::to_string(Convert<T>::to_value(s)) << std::endl;
+
+}
+
+int main()
+{
+  initlog();
+
+  std::map<int, double>                 t0;
+  std::unordered_map<int, double>       t1;
+  std::multimap<int, double>            t2;
+  std::unordered_multimap<int, double>  t3;
+
+  auto p = std::make_pair(5,1.0);
+  auto p2 = std::make_pair(5,2.0);
+  auto p3 = std::make_pair(1,3.0);
+
+  t0.insert(p);
+  t1.insert(p);
+  t2.insert(p);
+  t3.insert(p);
+
+  t0.insert(p3);
+  t1.insert(p3);
+  t2.insert(p3);
+  t3.insert(p3);
+
+  t2.insert(p2);
+  t3.insert(p2);
+
+  test(t0);
+  test(t1);
+  test(t2);
+  test(t3);
+
+  return 0;
+}
diff --git a/tests/base/patterns_tools_06.output b/tests/base/patterns_tools_06.output
new file mode 100644 (file)
index 0000000..f0eea70
--- /dev/null
@@ -0,0 +1,13 @@
+
+DEAL::Pattern  : [Map of <[Integer range -2147483648...2147483647 (inclusive)]>:<[Double -MAX_DOUBLE...MAX_DOUBLE (inclusive)]> of length 0...4294967295 (inclusive)]
+DEAL::To String: 1:3.000000, 5:1.000000
+DEAL::To value : 1:3.000000, 5:1.000000
+DEAL::Pattern  : [Map of <[Integer range -2147483648...2147483647 (inclusive)]>:<[Double -MAX_DOUBLE...MAX_DOUBLE (inclusive)]> of length 0...4294967295 (inclusive)]
+DEAL::To String: 5:1.000000, 1:3.000000
+DEAL::To value : 1:3.000000, 5:1.000000
+DEAL::Pattern  : [Map of <[Integer range -2147483648...2147483647 (inclusive)]>:<[Double -MAX_DOUBLE...MAX_DOUBLE (inclusive)]> of length 0...4294967295 (inclusive)]
+DEAL::To String: 1:3.000000, 5:1.000000, 5:2.000000
+DEAL::To value : 1:3.000000, 5:1.000000, 5:2.000000
+DEAL::Pattern  : [Map of <[Integer range -2147483648...2147483647 (inclusive)]>:<[Double -MAX_DOUBLE...MAX_DOUBLE (inclusive)]> of length 0...4294967295 (inclusive)]
+DEAL::To String: 1:3.000000, 5:1.000000, 5:2.000000
+DEAL::To value : 1:3.000000, 5:1.000000, 5:2.000000
diff --git a/tests/base/patterns_tools_07.cc b/tests/base/patterns_tools_07.cc
new file mode 100644 (file)
index 0000000..718f1a4
--- /dev/null
@@ -0,0 +1,80 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2005 - 2015 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.
+//
+// ---------------------------------------------------------------------
+
+
+#include "../tests.h"
+#include <deal.II/base/patterns_tools.h>
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/point.h>
+#include <deal.II/base/numbers.h>
+#include <boost/core/demangle.hpp>
+
+#include <memory>
+
+using namespace dealii;
+using namespace PatternTools;
+
+// Try conversion on arbitrary container types
+
+template<class T>
+void test(T t)
+{
+  auto p = Convert<T>::to_pattern();
+  deallog << "Pattern  : " << p->description() << std::endl;
+  auto s = Convert<T>::to_string(t);
+  deallog << "To String: " << s << std::endl;
+  deallog << "To value : " << Convert<T>::to_string(Convert<T>::to_value(s)) << std::endl;
+
+}
+
+int main()
+{
+  initlog();
+
+  std::vector            <unsigned int> t00;
+  std::deque             <unsigned int> t01;
+  std::list              <unsigned int> t02;
+  std::set               <unsigned int> t03;
+  std::multiset          <unsigned int> t04;
+  std::unordered_set     <unsigned int> t05;
+  std::unordered_multiset<unsigned int> t06;
+
+  t00.insert(t00.end(), 0);
+  t01.insert(t01.end(), 1);
+  t02.insert(t02.end(), 2);
+  t03.insert(t03.end(), 3);
+  t04.insert(t04.end(), 4);
+  t05.insert(t05.end(), 5);
+  t06.insert(t06.end(), 6);
+
+
+  t00.insert(t00.end(), 1);
+  t01.insert(t01.end(), 2);
+  t02.insert(t02.end(), 3);
+  t03.insert(t03.end(), 4);
+  t04.insert(t04.end(), 5);
+  t05.insert(t05.end(), 6);
+  t06.insert(t06.end(), 7);
+
+  test(t00);
+  test(t01);
+  test(t02);
+  test(t03);
+  test(t04);
+  test(t05);
+  test(t06);
+
+  return 0;
+}
diff --git a/tests/base/patterns_tools_07.output b/tests/base/patterns_tools_07.output
new file mode 100644 (file)
index 0000000..c8ca05d
--- /dev/null
@@ -0,0 +1,22 @@
+
+DEAL::Pattern  : [List of <[Integer]> of length 0...4294967295 (inclusive)]
+DEAL::To String: 0, 1
+DEAL::To value : 0, 1
+DEAL::Pattern  : [List of <[Integer]> of length 0...4294967295 (inclusive)]
+DEAL::To String: 1, 2
+DEAL::To value : 1, 2
+DEAL::Pattern  : [List of <[Integer]> of length 0...4294967295 (inclusive)]
+DEAL::To String: 2, 3
+DEAL::To value : 2, 3
+DEAL::Pattern  : [List of <[Integer]> of length 0...4294967295 (inclusive)]
+DEAL::To String: 3, 4
+DEAL::To value : 3, 4
+DEAL::Pattern  : [List of <[Integer]> of length 0...4294967295 (inclusive)]
+DEAL::To String: 4, 5
+DEAL::To value : 4, 5
+DEAL::Pattern  : [List of <[Integer]> of length 0...4294967295 (inclusive)]
+DEAL::To String: 5, 6
+DEAL::To value : 6, 5
+DEAL::Pattern  : [List of <[Integer]> of length 0...4294967295 (inclusive)]
+DEAL::To String: 6, 7
+DEAL::To value : 7, 6

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.