Flow123d  JS_before_hm-978-gfd4e3d8
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 #include <stdio.h> // for sprintf
23 #include <string> // for operator==, string
24 #include <vector> // for vector
25 #include <memory>
26 #include <armadillo>
27 #include "fields/field_algo_base.hh" // for FieldAlgorithmBase
28 #include "fields/field_values.hh" // for FieldValue<>::Enum, FieldValu...
29 #include "input/accessors.hh" // for ExcAccessorForNullStorage
30 #include "input/accessors_impl.hh" // for Record::val
31 #include "input/storage.hh" // for ExcStorageTypeMismatch
32 #include "input/type_record.hh" // for Record::ExcRecordKeyNotFound
33 #include "system/exceptions.hh" // for ExcAssertMsg::~ExcAssertMsg
34 #include "tools/time_governor.hh" // for TimeStep
35 
36 class FunctionParser;
37 template <int spacedim> class ElementAccessor;
38 class SurfaceDepth;
39 
40 using namespace std;
41 
42 /**
43  * Class representing fields given by runtime parsed formulas.
44  *
45  * Using library:
46  * http://warp.povusers.org/FunctionParser/
47  *
48  * TODO:
49  * correct support for discrete functions (use integer parser), actually we just convert double to int
50  *
51  */
52 template <int spacedim, class Value>
53 class FieldFormula : public FieldAlgorithmBase<spacedim, Value>
54 {
55 public:
58 
59  FieldFormula(unsigned int n_comp=0);
60 
61 
62  static const Input::Type::Record & get_input_type();
63 
64  virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData& init_data);
65 
66  /**
67  * For time dependent formulas returns always true. For time independent formulas returns true only for the first time.
68  */
69  bool set_time(const TimeStep &time) override;
70 
71  /**
72  * Create SurfaceDepth object if surface region is set.
73  *
74  * See also description of the FieldBase<...>::set_mesh.
75  */
76  void set_mesh(const Mesh *mesh, bool boundary_domain) override;
77 
78  /**
79  * Returns one value in one given point. ResultType can be used to avoid some costly calculation if the result is trivial.
80  */
81  virtual typename Value::return_type const &value(const Point &p, const ElementAccessor<spacedim> &elm);
82 
83  /**
84  * Returns std::vector of scalar values in several points at once.
85  */
86  virtual void value_list (const Armor::array &point_list, const ElementAccessor<spacedim> &elm,
88 
89 
90  virtual ~FieldFormula();
91 
92 private:
94 
95  /**
96  * Evaluate depth variable if it is contained in formula.
97  *
98  * Return arma vec of point coordinates extended by depth value (or zero if depth is not contained.
99  */
100  inline arma::vec eval_depth_var(const Point &p);
101 
102  // StringValue::return_type == StringTensor, which behaves like arma::mat<string>
104 
105  // Matrix of parsers corresponding to the formula matrix returned by formula_matrix_helper_
107 
108  /// Accessor to Input::Record
110 
111  /// Surface depth object calculate distance from surface.
112  std::shared_ptr<SurfaceDepth> surface_depth_;
113 
114  /// Flag indicates if depth variable 'd' is used in formula
116 
117  /// Flag indicates first call of set_time method, when FunctionParsers in parser_matrix_ must be initialized
119 
120  /// Registrar of class to factory
121  static const int registrar;
122 
123 
124 };
125 
126 
127 
128 #endif /* FIELD_FORMULA_HH_ */
ArmaVec< double, N > vec
Definition: armor.hh:861
Definition: mesh.h:78
std::shared_ptr< SurfaceDepth > surface_depth_
Surface depth object calculate distance from surface.
Helper struct stores data for initizalize descentants of FieldAlgorithmBase.
StringTensorInput< Value::NRows_, Value::NCols_ > STI
FieldAlgorithmBase< spacedim, Value >::Point Point
bool first_time_set_
Flag indicates first call of set_time method, when FunctionParsers in parser_matrix_ must be initiali...
StringTensor formula_matrix_
Basic time management class.
static constexpr bool value
Definition: json.hpp:87
bool has_depth_var_
Flag indicates if depth variable &#39;d&#39; is used in formula.
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
Space< spacedim >::Point Point
static const int registrar
Registrar of class to factory.
Input::Record in_rec_
Accessor to Input::Record.
FieldAlgorithmBase< spacedim, Value > FactoryBaseType
std::vector< std::vector< FunctionParser > > parser_matrix_
Record type proxy class.
Definition: type_record.hh:182
Representation of one time step..