libmesh解析
本工作只是尝试解析原libmesh的代码,供学习使用
 全部  命名空间 文件 函数 变量 类型定义 枚举 枚举值 友元 
parsed_fem_function.h
浏览该文件的文档.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2023 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 
31 #ifdef LIBMESH_HAVE_FPARSER
32 // FParser includes
33 #include "libmesh/fparser.hh"
34 #endif
35 
36 // C++ includes
37 #include <cctype>
38 #include <iomanip>
39 #include <sstream>
40 #include <string>
41 #include <vector>
42 #include <memory>
43 
44 namespace libMesh
45 {
46 
58 template <typename Output=Number>
59 class ParsedFEMFunction : public FEMFunctionBase<Output>
60 {
61 public:
62 
71  explicit
72  ParsedFEMFunction (const System & sys,
73  std::string expression,
74  const std::vector<std::string> * additional_vars=nullptr,
75  const std::vector<Output> * initial_vals=nullptr);
76 
87 
96 
102  void reparse (std::string expression);
103 
109  virtual void init_context (const FEMContext & c) override;
110 
116  virtual std::unique_ptr<FEMFunctionBase<Output>> clone () const override;
117 
127  virtual Output operator() (const FEMContext & c,
128  const Point & p,
129  const Real time = 0.) override;
130 
139  void operator() (const FEMContext & c,
140  const Point & p,
141  const Real time,
142  DenseVector<Output> & output) override;
143 
155  virtual Output component(const FEMContext & c,
156  unsigned int i,
157  const Point & p,
158  Real time=0.) override;
159 
165  const std::string & expression() { return _expression; }
166 
177  Output get_inline_value(std::string_view inline_var_name) const;
178 
190  void set_inline_value(std::string_view inline_var_name,
191  Output newval);
192 
193 protected:
199  void partial_reparse (std::string expression);
200 
209  std::size_t find_name (std::string_view varname,
210  std::string_view expr) const;
211 
219  void eval_args(const FEMContext & c,
220  const Point & p,
221  const Real time);
222 
223 #ifdef LIBMESH_HAVE_FPARSER
224 
233  inline Output eval(FunctionParserBase<Output> & parser,
234  std::string_view libmesh_dbg_var(function_name),
235  unsigned int libmesh_dbg_var(component_idx)) const;
236 #else // LIBMESH_HAVE_FPARSER
237 
246  inline Output eval(char & libmesh_dbg_var(parser),
247  std::string_view libmesh_dbg_var(function_name),
248  unsigned int libmesh_dbg_var(component_idx)) const;
249 #endif
250 
251 private:
252  const System & _sys; // 系统对象的常量引用。
253  std::string _expression; // 解析表达式字符串。
254  std::vector<std::string> _subexpressions; // 子表达式的字符串向量。
255  unsigned int _n_vars, // 变量数量。
256  _n_requested_vars, // 请求的变量数量。
257  _n_requested_grad_components, // 请求的梯度分量数量。
258  _n_requested_hess_components; // 请求的 Hessian 矩阵分量数量。
259  bool _requested_normals; // 是否请求法向量。
260 #ifdef LIBMESH_HAVE_FPARSER
261  std::vector<std::unique_ptr<FunctionParserBase<Output>>> parsers; // FunctionParserBase 对象的唯一指针向量。
262 #else
263  std::vector<char*> parsers; // 字符指针向量。
264 #endif
265  std::vector<Output> _spacetime; // 时空信息向量。
266 
267  // 标志,指示需要计算哪些变量
268 
269  // _need_var[v] 为真,当且仅当需要计算 value(v)
270  std::vector<bool> _need_var;
271 
272  // _need_var_grad[v*LIBMESH_DIM+i] 为真,当且仅当需要计算 grad(v, i)
273  std::vector<bool> _need_var_grad;
274 
275 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
276  // _need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+i*LIBMESH_DIM+j] 为真,
277  // 当且仅当需要计算 grad(v, i, j)
278  std::vector<bool> _need_var_hess;
279 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
280 
281  // 可以由函数解析器解析和处理的变量/值
282  std::string variables; // 可解析的变量字符串。
283  std::vector<std::string> _additional_vars; // 附加变量名称的字符串向量。
284  std::vector<Output> _initial_vals; // 初始值的输出向量。
285 };
286 
287 
288 /*----------------------- Inline functions ----------------------------------*/
289 
290 template <typename Output>
291 inline
293  std::string expression,
294  const std::vector<std::string> * additional_vars,
295  const std::vector<Output> * initial_vals) :
296  _sys(sys),
297  _expression (), // overridden by parse()
298  _n_vars(sys.n_vars()),
299  _n_requested_vars(0),
300  _n_requested_grad_components(0),
301  _n_requested_hess_components(0),
302  _requested_normals(false),
303  _need_var(_n_vars, false),
304  _need_var_grad(_n_vars*LIBMESH_DIM, false),
305 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
306  _need_var_hess(_n_vars*LIBMESH_DIM*LIBMESH_DIM, false),
307 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
308  _additional_vars (additional_vars ? *additional_vars : std::vector<std::string>()),
309  _initial_vals (initial_vals ? *initial_vals : std::vector<Output>())
310 {
311  this->reparse(std::move(expression));
312 }
313 
314 
315 template <typename Output>
316 inline
318  FEMFunctionBase<Output>(other),
319  _sys(other._sys),
320  _expression(other._expression),
321  _subexpressions(other._subexpressions),
322  _n_vars(other._n_vars),
323  _n_requested_vars(other._n_requested_vars),
324  _n_requested_grad_components(other._n_requested_grad_components),
325  _n_requested_hess_components(other._n_requested_hess_components),
326  _requested_normals(other._requested_normals),
327  _spacetime(other._spacetime),
328  _need_var(other._need_var),
329  _need_var_grad(other._need_var_grad),
330 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
331  _need_var_hess(other._need_var_hess),
332 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
333  variables(other.variables),
334  _additional_vars(other._additional_vars),
335  _initial_vals(other._initial_vals)
336 {
337  this->reparse(_expression);
338 }
339 
340 
341 template <typename Output>
342 inline
345 {
346  // We can only be assigned another ParsedFEMFunction defined on the same System
347  libmesh_assert(&_sys == &other._sys);
348 
349  // Use copy-and-swap idiom
350  ParsedFEMFunction<Output> tmp(other);
351  std::swap(tmp, *this);
352  return *this;
353 }
354 
355 
356 template <typename Output>
357 inline
358 void
359 ParsedFEMFunction<Output>::reparse (std::string expression)
360 {
361  variables = "x";
362 #if LIBMESH_DIM > 1
363  variables += ",y";
364 #endif
365 #if LIBMESH_DIM > 2
366  variables += ",z";
367 #endif
368  variables += ",t";
369 
370  for (unsigned int v=0; v != _n_vars; ++v)
371  {
372  const std::string & varname = _sys.variable_name(v);
373  std::size_t varname_i = find_name(varname, expression);
374 
375  // If we didn't find our variable name then let's go to the
376  // next.
377  if (varname_i == std::string::npos)
378  continue;
379 
380  _need_var[v] = true;
381  variables += ',';
382  variables += varname;
383  _n_requested_vars++;
384  }
385 
386  for (unsigned int v=0; v != _n_vars; ++v)
387  {
388  const std::string & varname = _sys.variable_name(v);
389 
390  for (unsigned int d=0; d != LIBMESH_DIM; ++d)
391  {
392  std::string gradname = std::string("grad_") +
393  "xyz"[d] + '_' + varname;
394  std::size_t gradname_i = find_name(gradname, expression);
395 
396  // If we didn't find that gradient component of our
397  // variable name then let's go to the next.
398  if (gradname_i == std::string::npos)
399  continue;
400 
401  _need_var_grad[v*LIBMESH_DIM+d] = true;
402  variables += ',';
403  variables += gradname;
404  _n_requested_grad_components++;
405  }
406  }
407 
408 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
409  for (unsigned int v=0; v != _n_vars; ++v)
410  {
411  const std::string & varname = _sys.variable_name(v);
412 
413  for (unsigned int d1=0; d1 != LIBMESH_DIM; ++d1)
414  for (unsigned int d2=0; d2 != LIBMESH_DIM; ++d2)
415  {
416  std::string hessname = std::string("hess_") +
417  "xyz"[d1] + "xyz"[d2] + '_' + varname;
418  std::size_t hessname_i = find_name(hessname, expression);
419 
420  // If we didn't find that hessian component of our
421  // variable name then let's go to the next.
422  if (hessname_i == std::string::npos)
423  continue;
424 
425  _need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+d1*LIBMESH_DIM+d2]
426  = true;
427  variables += ',';
428  variables += hessname;
429  _n_requested_hess_components++;
430  }
431  }
432 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
433 
434  {
435  std::size_t nx_i = find_name("n_x", expression);
436  std::size_t ny_i = find_name("n_y", expression);
437  std::size_t nz_i = find_name("n_z", expression);
438 
439  // If we found any requests for normal components then we'll
440  // compute normals
441  if (nx_i != std::string::npos ||
442  ny_i != std::string::npos ||
443  nz_i != std::string::npos)
444  {
445  _requested_normals = true;
446  variables += ',';
447  variables += "n_x";
448  if (LIBMESH_DIM > 1)
449  {
450  variables += ',';
451  variables += "n_y";
452  }
453  if (LIBMESH_DIM > 2)
454  {
455  variables += ',';
456  variables += "n_z";
457  }
458  }
459  }
460 
461  _spacetime.resize
462  (LIBMESH_DIM + 1 + _n_requested_vars +
463  _n_requested_grad_components + _n_requested_hess_components +
464  (_requested_normals ? LIBMESH_DIM : 0) +
465  _additional_vars.size());
466 
467  // If additional vars were passed, append them to the string
468  // that we send to the function parser. Also add them to the
469  // end of our spacetime vector
470  unsigned int offset = LIBMESH_DIM + 1 + _n_requested_vars +
471  _n_requested_grad_components + _n_requested_hess_components;
472 
473  for (auto i : index_range(_additional_vars))
474  {
475  variables += "," + _additional_vars[i];
476  // Initialize extra variables to the vector passed in or zero
477  // Note: The initial_vals vector can be shorter than the additional_vars vector
478  _spacetime[offset + i] =
479  (i < _initial_vals.size()) ? _initial_vals[i] : 0;
480  }
481 
482  this->partial_reparse(std::move(expression));
483 }
484 
485 template <typename Output>
486 inline
487 void
489 {
490  for (unsigned int v=0; v != _n_vars; ++v)
491  {
492  FEBase * elem_fe;
493  c.get_element_fe(v, elem_fe);
494  if (_n_requested_vars)
495  elem_fe->get_phi();
496  if (_n_requested_grad_components)
497  elem_fe->get_dphi();
498 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
499  if (_n_requested_hess_components)
500  elem_fe->get_d2phi();
501 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
502  }
503 
504  if (_requested_normals)
505  {
506  FEBase * side_fe;
507  c.get_side_fe(0, side_fe);
508 
509  side_fe->get_normals();
510 
511  // FIXME: this is a hack to support normals at quadrature
512  // points; we don't support normals elsewhere.
513  side_fe->get_xyz();
514  }
515 }
516 
517 template <typename Output>
518 inline
519 std::unique_ptr<FEMFunctionBase<Output>>
521 {
522  return std::make_unique<ParsedFEMFunction>
523  (_sys, _expression, &_additional_vars, &_initial_vals);
524 }
525 
526 template <typename Output>
527 inline
528 Output
530  const Point & p,
531  const Real time)
532 {
533  eval_args(c, p, time);
534 
535  return eval(*parsers[0], "f", 0);
536 }
537 
538 
539 
540 template <typename Output>
541 inline
542 void
544  const Point & p,
545  const Real time,
546  DenseVector<Output> & output)
547 {
548  eval_args(c, p, time);
549 
550  unsigned int size = output.size();
551 
552  libmesh_assert_equal_to (size, parsers.size());
553 
554  for (unsigned int i=0; i != size; ++i)
555  output(i) = eval(*parsers[i], "f", i);
556 }
557 
558 
559 template <typename Output>
560 inline
561 Output
563  unsigned int i,
564  const Point & p,
565  Real time)
566 {
567  eval_args(c, p, time);
568 
569  libmesh_assert_less (i, parsers.size());
570  return eval(*parsers[i], "f", i);
571 }
572 
573 template <typename Output>
574 inline
575 Output
576 ParsedFEMFunction<Output>::get_inline_value(std::string_view inline_var_name) const
577 {
578  libmesh_assert_greater (_subexpressions.size(), 0);
579 
580 #ifndef NDEBUG
581  bool found_var_name = false;
582 #endif
583  Output old_var_value(0.);
584 
585  for (const std::string & subexpression : _subexpressions)
586  {
587  const std::size_t varname_i =
588  find_name(inline_var_name, subexpression);
589  if (varname_i == std::string::npos)
590  continue;
591 
592  const std::size_t assignment_i =
593  subexpression.find(":", varname_i+1);
594 
595  libmesh_assert_not_equal_to(assignment_i, std::string::npos);
596 
597  libmesh_assert_equal_to(subexpression[assignment_i+1], '=');
598  for (std::size_t i = varname_i+1; i != assignment_i; ++i)
599  libmesh_assert_equal_to(subexpression[i], ' ');
600 
601  std::size_t end_assignment_i =
602  subexpression.find(";", assignment_i+1);
603 
604  libmesh_assert_not_equal_to(end_assignment_i, std::string::npos);
605 
606  std::string new_subexpression =
607  subexpression.substr(0, end_assignment_i+1).
608  append(inline_var_name);
609 
610 #ifdef LIBMESH_HAVE_FPARSER
611  // Parse and evaluate the new subexpression.
612  // Add the same constants as we used originally.
613  auto fp = std::make_unique<FunctionParserBase<Output>>();
614  fp->AddConstant("NaN", std::numeric_limits<Real>::quiet_NaN());
615  fp->AddConstant("pi", std::acos(Real(-1)));
616  fp->AddConstant("e", std::exp(Real(1)));
617  libmesh_error_msg_if
618  (fp->Parse(new_subexpression, variables) != -1, // -1 for success
619  "ERROR: FunctionParser is unable to parse modified expression: "
620  << new_subexpression << '\n' << fp->ErrorMsg());
621 
622  Output new_var_value = this->eval(*fp, new_subexpression, 0);
623 #ifdef NDEBUG
624  return new_var_value;
625 #else
626  if (found_var_name)
627  {
628  libmesh_assert_equal_to(old_var_value, new_var_value);
629  }
630  else
631  {
632  old_var_value = new_var_value;
633  found_var_name = true;
634  }
635 #endif
636 
637 #else
638  libmesh_error_msg("ERROR: This functionality requires fparser!");
639 #endif
640  }
641 
642  libmesh_assert(found_var_name);
643  return old_var_value;
644 }
645 
646 template <typename Output>
647 inline
648 void
649 ParsedFEMFunction<Output>::set_inline_value (std::string_view inline_var_name,
650  Output newval)
651 {
652  libmesh_assert_greater (_subexpressions.size(), 0);
653 
654 #ifndef NDEBUG
655  bool found_var_name = false;
656 #endif
657  for (auto s : index_range(_subexpressions))
658  {
659  const std::string & subexpression = _subexpressions[s];
660  const std::size_t varname_i =
661  find_name(inline_var_name, subexpression);
662  if (varname_i == std::string::npos)
663  continue;
664 
665 #ifndef NDEBUG
666  found_var_name = true;
667 #endif
668 
669  const std::size_t assignment_i =
670  subexpression.find(":", varname_i+1);
671 
672  libmesh_assert_not_equal_to(assignment_i, std::string::npos);
673 
674  libmesh_assert_equal_to(subexpression[assignment_i+1], '=');
675  for (std::size_t i = varname_i+1; i != assignment_i; ++i)
676  libmesh_assert_equal_to(subexpression[i], ' ');
677 
678  std::size_t end_assignment_i =
679  subexpression.find(";", assignment_i+1);
680 
681  libmesh_assert_not_equal_to(end_assignment_i, std::string::npos);
682 
683  std::ostringstream new_subexpression;
684  new_subexpression << subexpression.substr(0, assignment_i+2)
685  << std::setprecision(std::numeric_limits<Output>::digits10+2)
686 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
687  << '(' << newval.real() << '+'
688  << newval.imag() << 'i' << ')'
689 #else
690  << newval
691 #endif
692  << subexpression.substr(end_assignment_i,
693  std::string::npos);
694  _subexpressions[s] = new_subexpression.str();
695  }
696 
697  libmesh_assert(found_var_name);
698 
699  std::string new_expression;
700 
701  for (const auto & subexpression : _subexpressions)
702  {
703  new_expression += '{';
704  new_expression += subexpression;
705  new_expression += '}';
706  }
707 
708  this->partial_reparse(new_expression);
709 }
710 
711 
712 template <typename Output>
713 inline
714 void
716 {
717  _expression = std::move(expression);
718  _subexpressions.clear();
719  parsers.clear();
720 
721  size_t nextstart = 0, end = 0;
722 
723  while (end != std::string::npos)
724  {
725  // If we're past the end of the string, we can't make any more
726  // subparsers
727  if (nextstart >= _expression.size())
728  break;
729 
730  // If we're at the start of a brace delimited section, then we
731  // parse just that section:
732  if (_expression[nextstart] == '{')
733  {
734  nextstart++;
735  end = _expression.find('}', nextstart);
736  }
737  // otherwise we parse the whole thing
738  else
739  end = std::string::npos;
740 
741  // We either want the whole end of the string (end == npos) or
742  // a substring in the middle.
743  _subexpressions.push_back
744  (_expression.substr(nextstart, (end == std::string::npos) ?
745  std::string::npos : end - nextstart));
746 
747  // fparser can crash on empty expressions
748  libmesh_error_msg_if(_subexpressions.back().empty(),
749  "ERROR: FunctionParser is unable to parse empty expression.\n");
750 
751 
752 #ifdef LIBMESH_HAVE_FPARSER
753  // Parse (and optimize if possible) the subexpression.
754  // Add some basic constants, to Real precision.
755  auto fp = std::make_unique<FunctionParserBase<Output>>();
756  fp->AddConstant("NaN", std::numeric_limits<Real>::quiet_NaN());
757  fp->AddConstant("pi", std::acos(Real(-1)));
758  fp->AddConstant("e", std::exp(Real(1)));
759  libmesh_error_msg_if
760  (fp->Parse(_subexpressions.back(), variables) != -1, // -1 for success
761  "ERROR: FunctionParser is unable to parse expression: "
762  << _subexpressions.back() << '\n' << fp->ErrorMsg());
763  fp->Optimize();
764  parsers.push_back(std::move(fp));
765 #else
766  libmesh_error_msg("ERROR: This functionality requires fparser!");
767 #endif
768 
769  // If at end, use nextstart=maxSize. Else start at next
770  // character.
771  nextstart = (end == std::string::npos) ?
772  std::string::npos : end + 1;
773  }
774 }
775 
776 template <typename Output>
777 inline
778 std::size_t
779 ParsedFEMFunction<Output>::find_name (std::string_view varname,
780  std::string_view expr) const
781 {
782  const std::size_t namesize = varname.size();
783  std::size_t varname_i = expr.find(varname);
784 
785  while ((varname_i != std::string::npos) &&
786  (((varname_i > 0) &&
787  (std::isalnum(expr[varname_i-1]) ||
788  (expr[varname_i-1] == '_'))) ||
789  ((varname_i+namesize < expr.size()) &&
790  (std::isalnum(expr[varname_i+namesize]) ||
791  (expr[varname_i+namesize] == '_')))))
792  {
793  varname_i = expr.find(varname, varname_i+1);
794  }
795 
796  return varname_i;
797 }
798 
799 
800 
801 // Helper function for evaluating function arguments
802 template <typename Output>
803 inline
804 void
806  const Point & p,
807  const Real time)
808 {
809  _spacetime[0] = p(0);
810 #if LIBMESH_DIM > 1
811  _spacetime[1] = p(1);
812 #endif
813 #if LIBMESH_DIM > 2
814  _spacetime[2] = p(2);
815 #endif
816  _spacetime[LIBMESH_DIM] = time;
817 
818  unsigned int request_index = 0;
819  for (unsigned int v=0; v != _n_vars; ++v)
820  {
821  if (!_need_var[v])
822  continue;
823 
824  c.point_value(v, p, _spacetime[LIBMESH_DIM+1+request_index]);
825  request_index++;
826  }
827 
828  if (_n_requested_grad_components)
829  for (unsigned int v=0; v != _n_vars; ++v)
830  {
831  if (!_need_var_grad[v*LIBMESH_DIM]
832 #if LIBMESH_DIM > 1
833  && !_need_var_grad[v*LIBMESH_DIM+1]
834 #if LIBMESH_DIM > 2
835  && !_need_var_grad[v*LIBMESH_DIM+2]
836 #endif
837 #endif
838  )
839  continue;
840 
841  Gradient g;
842  c.point_gradient(v, p, g);
843 
844  for (unsigned int d=0; d != LIBMESH_DIM; ++d)
845  {
846  if (!_need_var_grad[v*LIBMESH_DIM+d])
847  continue;
848 
849  _spacetime[LIBMESH_DIM+1+request_index] = g(d);
850  request_index++;
851  }
852  }
853 
854 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
855  if (_n_requested_hess_components)
856  for (unsigned int v=0; v != _n_vars; ++v)
857  {
858  if (!_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM]
859 #if LIBMESH_DIM > 1
860  && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+1]
861  && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+2]
862  && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+3]
863 #if LIBMESH_DIM > 2
864  && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+4]
865  && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+5]
866  && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+6]
867  && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+7]
868  && !_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+8]
869 #endif
870 #endif
871  )
872  continue;
873 
874  Tensor h;
875  c.point_hessian(v, p, h);
876 
877  for (unsigned int d1=0; d1 != LIBMESH_DIM; ++d1)
878  for (unsigned int d2=0; d2 != LIBMESH_DIM; ++d2)
879  {
880  if (!_need_var_hess[v*LIBMESH_DIM*LIBMESH_DIM+d1*LIBMESH_DIM+d2])
881  continue;
882 
883  _spacetime[LIBMESH_DIM+1+request_index] = h(d1,d2);
884  request_index++;
885  }
886  }
887 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
888 
889  if (_requested_normals)
890  {
891  FEBase * side_fe;
892  c.get_side_fe(0, side_fe);
893 
894  const std::vector<Point> & normals = side_fe->get_normals();
895 
896  const std::vector<Point> & xyz = side_fe->get_xyz();
897 
898  libmesh_assert_equal_to(normals.size(), xyz.size());
899 
900  // We currently only support normals at quadrature points!
901 #ifndef NDEBUG
902  bool at_quadrature_point = false;
903 #endif
904  for (auto qp : index_range(normals))
905  {
906  if (p == xyz[qp])
907  {
908  const Point & n = normals[qp];
909  for (unsigned int d=0; d != LIBMESH_DIM; ++d)
910  {
911  _spacetime[LIBMESH_DIM+1+request_index] = n(d);
912  request_index++;
913  }
914 #ifndef NDEBUG
915  at_quadrature_point = true;
916 #endif
917  break;
918  }
919  }
920 
921  libmesh_assert(at_quadrature_point);
922  }
923 
924  // The remaining locations in _spacetime are currently set at construction
925  // but could potentially be made dynamic
926 }
927 
928 
929 // Evaluate the ith FunctionParser and check the result
930 #ifdef LIBMESH_HAVE_FPARSER
931 template <typename Output>
932 inline
933 Output
934 ParsedFEMFunction<Output>::eval (FunctionParserBase<Output> & parser,
935  std::string_view libmesh_dbg_var(function_name),
936  unsigned int libmesh_dbg_var(component_idx)) const
937 {
938 #ifndef NDEBUG
939  Output result = parser.Eval(_spacetime.data());
940  int error_code = parser.EvalError();
941  if (error_code)
942  {
943  libMesh::err << "ERROR: FunctionParser is unable to evaluate component "
944  << component_idx
945  << " for '"
946  << function_name;
947 
948  for (auto i : index_range(parsers))
949  if (parsers[i].get() == &parser)
950  libMesh::err << "' of expression '"
951  << _subexpressions[i];
952 
953  libMesh::err << "' with arguments:\n";
954  for (const auto & st : _spacetime)
955  libMesh::err << '\t' << st << '\n';
956  libMesh::err << '\n';
957 
958  // Currently no API to report error messages, we'll do it manually
959  std::string error_message = "Reason: ";
960 
961  switch (error_code)
962  {
963  case 1:
964  error_message += "Division by zero";
965  break;
966  case 2:
967  error_message += "Square Root error (negative value)";
968  break;
969  case 3:
970  error_message += "Log error (negative value)";
971  break;
972  case 4:
973  error_message += "Trigonometric error (asin or acos of illegal value)";
974  break;
975  case 5:
976  error_message += "Maximum recursion level reached";
977  break;
978  default:
979  error_message += "Unknown";
980  break;
981  }
982  libmesh_error_msg(error_message);
983  }
984 
985  return result;
986 #else
987  return parser.Eval(_spacetime.data());
988 #endif
989 }
990 #else // LIBMESH_HAVE_FPARSER
991 template <typename Output>
992 inline
993 Output
994 ParsedFEMFunction<Output>::eval (char & /*parser*/,
995  std::string_view /*function_name*/,
996  unsigned int /*component_idx*/) const
997 {
998  libmesh_error_msg("ERROR: This functionality requires fparser!");
999  return Output(0);
1000 }
1001 #endif
1002 
1003 
1004 } // namespace libMesh
1005 
1006 #endif // LIBMESH_PARSED_FEM_FUNCTION_H
Output get_inline_value(std::string_view inline_var_name) const
获取内联变量的值。
std::vector< bool > _need_var_grad
virtual void init_context(const FEMContext &c) override
初始化上下文。
std::vector< Output > _initial_vals
Output eval(FunctionParserBase< Output > &parser, std::string_view libmesh_dbg_var(function_name), unsigned int libmesh_dbg_var(component_idx)) const
评估第 i 个 FunctionParser 并检查结果。
std::vector< std::string > _additional_vars
此类定义了LIBMESH_DIM维的实数或复数空间中的向量。
Definition: tensor_tools.h:35
std::vector< std::string > _subexpressions
void partial_reparse(std::string expression)
用于重新解析表达式的辅助函数。
std::vector< Output > _spacetime
virtual std::unique_ptr< FEMFunctionBase< Output > > clone() const override
克隆函数。
virtual Output operator()(const FEMContext &c, const Point &p, const Real time=0.) override
调用运算符,用于计算解析函数的值。
ParsedFEMFunction 提供对 FEMSystem 中基于 FParser 的解析函数的支持。
void reparse(std::string expression)
重新解析使用新表达式。
std::vector< std::unique_ptr< FunctionParserBase< Output > > > parsers
ParsedFEMFunction(const System &sys, std::string expression, const std::vector< std::string > *additional_vars=nullptr, const std::vector< Output > *initial_vals=nullptr)
构造函数。
ParsedFEMFunction & operator=(const ParsedFEMFunction &)
赋值运算符。
std::vector< bool > _need_var_hess
void set_inline_value(std::string_view inline_var_name, Output newval)
更改内联变量的值。
void eval_args(const FEMContext &c, const Point &p, const Real time)
用于计算函数参数的辅助函数。
OStreamProxy err
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::size_t find_name(std::string_view varname, std::string_view expr) const
用于解析变量名称的辅助函数。
virtual Output component(const FEMContext &c, unsigned int i, const Point &p, Real time=0.) override
计算解析函数的第 i 个分量的值。
std::vector< char * > parsers
FEMFunctionBase是一个基类,用户可以从中派生出“函数样式”的对象,以在FEMSystem中使用。
virtual unsigned int size() const overridefinal
Definition: dense_vector.h:111
const std::string & expression()
获取解析表达式。
此类定义了LIBMESH_DIM维度的实数或复数空间中的张量。typedef RealTensorValue总是定义为实数值的张量, 而NumberTensorValue则根据库的配置定义为实数或复数值...
Definition: tensor_tools.h:37