Reference documentation for deal.II version 8.5.1
function_parser.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2005 - 2016 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
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>
21 #include <cmath>
22 #include <map>
23 
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
28 
29 #ifdef DEAL_II_WITH_MUPARSER
30 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
31 #include <muParser.h>
32 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
33 #else
34 
35 
36 
37 namespace fparser
38 {
39  class FunctionParser
40  {};
41 }
42 #endif
43 
44 DEAL_II_NAMESPACE_OPEN
45 
46 
47 
48 template <int dim>
50  const double initial_time,
51  const double h)
52  :
53  AutoDerivativeFunction<dim>(h, n_components, initial_time),
54  initialized (false),
55  n_vars (0)
56 {}
57 
58 
59 
60 template <int dim>
62 {}
63 
64 #ifdef DEAL_II_WITH_MUPARSER
65 
66 template <int dim>
67 void FunctionParser<dim>::initialize (const std::string &variables,
68  const std::vector<std::string> &expressions,
69  const std::map<std::string, double> &constants,
70  const bool time_dependent)
71 {
72  this->fp.clear(); // this will reset all thread-local objects
73 
74  this->constants = constants;
75  this->var_names = Utilities::split_string_list(variables, ',');
76  this->expressions = expressions;
77  AssertThrow(((time_dependent)?dim+1:dim) == var_names.size(),
78  ExcMessage("Wrong number of variables"));
79 
80  // We check that the number of
81  // components of this function
82  // matches the number of components
83  // passed in as a vector of
84  // strings.
85  AssertThrow(this->n_components == expressions.size(),
86  ExcInvalidExpressionSize(this->n_components,
87  expressions.size()) );
88 
89  // Now we define how many variables
90  // we expect to read in. We
91  // distinguish between two cases:
92  // Time dependent problems, and not
93  // time dependent problems. In the
94  // first case the number of
95  // variables is given by the
96  // dimension plus one. In the other
97  // case, the number of variables is
98  // equal to the dimension. Once we
99  // parsed the variables string, if
100  // none of this is the case, then
101  // an exception is thrown.
102  if (time_dependent)
103  n_vars = dim+1;
104  else
105  n_vars = dim;
106 
107  // create a parser object for the current thread we can then query
108  // in value() and vector_value(). this is not strictly necessary
109  // because a user may never call these functions on the current
110  // thread, but it gets us error messages about wrong formulas right
111  // away
112  init_muparser ();
113 
114  // finally set the initialization bit
115  initialized = true;
116 }
117 
118 
119 
120 namespace internal
121 {
122  // convert double into int
123  int mu_round(double val)
124  {
125  return static_cast<int>(val + ((val>=0.0) ? 0.5 : -0.5) );
126  }
127 
128  double mu_if(double condition, double thenvalue, double elsevalue)
129  {
130  if (mu_round(condition))
131  return thenvalue;
132  else
133  return elsevalue;
134  }
135 
136  double mu_or(double left, double right)
137  {
138  return (mu_round(left)) || (mu_round(right));
139  }
140 
141  double mu_and(double left, double right)
142  {
143  return (mu_round(left)) && (mu_round(right));
144  }
145 
146  double mu_int(double value)
147  {
148  return static_cast<double>(mu_round(value));
149  }
150 
151  double mu_ceil(double value)
152  {
153  return ceil(value);
154  }
155 
156  double mu_floor(double value)
157  {
158  return floor(value);
159  }
160 
161  double mu_cot(double value)
162  {
163  return 1.0/tan(value);
164  }
165 
166  double mu_csc(double value)
167  {
168  return 1.0/sin(value);
169  }
170 
171  double mu_sec(double value)
172  {
173  return 1.0/cos(value);
174  }
175 
176  double mu_log(double value)
177  {
178  return log(value);
179  }
180 
181  double mu_pow(double a, double b)
182  {
183  return std::pow(a, b);
184  }
185 
186  double mu_erfc(double value)
187  {
188  return boost::math::erfc(value);
189  }
190 
191  // returns a random value in the range [0,1] initializing the generator
192  // with the given seed
193  double mu_rand_seed(double seed)
194  {
195  static Threads::Mutex rand_mutex;
196  Threads::Mutex::ScopedLock lock(rand_mutex);
197 
198  static boost::random::uniform_real_distribution<> uniform_distribution(0,1);
199 
200  // for each seed an unique random number generator is created,
201  // which is initialized with the seed itself
202  static std::map<double, boost::random::mt19937> rng_map;
203 
204  if (rng_map.find(seed) == rng_map.end())
205  rng_map[seed] = boost::random::mt19937(static_cast<unsigned int>(seed));
206 
207  return uniform_distribution(rng_map[seed]);
208  }
209 
210  // returns a random value in the range [0,1]
211  double mu_rand()
212  {
213  static Threads::Mutex rand_mutex;
214  Threads::Mutex::ScopedLock lock(rand_mutex);
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);
218  }
219 
220 }
221 
222 
223 
224 template <int dim>
226 {
227  // check that we have not already initialized the parser on the
228  // current thread, i.e., that the current function is only called
229  // once per thread
230  Assert (fp.get().size()==0, ExcInternalError());
231 
232  // initialize the objects for the current thread (fp.get() and
233  // vars.get())
234  fp.get().resize(this->n_components);
235  vars.get().resize(var_names.size());
236  for (unsigned int component=0; component<this->n_components; ++component)
237  {
238  for (std::map< std::string, double >::const_iterator constant = constants.begin();
239  constant != constants.end(); ++constant)
240  {
241  fp.get()[component].DefineConst(constant->first.c_str(), constant->second);
242  }
243 
244  for (unsigned int iv=0; iv<var_names.size(); ++iv)
245  fp.get()[component].DefineVar(var_names[iv].c_str(), &vars.get()[iv]);
246 
247  // define some compatibility functions:
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);
262 
263  try
264  {
265  // muparser expects that functions have no
266  // space between the name of the function and the opening
267  // parenthesis. this is awkward because it is not backward
268  // compatible to the library we used to use before muparser
269  // (the fparser library) but also makes no real sense.
270  // consequently, in the expressions we set, remove any space
271  // we may find after function names
272  std::string transformed_expression = expressions[component];
273 
274  const char *function_names[] =
275  {
276  // functions predefined by muparser
277  "sin",
278  "cos",
279  "tan",
280  "asin",
281  "acos",
282  "atan",
283  "sinh",
284  "cosh",
285  "tanh",
286  "asinh",
287  "acosh",
288  "atanh",
289  "atan2",
290  "log2",
291  "log10",
292  "log",
293  "ln",
294  "exp",
295  "sqrt",
296  "sign",
297  "rint",
298  "abs",
299  "min",
300  "max",
301  "sum",
302  "avg",
303  // functions we define ourselves above
304  "if",
305  "int",
306  "ceil",
307  "cot",
308  "csc",
309  "floor",
310  "sec",
311  "pow",
312  "erfc",
313  "rand",
314  "rand_seed"
315  };
316  for (unsigned int f=0; f<sizeof(function_names)/sizeof(function_names[0]); ++f)
317  {
318  const std::string function_name = function_names[f];
319  const unsigned int function_name_length = function_name.size();
320 
321  std::string::size_type pos = 0;
322  while (true)
323  {
324  // try to find any occurrences of the function name
325  pos = transformed_expression.find (function_name, pos);
326  if (pos == std::string::npos)
327  break;
328 
329  // replace whitespace until there no longer is any
330  while ((pos+function_name_length<transformed_expression.size())
331  &&
332  ((transformed_expression[pos+function_name_length] == ' ')
333  ||
334  (transformed_expression[pos+function_name_length] == '\t')))
335  transformed_expression.erase (transformed_expression.begin()+pos+function_name_length);
336 
337  // move the current search position by the size of the
338  // actual function name
339  pos += function_name_length;
340  }
341  }
342 
343  // now use the transformed expression
344  fp.get()[component].SetExpr(transformed_expression);
345  }
346  catch (mu::ParserError &e)
347  {
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()));
354  }
355  }
356 }
357 
358 
359 
360 template <int dim>
361 void FunctionParser<dim>::initialize (const std::string &vars,
362  const std::string &expression,
363  const std::map<std::string, double> &constants,
364  const bool time_dependent)
365 {
366  initialize(vars, Utilities::split_string_list(expression, ';'),
367  constants, time_dependent);
368 }
369 
370 
371 
372 template <int dim>
374  const unsigned int component) const
375 {
376  Assert (initialized==true, ExcNotInitialized());
377  Assert (component < this->n_components,
378  ExcIndexRange(component, 0, this->n_components));
379 
380  // initialize the parser if that hasn't happened yet on the current thread
381  if (fp.get().size() == 0)
382  init_muparser();
383 
384  for (unsigned int i=0; i<dim; ++i)
385  vars.get()[i] = p(i);
386  if (dim != n_vars)
387  vars.get()[dim] = this->get_time();
388 
389  try
390  {
391  return fp.get()[component].Eval();
392  }
393  catch (mu::ParserError &e)
394  {
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()));
401  return 0.0;
402  }
403 }
404 
405 
406 
407 template <int dim>
409  Vector<double> &values) const
410 {
411  Assert (initialized==true, ExcNotInitialized());
412  Assert (values.size() == this->n_components,
413  ExcDimensionMismatch (values.size(), this->n_components));
414 
415 
416  // initialize the parser if that hasn't happened yet on the current thread
417  if (fp.get().size() == 0)
418  init_muparser();
419 
420  for (unsigned int i=0; i<dim; ++i)
421  vars.get()[i] = p(i);
422  if (dim != n_vars)
423  vars.get()[dim] = this->get_time();
424 
425  for (unsigned int component = 0; component < this->n_components;
426  ++component)
427  values(component) = fp.get()[component].Eval();
428 }
429 
430 #else
431 
432 
433 template <int dim>
434 void
435 FunctionParser<dim>::initialize(const std::string &,
436  const std::vector<std::string> &,
437  const std::map<std::string, double> &,
438  const bool)
439 {
440  Assert(false, ExcNeedsFunctionparser());
441 }
442 
443 template <int dim>
444 void
445 FunctionParser<dim>::initialize(const std::string &,
446  const std::string &,
447  const std::map<std::string, double> &,
448  const bool)
449 {
450  Assert(false, ExcNeedsFunctionparser());
451 }
452 
453 
454 
455 template <int dim>
457  const Point<dim> &, unsigned int) const
458 {
459  Assert(false, ExcNeedsFunctionparser());
460  return 0.;
461 }
462 
463 
464 template <int dim>
466  const Point<dim> &, Vector<double> &) const
467 {
468  Assert(false, ExcNeedsFunctionparser());
469 }
470 
471 
472 #endif
473 
474 // Explicit Instantiations.
475 
476 template class FunctionParser<1>;
477 template class FunctionParser<2>;
478 template class FunctionParser<3>;
479 
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)
Definition: exceptions.h:369
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:313
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=',')
Definition: utilities.cc:268
std::size_t size() const
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()