20 #ifndef LIBMESH_PARSED_FEM_FUNCTION_H
21 #define LIBMESH_PARSED_FEM_FUNCTION_H
23 #include "libmesh/libmesh_config.h"
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"
32 #ifdef LIBMESH_HAVE_FPARSER
34 #include "libmesh/fparser.hh"
56 template <
typename Output=Number>
67 const std::vector<std::string> * additional_vars=
nullptr,
68 const std::vector<Output> * initial_vals=
nullptr);
96 virtual std::unique_ptr<FEMFunctionBase<Output>>
clone ()
const override;
100 const Real time = 0.)
override;
110 Real time=0.)
override;
142 std::size_t
find_name (
const std::string & varname,
143 const std::string & expr)
const;
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;
170 #ifdef LIBMESH_HAVE_FPARSER
171 std::vector<FunctionParserBase<Output>>
parsers;
185 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
189 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
200 template <
typename Output>
203 const std::string & expression,
204 const std::vector<std::string> * additional_vars,
205 const std::vector<Output> * initial_vals) :
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),
218 _additional_vars (additional_vars ? *additional_vars :
std::vector<
std::string>()),
219 _initial_vals (initial_vals ? *initial_vals :
std::vector<Output>())
221 this->reparse(expression);
225 template <
typename Output>
239 for (
unsigned int v=0; v != _n_vars; ++v)
241 const std::string & varname = _sys.variable_name(v);
242 std::size_t varname_i = find_name(varname, expression);
246 if (varname_i == std::string::npos)
251 variables += varname;
255 for (
unsigned int v=0; v != _n_vars; ++v)
257 const std::string & varname = _sys.variable_name(v);
259 for (
unsigned int d=0; d != LIBMESH_DIM; ++d)
261 std::string gradname = std::string(
"grad_") +
262 "xyz"[d] +
'_' + varname;
263 std::size_t gradname_i = find_name(gradname, expression);
267 if (gradname_i == std::string::npos)
270 _need_var_grad[v*LIBMESH_DIM+d] =
true;
272 variables += gradname;
273 _n_requested_grad_components++;
277 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
278 for (
unsigned int v=0; v != _n_vars; ++v)
280 const std::string & varname = _sys.variable_name(v);
282 for (
unsigned int d1=0; d1 != LIBMESH_DIM; ++d1)
283 for (
unsigned int d2=0; d2 != LIBMESH_DIM; ++d2)
285 std::string hessname = std::string(
"hess_") +
286 "xyz"[d1] +
"xyz"[d2] +
'_' + varname;
287 std::size_t hessname_i = find_name(hessname, expression);
291 if (hessname_i == std::string::npos)
294 _need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+d1*LIBMESH_DIM+d2]
297 variables += hessname;
298 _n_requested_hess_components++;
301 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
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);
310 if (nx_i != std::string::npos ||
311 ny_i != std::string::npos ||
312 nz_i != std::string::npos)
314 _requested_normals =
true;
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());
339 unsigned int offset = LIBMESH_DIM + 1 + _n_requested_vars +
340 _n_requested_grad_components + _n_requested_hess_components;
344 variables +=
"," + _additional_vars[i];
347 _spacetime[offset + i] =
348 (i < _initial_vals.size()) ? _initial_vals[i] : 0;
351 this->partial_reparse(expression);
354 template <
typename Output>
359 for (
unsigned int v=0; v != _n_vars; ++v)
363 if (_n_requested_vars)
365 if (_n_requested_grad_components)
367 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
368 if (_n_requested_hess_components)
370 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
373 if (_requested_normals)
386 template <
typename Output>
388 std::unique_ptr<FEMFunctionBase<Output>>
391 return libmesh_make_unique<ParsedFEMFunction>
392 (_sys, _expression, &_additional_vars, &_initial_vals);
395 template <
typename Output>
402 eval_args(c, p, time);
404 return eval(parsers[0],
"f", 0);
409 template <
typename Output>
417 eval_args(c, p, time);
419 unsigned int size = output.
size();
421 libmesh_assert_equal_to (size, parsers.size());
423 for (
unsigned int i=0; i != size; ++i)
424 output(i) = eval(parsers[i],
"f", i);
428 template <
typename Output>
436 eval_args(c, p, time);
438 libmesh_assert_less (i, parsers.size());
439 return eval(parsers[i],
"f", i);
442 template <
typename Output>
447 libmesh_assert_greater (_subexpressions.size(), 0);
450 bool found_var_name =
false;
452 Output old_var_value(0.);
454 for (
const std::string & subexpression : _subexpressions)
456 const std::size_t varname_i =
457 find_name(inline_var_name, subexpression);
458 if (varname_i == std::string::npos)
461 const std::size_t assignment_i =
462 subexpression.find(
":", varname_i+1);
464 libmesh_assert_not_equal_to(assignment_i, std::string::npos);
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],
' ');
470 std::size_t end_assignment_i =
471 subexpression.find(
";", assignment_i+1);
473 libmesh_assert_not_equal_to(end_assignment_i, std::string::npos);
475 std::string new_subexpression =
476 subexpression.substr(0, end_assignment_i+1) +
479 #ifdef LIBMESH_HAVE_FPARSER
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)
488 (
"ERROR: FunctionParser is unable to parse modified expression: "
489 << new_subexpression <<
'\n' << fp.ErrorMsg());
491 Output new_var_value = this->eval(fp, new_subexpression, 0);
493 return new_var_value;
497 libmesh_assert_equal_to(old_var_value, new_var_value);
501 old_var_value = new_var_value;
502 found_var_name =
true;
507 libmesh_error_msg(
"ERROR: This functionality requires fparser!");
512 return old_var_value;
515 template <
typename Output>
521 libmesh_assert_greater (_subexpressions.size(), 0);
524 bool found_var_name =
false;
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)
535 found_var_name =
true;
538 const std::size_t assignment_i =
539 subexpression.find(
":", varname_i+1);
541 libmesh_assert_not_equal_to(assignment_i, std::string::npos);
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],
' ');
547 std::size_t end_assignment_i =
548 subexpression.find(
";", assignment_i+1);
550 libmesh_assert_not_equal_to(end_assignment_i, std::string::npos);
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' <<
')'
561 << subexpression.substr(end_assignment_i,
563 _subexpressions[s] = new_subexpression.str();
568 std::string new_expression;
570 for (
const auto & subexpression : _subexpressions)
572 new_expression +=
'{';
573 new_expression += subexpression;
574 new_expression +=
'}';
577 this->partial_reparse(new_expression);
581 template <
typename Output>
586 _expression = expression;
587 _subexpressions.clear();
590 size_t nextstart = 0,
end = 0;
592 while (
end != std::string::npos)
596 if (nextstart >= expression.size())
601 if (expression[nextstart] ==
'{')
604 end = expression.find(
'}', nextstart);
608 end = std::string::npos;
612 _subexpressions.push_back
613 (expression.substr(nextstart, (
end == std::string::npos) ?
614 std::string::npos :
end - nextstart));
617 if (_subexpressions.back().empty())
618 libmesh_error_msg(
"ERROR: FunctionParser is unable to parse empty expression.\n");
621 #ifdef LIBMESH_HAVE_FPARSER
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)
630 (
"ERROR: FunctionParser is unable to parse expression: "
631 << _subexpressions.back() <<
'\n' << fp.ErrorMsg());
633 parsers.push_back(fp);
635 libmesh_error_msg(
"ERROR: This functionality requires fparser!");
640 nextstart = (
end == std::string::npos) ?
641 std::string::npos :
end + 1;
645 template <
typename Output>
649 const std::string & expr)
const
651 const std::size_t namesize = varname.size();
652 std::size_t varname_i = expr.find(varname);
654 while ((varname_i != std::string::npos) &&
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] ==
'_')))))
662 varname_i = expr.find(varname, varname_i+1);
671 template <
typename Output>
678 _spacetime[0] = p(0);
680 _spacetime[1] = p(1);
683 _spacetime[2] = p(2);
685 _spacetime[LIBMESH_DIM] = time;
687 unsigned int request_index = 0;
688 for (
unsigned int v=0; v != _n_vars; ++v)
693 c.
point_value(v, p, _spacetime[LIBMESH_DIM+1+request_index]);
697 if (_n_requested_grad_components)
698 for (
unsigned int v=0; v != _n_vars; ++v)
700 if (!_need_var_grad[v*LIBMESH_DIM]
702 && !_need_var_grad[v*LIBMESH_DIM+1]
704 && !_need_var_grad[v*LIBMESH_DIM+2]
713 for (
unsigned int d=0; d != LIBMESH_DIM; ++d)
715 if (!_need_var_grad[v*LIBMESH_DIM+d])
718 _spacetime[LIBMESH_DIM+1+request_index] = g(d);
723 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
724 if (_n_requested_hess_components)
725 for (
unsigned int v=0; v != _n_vars; ++v)
727 if (!_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM]
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]
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]
746 for (
unsigned int d1=0; d1 != LIBMESH_DIM; ++d1)
747 for (
unsigned int d2=0; d2 != LIBMESH_DIM; ++d2)
749 if (!_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+d1*LIBMESH_DIM+d2])
752 _spacetime[LIBMESH_DIM+1+request_index] = h(d1,d2);
756 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
758 if (_requested_normals)
763 const std::vector<Point> & normals = side_fe->
get_normals();
765 const std::vector<Point> & xyz = side_fe->
get_xyz();
767 libmesh_assert_equal_to(normals.size(), xyz.size());
771 bool at_quadrature_point =
false;
777 const Point & n = normals[qp];
778 for (
unsigned int d=0; d != LIBMESH_DIM; ++d)
780 _spacetime[LIBMESH_DIM+1+request_index] = n(d);
784 at_quadrature_point =
true;
799 #ifdef LIBMESH_HAVE_FPARSER
800 template <
typename Output>
804 const std::string & libmesh_dbg_var(function_name),
805 unsigned int libmesh_dbg_var(component_idx))
const
808 Output result = parser.Eval(_spacetime.data());
809 int error_code = parser.EvalError();
812 libMesh::err <<
"ERROR: FunctionParser is unable to evaluate component "
814 <<
" of expression '"
816 <<
"' with arguments:\n";
817 for (
const auto & st : _spacetime)
822 std::string error_message =
"Reason: ";
827 error_message +=
"Division by zero";
830 error_message +=
"Square Root error (negative value)";
833 error_message +=
"Log error (negative value)";
836 error_message +=
"Trigonometric error (asin or acos of illegal value)";
839 error_message +=
"Maximum recursion level reached";
842 error_message +=
"Unknown";
845 libmesh_error_msg(error_message);
850 return parser.Eval(_spacetime.data());
853 #else // LIBMESH_HAVE_FPARSER
854 template <
typename Output>
858 const std::string & ,
861 libmesh_error_msg(
"ERROR: This functionality requires fparser!");
869 #endif // LIBMESH_PARSED_FEM_FUNCTION_H