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" 31 #ifdef LIBMESH_HAVE_FPARSER 33 #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);
91 virtual std::unique_ptr<FEMFunctionBase<Output>>
clone ()
const override;
95 const Real time = 0.)
override;
105 Real time=0.)
override;
137 std::size_t
find_name (std::string_view varname,
138 std::string_view expr)
const;
146 #ifdef LIBMESH_HAVE_FPARSER 147 inline Output
eval(FunctionParserBase<Output> & parser,
148 std::string_view libmesh_dbg_var(function_name),
149 unsigned int libmesh_dbg_var(component_idx))
const;
150 #else // LIBMESH_HAVE_FPARSER 151 inline Output
eval(
char & libmesh_dbg_var(parser),
152 std::string_view libmesh_dbg_var(function_name),
153 unsigned int libmesh_dbg_var(component_idx))
const;
165 #ifdef LIBMESH_HAVE_FPARSER 166 std::vector<std::unique_ptr<FunctionParserBase<Output>>>
parsers;
180 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 184 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 195 template <
typename Output>
198 std::string expression,
199 const std::vector<std::string> * additional_vars,
200 const std::vector<Output> * initial_vals) :
204 _n_requested_vars(0),
205 _n_requested_grad_components(0),
206 _n_requested_hess_components(0),
207 _requested_normals(false),
208 _need_var(_n_vars, false),
209 _need_var_grad(_n_vars*LIBMESH_DIM, false),
210 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
211 _need_var_hess(_n_vars*LIBMESH_DIM*LIBMESH_DIM, false),
213 _additional_vars (additional_vars ? *additional_vars :
std::vector<
std::string>()),
214 _initial_vals (initial_vals ? *initial_vals :
std::vector<Output>())
216 this->reparse(std::move(expression));
220 template <
typename Output>
225 _expression(other._expression),
226 _subexpressions(other._subexpressions),
227 _n_vars(other._n_vars),
228 _n_requested_vars(other._n_requested_vars),
229 _n_requested_grad_components(other._n_requested_grad_components),
230 _n_requested_hess_components(other._n_requested_hess_components),
231 _requested_normals(other._requested_normals),
232 _spacetime(other._spacetime),
233 _need_var(other._need_var),
234 _need_var_grad(other._need_var_grad),
235 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
236 _need_var_hess(other._need_var_hess),
238 variables(other.variables),
239 _additional_vars(other._additional_vars),
240 _initial_vals(other._initial_vals)
242 this->reparse(_expression);
246 template <
typename Output>
256 std::swap(tmp, *
this);
261 template <
typename Output>
275 for (
unsigned int v=0; v != _n_vars; ++v)
277 const std::string & varname = _sys.variable_name(v);
278 std::size_t varname_i = find_name(varname, expression);
282 if (varname_i == std::string::npos)
287 variables += varname;
291 for (
unsigned int v=0; v != _n_vars; ++v)
293 const std::string & varname = _sys.variable_name(v);
295 for (
unsigned int d=0; d != LIBMESH_DIM; ++d)
297 std::string gradname = std::string(
"grad_") +
298 "xyz"[d] +
'_' + varname;
299 std::size_t gradname_i = find_name(gradname, expression);
303 if (gradname_i == std::string::npos)
306 _need_var_grad[v*LIBMESH_DIM+d] =
true;
308 variables += gradname;
309 _n_requested_grad_components++;
313 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 314 for (
unsigned int v=0; v != _n_vars; ++v)
316 const std::string & varname = _sys.variable_name(v);
318 for (
unsigned int d1=0; d1 != LIBMESH_DIM; ++d1)
319 for (
unsigned int d2=0; d2 != LIBMESH_DIM; ++d2)
321 std::string hessname = std::string(
"hess_") +
322 "xyz"[d1] +
"xyz"[d2] +
'_' + varname;
323 std::size_t hessname_i = find_name(hessname, expression);
327 if (hessname_i == std::string::npos)
330 _need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+d1*LIBMESH_DIM+d2]
333 variables += hessname;
334 _n_requested_hess_components++;
337 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 340 std::size_t nx_i = find_name(
"n_x", expression);
341 std::size_t ny_i = find_name(
"n_y", expression);
342 std::size_t nz_i = find_name(
"n_z", expression);
346 if (nx_i != std::string::npos ||
347 ny_i != std::string::npos ||
348 nz_i != std::string::npos)
350 _requested_normals =
true;
367 (LIBMESH_DIM + 1 + _n_requested_vars +
368 _n_requested_grad_components + _n_requested_hess_components +
369 (_requested_normals ? LIBMESH_DIM : 0) +
370 _additional_vars.size());
375 unsigned int offset = LIBMESH_DIM + 1 + _n_requested_vars +
376 _n_requested_grad_components + _n_requested_hess_components;
380 variables +=
"," + _additional_vars[i];
383 _spacetime[offset + i] =
384 (i < _initial_vals.size()) ? _initial_vals[i] : 0;
387 this->partial_reparse(std::move(expression));
390 template <
typename Output>
395 for (
unsigned int v=0; v != _n_vars; ++v)
399 bool request_nothing =
true;
400 if (_n_requested_vars)
403 request_nothing =
false;
405 if (_n_requested_grad_components)
408 request_nothing =
false;
410 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 411 if (_n_requested_hess_components)
414 request_nothing =
false;
416 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 421 if (_requested_normals)
434 template <
typename Output>
436 std::unique_ptr<FEMFunctionBase<Output>>
439 return std::make_unique<ParsedFEMFunction>
440 (_sys, _expression, &_additional_vars, &_initial_vals);
443 template <
typename Output>
450 eval_args(c, p, time);
452 return eval(*parsers[0],
"f", 0);
457 template <
typename Output>
465 eval_args(c, p, time);
467 unsigned int size = output.
size();
469 libmesh_assert_equal_to (size, parsers.size());
471 for (
unsigned int i=0; i != size; ++i)
472 output(i) = eval(*parsers[i],
"f", i);
476 template <
typename Output>
484 eval_args(c, p, time);
486 libmesh_assert_less (i, parsers.size());
487 return eval(*parsers[i],
"f", i);
490 template <
typename Output>
495 libmesh_assert_greater (_subexpressions.size(), 0);
498 bool found_var_name =
false;
500 Output old_var_value(0.);
502 for (
const std::string & subexpression : _subexpressions)
504 const std::size_t varname_i =
505 find_name(inline_var_name, subexpression);
506 if (varname_i == std::string::npos)
509 const std::size_t assignment_i =
510 subexpression.find(
":", varname_i+1);
512 libmesh_assert_not_equal_to(assignment_i, std::string::npos);
514 libmesh_assert_equal_to(subexpression[assignment_i+1],
'=');
515 for (std::size_t i = varname_i+1; i != assignment_i; ++i)
516 libmesh_assert_equal_to(subexpression[i],
' ');
518 std::size_t end_assignment_i =
519 subexpression.find(
";", assignment_i+1);
521 libmesh_assert_not_equal_to(end_assignment_i, std::string::npos);
523 std::string new_subexpression =
524 subexpression.substr(0, end_assignment_i+1).
525 append(inline_var_name);
527 #ifdef LIBMESH_HAVE_FPARSER 530 auto fp = std::make_unique<FunctionParserBase<Output>>();
531 fp->AddConstant(
"NaN", std::numeric_limits<Real>::quiet_NaN());
532 fp->AddConstant(
"pi", std::acos(
Real(-1)));
533 fp->AddConstant(
"e", std::exp(
Real(1)));
535 (fp->Parse(new_subexpression, variables) != -1,
536 "ERROR: FunctionParser is unable to parse modified expression: " 537 << new_subexpression <<
'\n' << fp->ErrorMsg());
539 Output new_var_value = this->eval(*fp, new_subexpression, 0);
541 return new_var_value;
545 libmesh_assert_equal_to(old_var_value, new_var_value);
549 old_var_value = new_var_value;
550 found_var_name =
true;
555 libmesh_error_msg(
"ERROR: This functionality requires fparser!");
560 return old_var_value;
563 template <
typename Output>
569 libmesh_assert_greater (_subexpressions.size(), 0);
572 bool found_var_name =
false;
576 const std::string & subexpression = _subexpressions[s];
577 const std::size_t varname_i =
578 find_name(inline_var_name, subexpression);
579 if (varname_i == std::string::npos)
583 found_var_name =
true;
586 const std::size_t assignment_i =
587 subexpression.find(
":", varname_i+1);
589 libmesh_assert_not_equal_to(assignment_i, std::string::npos);
591 libmesh_assert_equal_to(subexpression[assignment_i+1],
'=');
592 for (std::size_t i = varname_i+1; i != assignment_i; ++i)
593 libmesh_assert_equal_to(subexpression[i],
' ');
595 std::size_t end_assignment_i =
596 subexpression.find(
";", assignment_i+1);
598 libmesh_assert_not_equal_to(end_assignment_i, std::string::npos);
600 std::ostringstream new_subexpression;
601 new_subexpression << subexpression.substr(0, assignment_i+2)
602 << std::setprecision(std::numeric_limits<Output>::digits10+2)
603 #ifdef LIBMESH_USE_COMPLEX_NUMBERS 604 <<
'(' << newval.real() <<
'+' 605 << newval.imag() <<
'i' <<
')' 609 << subexpression.substr(end_assignment_i,
611 _subexpressions[s] = new_subexpression.str();
616 std::string new_expression;
618 for (
const auto & subexpression : _subexpressions)
620 new_expression +=
'{';
621 new_expression += subexpression;
622 new_expression +=
'}';
625 this->partial_reparse(new_expression);
629 template <
typename Output>
634 _expression = std::move(expression);
635 _subexpressions.clear();
638 size_t nextstart = 0, end = 0;
640 while (end != std::string::npos)
644 if (nextstart >= _expression.size())
649 if (_expression[nextstart] ==
'{')
652 end = _expression.find(
'}', nextstart);
656 end = std::string::npos;
660 _subexpressions.push_back
661 (_expression.substr(nextstart, (end == std::string::npos) ?
662 std::string::npos : end - nextstart));
665 libmesh_error_msg_if(_subexpressions.back().empty(),
666 "ERROR: FunctionParser is unable to parse empty expression.\n");
669 #ifdef LIBMESH_HAVE_FPARSER 672 auto fp = std::make_unique<FunctionParserBase<Output>>();
673 fp->AddConstant(
"NaN", std::numeric_limits<Real>::quiet_NaN());
674 fp->AddConstant(
"pi", std::acos(
Real(-1)));
675 fp->AddConstant(
"e", std::exp(
Real(1)));
677 (fp->Parse(_subexpressions.back(), variables) != -1,
678 "ERROR: FunctionParser is unable to parse expression: " 679 << _subexpressions.back() <<
'\n' << fp->ErrorMsg());
681 parsers.push_back(std::move(fp));
683 libmesh_error_msg(
"ERROR: This functionality requires fparser!");
688 nextstart = (end == std::string::npos) ?
689 std::string::npos : end + 1;
693 template <
typename Output>
697 std::string_view expr)
const 699 const std::size_t namesize = varname.size();
700 std::size_t varname_i = expr.find(varname);
702 while ((varname_i != std::string::npos) &&
704 (std::isalnum(expr[varname_i-1]) ||
705 (expr[varname_i-1] ==
'_'))) ||
706 ((varname_i+namesize < expr.size()) &&
707 (std::isalnum(expr[varname_i+namesize]) ||
708 (expr[varname_i+namesize] ==
'_')))))
710 varname_i = expr.find(varname, varname_i+1);
719 template <
typename Output>
726 _spacetime[0] = p(0);
728 _spacetime[1] = p(1);
731 _spacetime[2] = p(2);
733 _spacetime[LIBMESH_DIM] = time;
735 unsigned int request_index = 0;
736 for (
unsigned int v=0; v != _n_vars; ++v)
741 c.
point_value(v, p, _spacetime[LIBMESH_DIM+1+request_index]);
745 if (_n_requested_grad_components)
746 for (
unsigned int v=0; v != _n_vars; ++v)
748 if (!_need_var_grad[v*LIBMESH_DIM]
750 && !_need_var_grad[v*LIBMESH_DIM+1]
752 && !_need_var_grad[v*LIBMESH_DIM+2]
761 for (
unsigned int d=0; d != LIBMESH_DIM; ++d)
763 if (!_need_var_grad[v*LIBMESH_DIM+d])
766 _spacetime[LIBMESH_DIM+1+request_index] = g(d);
771 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 772 if (_n_requested_hess_components)
773 for (
unsigned int v=0; v != _n_vars; ++v)
775 if (!_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM]
777 && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+1]
778 && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+2]
779 && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+3]
781 && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+4]
782 && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+5]
783 && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+6]
784 && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+7]
785 && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+8]
794 for (
unsigned int d1=0; d1 != LIBMESH_DIM; ++d1)
795 for (
unsigned int d2=0; d2 != LIBMESH_DIM; ++d2)
797 if (!_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+d1*LIBMESH_DIM+d2])
800 _spacetime[LIBMESH_DIM+1+request_index] = h(d1,d2);
804 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 806 if (_requested_normals)
811 const std::vector<Point> & normals = side_fe->
get_normals();
813 const std::vector<Point> & xyz = side_fe->
get_xyz();
815 libmesh_assert_equal_to(normals.size(), xyz.size());
819 bool at_quadrature_point =
false;
825 const Point & n = normals[qp];
826 for (
unsigned int d=0; d != LIBMESH_DIM; ++d)
828 _spacetime[LIBMESH_DIM+1+request_index] = n(d);
832 at_quadrature_point =
true;
847 #ifdef LIBMESH_HAVE_FPARSER 848 template <
typename Output>
852 std::string_view libmesh_dbg_var(function_name),
853 unsigned int libmesh_dbg_var(component_idx))
const 856 Output result = parser.Eval(_spacetime.data());
857 int error_code = parser.EvalError();
860 libMesh::err <<
"ERROR: FunctionParser is unable to evaluate component " 866 if (parsers[i].
get() == &parser)
868 << _subexpressions[i];
871 for (
const auto & st : _spacetime)
876 std::string error_message =
"Reason: ";
881 error_message +=
"Division by zero";
884 error_message +=
"Square Root error (negative value)";
887 error_message +=
"Log error (negative value)";
890 error_message +=
"Trigonometric error (asin or acos of illegal value)";
893 error_message +=
"Maximum recursion level reached";
896 error_message +=
"Unknown";
899 libmesh_error_msg(error_message);
904 return parser.Eval(_spacetime.data());
907 #else // LIBMESH_HAVE_FPARSER 908 template <
typename Output>
915 libmesh_error_msg(
"ERROR: This functionality requires fparser!");
923 #endif // LIBMESH_PARSED_FEM_FUNCTION_H
unsigned int _n_requested_vars
std::vector< bool > _need_var_grad
virtual void init_context(const FEMContext &c) override
Prepares a context object for use.
std::vector< Output > _initial_vals
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...
std::vector< std::string > _additional_vars
std::size_t find_name(std::string_view varname, std::string_view expr) const
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
std::vector< std::string > _subexpressions
void partial_reparse(std::string expression)
std::vector< Output > _spacetime
The libMesh namespace provides an interface to certain functionality in the library.
virtual std::unique_ptr< FEMFunctionBase< Output > > clone() const override
Number point_value(unsigned int var, const Point &p) 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
Output get_inline_value(std::string_view inline_var_name) const
ParsedFEMFunction provides support for FParser-based parsed functions in FEMSystem.
Manages consistently variables, degrees of freedom, and coefficient vectors.
void reparse(std::string expression)
Re-parse with new expression.
virtual ~ParsedFEMFunction()=default
std::vector< std::unique_ptr< FunctionParserBase< Output > > > parsers
Tensor point_hessian(unsigned int var, const Point &p) const
This class provides all data required for a physics package (e.g.
virtual_for_inffe const std::vector< Point > & get_xyz() const
const std::vector< std::vector< OutputTensor > > & get_d2phi() const
ParsedFEMFunction(const System &sys, std::string expression, const std::vector< std::string > *additional_vars=nullptr, const std::vector< Output > *initial_vals=nullptr)
Constructor.
ParsedFEMFunction & operator=(const ParsedFEMFunction &)
Constructors.
std::vector< bool > _need_var_hess
void set_inline_value(std::string_view inline_var_name, Output newval)
Changes the value of an inline variable.
void eval_args(const FEMContext &c, const Point &p, const Real time)
unsigned int _n_requested_hess_components
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::vector< bool > _need_var
virtual_for_inffe const std::vector< Point > & get_normals() const
virtual unsigned int size() const override final
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
Accessor for interior finite element object for variable var for the largest dimension in the mesh...
std::vector< char * > parsers
Gradient point_gradient(unsigned int var, const Point &p) const
FEMFunctionBase is a base class from which users can derive in order to define "function-like" object...
A Point defines a location in LIBMESH_DIM dimensional Real space.
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
This class forms the foundation from which generic finite elements may be derived.
const std::string & expression()
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
Output eval(FunctionParserBase< Output > &parser, std::string_view libmesh_dbg_var(function_name), unsigned int libmesh_dbg_var(component_idx)) const
const std::vector< std::vector< OutputShape > > & get_phi() const