Flow123d  build_with_4.0.3-24aa4ef
field_formula.hh
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.hh
15  * @brief
16  */
17 
18 #ifndef FIELD_FORMULA_HH_
19 #define FIELD_FORMULA_HH_
20 
21 
22 #define BOOST_MATH_DISABLE_FLOAT128
23 
24 
25 #include <stdio.h> // for sprintf
26 #include <string> // for operator==, string
27 #include <vector> // for vector
28 #include <memory>
29 #include <armadillo>
30 #include <map>
31 #include "fields/field_algo_base.hh" // for FieldAlgorithmBase
32 #include "fields/field_values.hh" // for FieldValue<>::Enum, FieldValu...
33 #include "fields/field_set.hh"
34 #include "input/accessors.hh" // for ExcAccessorForNullStorage
35 #include "input/accessors_impl.hh" // for Record::val
36 #include "input/storage.hh" // for ExcStorageTypeMismatch
37 #include "input/type_record.hh" // for Record::ExcRecordKeyNotFound
38 #include "input/input_exception.hh" // for ExcAssertMsg::~ExcAssertMsg
39 #include "system/exceptions.hh" // for ExcAssertMsg::~ExcAssertMsg
40 #include "tools/time_governor.hh" // for TimeStep
41 //#include "include/assert.hh" // bparser
42 #include "parser.hh" // bparser
43 
44 template <int spacedim> class ElementAccessor;
45 class SurfaceDepth;
46 
47 using namespace std;
48 
49 /**
50  * Class representing fields given by runtime parsed formulas.
51  *
52  * Using libraries:
53  * https://github.com/flow123d/bparser/
54  *
55  * Allows parsing:
56  * - base variables: coordinates (x,y,z), time (t), surface depth (d); constants: e, pi
57  * - standard functions: trigonometric functions, min and max function, exponential function, ternary operator etc.
58  * - expressions dependendent on other fields
59  */
60 template <int spacedim, class Value>
61 class FieldFormula : public FieldAlgorithmBase<spacedim, Value>
62 {
63 public:
66 
67  TYPEDEF_ERR_INFO(EI_Field, std::string);
68  DECLARE_INPUT_EXCEPTION(ExcNotDoubleField,
69  << "Can not use integer valued field " << EI_Field::qval << " in the formula: \n");
70 
71  TYPEDEF_ERR_INFO(EI_BParserMsg, std::string);
72  TYPEDEF_ERR_INFO(EI_Formula, std::string);
73  DECLARE_INPUT_EXCEPTION(ExcParserError,
74  << "Parsing in " << EI_BParserMsg::val << " in the formula: " << EI_Formula::qval << "\n");
75 
76  // Temporary exception of FParser. TODO remove at the same time as FParser
77  TYPEDEF_ERR_INFO(EI_FParserMsg, std::string);
78  TYPEDEF_ERR_INFO(EI_Row, unsigned int);
79  TYPEDEF_ERR_INFO(EI_Col, unsigned int);
80  DECLARE_INPUT_EXCEPTION(ExcFParserError,
81  << "ParserError: " << EI_FParserMsg::val << "\n in the FieldFormula[" << EI_Row::val
82  << "][" << EI_Row::val << "] == " << EI_Formula::qval << " \n");
83 
84  FieldFormula(unsigned int n_comp=0);
85 
86 
87  static const Input::Type::Record & get_input_type();
88 
89  virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData& init_data);
90 
91  /**
92  * For time dependent formulas returns always true. For time independent formulas returns true only for the first time.
93  */
94  bool set_time(const TimeStep &time) override;
95 
96  /**
97  * Create SurfaceDepth object if surface region is set.
98  *
99  * See also description of the FieldBase<...>::set_mesh.
100  */
101  void set_mesh(const Mesh *mesh) override;
102 
103  void cache_update(FieldValueCache<typename Value::element_type> &data_cache,
104  ElementCacheMap &cache_map, unsigned int region_patch_idx) override;
105 
106  /**
107  * Set reference of FieldSet.
108  */
109  std::vector<const FieldCommon *> set_dependency(FieldSet &field_set) override;
110 
111  /**
112  * Overload @p FieldAlgorithmBase::cache_reinit
113  *
114  * Reinit arena data member.
115  */
116  void cache_reinit(const ElementCacheMap &cache_map) override;
117 
118  virtual ~FieldFormula();
119 
120 private:
121  /**
122  * Evaluate depth variable if it is contained in formula.
123  *
124  * Return arma vec of point coordinates extended by depth value (or zero if depth is not contained.
125  */
126  inline arma::vec eval_depth_var(const Point &p);
127 
128  // formula expression, string is set to BParser
129  std::string formula_;
130 
131  // Parser evaluating scalar, vector or tensor expression returned by formula_ string.
132  bparser::Parser b_parser_;
133 
134  /// Accessor to Input::Record
136 
137  /// Surface depth object calculate distance from surface.
138  std::shared_ptr<SurfaceDepth> surface_depth_;
139 
140  /// Flag indicates if depth variable 'd' is used in formula - obsolete parameter of FParser
142 
143  /// Flag indicates if time variable 't' is used in formula - parameter of BParser
144  bool has_time_;
145 
146  /// Helper variable for construct of arena, holds sum of sizes (over shape) of all dependent fields.
148 
149  /// Arena object providing data arrays
150  bparser::ArenaAlloc * arena_alloc_;
151 
152  // BParser data arrays and variables
153  double *X_; ///< Coordinates vector
154  double *x_; ///< Coordinates x, part of previous array
155  double *y_; ///< Coordinates y, part of previous array
156  double *z_; ///< Coordinates z, part of previous array
157  double *d_; ///< Surface depth variable, used optionally if 'd' variable is set
158  double *res_; ///< Result vector of BParser
159  uint *subsets_; ///< Subsets indices in range 0 ... n-1
161 
162  /**
163  * Data of fields evaluated in expressions.
164  *
165  * Temporary data member, we need to copy data from FieldValueCaches to arrays allocated in arena.
166  */
167  std::unordered_map<const FieldCommon *, double *> eval_field_data_;
168 
169  /// Registrar of class to factory
170  static const int registrar;
171 
172 
173 };
174 
175 
176 #endif /* FIELD_FORMULA_HH_ */
Directing class of FieldValueCache.
Space< spacedim >::Point Point
TYPEDEF_ERR_INFO(EI_Field, std::string)
double * y_
Coordinates y, part of previous array.
std::vector< const FieldCommon * > required_fields_
double * x_
Coordinates x, part of previous array.
TYPEDEF_ERR_INFO(EI_Formula, std::string)
uint * subsets_
Subsets indices in range 0 ... n-1.
TYPEDEF_ERR_INFO(EI_Col, unsigned int)
DECLARE_INPUT_EXCEPTION(ExcNotDoubleField,<< "Can not use integer valued field "<< EI_Field::qval<< " in the formula: \n")
TYPEDEF_ERR_INFO(EI_FParserMsg, std::string)
double * res_
Result vector of BParser.
std::unordered_map< const FieldCommon *, double * > eval_field_data_
double * d_
Surface depth variable, used optionally if 'd' variable is set.
bool has_depth_var_
Flag indicates if depth variable 'd' is used in formula - obsolete parameter of FParser.
Input::Record in_rec_
Accessor to Input::Record.
std::string formula_
double * z_
Coordinates z, part of previous array.
std::shared_ptr< SurfaceDepth > surface_depth_
Surface depth object calculate distance from surface.
bool has_time_
Flag indicates if time variable 't' is used in formula - parameter of BParser.
bparser::ArenaAlloc * arena_alloc_
Arena object providing data arrays.
TYPEDEF_ERR_INFO(EI_Row, unsigned int)
FieldAlgorithmBase< spacedim, Value >::Point Point
bparser::Parser b_parser_
uint sum_shape_sizes_
Helper variable for construct of arena, holds sum of sizes (over shape) of all dependent fields.
DECLARE_INPUT_EXCEPTION(ExcFParserError,<< "ParserError: "<< EI_FParserMsg::val<< "\n in the FieldFormula["<< EI_Row::val<< "]["<< EI_Row::val<< "] == "<< EI_Formula::qval<< " \n")
double * X_
Coordinates vector.
static const int registrar
Registrar of class to factory.
TYPEDEF_ERR_INFO(EI_BParserMsg, std::string)
DECLARE_INPUT_EXCEPTION(ExcParserError,<< "Parsing in "<< EI_BParserMsg::val<< " in the formula: "<< EI_Formula::qval<< "\n")
FieldAlgorithmBase< spacedim, Value > FactoryBaseType
Container for various descendants of FieldCommonBase.
Definition: field_set.hh:159
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
Record type proxy class.
Definition: type_record.hh:182
Definition: mesh.h:362
Representation of one time step..
unsigned int uint
ArmaVec< double, N > vec
Definition: armor.hh:885
Helper struct stores data for initizalize descentants of FieldAlgorithmBase.
Basic time management class.