parsed_fem_function.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_PARSED_FEM_FUNCTION_H
21 #define LIBMESH_PARSED_FEM_FUNCTION_H
22 
23 #include "libmesh/libmesh_config.h"
24 
25 // Local Includes
27 #include "libmesh/point.h"
28 #include "libmesh/system.h"
29 
30 #ifdef LIBMESH_HAVE_FPARSER
31 // FParser includes
32 #include "libmesh/fparser.hh"
33 #endif
34 
35 // C++ includes
36 #include <cctype>
37 #include <iomanip>
38 #include <sstream>
39 #include <string>
40 #include <vector>
41 
42 namespace libMesh
43 {
44 
54 template <typename Output=Number>
55 class ParsedFEMFunction : public FEMFunctionBase<Output>
56 {
57 public:
58 
62  explicit
63  ParsedFEMFunction (const System & sys,
64  const std::string & expression,
65  const std::vector<std::string> * additional_vars=nullptr,
66  const std::vector<Output> * initial_vals=nullptr);
67 
73 
83  ParsedFEMFunction (const ParsedFEMFunction &) = default;
84  ParsedFEMFunction (ParsedFEMFunction &&) = default;
85  virtual ~ParsedFEMFunction () = default;
86 
90  void reparse (const std::string & expression);
91 
92  virtual void init_context (const FEMContext & c) override;
93 
94  virtual std::unique_ptr<FEMFunctionBase<Output>> clone () const override;
95 
96  virtual Output operator() (const FEMContext & c,
97  const Point & p,
98  const Real time = 0.) override;
99 
100  void operator() (const FEMContext & c,
101  const Point & p,
102  const Real time,
103  DenseVector<Output> & output) override;
104 
105  virtual Output component(const FEMContext & c,
106  unsigned int i,
107  const Point & p,
108  Real time=0.) override;
109 
110  const std::string & expression() { return _expression; }
111 
120  Output get_inline_value(const std::string & inline_var_name) const;
121 
132  void set_inline_value(const std::string & inline_var_name,
133  Output newval);
134 
135 protected:
136  // Helper function for reparsing minor changes to expression
137  void partial_reparse (const std::string & expression);
138 
139  // Helper function for parsing out variable names
140  std::size_t find_name (const std::string & varname,
141  const std::string & expr) const;
142 
143  // Helper function for evaluating function arguments
144  void eval_args(const FEMContext & c,
145  const Point & p,
146  const Real time);
147 
148  // Evaluate the ith FunctionParser and check the result
149 #ifdef LIBMESH_HAVE_FPARSER
150  inline Output eval(FunctionParserBase<Output> & parser,
151  const std::string & libmesh_dbg_var(function_name),
152  unsigned int libmesh_dbg_var(component_idx)) const;
153 #else // LIBMESH_HAVE_FPARSER
154  inline Output eval(char & libmesh_dbg_var(parser),
155  const std::string & libmesh_dbg_var(function_name),
156  unsigned int libmesh_dbg_var(component_idx)) const;
157 #endif
158 
159 private:
160  const System & _sys;
161  std::string _expression;
162  std::vector<std::string> _subexpressions;
163  unsigned int _n_vars,
168 #ifdef LIBMESH_HAVE_FPARSER
169  std::vector<FunctionParserBase<Output>> parsers;
170 #else
171  std::vector<char> parsers;
172 #endif
173  std::vector<Output> _spacetime;
174 
175  // Flags for which variables need to be computed
176 
177  // _need_var[v] is true iff value(v) is needed
178  std::vector<bool> _need_var;
179 
180  // _need_var_grad[v*LIBMESH_DIM+i] is true iff grad(v,i) is needed
181  std::vector<bool> _need_var_grad;
182 
183 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
184  // _need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+i*LIBMESH_DIM+j] is true
185  // iff grad(v,i,j) is needed
186  std::vector<bool> _need_var_hess;
187 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
188 
189  // Variables/values that can be parsed and handled by the function parser
190  std::string variables;
191  std::vector<std::string> _additional_vars;
192  std::vector<Output> _initial_vals;
193 };
194 
195 
196 /*----------------------- Inline functions ----------------------------------*/
197 
198 template <typename Output>
199 inline
201  const std::string & expression,
202  const std::vector<std::string> * additional_vars,
203  const std::vector<Output> * initial_vals) :
204  _sys(sys),
205  _expression (), // overridden by parse()
206  _n_vars(sys.n_vars()),
207  _n_requested_vars(0),
208  _n_requested_grad_components(0),
209  _n_requested_hess_components(0),
210  _requested_normals(false),
211  _need_var(_n_vars, false),
212  _need_var_grad(_n_vars*LIBMESH_DIM, false),
213 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
214  _need_var_hess(_n_vars*LIBMESH_DIM*LIBMESH_DIM, false),
215 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
216  _additional_vars (additional_vars ? *additional_vars : std::vector<std::string>()),
217  _initial_vals (initial_vals ? *initial_vals : std::vector<Output>())
218 {
219  this->reparse(expression);
220 }
221 
222 
223 template <typename Output>
224 inline
225 void
226 ParsedFEMFunction<Output>::reparse (const std::string & expression)
227 {
228  variables = "x";
229 #if LIBMESH_DIM > 1
230  variables += ",y";
231 #endif
232 #if LIBMESH_DIM > 2
233  variables += ",z";
234 #endif
235  variables += ",t";
236 
237  for (unsigned int v=0; v != _n_vars; ++v)
238  {
239  const std::string & varname = _sys.variable_name(v);
240  std::size_t varname_i = find_name(varname, expression);
241 
242  // If we didn't find our variable name then let's go to the
243  // next.
244  if (varname_i == std::string::npos)
245  continue;
246 
247  _need_var[v] = true;
248  variables += ',';
249  variables += varname;
250  _n_requested_vars++;
251  }
252 
253  for (unsigned int v=0; v != _n_vars; ++v)
254  {
255  const std::string & varname = _sys.variable_name(v);
256 
257  for (unsigned int d=0; d != LIBMESH_DIM; ++d)
258  {
259  std::string gradname = std::string("grad_") +
260  "xyz"[d] + '_' + varname;
261  std::size_t gradname_i = find_name(gradname, expression);
262 
263  // If we didn't find that gradient component of our
264  // variable name then let's go to the next.
265  if (gradname_i == std::string::npos)
266  continue;
267 
268  _need_var_grad[v*LIBMESH_DIM+d] = true;
269  variables += ',';
270  variables += gradname;
271  _n_requested_grad_components++;
272  }
273  }
274 
275 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
276  for (unsigned int v=0; v != _n_vars; ++v)
277  {
278  const std::string & varname = _sys.variable_name(v);
279 
280  for (unsigned int d1=0; d1 != LIBMESH_DIM; ++d1)
281  for (unsigned int d2=0; d2 != LIBMESH_DIM; ++d2)
282  {
283  std::string hessname = std::string("hess_") +
284  "xyz"[d1] + "xyz"[d2] + '_' + varname;
285  std::size_t hessname_i = find_name(hessname, expression);
286 
287  // If we didn't find that hessian component of our
288  // variable name then let's go to the next.
289  if (hessname_i == std::string::npos)
290  continue;
291 
292  _need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+d1*LIBMESH_DIM+d2]
293  = true;
294  variables += ',';
295  variables += hessname;
296  _n_requested_hess_components++;
297  }
298  }
299 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
300 
301  {
302  std::size_t nx_i = find_name("n_x", expression);
303  std::size_t ny_i = find_name("n_y", expression);
304  std::size_t nz_i = find_name("n_z", expression);
305 
306  // If we found any requests for normal components then we'll
307  // compute normals
308  if (nx_i != std::string::npos ||
309  ny_i != std::string::npos ||
310  nz_i != std::string::npos)
311  {
312  _requested_normals = true;
313  variables += ',';
314  variables += "n_x";
315  if (LIBMESH_DIM > 1)
316  {
317  variables += ',';
318  variables += "n_y";
319  }
320  if (LIBMESH_DIM > 2)
321  {
322  variables += ',';
323  variables += "n_z";
324  }
325  }
326  }
327 
328  _spacetime.resize
329  (LIBMESH_DIM + 1 + _n_requested_vars +
330  _n_requested_grad_components + _n_requested_hess_components +
331  (_requested_normals ? LIBMESH_DIM : 0) +
332  _additional_vars.size());
333 
334  // If additional vars were passed, append them to the string
335  // that we send to the function parser. Also add them to the
336  // end of our spacetime vector
337  unsigned int offset = LIBMESH_DIM + 1 + _n_requested_vars +
338  _n_requested_grad_components + _n_requested_hess_components;
339 
340  for (std::size_t i=0; i < _additional_vars.size(); ++i)
341  {
342  variables += "," + _additional_vars[i];
343  // Initialize extra variables to the vector passed in or zero
344  // Note: The initial_vals vector can be shorter than the additional_vars vector
345  _spacetime[offset + i] =
346  (i < _initial_vals.size()) ? _initial_vals[i] : 0;
347  }
348 
349  this->partial_reparse(expression);
350 }
351 
352 template <typename Output>
353 inline
354 void
356 {
357  for (unsigned int v=0; v != _n_vars; ++v)
358  {
359  FEBase * elem_fe;
360  c.get_element_fe(v, elem_fe);
361  if (_n_requested_vars)
362  elem_fe->get_phi();
363  if (_n_requested_grad_components)
364  elem_fe->get_dphi();
365 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
366  if (_n_requested_hess_components)
367  elem_fe->get_d2phi();
368 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
369  }
370 
371  if (_requested_normals)
372  {
373  FEBase * side_fe;
374  c.get_side_fe(0, side_fe);
375 
376  side_fe->get_normals();
377 
378  // FIXME: this is a hack to support normals at quadrature
379  // points; we don't support normals elsewhere.
380  side_fe->get_xyz();
381  }
382 }
383 
384 template <typename Output>
385 inline
386 std::unique_ptr<FEMFunctionBase<Output>>
388 {
389  return std::unique_ptr<FEMFunctionBase<Output>>
390  (new ParsedFEMFunction(_sys, _expression, &_additional_vars, &_initial_vals));
391 }
392 
393 template <typename Output>
394 inline
395 Output
397  const Point & p,
398  const Real time)
399 {
400  eval_args(c, p, time);
401 
402  return eval(parsers[0], "f", 0);
403 }
404 
405 
406 
407 template <typename Output>
408 inline
409 void
411  const Point & p,
412  const Real time,
413  DenseVector<Output> & output)
414 {
415  eval_args(c, p, time);
416 
417  unsigned int size = output.size();
418 
419  libmesh_assert_equal_to (size, parsers.size());
420 
421  for (unsigned int i=0; i != size; ++i)
422  output(i) = eval(parsers[i], "f", i);
423 }
424 
425 
426 template <typename Output>
427 inline
428 Output
430  unsigned int i,
431  const Point & p,
432  Real time)
433 {
434  eval_args(c, p, time);
435 
436  libmesh_assert_less (i, parsers.size());
437  return eval(parsers[i], "f", i);
438 }
439 
440 template <typename Output>
441 inline
442 Output
443 ParsedFEMFunction<Output>::get_inline_value(const std::string & inline_var_name) const
444 {
445  libmesh_assert_greater (_subexpressions.size(), 0);
446 
447 #ifndef NDEBUG
448  bool found_var_name = false;
449 #endif
450  Output old_var_value(0.);
451 
452  for (std::size_t s=0; s != _subexpressions.size(); ++s)
453  {
454  const std::string & subexpression = _subexpressions[s];
455  const std::size_t varname_i =
456  find_name(inline_var_name, subexpression);
457  if (varname_i == std::string::npos)
458  continue;
459 
460  const std::size_t assignment_i =
461  subexpression.find(":", varname_i+1);
462 
463  libmesh_assert_not_equal_to(assignment_i, std::string::npos);
464 
465  libmesh_assert_equal_to(subexpression[assignment_i+1], '=');
466  for (std::size_t i = varname_i+1; i != assignment_i; ++i)
467  libmesh_assert_equal_to(subexpression[i], ' ');
468 
469  std::size_t end_assignment_i =
470  subexpression.find(";", assignment_i+1);
471 
472  libmesh_assert_not_equal_to(end_assignment_i, std::string::npos);
473 
474  std::string new_subexpression =
475  subexpression.substr(0, end_assignment_i+1) +
476  inline_var_name;
477 
478 #ifdef LIBMESH_HAVE_FPARSER
479  // Parse and evaluate the new subexpression.
480  // Add the same constants as we used originally.
481  FunctionParserBase<Output> fp;
482  fp.AddConstant("NaN", std::numeric_limits<Real>::quiet_NaN());
483  fp.AddConstant("pi", std::acos(Real(-1)));
484  fp.AddConstant("e", std::exp(Real(1)));
485  if (fp.Parse(new_subexpression, variables) != -1) // -1 for success
486  libmesh_error_msg
487  ("ERROR: FunctionParser is unable to parse modified expression: "
488  << new_subexpression << '\n' << fp.ErrorMsg());
489 
490  Output new_var_value = this->eval(fp, new_subexpression, 0);
491 #ifdef NDEBUG
492  return new_var_value;
493 #else
494  if (found_var_name)
495  {
496  libmesh_assert_equal_to(old_var_value, new_var_value);
497  }
498  else
499  {
500  old_var_value = new_var_value;
501  found_var_name = true;
502  }
503 #endif
504 
505 #else
506  libmesh_error_msg("ERROR: This functionality requires fparser!");
507 #endif
508  }
509 
510  libmesh_assert(found_var_name);
511  return old_var_value;
512 }
513 
514 template <typename Output>
515 inline
516 void
517 ParsedFEMFunction<Output>::set_inline_value (const std::string & inline_var_name,
518  Output newval)
519 {
520  libmesh_assert_greater (_subexpressions.size(), 0);
521 
522 #ifndef NDEBUG
523  bool found_var_name = false;
524 #endif
525  for (std::size_t s=0; s != _subexpressions.size(); ++s)
526  {
527  const std::string & subexpression = _subexpressions[s];
528  const std::size_t varname_i =
529  find_name(inline_var_name, subexpression);
530  if (varname_i == std::string::npos)
531  continue;
532 
533 #ifndef NDEBUG
534  found_var_name = true;
535 #endif
536 
537  const std::size_t assignment_i =
538  subexpression.find(":", varname_i+1);
539 
540  libmesh_assert_not_equal_to(assignment_i, std::string::npos);
541 
542  libmesh_assert_equal_to(subexpression[assignment_i+1], '=');
543  for (std::size_t i = varname_i+1; i != assignment_i; ++i)
544  libmesh_assert_equal_to(subexpression[i], ' ');
545 
546  std::size_t end_assignment_i =
547  subexpression.find(";", assignment_i+1);
548 
549  libmesh_assert_not_equal_to(end_assignment_i, std::string::npos);
550 
551  std::ostringstream new_subexpression;
552  new_subexpression << subexpression.substr(0, assignment_i+2)
553  << std::setprecision(std::numeric_limits<Output>::digits10+2)
554 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
555  << '(' << newval.real() << '+'
556  << newval.imag() << 'i' << ')'
557 #else
558  << newval
559 #endif
560  << subexpression.substr(end_assignment_i,
561  std::string::npos);
562  _subexpressions[s] = new_subexpression.str();
563  }
564 
565  libmesh_assert(found_var_name);
566 
567  std::string new_expression;
568 
569  for (std::size_t s=0; s != _subexpressions.size(); ++s)
570  {
571  new_expression += '{';
572  new_expression += _subexpressions[s];
573  new_expression += '}';
574  }
575 
576  this->partial_reparse(new_expression);
577 }
578 
579 
580 template <typename Output>
581 inline
582 void
583 ParsedFEMFunction<Output>::partial_reparse (const std::string & expression)
584 {
585  _expression = expression;
586  _subexpressions.clear();
587  parsers.clear();
588 
589  size_t nextstart = 0, end = 0;
590 
591  while (end != std::string::npos)
592  {
593  // If we're past the end of the string, we can't make any more
594  // subparsers
595  if (nextstart >= expression.size())
596  break;
597 
598  // If we're at the start of a brace delimited section, then we
599  // parse just that section:
600  if (expression[nextstart] == '{')
601  {
602  nextstart++;
603  end = expression.find('}', nextstart);
604  }
605  // otherwise we parse the whole thing
606  else
607  end = std::string::npos;
608 
609  // We either want the whole end of the string (end == npos) or
610  // a substring in the middle.
611  _subexpressions.push_back
612  (expression.substr(nextstart, (end == std::string::npos) ?
613  std::string::npos : end - nextstart));
614 
615  // fparser can crash on empty expressions
616  if (_subexpressions.back().empty())
617  libmesh_error_msg("ERROR: FunctionParser is unable to parse empty expression.\n");
618 
619 
620 #ifdef LIBMESH_HAVE_FPARSER
621  // Parse (and optimize if possible) the subexpression.
622  // Add some basic constants, to Real precision.
623  FunctionParserBase<Output> fp;
624  fp.AddConstant("NaN", std::numeric_limits<Real>::quiet_NaN());
625  fp.AddConstant("pi", std::acos(Real(-1)));
626  fp.AddConstant("e", std::exp(Real(1)));
627  if (fp.Parse(_subexpressions.back(), variables) != -1) // -1 for success
628  libmesh_error_msg
629  ("ERROR: FunctionParser is unable to parse expression: "
630  << _subexpressions.back() << '\n' << fp.ErrorMsg());
631  fp.Optimize();
632  parsers.push_back(fp);
633 #else
634  libmesh_error_msg("ERROR: This functionality requires fparser!");
635 #endif
636 
637  // If at end, use nextstart=maxSize. Else start at next
638  // character.
639  nextstart = (end == std::string::npos) ?
640  std::string::npos : end + 1;
641  }
642 }
643 
644 template <typename Output>
645 inline
646 std::size_t
647 ParsedFEMFunction<Output>::find_name (const std::string & varname,
648  const std::string & expr) const
649 {
650  const std::size_t namesize = varname.size();
651  std::size_t varname_i = expr.find(varname);
652 
653  while ((varname_i != std::string::npos) &&
654  (((varname_i > 0) &&
655  (std::isalnum(expr[varname_i-1]) ||
656  (expr[varname_i-1] == '_'))) ||
657  ((varname_i+namesize < expr.size()) &&
658  (std::isalnum(expr[varname_i+namesize]) ||
659  (expr[varname_i+namesize] == '_')))))
660  {
661  varname_i = expr.find(varname, varname_i+1);
662  }
663 
664  return varname_i;
665 }
666 
667 
668 
669 // Helper function for evaluating function arguments
670 template <typename Output>
671 inline
672 void
674  const Point & p,
675  const Real time)
676 {
677  _spacetime[0] = p(0);
678 #if LIBMESH_DIM > 1
679  _spacetime[1] = p(1);
680 #endif
681 #if LIBMESH_DIM > 2
682  _spacetime[2] = p(2);
683 #endif
684  _spacetime[LIBMESH_DIM] = time;
685 
686  unsigned int request_index = 0;
687  for (unsigned int v=0; v != _n_vars; ++v)
688  {
689  if (!_need_var[v])
690  continue;
691 
692  c.point_value(v, p, _spacetime[LIBMESH_DIM+1+request_index]);
693  request_index++;
694  }
695 
696  if (_n_requested_grad_components)
697  for (unsigned int v=0; v != _n_vars; ++v)
698  {
699  if (!_need_var_grad[v*LIBMESH_DIM]
700 #if LIBMESH_DIM > 1
701  && !_need_var_grad[v*LIBMESH_DIM+1]
702 #if LIBMESH_DIM > 2
703  && !_need_var_grad[v*LIBMESH_DIM+2]
704 #endif
705 #endif
706  )
707  continue;
708 
709  Gradient g;
710  c.point_gradient(v, p, g);
711 
712  for (unsigned int d=0; d != LIBMESH_DIM; ++d)
713  {
714  if (!_need_var_grad[v*LIBMESH_DIM+d])
715  continue;
716 
717  _spacetime[LIBMESH_DIM+1+request_index] = g(d);
718  request_index++;
719  }
720  }
721 
722 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
723  if (_n_requested_hess_components)
724  for (unsigned int v=0; v != _n_vars; ++v)
725  {
726  if (!_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM]
727 #if LIBMESH_DIM > 1
728  && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+1]
729  && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+2]
730  && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+3]
731 #if LIBMESH_DIM > 2
732  && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+4]
733  && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+5]
734  && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+6]
735  && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+7]
736  && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+8]
737 #endif
738 #endif
739  )
740  continue;
741 
742  Tensor h;
743  c.point_hessian(v, p, h);
744 
745  for (unsigned int d1=0; d1 != LIBMESH_DIM; ++d1)
746  for (unsigned int d2=0; d2 != LIBMESH_DIM; ++d2)
747  {
748  if (!_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+d1*LIBMESH_DIM+d2])
749  continue;
750 
751  _spacetime[LIBMESH_DIM+1+request_index] = h(d1,d2);
752  request_index++;
753  }
754  }
755 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
756 
757  if (_requested_normals)
758  {
759  FEBase * side_fe;
760  c.get_side_fe(0, side_fe);
761 
762  const std::vector<Point> & normals = side_fe->get_normals();
763 
764  const std::vector<Point> & xyz = side_fe->get_xyz();
765 
766  libmesh_assert_equal_to(normals.size(), xyz.size());
767 
768  // We currently only support normals at quadrature points!
769 #ifndef NDEBUG
770  bool at_quadrature_point = false;
771 #endif
772  for (std::size_t qp = 0; qp != normals.size(); ++qp)
773  {
774  if (p == xyz[qp])
775  {
776  const Point & n = normals[qp];
777  for (unsigned int d=0; d != LIBMESH_DIM; ++d)
778  {
779  _spacetime[LIBMESH_DIM+1+request_index] = n(d);
780  request_index++;
781  }
782 #ifndef NDEBUG
783  at_quadrature_point = true;
784 #endif
785  break;
786  }
787  }
788 
789  libmesh_assert(at_quadrature_point);
790  }
791 
792  // The remaining locations in _spacetime are currently set at construction
793  // but could potentially be made dynamic
794 }
795 
796 
797 // Evaluate the ith FunctionParser and check the result
798 #ifdef LIBMESH_HAVE_FPARSER
799 template <typename Output>
800 inline
801 Output
802 ParsedFEMFunction<Output>::eval (FunctionParserBase<Output> & parser,
803  const std::string & libmesh_dbg_var(function_name),
804  unsigned int libmesh_dbg_var(component_idx)) const
805 {
806 #ifndef NDEBUG
807  Output result = parser.Eval(_spacetime.data());
808  int error_code = parser.EvalError();
809  if (error_code)
810  {
811  libMesh::err << "ERROR: FunctionParser is unable to evaluate component "
812  << component_idx
813  << " of expression '"
814  << function_name
815  << "' with arguments:\n";
816  for (std::size_t j=0; j<_spacetime.size(); ++j)
817  libMesh::err << '\t' << _spacetime[j] << '\n';
818  libMesh::err << '\n';
819 
820  // Currently no API to report error messages, we'll do it manually
821  std::string error_message = "Reason: ";
822 
823  switch (error_code)
824  {
825  case 1:
826  error_message += "Division by zero";
827  break;
828  case 2:
829  error_message += "Square Root error (negative value)";
830  break;
831  case 3:
832  error_message += "Log error (negative value)";
833  break;
834  case 4:
835  error_message += "Trigonometric error (asin or acos of illegal value)";
836  break;
837  case 5:
838  error_message += "Maximum recursion level reached";
839  break;
840  default:
841  error_message += "Unknown";
842  break;
843  }
844  libmesh_error_msg(error_message);
845  }
846 
847  return result;
848 #else
849  return parser.Eval(_spacetime.data());
850 #endif
851 }
852 #else // LIBMESH_HAVE_FPARSER
853 template <typename Output>
854 inline
855 Output
856 ParsedFEMFunction<Output>::eval (char & /*parser*/,
857  const std::string & /*function_name*/,
858  unsigned int /*component_idx*/) const
859 {
860  libmesh_error_msg("ERROR: This functionality requires fparser!");
861  return Output(0);
862 }
863 #endif
864 
865 
866 } // namespace libMesh
867 
868 #endif // LIBMESH_PARSED_FEM_FUNCTION_H
virtual unsigned int size() const override
Definition: dense_vector.h:92
void set_inline_value(const std::string &inline_var_name, Output newval)
std::vector< bool > _need_var_grad
virtual void init_context(const FEMContext &c) override
std::vector< Output > _initial_vals
const std::vector< Point > & get_xyz() const
Definition: fe_abstract.h:238
void get_side_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Definition: fem_context.h:299
void partial_reparse(const std::string &expression)
std::vector< std::string > _additional_vars
IterBase * end
std::vector< std::string > _subexpressions
std::vector< Output > _spacetime
const unsigned int n_vars
Definition: tecplot_io.C:69
virtual std::unique_ptr< FEMFunctionBase< Output > > clone() const override
ParsedFEMFunction & operator=(const ParsedFEMFunction &)=delete
Number point_value(unsigned int var, const Point &p) const
Definition: fem_context.C:768
std::size_t find_name(const std::string &varname, const std::string &expr) const
virtual Output operator()(const FEMContext &c, const Point &p, const Real time=0.) override
const std::vector< std::vector< OutputGradient > > & get_dphi() const
Definition: fe_base.h:215
const std::vector< Point > & get_normals() const
Definition: fe_abstract.h:385
void reparse(const std::string &expression)
Support for using parsed functions in FEMSystem.
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:92
Output eval(FunctionParserBase< Output > &parser, const std::string &libmesh_dbg_var(function_name), unsigned int libmesh_dbg_var(component_idx)) const
virtual ~ParsedFEMFunction()=default
Tensor point_hessian(unsigned int var, const Point &p) const
Definition: fem_context.C:866
Output get_inline_value(const std::string &inline_var_name) const
OStreamProxy err(std::cerr)
const std::vector< std::vector< OutputTensor > > & get_d2phi() const
Definition: fe_base.h:289
std::vector< FunctionParserBase< Output > > parsers
std::vector< bool > _need_var_hess
void eval_args(const FEMContext &c, const Point &p, const Real time)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Output component(const FEMContext &c, unsigned int i, const Point &p, Real time=0.) override
void get_element_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Definition: fem_context.h:262
ParsedFEMFunction(const System &sys, const std::string &expression, const std::vector< std::string > *additional_vars=nullptr, const std::vector< Output > *initial_vals=nullptr)
Gradient point_gradient(unsigned int var, const Point &p) const
Definition: fem_context.C:814
A geometric point in (x,y,z) space.
Definition: point.h:38
const std::string & expression()
const std::vector< std::vector< OutputShape > > & get_phi() const
Definition: fe_base.h:207