libMesh
Public Member Functions | Protected Member Functions | Private Attributes | List of all members
libMesh::ParsedFEMFunction< Output > Class Template Reference

ParsedFEMFunction provides support for FParser-based parsed functions in FEMSystem. More...

#include <parsed_fem_function.h>

Inheritance diagram for libMesh::ParsedFEMFunction< Output >:
[legend]

Public Member Functions

 ParsedFEMFunction (const System &sys, std::string expression, const std::vector< std::string > *additional_vars=nullptr, const std::vector< Output > *initial_vals=nullptr)
 Constructor. More...
 
ParsedFEMFunctionoperator= (const ParsedFEMFunction &)
 Constructors. More...
 
ParsedFEMFunctionoperator= (ParsedFEMFunction &&)=delete
 
 ParsedFEMFunction (const ParsedFEMFunction &)
 
 ParsedFEMFunction (ParsedFEMFunction &&)=default
 
virtual ~ParsedFEMFunction ()=default
 
void reparse (std::string expression)
 Re-parse with new expression. More...
 
virtual void init_context (const FEMContext &c) override
 Prepares a context object for use. More...
 
virtual std::unique_ptr< FEMFunctionBase< Output > > clone () const override
 
virtual Output operator() (const FEMContext &c, const Point &p, const Real time=0.) override
 
void operator() (const FEMContext &c, const Point &p, const Real time, DenseVector< Output > &output) override
 Evaluation function for time-dependent vector-valued functions. More...
 
virtual Output component (const FEMContext &c, unsigned int i, const Point &p, Real time=0.) override
 
const std::string & expression ()
 
Output get_inline_value (std::string_view inline_var_name) const
 
void set_inline_value (std::string_view inline_var_name, Output newval)
 Changes the value of an inline variable. More...
 
virtual void init ()
 Any post-construction initialization. More...
 
void operator() (const FEMContext &, const Point &p, DenseVector< Output > &output)
 Evaluation function for time-independent vector-valued functions. More...
 

Protected Member Functions

void partial_reparse (std::string expression)
 
std::size_t find_name (std::string_view varname, std::string_view expr) const
 
void eval_args (const FEMContext &c, const Point &p, const Real time)
 
Output eval (FunctionParserBase< Output > &parser, std::string_view libmesh_dbg_var(function_name), unsigned int libmesh_dbg_var(component_idx)) const
 
Output eval (char &libmesh_dbg_var(parser), std::string_view libmesh_dbg_var(function_name), unsigned int libmesh_dbg_var(component_idx)) const
 

Private Attributes

const System_sys
 
std::string _expression
 
std::vector< std::string > _subexpressions
 
unsigned int _n_vars
 
unsigned int _n_requested_vars
 
unsigned int _n_requested_grad_components
 
unsigned int _n_requested_hess_components
 
bool _requested_normals
 
std::vector< std::unique_ptr< FunctionParserBase< Output > > > parsers
 
std::vector< char * > parsers
 
std::vector< Output > _spacetime
 
std::vector< bool > _need_var
 
std::vector< bool > _need_var_grad
 
std::vector< bool > _need_var_hess
 
std::string variables
 
std::vector< std::string > _additional_vars
 
std::vector< Output > _initial_vals
 

Detailed Description

template<typename Output = Number>
class libMesh::ParsedFEMFunction< Output >

ParsedFEMFunction provides support for FParser-based parsed functions in FEMSystem.

All overridden virtual functions are documented in fem_function_base.h.

Author
Roy Stogner
Date
2014 Support for using parsed functions in FEMSystem.

Definition at line 57 of file parsed_fem_function.h.

Constructor & Destructor Documentation

◆ ParsedFEMFunction() [1/3]

template<typename Output>
libMesh::ParsedFEMFunction< Output >::ParsedFEMFunction ( const System sys,
std::string  expression,
const std::vector< std::string > *  additional_vars = nullptr,
const std::vector< Output > *  initial_vals = nullptr 
)
inlineexplicit

Constructor.

Definition at line 197 of file parsed_fem_function.h.

200  :
201  _sys(sys),
202  _expression (), // overridden by parse()
203  _n_vars(sys.n_vars()),
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),
212 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
213  _additional_vars (additional_vars ? *additional_vars : std::vector<std::string>()),
214  _initial_vals (initial_vals ? *initial_vals : std::vector<Output>())
215 {
216  this->reparse(std::move(expression));
217 }
std::vector< bool > _need_var_grad
std::vector< Output > _initial_vals
std::vector< std::string > _additional_vars
void reparse(std::string expression)
Re-parse with new expression.
std::vector< bool > _need_var_hess
const std::string & expression()

◆ ParsedFEMFunction() [2/3]

template<typename Output>
libMesh::ParsedFEMFunction< Output >::ParsedFEMFunction ( const ParsedFEMFunction< Output > &  other)
inline

Definition at line 222 of file parsed_fem_function.h.

222  :
223  FEMFunctionBase<Output>(other),
224  _sys(other._sys),
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),
237 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
238  variables(other.variables),
239  _additional_vars(other._additional_vars),
240  _initial_vals(other._initial_vals)
241 {
242  this->reparse(_expression);
243 }
std::vector< bool > _need_var_grad
std::vector< Output > _initial_vals
std::vector< std::string > _additional_vars
std::vector< std::string > _subexpressions
std::vector< Output > _spacetime
void reparse(std::string expression)
Re-parse with new expression.
std::vector< bool > _need_var_hess

◆ ParsedFEMFunction() [3/3]

template<typename Output = Number>
libMesh::ParsedFEMFunction< Output >::ParsedFEMFunction ( ParsedFEMFunction< Output > &&  )
default

◆ ~ParsedFEMFunction()

template<typename Output = Number>
virtual libMesh::ParsedFEMFunction< Output >::~ParsedFEMFunction ( )
virtualdefault

Member Function Documentation

◆ clone()

template<typename Output >
std::unique_ptr< FEMFunctionBase< Output > > libMesh::ParsedFEMFunction< Output >::clone ( ) const
inlineoverridevirtual
Returns
A new copy of the function.

The new copy should be as "deep" as necessary to allow independent destruction and simultaneous evaluations of the copies in different threads.

Implements libMesh::FEMFunctionBase< Output >.

Definition at line 437 of file parsed_fem_function.h.

438 {
439  return std::make_unique<ParsedFEMFunction>
441 }
std::vector< Output > _initial_vals
std::vector< std::string > _additional_vars

◆ component()

template<typename Output >
Output libMesh::ParsedFEMFunction< Output >::component ( const FEMContext context,
unsigned int  i,
const Point p,
Real  time = 0. 
)
inlineoverridevirtual
Returns
The vector component i at coordinate p and time time.
Note
Subclasses aren't required to override this, since the default implementation is based on the full vector evaluation, which is often correct.
Subclasses are recommended to override this, since the default implementation is based on a vector evaluation, which is usually unnecessarily inefficient.

Reimplemented from libMesh::FEMFunctionBase< Output >.

Definition at line 479 of file parsed_fem_function.h.

483 {
484  eval_args(c, p, time);
485 
486  libmesh_assert_less (i, parsers.size());
487  return eval(*parsers[i], "f", i);
488 }
std::vector< std::unique_ptr< FunctionParserBase< Output > > > parsers
void eval_args(const FEMContext &c, const Point &p, const Real time)
Output eval(FunctionParserBase< Output > &parser, std::string_view libmesh_dbg_var(function_name), unsigned int libmesh_dbg_var(component_idx)) const

◆ eval() [1/2]

template<typename Output>
Output libMesh::ParsedFEMFunction< Output >::eval ( FunctionParserBase< Output > &  parser,
std::string_view   libmesh_dbg_varfunction_name,
unsigned int   libmesh_dbg_varcomponent_idx 
) const
inlineprotected

Definition at line 851 of file parsed_fem_function.h.

854 {
855 #ifndef NDEBUG
856  Output result = parser.Eval(_spacetime.data());
857  int error_code = parser.EvalError();
858  if (error_code)
859  {
860  libMesh::err << "ERROR: FunctionParser is unable to evaluate component "
861  << component_idx
862  << " for '"
863  << function_name;
864 
865  for (auto i : index_range(parsers))
866  if (parsers[i].get() == &parser)
867  libMesh::err << "' of expression '"
868  << _subexpressions[i];
869 
870  libMesh::err << "' with arguments:\n";
871  for (const auto & st : _spacetime)
872  libMesh::err << '\t' << st << '\n';
873  libMesh::err << '\n';
874 
875  // Currently no API to report error messages, we'll do it manually
876  std::string error_message = "Reason: ";
877 
878  switch (error_code)
879  {
880  case 1:
881  error_message += "Division by zero";
882  break;
883  case 2:
884  error_message += "Square Root error (negative value)";
885  break;
886  case 3:
887  error_message += "Log error (negative value)";
888  break;
889  case 4:
890  error_message += "Trigonometric error (asin or acos of illegal value)";
891  break;
892  case 5:
893  error_message += "Maximum recursion level reached";
894  break;
895  default:
896  error_message += "Unknown";
897  break;
898  }
899  libmesh_error_msg(error_message);
900  }
901 
902  return result;
903 #else
904  return parser.Eval(_spacetime.data());
905 #endif
906 }
OStreamProxy err
std::vector< std::string > _subexpressions
std::vector< Output > _spacetime
std::vector< std::unique_ptr< FunctionParserBase< Output > > > parsers
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:117

◆ eval() [2/2]

template<typename Output = Number>
Output libMesh::ParsedFEMFunction< Output >::eval ( char &  libmesh_dbg_varparser,
std::string_view   libmesh_dbg_varfunction_name,
unsigned int   libmesh_dbg_varcomponent_idx 
) const
inlineprotected

◆ eval_args()

template<typename Output >
void libMesh::ParsedFEMFunction< Output >::eval_args ( const FEMContext c,
const Point p,
const Real  time 
)
inlineprotected

Definition at line 722 of file parsed_fem_function.h.

725 {
726  _spacetime[0] = p(0);
727 #if LIBMESH_DIM > 1
728  _spacetime[1] = p(1);
729 #endif
730 #if LIBMESH_DIM > 2
731  _spacetime[2] = p(2);
732 #endif
733  _spacetime[LIBMESH_DIM] = time;
734 
735  unsigned int request_index = 0;
736  for (unsigned int v=0; v != _n_vars; ++v)
737  {
738  if (!_need_var[v])
739  continue;
740 
741  c.point_value(v, p, _spacetime[LIBMESH_DIM+1+request_index]);
742  request_index++;
743  }
744 
746  for (unsigned int v=0; v != _n_vars; ++v)
747  {
748  if (!_need_var_grad[v*LIBMESH_DIM]
749 #if LIBMESH_DIM > 1
750  && !_need_var_grad[v*LIBMESH_DIM+1]
751 #if LIBMESH_DIM > 2
752  && !_need_var_grad[v*LIBMESH_DIM+2]
753 #endif
754 #endif
755  )
756  continue;
757 
758  Gradient g;
759  c.point_gradient(v, p, g);
760 
761  for (unsigned int d=0; d != LIBMESH_DIM; ++d)
762  {
763  if (!_need_var_grad[v*LIBMESH_DIM+d])
764  continue;
765 
766  _spacetime[LIBMESH_DIM+1+request_index] = g(d);
767  request_index++;
768  }
769  }
770 
771 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
773  for (unsigned int v=0; v != _n_vars; ++v)
774  {
775  if (!_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM]
776 #if LIBMESH_DIM > 1
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]
780 #if LIBMESH_DIM > 2
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]
786 #endif
787 #endif
788  )
789  continue;
790 
791  Tensor h;
792  c.point_hessian(v, p, h);
793 
794  for (unsigned int d1=0; d1 != LIBMESH_DIM; ++d1)
795  for (unsigned int d2=0; d2 != LIBMESH_DIM; ++d2)
796  {
797  if (!_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+d1*LIBMESH_DIM+d2])
798  continue;
799 
800  _spacetime[LIBMESH_DIM+1+request_index] = h(d1,d2);
801  request_index++;
802  }
803  }
804 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
805 
806  if (_requested_normals)
807  {
808  FEBase * side_fe;
809  c.get_side_fe(0, side_fe);
810 
811  const std::vector<Point> & normals = side_fe->get_normals();
812 
813  const std::vector<Point> & xyz = side_fe->get_xyz();
814 
815  libmesh_assert_equal_to(normals.size(), xyz.size());
816 
817  // We currently only support normals at quadrature points!
818 #ifndef NDEBUG
819  bool at_quadrature_point = false;
820 #endif
821  for (auto qp : index_range(normals))
822  {
823  if (p == xyz[qp])
824  {
825  const Point & n = normals[qp];
826  for (unsigned int d=0; d != LIBMESH_DIM; ++d)
827  {
828  _spacetime[LIBMESH_DIM+1+request_index] = n(d);
829  request_index++;
830  }
831 #ifndef NDEBUG
832  at_quadrature_point = true;
833 #endif
834  break;
835  }
836  }
837 
838  libmesh_assert(at_quadrature_point);
839  }
840 
841  // The remaining locations in _spacetime are currently set at construction
842  // but could potentially be made dynamic
843 }
std::vector< bool > _need_var_grad
std::vector< Output > _spacetime
NumberVectorValue Gradient
libmesh_assert(ctx)
FEGenericBase< Real > FEBase
std::vector< bool > _need_var_hess
NumberTensorValue Tensor
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:117

◆ expression()

template<typename Output = Number>
const std::string& libMesh::ParsedFEMFunction< Output >::expression ( )
inline

Definition at line 107 of file parsed_fem_function.h.

107 { return _expression; }

◆ find_name()

template<typename Output >
std::size_t libMesh::ParsedFEMFunction< Output >::find_name ( std::string_view  varname,
std::string_view  expr 
) const
inlineprotected

Definition at line 696 of file parsed_fem_function.h.

698 {
699  const std::size_t namesize = varname.size();
700  std::size_t varname_i = expr.find(varname);
701 
702  while ((varname_i != std::string::npos) &&
703  (((varname_i > 0) &&
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] == '_')))))
709  {
710  varname_i = expr.find(varname, varname_i+1);
711  }
712 
713  return varname_i;
714 }

◆ get_inline_value()

template<typename Output >
Output libMesh::ParsedFEMFunction< Output >::get_inline_value ( std::string_view  inline_var_name) const
inline
Returns
The value of an inline variable.
Note
Will only be correct if the inline variable value is independent of input variables, if the inline variable is not redefined within any subexpression, and if the inline variable takes the same value within any subexpressions where it appears.

Definition at line 493 of file parsed_fem_function.h.

Referenced by libMesh::ParsedFEMFunctionParameter< T >::get(), ParsedFEMFunctionTest::testInlineGetter(), and ParsedFEMFunctionTest::testInlineSetter().

494 {
495  libmesh_assert_greater (_subexpressions.size(), 0);
496 
497 #ifndef NDEBUG
498  bool found_var_name = false;
499 #endif
500  Output old_var_value(0.);
501 
502  for (const std::string & subexpression : _subexpressions)
503  {
504  const std::size_t varname_i =
505  find_name(inline_var_name, subexpression);
506  if (varname_i == std::string::npos)
507  continue;
508 
509  const std::size_t assignment_i =
510  subexpression.find(":", varname_i+1);
511 
512  libmesh_assert_not_equal_to(assignment_i, std::string::npos);
513 
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], ' ');
517 
518  std::size_t end_assignment_i =
519  subexpression.find(";", assignment_i+1);
520 
521  libmesh_assert_not_equal_to(end_assignment_i, std::string::npos);
522 
523  std::string new_subexpression =
524  subexpression.substr(0, end_assignment_i+1).
525  append(inline_var_name);
526 
527 #ifdef LIBMESH_HAVE_FPARSER
528  // Parse and evaluate the new subexpression.
529  // Add the same constants as we used originally.
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)));
534  libmesh_error_msg_if
535  (fp->Parse(new_subexpression, variables) != -1, // -1 for success
536  "ERROR: FunctionParser is unable to parse modified expression: "
537  << new_subexpression << '\n' << fp->ErrorMsg());
538 
539  Output new_var_value = this->eval(*fp, new_subexpression, 0);
540 #ifdef NDEBUG
541  return new_var_value;
542 #else
543  if (found_var_name)
544  {
545  libmesh_assert_equal_to(old_var_value, new_var_value);
546  }
547  else
548  {
549  old_var_value = new_var_value;
550  found_var_name = true;
551  }
552 #endif
553 
554 #else
555  libmesh_error_msg("ERROR: This functionality requires fparser!");
556 #endif
557  }
558 
559  libmesh_assert(found_var_name);
560  return old_var_value;
561 }
std::size_t find_name(std::string_view varname, std::string_view expr) const
std::vector< std::string > _subexpressions
libmesh_assert(ctx)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Output eval(FunctionParserBase< Output > &parser, std::string_view libmesh_dbg_var(function_name), unsigned int libmesh_dbg_var(component_idx)) const

◆ init()

template<typename Output = Number>
virtual void libMesh::FEMFunctionBase< Output >::init ( )
inlinevirtualinherited

Any post-construction initialization.

Reimplemented in libMesh::WrappedFunctor< Output >.

Definition at line 69 of file fem_function_base.h.

69 {}

◆ init_context()

template<typename Output >
void libMesh::ParsedFEMFunction< Output >::init_context ( const FEMContext )
inlineoverridevirtual

Prepares a context object for use.

Most problems will want to reimplement this for efficiency, in order to call FE::get_*() as their particular function requires.

Reimplemented from libMesh::FEMFunctionBase< Output >.

Definition at line 393 of file parsed_fem_function.h.

394 {
395  for (unsigned int v=0; v != _n_vars; ++v)
396  {
397  FEBase * elem_fe;
398  c.get_element_fe(v, elem_fe);
399  bool request_nothing = true;
400  if (_n_requested_vars)
401  {
402  elem_fe->get_phi();
403  request_nothing = false;
404  }
406  {
407  elem_fe->get_dphi();
408  request_nothing = false;
409  }
410 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
412  {
413  elem_fe->get_d2phi();
414  request_nothing = false;
415  }
416 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
417  if (request_nothing)
418  elem_fe->get_nothing();
419  }
420 
421  if (_requested_normals)
422  {
423  FEBase * side_fe;
424  c.get_side_fe(0, side_fe);
425 
426  side_fe->get_normals();
427 
428  // FIXME: this is a hack to support normals at quadrature
429  // points; we don't support normals elsewhere.
430  side_fe->get_xyz();
431  }
432 }
FEGenericBase< Real > FEBase

◆ operator()() [1/3]

template<typename Output >
Output libMesh::ParsedFEMFunction< Output >::operator() ( const FEMContext ,
const Point p,
const Real  time = 0. 
)
inlineoverridevirtual
Returns
The scalar function value at coordinate p and time time, which defaults to zero.

Pure virtual, so you have to override it.

Implements libMesh::FEMFunctionBase< Output >.

Definition at line 446 of file parsed_fem_function.h.

449 {
450  eval_args(c, p, time);
451 
452  return eval(*parsers[0], "f", 0);
453 }
std::vector< std::unique_ptr< FunctionParserBase< Output > > > parsers
void eval_args(const FEMContext &c, const Point &p, const Real time)
Output eval(FunctionParserBase< Output > &parser, std::string_view libmesh_dbg_var(function_name), unsigned int libmesh_dbg_var(component_idx)) const

◆ operator()() [2/3]

template<typename Output>
void libMesh::ParsedFEMFunction< Output >::operator() ( const FEMContext ,
const Point p,
const Real  time,
DenseVector< Output > &  output 
)
inlineoverridevirtual

Evaluation function for time-dependent vector-valued functions.

Sets output values in the passed-in output DenseVector.

Pure virtual, so you have to override it.

Implements libMesh::FEMFunctionBase< Output >.

Definition at line 460 of file parsed_fem_function.h.

464 {
465  eval_args(c, p, time);
466 
467  unsigned int size = output.size();
468 
469  libmesh_assert_equal_to (size, parsers.size());
470 
471  for (unsigned int i=0; i != size; ++i)
472  output(i) = eval(*parsers[i], "f", i);
473 }
std::vector< std::unique_ptr< FunctionParserBase< Output > > > parsers
void eval_args(const FEMContext &c, const Point &p, const Real time)
virtual unsigned int size() const override final
Definition: dense_vector.h:104
Output eval(FunctionParserBase< Output > &parser, std::string_view libmesh_dbg_var(function_name), unsigned int libmesh_dbg_var(component_idx)) const

◆ operator()() [3/3]

template<typename Output >
void libMesh::FEMFunctionBase< Output >::operator() ( const FEMContext context,
const Point p,
DenseVector< Output > &  output 
)
inlineinherited

Evaluation function for time-independent vector-valued functions.

Sets output values in the passed-in output DenseVector.

Definition at line 149 of file fem_function_base.h.

152 {
153  // Call the time-dependent function with t=0.
154  this->operator()(context, p, 0., output);
155 }
virtual Output operator()(const FEMContext &, const Point &p, const Real time=0.)=0

◆ operator=() [1/2]

template<typename Output >
ParsedFEMFunction< Output > & libMesh::ParsedFEMFunction< Output >::operator= ( const ParsedFEMFunction< Output > &  other)
inline

Constructors.

  • This class contains a const reference so it can't be default copy or move-assigned. We manually implement the former
  • This class contains unique_ptrs so it can't be default copy constructed or assigned.
  • It can still be default moved and deleted.

Definition at line 249 of file parsed_fem_function.h.

250 {
251  // We can only be assigned another ParsedFEMFunction defined on the same System
252  libmesh_assert(&_sys == &other._sys);
253 
254  // Use copy-and-swap idiom
255  ParsedFEMFunction<Output> tmp(other);
256  std::swap(tmp, *this);
257  return *this;
258 }
libmesh_assert(ctx)

◆ operator=() [2/2]

template<typename Output = Number>
ParsedFEMFunction& libMesh::ParsedFEMFunction< Output >::operator= ( ParsedFEMFunction< Output > &&  )
delete

◆ partial_reparse()

template<typename Output >
void libMesh::ParsedFEMFunction< Output >::partial_reparse ( std::string  expression)
inlineprotected

Definition at line 632 of file parsed_fem_function.h.

633 {
634  _expression = std::move(expression);
635  _subexpressions.clear();
636  parsers.clear();
637 
638  size_t nextstart = 0, end = 0;
639 
640  while (end != std::string::npos)
641  {
642  // If we're past the end of the string, we can't make any more
643  // subparsers
644  if (nextstart >= _expression.size())
645  break;
646 
647  // If we're at the start of a brace delimited section, then we
648  // parse just that section:
649  if (_expression[nextstart] == '{')
650  {
651  nextstart++;
652  end = _expression.find('}', nextstart);
653  }
654  // otherwise we parse the whole thing
655  else
656  end = std::string::npos;
657 
658  // We either want the whole end of the string (end == npos) or
659  // a substring in the middle.
660  _subexpressions.push_back
661  (_expression.substr(nextstart, (end == std::string::npos) ?
662  std::string::npos : end - nextstart));
663 
664  // fparser can crash on empty expressions
665  libmesh_error_msg_if(_subexpressions.back().empty(),
666  "ERROR: FunctionParser is unable to parse empty expression.\n");
667 
668 
669 #ifdef LIBMESH_HAVE_FPARSER
670  // Parse (and optimize if possible) the subexpression.
671  // Add some basic constants, to Real precision.
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)));
676  libmesh_error_msg_if
677  (fp->Parse(_subexpressions.back(), variables) != -1, // -1 for success
678  "ERROR: FunctionParser is unable to parse expression: "
679  << _subexpressions.back() << '\n' << fp->ErrorMsg());
680  fp->Optimize();
681  parsers.push_back(std::move(fp));
682 #else
683  libmesh_error_msg("ERROR: This functionality requires fparser!");
684 #endif
685 
686  // If at end, use nextstart=maxSize. Else start at next
687  // character.
688  nextstart = (end == std::string::npos) ?
689  std::string::npos : end + 1;
690  }
691 }
std::vector< std::string > _subexpressions
std::vector< std::unique_ptr< FunctionParserBase< Output > > > parsers
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const std::string & expression()

◆ reparse()

template<typename Output >
void libMesh::ParsedFEMFunction< Output >::reparse ( std::string  expression)
inline

Re-parse with new expression.

Definition at line 264 of file parsed_fem_function.h.

265 {
266  variables = "x";
267 #if LIBMESH_DIM > 1
268  variables += ",y";
269 #endif
270 #if LIBMESH_DIM > 2
271  variables += ",z";
272 #endif
273  variables += ",t";
274 
275  for (unsigned int v=0; v != _n_vars; ++v)
276  {
277  const std::string & varname = _sys.variable_name(v);
278  std::size_t varname_i = find_name(varname, expression);
279 
280  // If we didn't find our variable name then let's go to the
281  // next.
282  if (varname_i == std::string::npos)
283  continue;
284 
285  _need_var[v] = true;
286  variables += ',';
287  variables += varname;
289  }
290 
291  for (unsigned int v=0; v != _n_vars; ++v)
292  {
293  const std::string & varname = _sys.variable_name(v);
294 
295  for (unsigned int d=0; d != LIBMESH_DIM; ++d)
296  {
297  std::string gradname = std::string("grad_") +
298  "xyz"[d] + '_' + varname;
299  std::size_t gradname_i = find_name(gradname, expression);
300 
301  // If we didn't find that gradient component of our
302  // variable name then let's go to the next.
303  if (gradname_i == std::string::npos)
304  continue;
305 
306  _need_var_grad[v*LIBMESH_DIM+d] = true;
307  variables += ',';
308  variables += gradname;
310  }
311  }
312 
313 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
314  for (unsigned int v=0; v != _n_vars; ++v)
315  {
316  const std::string & varname = _sys.variable_name(v);
317 
318  for (unsigned int d1=0; d1 != LIBMESH_DIM; ++d1)
319  for (unsigned int d2=0; d2 != LIBMESH_DIM; ++d2)
320  {
321  std::string hessname = std::string("hess_") +
322  "xyz"[d1] + "xyz"[d2] + '_' + varname;
323  std::size_t hessname_i = find_name(hessname, expression);
324 
325  // If we didn't find that hessian component of our
326  // variable name then let's go to the next.
327  if (hessname_i == std::string::npos)
328  continue;
329 
330  _need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+d1*LIBMESH_DIM+d2]
331  = true;
332  variables += ',';
333  variables += hessname;
335  }
336  }
337 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
338 
339  {
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);
343 
344  // If we found any requests for normal components then we'll
345  // compute normals
346  if (nx_i != std::string::npos ||
347  ny_i != std::string::npos ||
348  nz_i != std::string::npos)
349  {
350  _requested_normals = true;
351  variables += ',';
352  variables += "n_x";
353  if (LIBMESH_DIM > 1)
354  {
355  variables += ',';
356  variables += "n_y";
357  }
358  if (LIBMESH_DIM > 2)
359  {
360  variables += ',';
361  variables += "n_z";
362  }
363  }
364  }
365 
366  _spacetime.resize
367  (LIBMESH_DIM + 1 + _n_requested_vars +
369  (_requested_normals ? LIBMESH_DIM : 0) +
370  _additional_vars.size());
371 
372  // If additional vars were passed, append them to the string
373  // that we send to the function parser. Also add them to the
374  // end of our spacetime vector
375  unsigned int offset = LIBMESH_DIM + 1 + _n_requested_vars +
377 
378  for (auto i : index_range(_additional_vars))
379  {
380  variables += "," + _additional_vars[i];
381  // Initialize extra variables to the vector passed in or zero
382  // Note: The initial_vals vector can be shorter than the additional_vars vector
383  _spacetime[offset + i] =
384  (i < _initial_vals.size()) ? _initial_vals[i] : 0;
385  }
386 
387  this->partial_reparse(std::move(expression));
388 }
std::vector< bool > _need_var_grad
std::vector< Output > _initial_vals
std::vector< std::string > _additional_vars
std::size_t find_name(std::string_view varname, std::string_view expr) const
void partial_reparse(std::string expression)
std::vector< Output > _spacetime
const std::string & variable_name(const unsigned int i) const
Definition: system.h:2478
std::vector< bool > _need_var_hess
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:117
const std::string & expression()

◆ set_inline_value()

template<typename Output>
void libMesh::ParsedFEMFunction< Output >::set_inline_value ( std::string_view  inline_var_name,
Output  newval 
)
inline

Changes the value of an inline variable.

Note
Forever after, the variable value will take the given constant, independent of input variables, in every subexpression where it is already defined.
Currently only works if the inline variable is not redefined within any one subexpression.

Definition at line 566 of file parsed_fem_function.h.

Referenced by libMesh::ParsedFEMFunctionParameter< T >::set(), and ParsedFEMFunctionTest::testInlineSetter().

568 {
569  libmesh_assert_greater (_subexpressions.size(), 0);
570 
571 #ifndef NDEBUG
572  bool found_var_name = false;
573 #endif
574  for (auto s : index_range(_subexpressions))
575  {
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)
580  continue;
581 
582 #ifndef NDEBUG
583  found_var_name = true;
584 #endif
585 
586  const std::size_t assignment_i =
587  subexpression.find(":", varname_i+1);
588 
589  libmesh_assert_not_equal_to(assignment_i, std::string::npos);
590 
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], ' ');
594 
595  std::size_t end_assignment_i =
596  subexpression.find(";", assignment_i+1);
597 
598  libmesh_assert_not_equal_to(end_assignment_i, std::string::npos);
599 
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' << ')'
606 #else
607  << newval
608 #endif
609  << subexpression.substr(end_assignment_i,
610  std::string::npos);
611  _subexpressions[s] = new_subexpression.str();
612  }
613 
614  libmesh_assert(found_var_name);
615 
616  std::string new_expression;
617 
618  for (const auto & subexpression : _subexpressions)
619  {
620  new_expression += '{';
621  new_expression += subexpression;
622  new_expression += '}';
623  }
624 
625  this->partial_reparse(new_expression);
626 }
std::size_t find_name(std::string_view varname, std::string_view expr) const
std::vector< std::string > _subexpressions
void partial_reparse(std::string expression)
libmesh_assert(ctx)
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:117

Member Data Documentation

◆ _additional_vars

template<typename Output = Number>
std::vector<std::string> libMesh::ParsedFEMFunction< Output >::_additional_vars
private

Definition at line 188 of file parsed_fem_function.h.

◆ _expression

template<typename Output = Number>
std::string libMesh::ParsedFEMFunction< Output >::_expression
private

Definition at line 158 of file parsed_fem_function.h.

Referenced by libMesh::ParsedFEMFunction< T >::expression().

◆ _initial_vals

template<typename Output = Number>
std::vector<Output> libMesh::ParsedFEMFunction< Output >::_initial_vals
private

Definition at line 189 of file parsed_fem_function.h.

◆ _n_requested_grad_components

template<typename Output = Number>
unsigned int libMesh::ParsedFEMFunction< Output >::_n_requested_grad_components
private

Definition at line 160 of file parsed_fem_function.h.

◆ _n_requested_hess_components

template<typename Output = Number>
unsigned int libMesh::ParsedFEMFunction< Output >::_n_requested_hess_components
private

Definition at line 160 of file parsed_fem_function.h.

◆ _n_requested_vars

template<typename Output = Number>
unsigned int libMesh::ParsedFEMFunction< Output >::_n_requested_vars
private

Definition at line 160 of file parsed_fem_function.h.

◆ _n_vars

template<typename Output = Number>
unsigned int libMesh::ParsedFEMFunction< Output >::_n_vars
private

Definition at line 160 of file parsed_fem_function.h.

◆ _need_var

template<typename Output = Number>
std::vector<bool> libMesh::ParsedFEMFunction< Output >::_need_var
private

Definition at line 175 of file parsed_fem_function.h.

◆ _need_var_grad

template<typename Output = Number>
std::vector<bool> libMesh::ParsedFEMFunction< Output >::_need_var_grad
private

Definition at line 178 of file parsed_fem_function.h.

◆ _need_var_hess

template<typename Output = Number>
std::vector<bool> libMesh::ParsedFEMFunction< Output >::_need_var_hess
private

Definition at line 183 of file parsed_fem_function.h.

◆ _requested_normals

template<typename Output = Number>
bool libMesh::ParsedFEMFunction< Output >::_requested_normals
private

Definition at line 164 of file parsed_fem_function.h.

◆ _spacetime

template<typename Output = Number>
std::vector<Output> libMesh::ParsedFEMFunction< Output >::_spacetime
private

Definition at line 170 of file parsed_fem_function.h.

◆ _subexpressions

template<typename Output = Number>
std::vector<std::string> libMesh::ParsedFEMFunction< Output >::_subexpressions
private

Definition at line 159 of file parsed_fem_function.h.

◆ _sys

template<typename Output = Number>
const System& libMesh::ParsedFEMFunction< Output >::_sys
private

Definition at line 157 of file parsed_fem_function.h.

Referenced by libMesh::ParsedFEMFunction< T >::operator=().

◆ parsers [1/2]

template<typename Output = Number>
std::vector<std::unique_ptr<FunctionParserBase<Output> > > libMesh::ParsedFEMFunction< Output >::parsers
private

Definition at line 166 of file parsed_fem_function.h.

◆ parsers [2/2]

template<typename Output = Number>
std::vector<char*> libMesh::ParsedFEMFunction< Output >::parsers
private

Definition at line 168 of file parsed_fem_function.h.

◆ variables

template<typename Output = Number>
std::string libMesh::ParsedFEMFunction< Output >::variables
private

Definition at line 187 of file parsed_fem_function.h.


The documentation for this class was generated from the following file: