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