Flow123d  JS_before_hm-2208-gb971e62bf
field_formula.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file field_formula.cc
15  * @brief
16  */
17 
18 
19 
20 #include "fields/field_formula.hh"
21 #include "fields/field_instances.hh" // for instantiation macros
22 #include "fields/surface_depth.hh"
23 #include "fparser.hh"
24 #include "input/input_type.hh"
25 #include "include/arena_alloc.hh" // bparser
26 #include <boost/algorithm/string/replace.hpp>
27 
28 
29 /// Implementation.
30 
31 namespace it = Input::Type;
32 
33 FLOW123D_FORCE_LINK_IN_CHILD(field_formula)
34 
35 
36 template <int spacedim, class Value>
38 {
39 
40  return it::Record("FieldFormula", FieldAlgorithmBase<spacedim,Value>::template_name()+" Field given by runtime interpreted formula.")
43  .declare_key("value", STI::get_input_type() , it::Default::obligatory(),
44  "String, array of strings, or matrix of strings with formulas for individual "
45  "entries of scalar, vector, or tensor value respectively.\n"
46  "For vector values, you can use just one string to enter homogeneous vector.\n"
47  "For square (($N\\times N$))-matrix values, you can use:\n\n"
48  " - array of strings of size (($N$)) to enter diagonal matrix\n"
49  " - array of strings of size (($\\frac12N(N+1)$)) to enter symmetric matrix (upper triangle, row by row)\n"
50  " - just one string to enter (spatially variable) multiple of the unit matrix.\n"
51  "Formula can contain variables ```x,y,z,t,d``` and usual operators and functions." )
52  //.declare_key("unit", FieldAlgorithmBase<spacedim, Value>::get_input_type_unit_si(), it::Default::optional(),
53  // "Definition of unit.")
54  .declare_key("surface_direction", it::String(), it::Default("\"0 0 1\""),
55  "The vector used to project evaluation point onto the surface.")
56  .declare_key("surface_region", it::String(), it::Default::optional(),
57  "The name of region set considered as the surface. You have to set surface region if you "
58  "want to use formula variable ```d```.")
59  .allow_auto_conversion("value")
60  .close();
61 }
62 
63 
64 
65 template <int spacedim, class Value>
67  Input::register_class< FieldFormula<spacedim, Value>, unsigned int >("FieldFormula") +
69 
70 
71 
72 template <int spacedim, class Value>
74 : FieldAlgorithmBase<spacedim, Value>(n_comp),
75  formula_matrix_(this->value_.n_rows(), this->value_.n_cols()),
76  first_time_set_(true), arena_alloc_(nullptr)
77 {
78  this->is_constant_in_space_ = false;
79  parser_matrix_.resize(this->value_.n_rows());
80  for(unsigned int row=0; row < this->value_.n_rows(); row++) {
81  parser_matrix_[row].resize(this->value_.n_cols());
82  }
83  b_parser_.reserve(this->value_.n_rows()*this->value_.n_cols());
84  for(unsigned int i=0; i < this->value_.n_rows()*this->value_.n_cols(); i++) {
85  b_parser_.emplace_back( 1.1 * CacheMapElementNumber::get() );
86  }
87 }
88 
89 
90 
91 template <int spacedim, class Value>
93  this->init_unit_conversion_coefficient(rec, init_data);
94 
95  // read formulas form input
96  STI::init_from_input( formula_matrix_, rec.val<typename STI::AccessType>("value") );
97  in_rec_ = rec;
98 }
99 
100 
101 template <int spacedim, class Value>
103 
104  if (!time.use_fparser_) {
105  this->time_=time;
106  this->is_constant_in_space_ = false;
107  return true;
108  }
109 
110  /* OLD FPARSER CODE */
111  bool any_parser_changed = false;
112  std::string value_input_address = in_rec_.address_string();
113  has_depth_var_ = false;
114  this->is_constant_in_space_ = true; // set flag to true, then if found 'x', 'y', 'z' or 'd' reset to false
115 
116 
117  std::string vars = string("x,y,z").substr(0, 2*spacedim-1);
118  std::vector<bool> time_dependent(this->value_.n_rows() * this->value_.n_cols(), false);
119  // update parsers
120  for(unsigned int row=0; row < this->value_.n_rows(); row++)
121  for(unsigned int col=0; col < this->value_.n_cols(); col++) {
122  // get all variable names from the formula
123  std::vector<std::string> var_list;
124 
125  FunctionParser tmp_parser;
126  tmp_parser.AddConstant("Pi", 3.14159265358979323846);
127  tmp_parser.AddConstant("E", 2.71828182845904523536);
128 
129 
130 #pragma GCC diagnostic push
131 #pragma GCC diagnostic ignored "-Wunused-variable"
132  {
133  int err=tmp_parser.ParseAndDeduceVariables(formula_matrix_.at(row,col), var_list);
134  ASSERT(err != FunctionParser::FP_NO_ERROR)(tmp_parser.ErrorMsg()).error("ParseAndDeduceVariables error\n");
135  }
136 #pragma GCC diagnostic pop
137 
138  for(std::string &var_name : var_list ) {
139  if (var_name == std::string("t") ) time_dependent[row*this->value_.n_rows()+col]=true;
140  else if (var_name == std::string("d") ) {
141  this->is_constant_in_space_ = false;
142  if (surface_depth_)
143  has_depth_var_=true;
144  else
145  WarningOut().fmt("Unset surface region. Variable '{}' in the FieldFormula[{}][{}] == '{}' will be set to zero\n at the input address:\n {} \n",
146  var_name, row, col, formula_matrix_.at(row,col), value_input_address );
147  }
148  else if (var_name == "x" || var_name == "y" || var_name == "z") {
149  this->is_constant_in_space_ = false;
150  continue;
151  }
152  else
153  WarningOut().fmt("Unknown variable '{}' in the FieldFormula[{}][{}] == '{}'\n at the input address:\n {} \n",
154  var_name, row, col, formula_matrix_.at(row,col), value_input_address );
155  }
156 
157  // Seems that we can not just add 't' constant to tmp_parser, since it was already Parsed.
158  parser_matrix_[row][col].AddConstant("Pi", 3.14159265358979323846);
159  parser_matrix_[row][col].AddConstant("E", 2.71828182845904523536);
160  if (time_dependent[row*this->value_.n_rows()+col]) {
161  parser_matrix_[row][col].AddConstant("t", time.end());
162  }
163  }
164 
165  if (has_depth_var_)
166  vars += string(",d");
167 
168  // update parsers
169  for(unsigned int row=0; row < this->value_.n_rows(); row++)
170  for(unsigned int col=0; col < this->value_.n_cols(); col++) {
171  // TODO:
172  // - possibly add user defined constants and units here ...
173  // - optimization; possibly parse only if time_dependent || formula_matrix[][] has changed ...
174  //parser_matrix_[row][col] = tmp_parser;
175  if (time_dependent[row*this->value_.n_rows()+col] || first_time_set_ ) {
176  parser_matrix_[row][col].Parse(formula_matrix_.at(row,col), vars);
177 
178  if ( parser_matrix_[row][col].GetParseErrorType() != FunctionParser::FP_NO_ERROR ) {
179  THROW( ExcFParserError() << EI_FParserMsg(parser_matrix_[row][col].ErrorMsg()) << EI_Row(row)
180  << EI_Col(col) << EI_Formula(formula_matrix_.at(row,col)) );
181  }
182 
183  parser_matrix_[row][col].Optimize();
184  any_parser_changed = true;
185  }
186 
187 
188  }
189 
190  first_time_set_ = false;
191  this->time_=time;
192  return any_parser_changed;
193 }
194 
195 
196 template <int spacedim, class Value>
197 void FieldFormula<spacedim, Value>::set_mesh(const Mesh *mesh, FMT_UNUSED bool boundary_domain) {
198  // create SurfaceDepth object if surface region is set
199  std::string surface_region;
200  if ( in_rec_.opt_val("surface_region", surface_region) ) {
201  surface_depth_ = std::make_shared<SurfaceDepth>(mesh, surface_region, in_rec_.val<std::string>("surface_direction"));
202  }
203 }
204 
205 
206 /**
207  * Returns one value in one given point. ResultType can be used to avoid some costly calculation if the result is trivial.
208  */
209 template <int spacedim, class Value>
210 typename Value::return_type const & FieldFormula<spacedim, Value>::value(const Point &p, FMT_UNUSED const ElementAccessor<spacedim> &elm)
211 {
212 
213  auto p_depth = this->eval_depth_var(p);
214  for(unsigned int row=0; row < this->value_.n_rows(); row++)
215  for(unsigned int col=0; col < this->value_.n_cols(); col++) {
216  this->value_(row,col) = this->unit_conversion_coefficient_ * parser_matrix_[row][col].Eval(p_depth.memptr());
217  }
218  return this->r_value_;
219 }
220 
221 
222 /**
223  * Returns std::vector of scalar values in several points at once.
224  */
225 template <int spacedim, class Value>
228 {
229  ASSERT_EQ( point_list.size(), value_list.size() );
230  ASSERT_DBG( point_list.n_rows() == spacedim && point_list.n_cols() == 1).error("Invalid point size.\n");
231  for(unsigned int i=0; i< point_list.size(); i++) {
232  Value envelope(value_list[i]);
233  ASSERT_EQ( envelope.n_rows(), this->value_.n_rows() )(i)(envelope.n_rows())(this->value_.n_rows())
234  .error("value_list['i'] has wrong number of rows\n");
235  auto p_depth = this->eval_depth_var(point_list.vec<spacedim>(i));
236 
237  for(unsigned int row=0; row < this->value_.n_rows(); row++)
238  for(unsigned int col=0; col < this->value_.n_cols(); col++) {
239  envelope(row,col) = this->unit_conversion_coefficient_ * parser_matrix_[row][col].Eval(p_depth.memptr());
240  }
241  }
242 }
243 
244 
245 template <int spacedim, class Value>
247  ElementCacheMap &cache_map, unsigned int region_patch_idx)
248 {
249  unsigned int reg_chunk_begin = cache_map.region_chunk_begin(region_patch_idx);
250  unsigned int reg_chunk_end = cache_map.region_chunk_end(region_patch_idx);
251 
252  for (unsigned int i=reg_chunk_begin; i<reg_chunk_end; ++i) {
253  res_[i] = 0.0;
254  }
255  for (auto it : eval_field_data_) {
256  // Copy data from dependent fields to arena. Temporary solution.
257  // TODO hold field data caches in arena, remove this step
258  auto value_cache = it.first->value_cache();
259  for (unsigned int i=reg_chunk_begin; i<reg_chunk_end; ++i) {
260  if (it.first->name() == "X") {
261  x_[i] = value_cache->template vec<3>(i)(0);
262  y_[i] = value_cache->template vec<3>(i)(1);
263  z_[i] = value_cache->template vec<3>(i)(2);
264  } else
265  it.second[i] = value_cache->data_[i];
266  }
267  }
268 
269  // Get vector of subsets as subarray
270  uint subsets_begin = reg_chunk_begin / ElementCacheMap::simd_size_double;
271  uint subsets_end = reg_chunk_end / ElementCacheMap::simd_size_double;
272  std::vector<uint> subset_vec;
273  subset_vec.assign(subsets_ + subsets_begin, subsets_ + subsets_end);
274 
275  for(unsigned int row=0; row < this->value_.n_rows(); row++)
276  for(unsigned int col=0; col < this->value_.n_cols(); col++) {
277  b_parser_[row*this->value_.n_cols()+col].set_subset(subset_vec);
278  b_parser_[row*this->value_.n_cols()+col].run();
279  for (unsigned int i=reg_chunk_begin; i<reg_chunk_end; ++i) {
280  auto cache_val = data_cache.template mat<Value::NRows_, Value::NCols_>(i);
281  cache_val(row, col) = res_[i];
282  data_cache.set(i) = cache_val;
283  }
284  }
285 }
286 
287 
288 template <int spacedim, class Value>
290 {
291  if (surface_depth_ && has_depth_var_) {
292  // add value of depth
293  arma::vec p_depth(spacedim+1);
294  p_depth.subvec(0,spacedim-1) = p;
295  try {
296  p_depth(spacedim) = surface_depth_->compute_distance(p);
297  } catch (SurfaceDepth::ExcTooLargeSnapDistance &e) {
298  e << SurfaceDepth::EI_FieldTime(this->time_.end());
299  e << in_rec_.ei_address();
300  throw;
301  }
302  return p_depth;
303  } else {
304  return p;
305  }
306 }
307 
309  uint r = 1;
310  for (auto i : shape) r *= i;
311  return r;
312 }
313 
314 
315 template <int spacedim, class Value>
317  required_fields_.clear(); // returned value
318 
319  std::vector<std::string> variables;
320  for(unsigned int row=0; row < this->value_.n_rows(); row++)
321  for(unsigned int col=0; col < this->value_.n_cols(); col++) {
322  // set expression and data to BParser
323  unsigned int i_p = row*this->value_.n_cols()+col;
324  //b_parser_[i_p].parse(formula_matrix_.at(row,col));
325  std::string expr = formula_matrix_.at(row,col); // Need replace some operations to make them compatible with BParser.
326  // It will be solved by conversion script after remove fparser, but
327  // we mix using of BParser and fparser and need this solution now.
328  boost::replace_all(expr, "^", "**"); // power function
329  boost::replace_all(expr, "max(", "maximum("); // max function
330  boost::replace_all(expr, "min(", "minimum("); // min function
331  boost::replace_all(expr, "Pi", "pi"); // Math.pi
332  boost::replace_all(expr, "E", "e"); // Math.e
333  { // ternary operator
334  std::string pref("if(");
335  auto res = std::mismatch(pref.begin(), pref.end(), expr.begin());
336  if ( (res.first == pref.end()) && (expr.back() == ')') ) {
337  std::string subexpr = expr.substr(3, expr.size()-4);
338  std::string delimiter = ",";
339  std::string cond = subexpr.substr(0, subexpr.find(delimiter));
340  subexpr.erase(0, cond.size()+1);
341  std::string if_case = subexpr.substr(0, subexpr.find(delimiter));
342  std::string else_case = subexpr.substr(if_case.size()+1);
343  expr = "(" + if_case + " if " + cond + " else " + else_case +")";
344  }
345  }
346  try {
347  b_parser_[i_p].parse( expr );
348  } catch (std::exception const& e) {
349  if (typeid(e) == typeid(bparser::Exception))
350  THROW( ExcParserError() << EI_BParserMsg(e.what()) << EI_Formula(expr) << Input::EI_Address( in_rec_.address_string() ) );
351  else throw;
352  }
353  auto free_symbols = b_parser_[i_p].free_symbols();
354  variables.insert(variables.end(), free_symbols.begin(), free_symbols.end());
355  }
356 
357  std::sort( variables.begin(), variables.end() );
358  variables.erase( std::unique( variables.begin(), variables.end() ), variables.end() );
359  has_time_=false;
360  sum_shape_sizes_=0; // scecifies size of arena
361  for (auto var : variables) {
362  if (var == "x" || var == "y" || var == "z") {
363  required_fields_.push_back( field_set.field("X") );
364  sum_shape_sizes_ += spacedim;
365  }
366  else if (var == "t") has_time_ = true;
367  else {
368  auto field_ptr = field_set.field(var);
369  if (field_ptr != nullptr) required_fields_.push_back( field_ptr );
370  else {
371  field_ptr = field_set.user_field(var, this->time_);
372  if (field_ptr != nullptr) required_fields_.push_back( field_ptr );
373  else THROW( ExcUnknownField() << EI_Field(var) << Input::EI_Address( in_rec_.address_string() ) );
374  }
375  // TODO: Test the exception, report input line of the formula.
376  if (field_ptr->value_cache() == nullptr) THROW( ExcNotDoubleField() << EI_Field(var) << Input::EI_Address( in_rec_.address_string() ) );
377  // TODO: Test the exception, report input line of the formula.
378 
379  sum_shape_sizes_ += n_shape( field_ptr->shape_ );
380  if (var == "d") {
381  field_set.set_surface_depth(this->surface_depth_);
382  }
383  }
384  }
385 
386  return required_fields_;
387 }
388 
389 
390 template <int spacedim, class Value>
392 {
393  if (arena_alloc_!=nullptr) {
394  delete arena_alloc_;
395  }
396  eval_field_data_.clear();
397  uint vec_size = 1.1 * CacheMapElementNumber::get();
398  while (vec_size%ElementCacheMap::simd_size_double > 0) vec_size++; // alignment of block size
399  // number of subset alignment to block size
401  uint n_vectors = sum_shape_sizes_ + 1; // needs add space of result vector
402  arena_alloc_ = new bparser::ArenaAlloc(ElementCacheMap::simd_size_double, n_vectors * vec_size * sizeof(double) + n_subsets * sizeof(uint));
403  res_ = arena_alloc_->create_array<double>(vec_size);
404  for (auto field : required_fields_) {
405  std::string field_name = field->name();
406  eval_field_data_[field] = arena_alloc_->create_array<double>(n_shape( field->shape_ ) * vec_size);
407  if (field_name == "X") {
408  x_ = eval_field_data_[field] + 0;
409  y_ = eval_field_data_[field] + vec_size;
410  z_ = eval_field_data_[field] + 2*vec_size;
411  }
412  }
413  subsets_ = arena_alloc_->create_array<uint>(n_subsets);
414 
415  for(unsigned int row=0; row < this->value_.n_rows(); row++)
416  for(unsigned int col=0; col < this->value_.n_cols(); col++) {
417  // set expression and data to BParser
418  unsigned int i_p = row*this->value_.n_cols()+col;
419  if (has_time_) {
420  b_parser_[i_p].set_constant("t", {}, {this->time_.end()});
421  }
422  for (auto field : required_fields_) {
423  std::string field_name = field->name();
424  if (field_name == "X") {
425  b_parser_[i_p].set_variable("x", {}, x_);
426  b_parser_[i_p].set_variable("y", {}, y_);
427  b_parser_[i_p].set_variable("z", {}, z_);
428  } else
429  b_parser_[i_p].set_variable(field_name, {}, eval_field_data_[field]);
430  }
431  b_parser_[i_p].set_variable("_result_", {}, res_);
432  b_parser_[i_p].compile();
433  }
434  for (uint i=0; i<n_subsets; ++i)
435  subsets_[i] = i;
436 }
437 
438 
439 template <int spacedim, class Value>
441 
442 
443 // Instantiations of FieldFormula
surface_depth.hh
CacheMapElementNumber::get
static unsigned int get()
Return number of stored elements.
Definition: field_value_cache.hh:93
FieldFormula::eval_depth_var
arma::vec eval_depth_var(const Point &p)
Definition: field_formula.cc:289
Armor::vec
ArmaVec< double, N > vec
Definition: armor.hh:885
FieldFormula::b_parser_
std::vector< bparser::Parser > b_parser_
Definition: field_formula.hh:151
ASSERT
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
Definition: asserts.hh:347
ElementCacheMap::simd_size_double
static const unsigned int simd_size_double
Definition: field_value_cache.hh:158
ElementCacheMap
Directing class of FieldValueCache.
Definition: field_value_cache.hh:151
Input::Record::val
const Ret val(const string &key) const
Definition: accessors_impl.hh:31
Armor::Array::vec
ArmaVec< Type, nr > vec(uint mat_index) const
Definition: armor.hh:821
FieldFormula::cache_update
void cache_update(FieldValueCache< typename Value::element_type > &data_cache, ElementCacheMap &cache_map, unsigned int region_patch_idx) override
Definition: field_formula.cc:246
FieldSet::user_field
FieldCommon * user_field(const std::string &field_name, const TimeStep &time)
Definition: field_set.cc:177
ASSERT_DBG
#define ASSERT_DBG(expr)
Definition: include_fadbad.hh:28
FLOW123D_FORCE_LINK_IN_CHILD
#define FLOW123D_FORCE_LINK_IN_CHILD(x)
Definition: global_defs.h:157
THROW
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
Definition: exceptions.hh:53
FieldAlgorithmBase::is_constant_in_space_
bool is_constant_in_space_
Flag detects that field is only dependent on time.
Definition: field_algo_base.hh:289
std::vector< bool >
ElementAccessor
Definition: dh_cell_accessor.hh:32
FieldSet::set_surface_depth
void set_surface_depth(std::shared_ptr< SurfaceDepth > surface_depth)
Set surface depth object to "d" field.
Definition: field_set.hh:374
FieldFormula::set_mesh
void set_mesh(const Mesh *mesh, bool boundary_domain) override
Definition: field_formula.cc:197
uint
unsigned int uint
Definition: mh_dofhandler.hh:101
TimeStep::use_fparser_
bool use_fparser_
Definition: time_governor.hh:235
FieldAlgorithmBase::Point
Space< spacedim >::Point Point
Definition: field_algo_base.hh:115
FieldFormula::~FieldFormula
virtual ~FieldFormula()
Definition: field_formula.cc:440
n_shape
uint n_shape(std::vector< uint > shape)
Definition: field_formula.cc:308
Armor::Array::n_rows
uint n_rows() const
Definition: armor.hh:715
ElementCacheMap::region_chunk_end
unsigned int region_chunk_end(unsigned int region_patch_idx) const
Return end position of region chunk in FieldValueCache.
Definition: field_value_cache.hh:257
FieldFormula::FieldFormula
FieldFormula(unsigned int n_comp=0)
Definition: field_formula.cc:73
Input::Type::Default
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
Input::Type::Record::derive_from
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:196
TimeStep::end
double end() const
Definition: time_governor.hh:161
FieldFormula::get_input_type
static const Input::Type::Record & get_input_type()
Implementation.
Definition: field_formula.cc:37
Input::Record
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
FieldAlgoBaseInitData
Helper struct stores data for initizalize descentants of FieldAlgorithmBase.
Definition: field_algo_base.hh:81
Input::Type::Record::allow_auto_conversion
virtual Record & allow_auto_conversion(const string &from_key)
Allows shorter input of the Record providing only value of the from_key given as the parameter.
Definition: type_record.cc:133
TimeStep
Representation of one time step..
Definition: time_governor.hh:123
FieldAlgorithmBase::value_
Value value_
Last value, prevents passing large values (vectors) by value.
Definition: field_algo_base.hh:280
FieldFormula::set_time
bool set_time(const TimeStep &time) override
Definition: field_formula.cc:102
Input::Type::Default::obligatory
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
ASSERT_EQ
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual)
Definition: asserts.hh:328
FieldFormula::init_from_input
virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData &init_data)
Definition: field_formula.cc:92
Input::Type::Record::declare_key
Record & declare_key(const string &key, std::shared_ptr< TypeBase > type, const Default &default_value, const string &description, TypeBase::attribute_map key_attributes=TypeBase::attribute_map())
Declares a new key of the Record.
Definition: type_record.cc:503
FieldSet
Container for various descendants of FieldCommonBase.
Definition: field_set.hh:159
FieldFormula::cache_reinit
void cache_reinit(const ElementCacheMap &cache_map) override
Definition: field_formula.cc:391
Armor::Array::n_cols
uint n_cols() const
Definition: armor.hh:720
Armor::Array::set
ArrayMatSet set(uint index)
Definition: armor.hh:838
Input::Type::Record::close
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:304
Armor::Array::size
unsigned int size() const
Definition: armor.hh:728
Input::Type
Definition: balance.hh:41
Input::Type::Record
Record type proxy class.
Definition: type_record.hh:182
Value
@ Value
Definition: finite_element.hh:43
input_type.hh
ElementCacheMap::region_chunk_begin
unsigned int region_chunk_begin(unsigned int region_patch_idx) const
Return begin position of region chunk in FieldValueCache.
Definition: field_value_cache.hh:251
Mesh
Definition: mesh.h:355
field_formula.hh
Input::Type::Record::copy_keys
Record & copy_keys(const Record &other)
Copy keys from other record.
Definition: type_record.cc:216
FieldFormula::parser_matrix_
std::vector< std::vector< FunctionParser > > parser_matrix_
Definition: field_formula.hh:150
FieldFormula::set_dependency
std::vector< const FieldCommon * > set_dependency(FieldSet &field_set) override
Definition: field_formula.cc:316
FieldAlgorithmBase
Definition: field_algo_base.hh:112
Input::Type::String
Class for declaration of the input data that are in string format.
Definition: type_base.hh:582
Input::Array
Accessor to input data conforming to declared Array.
Definition: accessors.hh:566
WarningOut
#define WarningOut()
Macro defining 'warning' record of log.
Definition: logger.hh:278
INSTANCE_ALL
#define INSTANCE_ALL(field)
Definition: field_instances.hh:43
field_instances.hh
Armor::Array< double >
FieldFormula::value_list
virtual void value_list(const Armor::array &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list)
Definition: field_formula.cc:226
Input::Type::Default::optional
static Default optional()
The factory function to make an empty default value which is optional.
Definition: type_record.hh:124
FieldFormula::value
virtual const Value::return_type & value(const Point &p, const ElementAccessor< spacedim > &elm)
Definition: field_formula.cc:210
FieldSet::field
FieldCommon * field(const std::string &field_name) const
Definition: field_set.cc:169
FieldFormula
Definition: field_formula.hh:63
FMT_UNUSED
#define FMT_UNUSED
Definition: posix.h:75