]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add Patterns::Map. Also fix a number of minor issues with Patterns::List.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Wed, 1 Aug 2012 16:15:38 +0000 (16:15 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Wed, 1 Aug 2012 16:15:38 +0000 (16:15 +0000)
git-svn-id: https://svn.dealii.org/trunk@25745 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/include/deal.II/base/parameter_handler.h
deal.II/source/base/parameter_handler.cc
tests/bits/parameter_handler_13.cc [new file with mode: 0644]
tests/bits/parameter_handler_13/cmp/generic [new file with mode: 0644]
tests/bits/parameter_handler_13/prm [new file with mode: 0644]

index b6d19a3febd570b4fc20a1578bbe8142715e9283..0cbea7afb43fa220d42f441faaadea0a4ad3a549 100644 (file)
@@ -237,6 +237,12 @@ enabled due to a missing include file in file
 <a name="specific"></a>
 <h3>Specific improvements</h3>
 
+<ol>
+<li> New: The Patterns::Map pattern allows to describe maps from keys
+to values in input files.
+<br>
+(Wolfgang Bangerth, 2012/08/01)
+
 <ol>
 <li> Fixed: DoFTools::make_zero_boundary_constraints now also works for parallel distributed triangulations.
 <br>
index 66754146a26072d3578cc614121adc4f6624b323..e7643d852da77d55e5249467e5af13e1782e0277 100644 (file)
@@ -650,7 +650,7 @@ namespace Patterns
       const unsigned int min_elements;
 
                                        /**
-                                        * Minimum number of elements
+                                        * Maximum number of elements
                                         * the list must have.
                                         */
       const unsigned int max_elements;
@@ -661,6 +661,143 @@ namespace Patterns
       static const char* description_init;
   };
 
+
+  /**
+   * This pattern matches a list of
+   * comma-separated values each of which
+   * denotes a pair of key and value. Both key and value
+   * have to match a pattern given to the
+   * constructor. For each entry of the map,
+   * parameters have to be entered in the form
+   * <code>key: value</code>. In other words, a
+   * map is described in the form
+   * <code>key1: value1, key2: value2, key3: value3, ...</code>.
+   *
+   * With two additional
+   * parameters, the number of elements this
+   * list has to have can be specified. If
+   * none is specified, the map may have
+   * zero or more entries.
+   */
+  class Map : public PatternBase
+  {
+  public:
+    /**
+     * Maximal integer value. If
+     * the numeric_limits class
+     * is available use this
+     * information to obtain the
+     * extremal values, otherwise
+     * set it so that this class
+     * understands that all values
+     * are allowed.
+     */
+    static const unsigned int max_int_value;
+
+    /**
+     * Constructor. Take the
+     * given parameter as the
+     * specification of valid
+     * elements of the list.
+     *
+     * The two other arguments can
+     * be used to denote minimal
+     * and maximal allowable
+     * lengths of the list.
+     */
+    Map (const PatternBase  &key_pattern,
+         const PatternBase  &value_pattern,
+         const unsigned int  min_elements = 0,
+         const unsigned int  max_elements = max_int_value);
+
+    /**
+     * Destructor.
+     */
+    virtual ~Map ();
+
+    /**
+     * Return <tt>true</tt> if the
+     * string is a comma-separated
+     * list of strings each of
+     * which match the pattern
+     * given to the constructor.
+     */
+    virtual bool match (const std::string &test_string) const;
+
+    /**
+     * Return a description of
+     * the pattern that valid
+     * strings are expected to
+     * match.
+     */
+    virtual std::string description () const;
+
+    /**
+     * Return a copy of the
+     * present object, which is
+     * newly allocated on the
+     * heap. Ownership of that
+     * object is transferred to
+     * the caller of this
+     * function.
+     */
+    virtual PatternBase * clone () const;
+
+    /**
+     * Creates new object if the start of
+     * description matches
+     * description_init.  Ownership of that
+     * object is transferred to the caller
+     * of this function.
+     */
+    static Map* create (const std::string& description);
+
+    /**
+     * Determine an estimate for
+     * the memory consumption (in
+     * bytes) of this object.
+     */
+    std::size_t memory_consumption () const;
+
+    /** @addtogroup Exceptions
+     * @{ */
+
+    /**
+     * Exception.
+     */
+    DeclException2 (ExcInvalidRange,
+        int, int,
+        << "The values " << arg1 << " and " << arg2
+        << " do not form a valid range.");
+    //@}
+  private:
+    /**
+     * Copy of the patterns that
+     * each key and each value of the map has
+     * to satisfy.
+     */
+    PatternBase *key_pattern;
+    PatternBase *value_pattern;
+
+    /**
+     * Minimum number of elements
+     * the list must have.
+     */
+    const unsigned int min_elements;
+
+    /**
+     * Maximum number of elements
+     * the list must have.
+     */
+    const unsigned int max_elements;
+
+    /**
+     * Initial part of description
+     */
+    static const char* description_init;
+  };
+
+
                                    /**
                                     * This class is much like the
                                     * Selection class, but it
index bfa6be40712e6aa0916f9d22b637beab348b4e1c..1a628bfea40a8b066aae5f6ac6b69a63b733e1e0 100644 (file)
@@ -559,7 +559,7 @@ namespace Patterns
   std::size_t
   List::memory_consumption () const
   {
-    return (sizeof(PatternBase) +
+    return (sizeof(*this) +
             MemoryConsumption::memory_consumption(*pattern));
   }
 
@@ -577,7 +577,7 @@ namespace Patterns
       std::string str;
       std::getline(is, str, '>');
 
-      std::auto_ptr<PatternBase> base_pattern (pattern_factory(str));
+      std_cxx1x::shared_ptr<PatternBase> base_pattern (pattern_factory(str));
 
       is.ignore(strlen(" of length "));
       if(!(is >> min_elements))
@@ -585,7 +585,7 @@ namespace Patterns
 
       is.ignore(strlen("..."));
       if(!(is >> max_elements))
-        return new List(*base_pattern);
+        return new List(*base_pattern, min_elements);
 
       return new List(*base_pattern, min_elements, max_elements);
     }
@@ -595,6 +595,183 @@ namespace Patterns
 
 
 
+  const unsigned int Map::max_int_value =
+#ifdef HAVE_STD_NUMERIC_LIMITS
+          std::numeric_limits<unsigned int>::max();
+#else
+          numbers::invalid_unsigned_int;
+#endif
+
+
+  const char* Map::description_init = "[Map";
+
+
+  Map::Map (const PatternBase  &p_key,
+            const PatternBase  &p_value,
+            const unsigned int  min_elements,
+            const unsigned int  max_elements)
+                  :
+                  key_pattern (p_key.clone()),
+                  value_pattern (p_value.clone()),
+                  min_elements (min_elements),
+                  max_elements (max_elements)
+  {
+    Assert (min_elements <= max_elements,
+            ExcInvalidRange (min_elements, max_elements));
+  }
+
+
+
+  Map::~Map ()
+  {
+    delete key_pattern;
+    key_pattern = 0;
+
+    delete value_pattern;
+    value_pattern = 0;
+  }
+
+
+
+  bool Map::match (const std::string &test_string_list) const
+  {
+    std::string tmp = test_string_list;
+    std::vector<std::string> split_list;
+    split_list.reserve (std::count (tmp.begin(), tmp.end(), ',')+1);
+
+                                     // first split the input list at comma sites
+    while (tmp.length() != 0)
+      {
+        std::string map_entry;
+        map_entry = tmp;
+
+        if (map_entry.find(",") != std::string::npos)
+          {
+            map_entry.erase (map_entry.find(","), std::string::npos);
+            tmp.erase (0, tmp.find(",")+1);
+          }
+        else
+          tmp = "";
+
+        while ((map_entry.length() != 0) &&
+               (std::isspace (map_entry[0])))
+          map_entry.erase (0,1);
+
+        while (std::isspace (map_entry[map_entry.length()-1]))
+          map_entry.erase (map_entry.length()-1, 1);
+
+        split_list.push_back (map_entry);
+      };
+
+    if ((split_list.size() < min_elements) ||
+        (split_list.size() > max_elements))
+      return false;
+
+                                     // check the different possibilities
+    for (std::vector<std::string>::const_iterator
+           test_string = split_list.begin();
+         test_string != split_list.end(); ++test_string)
+      {
+        // separate key and value from the test_string
+        if (test_string->find(":") == std::string::npos)
+          return false;
+
+        // we know now that there is a ':', so split the string there
+        // and trim spaces
+        std::string key = *test_string;
+        key.erase (key.find(":"), std::string::npos);
+        while ((key.length() > 0) && (std::isspace (key[key.length()-1])))
+          key.erase (key.length()-1, 1);
+
+        std::string value = *test_string;
+        value.erase (0, value.find(":")+1);
+        while ((value.length() > 0) && (std::isspace (value[0])))
+          value.erase (0, 1);
+
+        // then verify that the patterns are satisfied
+        if (key_pattern->match (key) == false)
+          return false;
+        if (value_pattern->match (value) == false)
+          return false;
+      }
+
+    return true;
+  }
+
+
+
+  std::string Map::description () const
+  {
+    std::ostringstream description;
+
+    description << description_init
+                << " map of <"
+                << key_pattern->description() << ":"
+                << value_pattern->description() << ">"
+                << " of length " << min_elements << "..." << max_elements
+                << " (inclusive)"
+                << "]";
+
+    return description.str();
+  }
+
+
+
+  PatternBase *
+  Map::clone () const
+  {
+    return new Map(*key_pattern, *value_pattern, min_elements, max_elements);
+  }
+
+
+  std::size_t
+  Map::memory_consumption () const
+  {
+    return (sizeof(*this) +
+            MemoryConsumption::memory_consumption (*key_pattern) +
+            MemoryConsumption::memory_consumption (*value_pattern));
+  }
+
+
+
+  Map* Map::create (const std::string& description)
+  {
+    if(description.compare(0, std::strlen(description_init), description_init) == 0)
+    {
+      int min_elements, max_elements;
+
+      std::istringstream is(description);
+      is.ignore(strlen(description_init) + strlen(" map of <"));
+
+      std::string str;
+      std::getline(is, str, '>');
+
+      // split 'str' into key and value
+      std::string key = str;
+      key.erase (key.find(":"), std::string::npos);
+
+      std::string value = str;
+      value.erase (0, value.find(":")+1);
+
+      std_cxx1x::shared_ptr<PatternBase> key_pattern (pattern_factory(key));
+      std_cxx1x::shared_ptr<PatternBase> value_pattern (pattern_factory(value));
+
+      is.ignore(strlen(" of length "));
+      if(!(is >> min_elements))
+        return new Map(*key_pattern, *value_pattern);
+
+      is.ignore(strlen("..."));
+      if(!(is >> max_elements))
+        return new Map(*key_pattern, *value_pattern, min_elements);
+
+      return new Map(*key_pattern, *value_pattern, min_elements, max_elements);
+    }
+    else
+      return 0;
+  }
+
+
+
   const char* MultipleSelection::description_init = "[MultipleSelection";
 
 
diff --git a/tests/bits/parameter_handler_13.cc b/tests/bits/parameter_handler_13.cc
new file mode 100644 (file)
index 0000000..9a9b258
--- /dev/null
@@ -0,0 +1,47 @@
+//----------------------------  parameter_handler_13.cc  ---------------------------
+//    $Id$
+//    Version: $Name$
+//
+//    Copyright (C) 2002, 2003, 2004, 2005, 2010, 2012 by the deal.II authors
+//
+//    This file is subject to QPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//----------------------------  parameter_handler_13.cc  ---------------------------
+
+
+// check the Patterns::Map pattern
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/parameter_handler.h>
+#include <fstream>
+
+void check (const char *p)
+{
+  ParameterHandler prm;
+  prm.declare_entry ("test_13", "-1:a, 0:b, 1:c",
+                     Patterns::Map(Patterns::Integer(-1,1),
+                                  Patterns::Selection("a|b|c"),
+                                  2,3));
+
+  std::ifstream in(p);
+  prm.read_input (in);
+
+  deallog << "test_13=" << prm.get ("test_13") << std::endl;
+}
+
+
+int main ()
+{
+  std::ofstream logfile("parameter_handler_13/output");
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+  deallog.threshold_double(1.e-10);
+
+  check ("parameter_handler_13/prm");
+
+  return 0;
+}
diff --git a/tests/bits/parameter_handler_13/cmp/generic b/tests/bits/parameter_handler_13/cmp/generic
new file mode 100644 (file)
index 0000000..3602280
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL::test_13=1 : b, 0: a, 1 : c
diff --git a/tests/bits/parameter_handler_13/prm b/tests/bits/parameter_handler_13/prm
new file mode 100644 (file)
index 0000000..238e1eb
--- /dev/null
@@ -0,0 +1 @@
+set test_13 = 1 : b, 0: a, 1 : c

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.