17 #include <deal.II/base/function_parser.h> 18 #include <deal.II/base/utilities.h> 19 #include <deal.II/base/thread_management.h> 20 #include <deal.II/lac/vector.h> 24 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
25 #include <boost/random.hpp> 26 #include <boost/math/special_functions/erf.hpp> 27 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
29 #ifdef DEAL_II_WITH_MUPARSER 30 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
32 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
44 DEAL_II_NAMESPACE_OPEN
50 const double initial_time,
64 #ifdef DEAL_II_WITH_MUPARSER 68 const std::vector<std::string> &expressions,
69 const std::map<std::string, double> &constants,
70 const bool time_dependent)
74 this->constants = constants;
76 this->expressions = expressions;
77 AssertThrow(((time_dependent)?dim+1:dim) == var_names.size(),
87 expressions.size()) );
123 int mu_round(
double val)
125 return static_cast<int>(val + ((val>=0.0) ? 0.5 : -0.5) );
128 double mu_if(
double condition,
double thenvalue,
double elsevalue)
130 if (mu_round(condition))
136 double mu_or(
double left,
double right)
138 return (mu_round(left)) || (mu_round(right));
141 double mu_and(
double left,
double right)
143 return (mu_round(left)) && (mu_round(right));
146 double mu_int(
double value)
148 return static_cast<double>(mu_round(value));
151 double mu_ceil(
double value)
156 double mu_floor(
double value)
161 double mu_cot(
double value)
163 return 1.0/tan(value);
166 double mu_csc(
double value)
168 return 1.0/sin(value);
171 double mu_sec(
double value)
173 return 1.0/cos(value);
176 double mu_log(
double value)
181 double mu_pow(
double a,
double b)
183 return std::pow(a, b);
186 double mu_erfc(
double value)
188 return boost::math::erfc(value);
193 double mu_rand_seed(
double seed)
198 static boost::random::uniform_real_distribution<> uniform_distribution(0,1);
202 static std::map<double, boost::random::mt19937> rng_map;
204 if (rng_map.find(seed) == rng_map.end())
205 rng_map[seed] = boost::random::mt19937(static_cast<unsigned int>(seed));
207 return uniform_distribution(rng_map[seed]);
215 static boost::random::uniform_real_distribution<> uniform_distribution(0,1);
216 static boost::random::mt19937 rng(static_cast<unsigned long>(std::time(0)));
217 return uniform_distribution(rng);
235 vars.get().resize(var_names.size());
236 for (
unsigned int component=0; component<this->
n_components; ++component)
238 for (std::map< std::string, double >::const_iterator constant = constants.begin();
239 constant != constants.end(); ++constant)
241 fp.get()[component].DefineConst(constant->first.c_str(), constant->second);
244 for (
unsigned int iv=0; iv<var_names.size(); ++iv)
245 fp.get()[component].DefineVar(var_names[iv].c_str(), &vars.get()[iv]);
248 fp.get()[component].DefineFun(
"if",internal::mu_if,
true);
249 fp.get()[component].DefineOprt(
"|", internal::mu_or, 1);
250 fp.get()[component].DefineOprt(
"&", internal::mu_and, 2);
251 fp.get()[component].DefineFun(
"int", internal::mu_int,
true);
252 fp.get()[component].DefineFun(
"ceil", internal::mu_ceil,
true);
253 fp.get()[component].DefineFun(
"cot", internal::mu_cot,
true);
254 fp.get()[component].DefineFun(
"csc", internal::mu_csc,
true);
255 fp.get()[component].DefineFun(
"floor", internal::mu_floor,
true);
256 fp.get()[component].DefineFun(
"sec", internal::mu_sec,
true);
257 fp.get()[component].DefineFun(
"log", internal::mu_log,
true);
258 fp.get()[component].DefineFun(
"pow", internal::mu_pow,
true);
259 fp.get()[component].DefineFun(
"erfc", internal::mu_erfc,
true);
260 fp.get()[component].DefineFun(
"rand_seed", internal::mu_rand_seed,
true);
261 fp.get()[component].DefineFun(
"rand", internal::mu_rand,
true);
272 std::string transformed_expression = expressions[component];
274 const char *function_names[] =
316 for (
unsigned int f=0; f<
sizeof(function_names)/
sizeof(function_names[0]); ++f)
318 const std::string function_name = function_names[f];
319 const unsigned int function_name_length = function_name.size();
321 std::string::size_type pos = 0;
325 pos = transformed_expression.find (function_name, pos);
326 if (pos == std::string::npos)
330 while ((pos+function_name_length<transformed_expression.size())
332 ((transformed_expression[pos+function_name_length] ==
' ')
334 (transformed_expression[pos+function_name_length] ==
'\t')))
335 transformed_expression.erase (transformed_expression.begin()+pos+function_name_length);
339 pos += function_name_length;
344 fp.get()[component].SetExpr(transformed_expression);
346 catch (mu::ParserError &e)
348 std::cerr <<
"Message: <" << e.GetMsg() <<
">\n";
349 std::cerr <<
"Formula: <" << e.GetExpr() <<
">\n";
350 std::cerr <<
"Token: <" << e.GetToken() <<
">\n";
351 std::cerr <<
"Position: <" << e.GetPos() <<
">\n";
352 std::cerr <<
"Errc: <" << e.GetCode() <<
">" << std::endl;
353 AssertThrow(
false, ExcParseError(e.GetCode(), e.GetMsg().c_str()));
362 const std::string &expression,
363 const std::map<std::string, double> &constants,
364 const bool time_dependent)
367 constants, time_dependent);
374 const unsigned int component)
const 381 if (fp.get().size() == 0)
384 for (
unsigned int i=0; i<dim; ++i)
385 vars.get()[i] = p(i);
387 vars.get()[dim] = this->get_time();
391 return fp.get()[component].Eval();
393 catch (mu::ParserError &e)
395 std::cerr <<
"Message: <" << e.GetMsg() <<
">\n";
396 std::cerr <<
"Formula: <" << e.GetExpr() <<
">\n";
397 std::cerr <<
"Token: <" << e.GetToken() <<
">\n";
398 std::cerr <<
"Position: <" << e.GetPos() <<
">\n";
399 std::cerr <<
"Errc: <" << e.GetCode() <<
">" << std::endl;
400 AssertThrow(
false, ExcParseError(e.GetCode(), e.GetMsg().c_str()));
417 if (fp.get().size() == 0)
420 for (
unsigned int i=0; i<dim; ++i)
421 vars.get()[i] = p(i);
423 vars.get()[dim] = this->get_time();
425 for (
unsigned int component = 0; component < this->
n_components;
427 values(component) = fp.get()[component].Eval();
436 const std::vector<std::string> &,
437 const std::map<std::string, double> &,
447 const std::map<std::string, double> &,
480 DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNeedsFunctionparser()
FunctionParser(const unsigned int n_components=1, const double initial_time=0.0, const double h=1e-8)
virtual double value(const Point< dim > &p, const unsigned int component=0) const
static ::ExceptionBase & ExcNotInitialized()
#define AssertThrow(cond, exc)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::vector< std::string > split_string_list(const std::string &s, const char delimiter=',')
void init_muparser() const
virtual void vector_value(const Point< dim > &p, Vector< double > &values) const
void initialize(const std::string &vars, const std::vector< std::string > &expressions, const ConstMap &constants, const bool time_dependent=false)
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
static ::ExceptionBase & ExcInternalError()