]> https://gitweb.dealii.org/ - dealii.git/commitdiff
extend Patterns::Tools::Convert to understand ComponentMask and FunctionParser
authorDenis Davydov <davydden@gmail.com>
Wed, 19 Sep 2018 21:11:58 +0000 (23:11 +0200)
committerDenis Davydov <davydden@gmail.com>
Wed, 19 Sep 2018 21:11:58 +0000 (23:11 +0200)
include/deal.II/base/function_parser.h
include/deal.II/base/patterns.h
source/base/function_parser.cc
tests/parameter_handler/patterns_14.cc [new file with mode: 0644]
tests/parameter_handler/patterns_14.output [new file with mode: 0644]
tests/parameter_handler/patterns_15.cc [new file with mode: 0644]
tests/parameter_handler/patterns_15.output [new file with mode: 0644]
tests/parameter_handler/patterns_16.cc [new file with mode: 0644]
tests/parameter_handler/patterns_16.output [new file with mode: 0644]

index 31de68488f6beaf03090fe36dcf103c592dfb967..b1d4c1c574a9d6cc47116ae872fcda58b2e0110b 100644 (file)
@@ -341,6 +341,13 @@ public:
   virtual void
   vector_value(const Point<dim> &p, Vector<double> &values) const override;
 
+  /**
+   * Return an array of function expressions (one per component), used to
+   * initialize this function.
+   */
+  const std::vector<std::string> &
+  get_expressions() const;
+
   /**
    * @addtogroup Exceptions
    * @{
index 9414567e8a4bbe1b872e153cba03f3ed0290b727..3b30fcfa3bb3662be9967a1c7b1749aee60d5537 100644 (file)
@@ -28,6 +28,8 @@
 #include <deal.II/base/template_constraints.h>
 #include <deal.II/base/utilities.h>
 
+#include <deal.II/fe/component_mask.h>
+
 #include <boost/archive/basic_archive.hpp>
 #include <boost/core/demangle.hpp>
 #include <boost/property_tree/ptree_fwd.hpp>
@@ -52,6 +54,8 @@ DEAL_II_NAMESPACE_OPEN
 // forward declarations for interfaces and friendship
 class LogStream;
 class MultipleParameterLoop;
+template <int dim>
+class FunctionParser;
 
 /**
  * Namespace for a few classes that act as patterns for the ParameterHandler
@@ -1715,6 +1719,23 @@ namespace Patterns
         static constexpr int map_rank  = RankInfo<Number>::map_rank;
       };
 
+      // Rank of FunctionParser
+      template <int dim>
+      struct RankInfo<std::unique_ptr<FunctionParser<dim>>>
+      {
+        static constexpr int list_rank = 1;
+        static constexpr int map_rank  = 0;
+      };
+
+      // Rank of ComponentMask
+      template <>
+      struct RankInfo<ComponentMask>
+      {
+        static constexpr int list_rank = 1;
+        static constexpr int map_rank  = 0;
+      };
+
+      // Rank of std::pair
       template <class Key, class Value>
       struct RankInfo<std::pair<Key, Value>>
       {
@@ -1982,6 +2003,106 @@ namespace Patterns
       }
     };
 
+    // Functions::FunctionParser
+    template <int dim>
+    struct Convert<std::unique_ptr<FunctionParser<dim>>>
+    {
+      using T = std::unique_ptr<FunctionParser<dim>>;
+
+      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>(
+          Patterns::Anything(),
+          1,
+          Patterns::List::max_int_value,
+          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())
+      {
+        auto p = dynamic_cast<const Patterns::List *>(pattern.get());
+        AssertThrow(p,
+                    ExcMessage("I need a List pattern to convert a string "
+                               "to a List compatbile type."));
+
+        const auto &expressions = t->get_expressions();
+        if (expressions.size() == 0)
+          return std::string();
+
+        std::string s = expressions[0];
+        for (unsigned int i = 1; i < expressions.size(); ++i)
+          s = s + p->get_separator() + expressions[i];
+
+        AssertThrow(pattern->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())
+      {
+        AssertThrow(pattern->match(s), ExcNoMatch(s, pattern.get()));
+
+        auto p = dynamic_cast<const Patterns::List *>(pattern.get());
+        AssertThrow(p,
+                    ExcMessage("I need a List pattern to convert a string "
+                               "to a List compatbile type."));
+
+        const auto expressions =
+          Utilities::split_string_list(s, p->get_separator());
+
+        T t = std_cxx14::make_unique<FunctionParser<dim>>(expressions.size());
+        const std::string var =
+          FunctionParser<dim>::default_variable_names() + ",t";
+        const typename FunctionParser<dim>::ConstMap constants;
+        t->initialize(var, expressions, constants, true);
+        return t;
+      }
+    };
+
+    // ComponentMask
+    template <>
+    struct Convert<ComponentMask>
+    {
+      using T = ComponentMask;
+
+      static std::unique_ptr<Patterns::PatternBase>
+      to_pattern()
+      {
+        return Convert<std::vector<bool>>::to_pattern();
+      }
+
+      static std::string
+      to_string(const T &                                     t,
+                const std::unique_ptr<Patterns::PatternBase> &pattern =
+                  Convert<T>::to_pattern())
+      {
+        std::vector<bool> mask(t.size());
+        for (unsigned int i = 0; i < t.size(); ++i)
+          mask[i] = t[i];
+
+        return Convert<std::vector<bool>>::to_string(mask, pattern);
+      }
+
+      static T
+      to_value(const std::string &                           s,
+               const std::unique_ptr<Patterns::PatternBase> &pattern =
+                 Convert<T>::to_pattern())
+      {
+        const auto mask = Convert<std::vector<bool>>::to_value(s, pattern);
+        return ComponentMask(mask);
+      }
+    };
+
     // Complex numbers
     template <class Number>
     struct Convert<std::complex<Number>>
index 754887ba1d3e35a92a886b0aacd39564bffc8616..24a63b19d5cd25eed729be626d3fe5e05143541d 100644 (file)
@@ -42,6 +42,14 @@ namespace fparser
 DEAL_II_NAMESPACE_OPEN
 
 
+template <int dim>
+const std::vector<std::string> &
+FunctionParser<dim>::get_expressions() const
+{
+  return expressions;
+}
+
+
 
 template <int dim>
 FunctionParser<dim>::FunctionParser(const unsigned int n_components,
diff --git a/tests/parameter_handler/patterns_14.cc b/tests/parameter_handler/patterns_14.cc
new file mode 100644 (file)
index 0000000..20e4e68
--- /dev/null
@@ -0,0 +1,56 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2005 - 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.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+// test add_parameters with ParsedFunction.
+
+#include <deal.II/base/function_parser.h>
+#include <deal.II/base/parameter_handler.h>
+#include <deal.II/base/std_cxx14/memory.h>
+
+#include <memory>
+
+#include "../tests.h"
+
+using namespace Patterns;
+using namespace Patterns::Tools;
+
+int
+main()
+{
+  initlog();
+
+  typedef std::unique_ptr<FunctionParser<3>> T;
+
+  T a;
+  a = Convert<T>::to_value("x*y,y-t");
+
+  ParameterHandler prm;
+  prm.add_parameter("A function", a);
+
+  prm.log_parameters(deallog);
+
+  prm.set("A function", "x*4,y*y*x+t*24");
+
+  deallog << "After ParameterHandler::set =========================="
+          << std::endl
+          << std::endl;
+  prm.log_parameters(deallog);
+
+  deallog << "Actual variables            =========================="
+          << std::endl
+          << std::endl;
+
+  deallog << Convert<T>::to_string(a) << std::endl;
+}
diff --git a/tests/parameter_handler/patterns_14.output b/tests/parameter_handler/patterns_14.output
new file mode 100644 (file)
index 0000000..67a3766
--- /dev/null
@@ -0,0 +1,8 @@
+
+DEAL:parameters::A function: x*y,y-t
+DEAL::After ParameterHandler::set ==========================
+DEAL::
+DEAL:parameters::A function: x*4,y*y*x+t*24
+DEAL::Actual variables            ==========================
+DEAL::
+DEAL::x*4,y*y*x+t*24
diff --git a/tests/parameter_handler/patterns_15.cc b/tests/parameter_handler/patterns_15.cc
new file mode 100644 (file)
index 0000000..ff52373
--- /dev/null
@@ -0,0 +1,57 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2005 - 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.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+// test add_parameters with ComponentMask.
+
+#include <deal.II/base/parameter_handler.h>
+#include <deal.II/base/std_cxx14/memory.h>
+
+#include <deal.II/fe/component_mask.h>
+
+#include <memory>
+
+#include "../tests.h"
+
+using namespace Patterns;
+using namespace Patterns::Tools;
+
+int
+main()
+{
+  initlog();
+
+  typedef ComponentMask T;
+
+  T a;
+  a = Convert<T>::to_value("true,false,true");
+
+  ParameterHandler prm;
+  prm.add_parameter("A mask", a);
+
+  prm.log_parameters(deallog);
+
+  prm.set("A mask", "false,false,true");
+
+  deallog << "After ParameterHandler::set =========================="
+          << std::endl
+          << std::endl;
+  prm.log_parameters(deallog);
+
+  deallog << "Actual variables            =========================="
+          << std::endl
+          << std::endl;
+
+  deallog << Convert<T>::to_string(a) << std::endl;
+}
diff --git a/tests/parameter_handler/patterns_15.output b/tests/parameter_handler/patterns_15.output
new file mode 100644 (file)
index 0000000..91f8e92
--- /dev/null
@@ -0,0 +1,8 @@
+
+DEAL:parameters::A mask: true, false, true
+DEAL::After ParameterHandler::set ==========================
+DEAL::
+DEAL:parameters::A mask: false,false,true
+DEAL::Actual variables            ==========================
+DEAL::
+DEAL::false, false, true
diff --git a/tests/parameter_handler/patterns_16.cc b/tests/parameter_handler/patterns_16.cc
new file mode 100644 (file)
index 0000000..f89e409
--- /dev/null
@@ -0,0 +1,58 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2005 - 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.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+// test add_parameters with std::map<unsigned int,ParsedFunction>
+
+#include <deal.II/base/function_parser.h>
+#include <deal.II/base/parameter_handler.h>
+#include <deal.II/base/std_cxx14/memory.h>
+
+#include <deal.II/fe/component_mask.h>
+
+#include <memory>
+
+#include "../tests.h"
+
+using namespace Patterns;
+using namespace Patterns::Tools;
+
+int
+main()
+{
+  initlog();
+
+  typedef std::map<types::boundary_id, std::unique_ptr<FunctionParser<3>>> T;
+
+  T a;
+  a = Convert<T>::to_value("0:x,y,z*t");
+
+  ParameterHandler prm;
+  prm.add_parameter("Boundary conditions", a);
+
+  prm.log_parameters(deallog);
+
+  prm.set("Boundary conditions", "0:x*x*x,y*y*y,z*z*z*t;1:0,x*y,t");
+
+  deallog << "After ParameterHandler::set =========================="
+          << std::endl
+          << std::endl;
+  prm.log_parameters(deallog);
+
+  deallog << "Actual variables            =========================="
+          << std::endl
+          << std::endl;
+
+  deallog << Convert<T>::to_string(a) << std::endl;
+}
diff --git a/tests/parameter_handler/patterns_16.output b/tests/parameter_handler/patterns_16.output
new file mode 100644 (file)
index 0000000..9a822a8
--- /dev/null
@@ -0,0 +1,8 @@
+
+DEAL:parameters::Boundary conditions: 0:x,y,z*t
+DEAL::After ParameterHandler::set ==========================
+DEAL::
+DEAL:parameters::Boundary conditions: 0:x*x*x,y*y*y,z*z*z*t;1:0,x*y,t
+DEAL::Actual variables            ==========================
+DEAL::
+DEAL::0:x*x*x,y*y*y,z*z*z*t; 1:0,x*y,t

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.