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
* @{
#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>
// 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
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>>
{
}
};
+ // 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>>
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,
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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;
+}
--- /dev/null
+
+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
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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;
+}
--- /dev/null
+
+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
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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;
+}
--- /dev/null
+
+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