14#include <range/v3/range/conversion.hpp>
15#include <range/v3/view/transform.hpp>
18#include <unordered_set>
30 : exprtk::ifunction<double>(1),
_curve(curve)
32 exprtk::disable_has_side_effects(*
this);
44 exprtk::symbol_table<T>& symbol_table,
45 std::vector<std::string>
const& string_expressions)
47 exprtk::parser<T> parser;
49 std::vector<exprtk::expression<T>> expressions(string_expressions.size());
50 for (
unsigned i = 0; i < string_expressions.size(); ++i)
52 expressions[i].register_symbol_table(symbol_table);
53 if (!parser.compile(string_expressions[i], expressions[i]))
55 OGS_FATAL(
"Error: {:s}\tExpression: {:s}\n",
57 string_expressions[i]);
64class Function::Implementation
71 std::vector<std::string>
const& variables,
72 std::vector<std::string>
const& value_string_expressions,
73 std::vector<std::pair<std::string, std::vector<std::string>>>
const&
74 dvalue_string_expressions,
76 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>>
const&
83 std::vector<std::string>
const& variables,
85 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>>
const&
95 std::vector<std::pair<Variable, std::vector<Expression>>>
109 std::vector<std::string>
const& variables,
110 std::map<std::string,
111 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>>
const&
114 exprtk::symbol_table<double> symbol_table;
115 symbol_table.add_constants();
117 std::unordered_set<std::string> curve_names;
118 for (
auto const& curve : curves)
120 curve_names.insert(curve.first);
123 std::unordered_set<std::string> used_curves;
125 symbol_table.create_variable(
"t");
127 for (
auto const& v : variables)
133 else if (v ==
"x" || v ==
"y" || v ==
"z")
135 symbol_table.create_variable(
"x");
136 symbol_table.create_variable(
"y");
137 symbol_table.create_variable(
"z");
138 spatial_position_is_required =
true;
140 else if (curve_names.contains(v))
142 used_curves.insert(v);
146 auto add_scalar = [&v, &symbol_table](
double& value)
147 { symbol_table.add_variable(v, value); };
150 [&v, &symbol_table](
double* ptr, std::size_t
const size)
151 { symbol_table.add_vector(v, ptr, size); };
155 { add_scalar(*address); },
158 auto constexpr size =
160 auto& result = address->template emplace<
161 Eigen::Matrix<double, size, 1>>();
162 add_vector(result.data(), size);
167 auto& result = address->template emplace<
168 Eigen::Matrix<double, size, 1>>();
169 add_vector(result.data(), size);
173 variable_array.visitVariable(add_any_variable, variable);
177 for (
const auto&
name : used_curves)
179 const auto& curve_ptr = curves.at(
name);
182 for (
auto& [
name, wrapper] : _curve_wrappers)
184 symbol_table.add_function(
name, wrapper);
191 std::vector<std::string>
const& variables,
192 std::vector<std::string>
const& value_string_expressions,
193 std::vector<std::pair<std::string, std::vector<std::string>>>
const&
194 dvalue_string_expressions,
195 std::map<std::string,
196 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>>
const&
199 auto symbol_table = createSymbolTable(variables, curves);
206 for (
auto const& [variable_name, string_expressions] :
207 dvalue_string_expressions)
209 dvalue_expressions.emplace_back(
219 for (
auto const& variable : variables)
224 double const value = *std::get<VariableArray::Scalar const*>(
227 if (std::isnan(value))
230 "Function property: Scalar variable '{:s}' is not "
237 auto assign_kelvin_vector = [&variable, &new_variable_array](
240 auto assign_value = [&destination = *address,
241 &variable]<
typename S>(S
const& source)
243 if constexpr (std::is_same_v<S, std::monostate>)
246 "Function property: Kelvin vector variable '{:s}' is "
252 if (std::holds_alternative<S>(destination))
260 "Function property: Mismatch of Kelvin vector "
261 "sizes for variable {:s}.",
267 std::visit(assign_value,
268 *std::get<VariableArray::KelvinVector const*>(
271 auto assign_deformation_gradient =
275 auto assign_value = [&destination = *address,
276 &variable]<
typename S>(S
const& source)
278 if constexpr (std::is_same_v<S, std::monostate>)
281 "Function property: Vectorized tensor variable '{:s}' "
282 "is not initialized.",
287 if (std::holds_alternative<S>(destination))
289 std::get<S>(destination) = source;
294 "Function property: Mismatch of vectorized tensor "
295 "sizes for variable {:s}.",
301 std::visit(assign_value,
302 *std::get<VariableArray::DeformationGradient const*>(
308 assign_deformation_gradient},
314 std::vector<Variable>
const& variables,
317 std::vector<exprtk::expression<double>>
const& expressions,
319 bool const spatial_position_is_required)
321 std::vector<double> result(expressions.size());
324 std::lock_guard lock_guard(mutex);
329 begin(expressions), end(expressions), begin(result),
330 [t, &pos, spatial_position_is_required](
auto const& e)
332 auto& symbol_table = e.get_symbol_table();
333 symbol_table.get_variable(
"t")->ref() = t;
335 if (spatial_position_is_required)
340 "FunctionParameter: The spatial position "
345 symbol_table.get_variable(
"x")->ref() = coords[0];
346 symbol_table.get_variable(
"y")->ref() = coords[1];
347 symbol_table.get_variable(
"z")->ref() = coords[2];
354 switch (result.size())
362 return Eigen::Vector2d{result[0], result[1]};
366 return Eigen::Vector3d{result[0], result[1], result[2]};
370 Eigen::Matrix<double, 2, 2> m;
371 m = Eigen::Map<Eigen::Matrix<double, 2, 2>
const>(result.data(), 2,
377 Eigen::Matrix<double, 3, 3> m;
378 m = Eigen::Map<Eigen::Matrix<double, 3, 3>
const>(result.data(), 3,
383 OGS_FATAL(
"Cannot convert a vector of size {} to a PropertyDataType",
388 std::vector<std::string>
const& value_string_expressions,
389 std::vector<std::pair<std::string, std::vector<std::string>>>
const&
390 dvalue_string_expressions)
392 std::vector<std::string> variables;
394 auto collect_variables = [&](
auto string_expressions)
396 for (
auto const& string_expression : string_expressions)
398 if (!exprtk::collect_variables(string_expression, variables))
401 "Collecting variables from expression '{}' didn't work.",
407 collect_variables(value_string_expressions);
408 for (
auto const& var_name_string_expression : dvalue_string_expressions)
410 collect_variables(var_name_string_expression.second);
419 std::vector<std::string>
const& value_string_expressions,
420 std::vector<std::pair<std::string, std::vector<std::string>>>
const&
421 dvalue_string_expressions,
422 std::map<std::string,
423 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>>
const&
428 auto const variables =
432 static std::unordered_set<std::string> filter_not_variables = {
"t",
"x",
434 for (
auto const& curve : curves)
436 filter_not_variables.insert(curve.first);
441 std::views::filter([](
const std::string& s)
442 {
return !filter_not_variables.contains(s); }) |
443 std::views::transform([](std::string
const& s)
445 ranges::to<std::vector>;
447 impl2_ = std::make_unique<Implementation<2>>(
448 variables, value_string_expressions, dvalue_string_expressions, curves);
449 impl3_ = std::make_unique<Implementation<3>>(
450 variables, value_string_expressions, dvalue_string_expressions, curves);
457 if (variable_array.
is2D())
461 if (variable_array.
is3D())
467 "Variable array has vectors for 2 and 3 dimensions simultaneously. "
468 "Mixed dimensions cannot be dealt within Function evaluation.");
473 double const t,
double const )
const
479 impl_ptr->value_expressions,
480 impl_ptr->variable_array,
mutex_,
481 impl_ptr->spatial_position_is_required);
489 double const t,
double const )
const
494 auto const it = std::find_if(begin(impl_ptr->dvalue_expressions),
495 end(impl_ptr->dvalue_expressions),
496 [&variable](
auto const& v)
497 { return v.first == variable; });
499 if (it == end(impl_ptr->dvalue_expressions))
502 "Requested derivative with respect to the variable {:s} "
504 "provided for Function-type property {:s}.",
509 it->second, impl_ptr->variable_array,
511 impl_ptr->spatial_position_is_required);
Definition of the PiecewiseLinearInterpolation class.
exprtk::symbol_table< double > createSymbolTable(std::vector< std::string > const &variables, std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation > > const &curves)
std::vector< Expression > value_expressions
VariableArray variable_array
bool spatial_position_is_required
Implementation(std::vector< std::string > const &variables, std::vector< std::string > const &value_string_expressions, std::vector< std::pair< std::string, std::vector< std::string > > > const &dvalue_string_expressions, std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation > > const &curves)
std::vector< std::pair< Variable, std::vector< Expression > > > dvalue_expressions
exprtk::expression< double > Expression
std::map< std::string, CurveWrapper > _curve_wrappers
std::unique_ptr< Implementation< 3 > > impl3_
std::vector< Variable > variables_
Variables used in the exprtk expressions.
std::variant< Function::Implementation< 2 > *, Function::Implementation< 3 > * > getImplementationForDimensionOfVariableArray(VariableArray const &variable_array) const
Function(std::string name, std::vector< std::string > const &value_string_expressions, std::vector< std::pair< std::string, std::vector< std::string > > > const &dvalue_string_expressions, std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation > > const &curves)
PropertyDataType dValue(VariableArray const &variable_array, Variable const variable, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const override
std::unique_ptr< Implementation< 2 > > impl2_
virtual PropertyDataType value() const
std::variant< std::monostate, Eigen::Vector< double, 5 >, Eigen::Vector< double, 9 > > DeformationGradient
std::variant< std::monostate, Eigen::Vector< double, 4 >, Eigen::Vector< double, 6 > > KelvinVector
VariablePointerConst address_of(Variable const v) const
auto visitVariable(Visitor &&visitor, Variable const variable)
double getValue(double pnt_to_interpolate) const
Calculates the interpolation value.
std::optional< MathLib::Point3d > const getCoordinates() const
void makeVectorUnique(std::vector< T > &v)
static std::vector< exprtk::expression< T > > compileExpressions(exprtk::symbol_table< T > &symbol_table, std::vector< std::string > const &string_expressions)
static std::vector< std::string > collectVariables(std::vector< std::string > const &value_string_expressions, std::vector< std::pair< std::string, std::vector< std::string > > > const &dvalue_string_expressions)
static void updateVariableArrayValues(std::vector< Variable > const &variables, VariableArray const &new_variable_array, VariableArray &variable_array)
static PropertyDataType evaluateExpressions(std::vector< Variable > const &variables, VariableArray const &new_variable_array, ParameterLib::SpatialPosition const &pos, double const t, std::vector< exprtk::expression< double > > const &expressions, VariableArray &variable_array, std::mutex &mutex, bool const spatial_position_is_required)
static const std::array< std::string, static_cast< int >(Variable::number_of_variables)> variable_enum_to_string
Variable convertStringToVariable(std::string const &string)
std::variant< double, Eigen::Matrix< double, 2, 1 >, Eigen::Matrix< double, 3, 1 >, Eigen::Matrix< double, 2, 2 >, Eigen::Matrix< double, 3, 3 >, Eigen::Matrix< double, 4, 1 >, Eigen::Matrix< double, 6, 1 >, Eigen::MatrixXd > PropertyDataType
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
constexpr int size(int const displacement_dim)
Vectorized tensor size for given displacement dimension.
double operator()(const double &t) override
CurveWrapper(const MathLib::PiecewiseLinearInterpolation &curve)
const MathLib::PiecewiseLinearInterpolation & _curve