20 #ifndef LIBMESH_PARSED_FEM_FUNCTION_H 21 #define LIBMESH_PARSED_FEM_FUNCTION_H 30 #ifdef LIBMESH_HAVE_FPARSER 32 #include "libmesh/fparser.hh" 54 template <
typename Output=Number>
65 const std::vector<std::string> * additional_vars=
nullptr,
66 const std::vector<Output> * initial_vals=
nullptr);
94 virtual std::unique_ptr<FEMFunctionBase<Output>>
clone ()
const override;
98 const Real time = 0.)
override;
108 Real time=0.)
override;
140 std::size_t
find_name (
const std::string & varname,
141 const std::string & expr)
const;
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;
168 #ifdef LIBMESH_HAVE_FPARSER 169 std::vector<FunctionParserBase<Output>>
parsers;
183 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 187 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 198 template <
typename Output>
201 const std::string & expression,
202 const std::vector<std::string> * additional_vars,
203 const std::vector<Output> * initial_vals) :
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),
216 _additional_vars (additional_vars ? *additional_vars :
std::vector<
std::string>()),
217 _initial_vals (initial_vals ? *initial_vals :
std::vector<Output>())
219 this->reparse(expression);
223 template <
typename Output>
237 for (
unsigned int v=0; v != _n_vars; ++v)
239 const std::string & varname = _sys.variable_name(v);
240 std::size_t varname_i = find_name(varname, expression);
244 if (varname_i == std::string::npos)
249 variables += varname;
253 for (
unsigned int v=0; v != _n_vars; ++v)
255 const std::string & varname = _sys.variable_name(v);
257 for (
unsigned int d=0; d != LIBMESH_DIM; ++d)
259 std::string gradname = std::string(
"grad_") +
260 "xyz"[d] +
'_' + varname;
261 std::size_t gradname_i = find_name(gradname, expression);
265 if (gradname_i == std::string::npos)
268 _need_var_grad[v*LIBMESH_DIM+d] =
true;
270 variables += gradname;
271 _n_requested_grad_components++;
275 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 276 for (
unsigned int v=0; v != _n_vars; ++v)
278 const std::string & varname = _sys.variable_name(v);
280 for (
unsigned int d1=0; d1 != LIBMESH_DIM; ++d1)
281 for (
unsigned int d2=0; d2 != LIBMESH_DIM; ++d2)
283 std::string hessname = std::string(
"hess_") +
284 "xyz"[d1] +
"xyz"[d2] +
'_' + varname;
285 std::size_t hessname_i = find_name(hessname, expression);
289 if (hessname_i == std::string::npos)
292 _need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+d1*LIBMESH_DIM+d2]
295 variables += hessname;
296 _n_requested_hess_components++;
299 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 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);
308 if (nx_i != std::string::npos ||
309 ny_i != std::string::npos ||
310 nz_i != std::string::npos)
312 _requested_normals =
true;
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());
337 unsigned int offset = LIBMESH_DIM + 1 + _n_requested_vars +
338 _n_requested_grad_components + _n_requested_hess_components;
340 for (std::size_t i=0; i < _additional_vars.size(); ++i)
342 variables +=
"," + _additional_vars[i];
345 _spacetime[offset + i] =
346 (i < _initial_vals.size()) ? _initial_vals[i] : 0;
349 this->partial_reparse(expression);
352 template <
typename Output>
357 for (
unsigned int v=0; v != _n_vars; ++v)
361 if (_n_requested_vars)
363 if (_n_requested_grad_components)
365 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 366 if (_n_requested_hess_components)
368 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 371 if (_requested_normals)
384 template <
typename Output>
386 std::unique_ptr<FEMFunctionBase<Output>>
389 return std::unique_ptr<FEMFunctionBase<Output>>
393 template <
typename Output>
400 eval_args(c, p, time);
402 return eval(parsers[0],
"f", 0);
407 template <
typename Output>
415 eval_args(c, p, time);
417 unsigned int size = output.
size();
419 libmesh_assert_equal_to (size, parsers.size());
421 for (
unsigned int i=0; i != size; ++i)
422 output(i) = eval(parsers[i],
"f", i);
426 template <
typename Output>
434 eval_args(c, p, time);
436 libmesh_assert_less (i, parsers.size());
437 return eval(parsers[i],
"f", i);
440 template <
typename Output>
445 libmesh_assert_greater (_subexpressions.size(), 0);
448 bool found_var_name =
false;
450 Output old_var_value(0.);
452 for (std::size_t s=0; s != _subexpressions.size(); ++s)
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)
460 const std::size_t assignment_i =
461 subexpression.find(
":", varname_i+1);
463 libmesh_assert_not_equal_to(assignment_i, std::string::npos);
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],
' ');
469 std::size_t end_assignment_i =
470 subexpression.find(
";", assignment_i+1);
472 libmesh_assert_not_equal_to(end_assignment_i, std::string::npos);
474 std::string new_subexpression =
475 subexpression.substr(0, end_assignment_i+1) +
478 #ifdef LIBMESH_HAVE_FPARSER 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)
487 (
"ERROR: FunctionParser is unable to parse modified expression: " 488 << new_subexpression <<
'\n' << fp.ErrorMsg());
490 Output new_var_value = this->eval(fp, new_subexpression, 0);
492 return new_var_value;
496 libmesh_assert_equal_to(old_var_value, new_var_value);
500 old_var_value = new_var_value;
501 found_var_name =
true;
506 libmesh_error_msg(
"ERROR: This functionality requires fparser!");
510 libmesh_assert(found_var_name);
511 return old_var_value;
514 template <
typename Output>
520 libmesh_assert_greater (_subexpressions.size(), 0);
523 bool found_var_name =
false;
525 for (std::size_t s=0; s != _subexpressions.size(); ++s)
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)
534 found_var_name =
true;
537 const std::size_t assignment_i =
538 subexpression.find(
":", varname_i+1);
540 libmesh_assert_not_equal_to(assignment_i, std::string::npos);
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],
' ');
546 std::size_t end_assignment_i =
547 subexpression.find(
";", assignment_i+1);
549 libmesh_assert_not_equal_to(end_assignment_i, std::string::npos);
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' <<
')' 560 << subexpression.substr(end_assignment_i,
562 _subexpressions[s] = new_subexpression.str();
565 libmesh_assert(found_var_name);
567 std::string new_expression;
569 for (std::size_t s=0; s != _subexpressions.size(); ++s)
571 new_expression +=
'{';
572 new_expression += _subexpressions[s];
573 new_expression +=
'}';
576 this->partial_reparse(new_expression);
580 template <
typename Output>
585 _expression = expression;
586 _subexpressions.clear();
589 size_t nextstart = 0,
end = 0;
591 while (
end != std::string::npos)
595 if (nextstart >= expression.size())
600 if (expression[nextstart] ==
'{')
603 end = expression.find(
'}', nextstart);
607 end = std::string::npos;
611 _subexpressions.push_back
612 (expression.substr(nextstart, (
end == std::string::npos) ?
613 std::string::npos :
end - nextstart));
616 if (_subexpressions.back().empty())
617 libmesh_error_msg(
"ERROR: FunctionParser is unable to parse empty expression.\n");
620 #ifdef LIBMESH_HAVE_FPARSER 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)
629 (
"ERROR: FunctionParser is unable to parse expression: " 630 << _subexpressions.back() <<
'\n' << fp.ErrorMsg());
632 parsers.push_back(fp);
634 libmesh_error_msg(
"ERROR: This functionality requires fparser!");
639 nextstart = (
end == std::string::npos) ?
640 std::string::npos :
end + 1;
644 template <
typename Output>
648 const std::string & expr)
const 650 const std::size_t namesize = varname.size();
651 std::size_t varname_i = expr.find(varname);
653 while ((varname_i != std::string::npos) &&
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] ==
'_')))))
661 varname_i = expr.find(varname, varname_i+1);
670 template <
typename Output>
677 _spacetime[0] = p(0);
679 _spacetime[1] = p(1);
682 _spacetime[2] = p(2);
684 _spacetime[LIBMESH_DIM] = time;
686 unsigned int request_index = 0;
687 for (
unsigned int v=0; v != _n_vars; ++v)
692 c.
point_value(v, p, _spacetime[LIBMESH_DIM+1+request_index]);
696 if (_n_requested_grad_components)
697 for (
unsigned int v=0; v != _n_vars; ++v)
699 if (!_need_var_grad[v*LIBMESH_DIM]
701 && !_need_var_grad[v*LIBMESH_DIM+1]
703 && !_need_var_grad[v*LIBMESH_DIM+2]
712 for (
unsigned int d=0; d != LIBMESH_DIM; ++d)
714 if (!_need_var_grad[v*LIBMESH_DIM+d])
717 _spacetime[LIBMESH_DIM+1+request_index] = g(d);
722 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 723 if (_n_requested_hess_components)
724 for (
unsigned int v=0; v != _n_vars; ++v)
726 if (!_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM]
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]
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]
745 for (
unsigned int d1=0; d1 != LIBMESH_DIM; ++d1)
746 for (
unsigned int d2=0; d2 != LIBMESH_DIM; ++d2)
748 if (!_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+d1*LIBMESH_DIM+d2])
751 _spacetime[LIBMESH_DIM+1+request_index] = h(d1,d2);
755 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 757 if (_requested_normals)
762 const std::vector<Point> & normals = side_fe->
get_normals();
764 const std::vector<Point> & xyz = side_fe->
get_xyz();
766 libmesh_assert_equal_to(normals.size(), xyz.size());
770 bool at_quadrature_point =
false;
772 for (std::size_t qp = 0; qp != normals.size(); ++qp)
776 const Point & n = normals[qp];
777 for (
unsigned int d=0; d != LIBMESH_DIM; ++d)
779 _spacetime[LIBMESH_DIM+1+request_index] = n(d);
783 at_quadrature_point =
true;
789 libmesh_assert(at_quadrature_point);
798 #ifdef LIBMESH_HAVE_FPARSER 799 template <
typename Output>
803 const std::string & libmesh_dbg_var(function_name),
804 unsigned int libmesh_dbg_var(component_idx))
const 807 Output result = parser.Eval(_spacetime.data());
808 int error_code = parser.EvalError();
811 libMesh::err <<
"ERROR: FunctionParser is unable to evaluate component " 813 <<
" of expression '" 815 <<
"' with arguments:\n";
816 for (std::size_t j=0; j<_spacetime.size(); ++j)
821 std::string error_message =
"Reason: ";
826 error_message +=
"Division by zero";
829 error_message +=
"Square Root error (negative value)";
832 error_message +=
"Log error (negative value)";
835 error_message +=
"Trigonometric error (asin or acos of illegal value)";
838 error_message +=
"Maximum recursion level reached";
841 error_message +=
"Unknown";
844 libmesh_error_msg(error_message);
849 return parser.Eval(_spacetime.data());
852 #else // LIBMESH_HAVE_FPARSER 853 template <
typename Output>
857 const std::string & ,
860 libmesh_error_msg(
"ERROR: This functionality requires fparser!");
868 #endif // LIBMESH_PARSED_FEM_FUNCTION_H virtual unsigned int size() const override
void set_inline_value(const std::string &inline_var_name, Output newval)
unsigned int _n_requested_vars
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
void get_side_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
void partial_reparse(const std::string &expression)
std::vector< std::string > _additional_vars
std::vector< std::string > _subexpressions
std::vector< Output > _spacetime
const unsigned int n_vars
virtual std::unique_ptr< FEMFunctionBase< Output > > clone() const override
ParsedFEMFunction & operator=(const ParsedFEMFunction &)=delete
Number point_value(unsigned int var, const Point &p) const
std::size_t find_name(const std::string &varname, const std::string &expr) const
unsigned int _n_requested_grad_components
virtual Output operator()(const FEMContext &c, const Point &p, const Real time=0.) override
const std::vector< std::vector< OutputGradient > > & get_dphi() const
const std::vector< Point > & get_normals() const
void reparse(const std::string &expression)
Support for using parsed functions in FEMSystem.
Manages consistently variables, degrees of freedom, and coefficient vectors.
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
Output get_inline_value(const std::string &inline_var_name) const
OStreamProxy err(std::cerr)
const std::vector< std::vector< OutputTensor > > & get_d2phi() const
std::vector< FunctionParserBase< Output > > parsers
std::vector< bool > _need_var_hess
void eval_args(const FEMContext &c, const Point &p, const Real time)
unsigned int _n_requested_hess_components
std::vector< char > parsers
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::vector< bool > _need_var
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
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
A geometric point in (x,y,z) space.
const std::string & expression()
const std::vector< std::vector< OutputShape > > & get_phi() const