18#include <range/v3/range/conversion.hpp>
19#include <range/v3/view/transform.hpp>
22#include <unordered_set>
35 : exprtk::ifunction<double>(1),
_curve(curve)
37 exprtk::disable_has_side_effects(*
this);
49 exprtk::symbol_table<T>& symbol_table,
50 std::vector<std::string>
const& string_expressions)
52 exprtk::parser<T> parser;
54 std::vector<exprtk::expression<T>> expressions(string_expressions.size());
55 for (
unsigned i = 0; i < string_expressions.size(); ++i)
57 expressions[i].register_symbol_table(symbol_table);
58 if (!parser.compile(string_expressions[i], expressions[i]))
60 OGS_FATAL(
"Error: {:s}\tExpression: {:s}\n",
62 string_expressions[i]);
69class Function::Implementation
77 std::vector<std::string>
const& variables,
78 std::vector<std::string>
const& value_string_expressions,
79 std::vector<std::pair<std::string, std::vector<std::string>>>
const&
80 dvalue_string_expressions,
82 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>>
const&
89 std::vector<std::string>
const& variables,
91 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>>
const&
98 std::vector<std::string>
const& variables,
100 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>>
const&
110 std::vector<std::vector<std::pair<Variable, std::vector<Expression>>>>
124 std::vector<std::string>
const& variables,
125 std::map<std::string,
126 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>>
const&
130 exprtk::symbol_table<double> symbol_table;
131 symbol_table.add_constants();
133 std::unordered_set<std::string> curve_names;
134 for (
auto const& curve : curves)
136 curve_names.insert(curve.first);
139 std::unordered_set<std::string> used_curves;
141 symbol_table.create_variable(
"t");
143 for (
auto const& v : variables)
149 else if (v ==
"x" || v ==
"y" || v ==
"z")
151 symbol_table.create_variable(
"x");
152 symbol_table.create_variable(
"y");
153 symbol_table.create_variable(
"z");
154 spatial_position_is_required =
true;
156 else if (curve_names.contains(v))
158 used_curves.insert(v);
162 auto add_scalar = [&v, &symbol_table](
double& value)
163 { symbol_table.add_variable(v, value); };
166 [&v, &symbol_table](
double* ptr, std::size_t
const size)
167 { symbol_table.add_vector(v, ptr, size); };
171 { add_scalar(*address); },
174 auto constexpr size =
176 auto& result = address->template emplace<
177 Eigen::Matrix<double, size, 1>>();
178 add_vector(result.data(), size);
183 auto& result = address->template emplace<
184 Eigen::Matrix<double, size, 1>>();
185 add_vector(result.data(), size);
193 for (
const auto&
name : used_curves)
195 const auto& curve_ptr = curves.at(
name);
198 for (
auto& [
name, wrapper] : _curve_wrappers)
200 symbol_table.add_function(
name, wrapper);
205std::vector<exprtk::symbol_table<double>>
208 std::vector<std::string>
const& variables,
209 std::map<std::string,
210 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>>
const&
213 std::vector<exprtk::symbol_table<double>> symbol_tables(num_threads);
215 for (
int thread_id = 0; thread_id < num_threads; ++thread_id)
217 symbol_tables[thread_id] =
218 createSymbolTable(variables, curves, variable_arrays[thread_id]);
221 return symbol_tables;
227 std::vector<std::string>
const& variables,
228 std::vector<std::string>
const& value_string_expressions,
229 std::vector<std::pair<std::string, std::vector<std::string>>>
const&
230 dvalue_string_expressions,
231 std::map<std::string,
232 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>>
const&
235 variable_arrays.resize(num_threads);
236 value_expressions.resize(num_threads);
237 dvalue_expressions.resize(num_threads);
239 auto symbol_tables = createSymbolTables(num_threads, variables, curves);
242 for (
int thread_id = 0; thread_id < num_threads; ++thread_id)
245 symbol_tables[thread_id], value_string_expressions);
249 for (
int thread_id = 0; thread_id < num_threads; ++thread_id)
251 for (
auto const& [variable_name, string_expressions] :
252 dvalue_string_expressions)
254 dvalue_expressions[thread_id].emplace_back(
257 string_expressions));
266 for (
auto const& variable : variables)
271 double const value = *std::get<VariableArray::Scalar const*>(
274 if (std::isnan(value))
277 "Function property: Scalar variable '{:s}' is not "
284 auto assign_kelvin_vector = [&variable, &new_variable_array](
287 auto assign_value = [&destination = *address,
288 &variable]<
typename S>(S
const& source)
290 if constexpr (std::is_same_v<S, std::monostate>)
293 "Function property: Kelvin vector variable '{:s}' is "
299 if (std::holds_alternative<S>(destination))
307 "Function property: Mismatch of Kelvin vector "
308 "sizes for variable {:s}.",
314 std::visit(assign_value,
315 *std::get<VariableArray::KelvinVector const*>(
318 auto assign_deformation_gradient =
322 auto assign_value = [&destination = *address,
323 &variable]<
typename S>(S
const& source)
325 if constexpr (std::is_same_v<S, std::monostate>)
328 "Function property: Vectorized tensor variable '{:s}' "
329 "is not initialized.",
334 if (std::holds_alternative<S>(destination))
336 std::get<S>(destination) = source;
341 "Function property: Mismatch of vectorized tensor "
342 "sizes for variable {:s}.",
348 std::visit(assign_value,
349 *std::get<VariableArray::DeformationGradient const*>(
355 assign_deformation_gradient},
361 std::vector<Variable>
const& variables,
364 std::vector<exprtk::expression<double>>
const& expressions,
365 VariableArray& variable_array,
bool const spatial_position_is_required)
367 std::vector<double> result(expressions.size());
371 std::transform(begin(expressions), end(expressions), begin(result),
372 [t, &pos, spatial_position_is_required](
auto const& e)
374 auto& symbol_table = e.get_symbol_table();
375 symbol_table.get_variable(
"t")->ref() = t;
377 if (spatial_position_is_required)
382 "FunctionParameter: The spatial position "
387 symbol_table.get_variable(
"x")->ref() = coords[0];
388 symbol_table.get_variable(
"y")->ref() = coords[1];
389 symbol_table.get_variable(
"z")->ref() = coords[2];
395 switch (result.size())
403 return Eigen::Vector2d{result[0], result[1]};
407 return Eigen::Vector3d{result[0], result[1], result[2]};
411 Eigen::Matrix<double, 2, 2> m;
412 m = Eigen::Map<Eigen::Matrix<double, 2, 2>
const>(result.data(), 2,
418 Eigen::Matrix<double, 3, 3> m;
419 m = Eigen::Map<Eigen::Matrix<double, 3, 3>
const>(result.data(), 3,
424 OGS_FATAL(
"Cannot convert a vector of size {} to a PropertyDataType",
429 std::vector<std::string>
const& value_string_expressions,
430 std::vector<std::pair<std::string, std::vector<std::string>>>
const&
431 dvalue_string_expressions)
433 std::vector<std::string> variables;
435 auto collect_variables = [&](
auto string_expressions)
437 for (
auto const& string_expression : string_expressions)
439 if (!exprtk::collect_variables(string_expression, variables))
442 "Collecting variables from expression '{}' didn't work.",
448 collect_variables(value_string_expressions);
449 for (
auto const& var_name_string_expression : dvalue_string_expressions)
451 collect_variables(var_name_string_expression.second);
460 std::vector<std::string>
const& value_string_expressions,
461 std::vector<std::pair<std::string, std::vector<std::string>>>
const&
462 dvalue_string_expressions,
463 std::map<std::string,
464 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>>
const&
469 auto const variables =
473 static std::unordered_set<std::string> filter_not_variables = {
"t",
"x",
475 for (
auto const& curve : curves)
477 filter_not_variables.insert(curve.first);
482 std::views::filter([](
const std::string& s)
483 {
return !filter_not_variables.contains(s); }) |
484 std::views::transform([](std::string
const& s)
486 ranges::to<std::vector>;
488 auto const get_number_omp_threads = []()
491 return omp_get_max_threads();
498 get_number_omp_threads());
500 impl2_ = std::make_unique<Implementation<2>>(
501 num_threads, variables, value_string_expressions,
502 dvalue_string_expressions, curves);
503 impl3_ = std::make_unique<Implementation<3>>(
504 num_threads, variables, value_string_expressions,
505 dvalue_string_expressions, curves);
512 if (variable_array.
is2D())
516 if (variable_array.
is3D())
522 "Variable array has vectors for 2 and 3 dimensions simultaneously. "
523 "Mixed dimensions cannot be dealt within Function evaluation.");
528 double const t,
double const )
const
531 int const thread_id = omp_get_thread_num();
533 int const thread_id = 0;
539 static_cast<int>(impl_ptr->value_expressions.size()))
542 "In Function-type property '{:s}' evaluation the "
543 "OMP-thread with id {:d} exceeds the number of allocated "
545 name_, thread_id, impl_ptr->value_expressions.size());
548 impl_ptr->value_expressions[thread_id],
549 impl_ptr->variable_arrays[thread_id],
550 impl_ptr->spatial_position_is_required);
558 double const t,
double const )
const
561 int const thread_id = omp_get_thread_num();
563 int const thread_id = 0;
569 static_cast<int>(impl_ptr->dvalue_expressions.size()))
572 "In Function-type property '{:s}' evaluation the "
573 "OMP-thread with id {:d} exceeds the number of allocated "
575 name_, thread_id, impl_ptr->value_expressions.size());
577 auto const it = std::find_if(
578 begin(impl_ptr->dvalue_expressions[thread_id]),
579 end(impl_ptr->dvalue_expressions[thread_id]),
580 [&variable](
auto const& v) { return v.first == variable; });
582 if (it == end(impl_ptr->dvalue_expressions[thread_id]))
585 "Requested derivative with respect to the variable {:s} "
587 "provided for Function-type property {:s}.",
593 impl_ptr->variable_arrays[thread_id],
594 impl_ptr->spatial_position_is_required);
Definition of the PiecewiseLinearInterpolation class.
std::vector< VariableArray > variable_arrays
bool spatial_position_is_required
std::vector< exprtk::symbol_table< double > > createSymbolTables(int num_threads, std::vector< std::string > const &variables, std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation > > const &curves)
Create symbol tables for all threads.
Implementation(int num_threads, 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)
exprtk::expression< double > Expression
std::map< std::string, CurveWrapper > _curve_wrappers
std::vector< std::vector< std::pair< Variable, std::vector< Expression > > > > dvalue_expressions
exprtk::symbol_table< double > createSymbolTable(std::vector< std::string > const &variables, std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation > > const &curves, VariableArray &variable_array)
std::vector< std::vector< Expression > > value_expressions
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
int getNumberOfAssemblyThreads()
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 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, bool const spatial_position_is_required)
static void updateVariableArrayValues(std::vector< Variable > const &variables, VariableArray const &new_variable_array, VariableArray &variable_array)
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