Loading [MathJax]/extensions/tex2jax.js
OGS
Function.cpp
Go to the documentation of this file.
1
11
12#include <exprtk.hpp>
13#include <numeric>
14#include <range/v3/range/conversion.hpp>
15#include <range/v3/view/transform.hpp>
16#include <ranges>
17#include <typeinfo>
18#include <unordered_set>
19
20#include "BaseLib/Algorithm.h"
24
25namespace MaterialPropertyLib
26{
27struct CurveWrapper : public exprtk::ifunction<double>
28{
30 : exprtk::ifunction<double>(1), _curve(curve)
31 {
32 exprtk::disable_has_side_effects(*this);
33 }
34
35 double operator()(const double& t) override { return _curve.getValue(t); }
36
37private:
39};
40// Passing symbol table by reference as required by register_symbol_table()
41// call.
42template <typename T>
43static std::vector<exprtk::expression<T>> compileExpressions(
44 exprtk::symbol_table<T>& symbol_table,
45 std::vector<std::string> const& string_expressions)
46{
47 exprtk::parser<T> parser;
48
49 std::vector<exprtk::expression<T>> expressions(string_expressions.size());
50 for (unsigned i = 0; i < string_expressions.size(); ++i)
51 {
52 expressions[i].register_symbol_table(symbol_table);
53 if (!parser.compile(string_expressions[i], expressions[i]))
54 {
55 OGS_FATAL("Error: {:s}\tExpression: {:s}\n",
56 parser.error(),
57 string_expressions[i]);
58 }
59 }
60 return expressions;
61}
62
63template <int D>
64class Function::Implementation
65{
66public:
67 using Expression = exprtk::expression<double>;
68
69public:
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,
75 std::map<std::string,
76 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
77 curves);
78
79private:
82 exprtk::symbol_table<double> createSymbolTable(
83 std::vector<std::string> const& variables,
84 std::map<std::string,
85 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
86 curves);
87
88public:
91 std::vector<Expression> value_expressions;
92
95 std::vector<std::pair<Variable, std::vector<Expression>>>
97
101
102 std::map<std::string, CurveWrapper> _curve_wrappers;
103
105};
106
107template <int D>
109 std::vector<std::string> const& variables,
110 std::map<std::string,
111 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
112 curves)
113{
114 exprtk::symbol_table<double> symbol_table;
115 symbol_table.add_constants();
116
117 std::unordered_set<std::string> curve_names;
118 for (auto const& curve : curves)
119 {
120 curve_names.insert(curve.first);
121 }
122
123 std::unordered_set<std::string> used_curves;
124
125 symbol_table.create_variable("t");
126
127 for (auto const& v : variables)
128 {
129 if (v == "t")
130 {
131 continue;
132 }
133 else if (v == "x" || v == "y" || v == "z")
134 {
135 symbol_table.create_variable("x");
136 symbol_table.create_variable("y");
137 symbol_table.create_variable("z");
138 spatial_position_is_required = true;
139 }
140 else if (curve_names.contains(v))
141 {
142 used_curves.insert(v);
143 }
144 else
145 {
146 auto add_scalar = [&v, &symbol_table](double& value)
147 { symbol_table.add_variable(v, value); };
148
149 auto add_vector =
150 [&v, &symbol_table](double* ptr, std::size_t const size)
151 { symbol_table.add_vector(v, ptr, size); };
152
153 auto add_any_variable = BaseLib::Overloaded{
154 [&add_scalar](VariableArray::Scalar* address)
155 { add_scalar(*address); },
156 [&add_vector](VariableArray::KelvinVector* address)
157 {
158 auto constexpr size =
160 auto& result = address->template emplace<
161 Eigen::Matrix<double, size, 1>>();
162 add_vector(result.data(), size);
163 },
164 [&add_vector](VariableArray::DeformationGradient* address)
165 {
166 auto constexpr size = MathLib::VectorizedTensor::size(D);
167 auto& result = address->template emplace<
168 Eigen::Matrix<double, size, 1>>();
169 add_vector(result.data(), size);
170 }};
171
172 Variable const variable = convertStringToVariable(v);
173 variable_array.visitVariable(add_any_variable, variable);
174 }
175 }
176
177 for (const auto& name : used_curves)
178 {
179 const auto& curve_ptr = curves.at(name);
180 _curve_wrappers.emplace(name, CurveWrapper(*curve_ptr));
181 }
182 for (auto& [name, wrapper] : _curve_wrappers)
183 {
184 symbol_table.add_function(name, wrapper);
185 }
186 return symbol_table;
187}
188
189template <int D>
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&
197 curves)
198{
199 auto symbol_table = createSymbolTable(variables, curves);
200
201 // value expressions.
202 value_expressions =
203 compileExpressions(symbol_table, value_string_expressions);
204
205 // dValue expressions.
206 for (auto const& [variable_name, string_expressions] :
207 dvalue_string_expressions)
208 {
209 dvalue_expressions.emplace_back(
210 convertStringToVariable(variable_name),
211 compileExpressions(symbol_table, string_expressions));
212 }
213}
214
215static void updateVariableArrayValues(std::vector<Variable> const& variables,
216 VariableArray const& new_variable_array,
217 VariableArray& variable_array)
218{
219 for (auto const& variable : variables)
220 {
221 auto assign_scalar =
222 [&variable, &new_variable_array](VariableArray::Scalar* address)
223 {
224 double const value = *std::get<VariableArray::Scalar const*>(
225 new_variable_array.address_of(variable));
226
227 if (std::isnan(value))
228 {
229 OGS_FATAL(
230 "Function property: Scalar variable '{:s}' is not "
231 "initialized.",
232 variable_enum_to_string[static_cast<int>(variable)]);
233 }
234
235 *address = value;
236 };
237 auto assign_kelvin_vector = [&variable, &new_variable_array](
239 {
240 auto assign_value = [&destination = *address,
241 &variable]<typename S>(S const& source)
242 {
243 if constexpr (std::is_same_v<S, std::monostate>)
244 {
245 OGS_FATAL(
246 "Function property: Kelvin vector variable '{:s}' is "
247 "not initialized.",
248 variable_enum_to_string[static_cast<int>(variable)]);
249 }
250 else
251 {
252 if (std::holds_alternative<S>(destination))
253 {
254 std::get<S>(destination) = MathLib::KelvinVector::
256 }
257 else
258 {
259 OGS_FATAL(
260 "Function property: Mismatch of Kelvin vector "
261 "sizes for variable {:s}.",
262 variable_enum_to_string[static_cast<int>(
263 variable)]);
264 }
265 }
266 };
267 std::visit(assign_value,
268 *std::get<VariableArray::KelvinVector const*>(
269 new_variable_array.address_of(variable)));
270 };
271 auto assign_deformation_gradient =
272 [&variable,
273 &new_variable_array](VariableArray::DeformationGradient* address)
274 {
275 auto assign_value = [&destination = *address,
276 &variable]<typename S>(S const& source)
277 {
278 if constexpr (std::is_same_v<S, std::monostate>)
279 {
280 OGS_FATAL(
281 "Function property: Vectorized tensor variable '{:s}' "
282 "is not initialized.",
283 variable_enum_to_string[static_cast<int>(variable)]);
284 }
285 else
286 {
287 if (std::holds_alternative<S>(destination))
288 {
289 std::get<S>(destination) = source;
290 }
291 else
292 {
293 OGS_FATAL(
294 "Function property: Mismatch of vectorized tensor "
295 "sizes for variable {:s}.",
296 variable_enum_to_string[static_cast<int>(
297 variable)]);
298 }
299 }
300 };
301 std::visit(assign_value,
302 *std::get<VariableArray::DeformationGradient const*>(
303 new_variable_array.address_of(variable)));
304 };
305
306 variable_array.visitVariable(
307 BaseLib::Overloaded{assign_scalar, assign_kelvin_vector,
308 assign_deformation_gradient},
309 variable);
310 }
311}
312
314 std::vector<Variable> const& variables,
315 VariableArray const& new_variable_array,
316 ParameterLib::SpatialPosition const& pos, double const t,
317 std::vector<exprtk::expression<double>> const& expressions,
318 VariableArray& variable_array, std::mutex& mutex,
319 bool const spatial_position_is_required)
320{
321 std::vector<double> result(expressions.size());
322
323 {
324 std::lock_guard lock_guard(mutex);
325 updateVariableArrayValues(variables, new_variable_array,
326 variable_array);
327
328 std::transform(
329 begin(expressions), end(expressions), begin(result),
330 [t, &pos, spatial_position_is_required](auto const& e)
331 {
332 auto& symbol_table = e.get_symbol_table();
333 symbol_table.get_variable("t")->ref() = t;
334
335 if (spatial_position_is_required)
336 {
337 if (!pos.getCoordinates())
338 {
339 OGS_FATAL(
340 "FunctionParameter: The spatial position "
341 "has to be set by "
342 "coordinates.");
343 }
344 auto const coords = pos.getCoordinates().value();
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];
348 }
349
350 return e.value();
351 });
352 }
353
354 switch (result.size())
355 {
356 case 1:
357 {
358 return result[0];
359 }
360 case 2:
361 {
362 return Eigen::Vector2d{result[0], result[1]};
363 }
364 case 3:
365 {
366 return Eigen::Vector3d{result[0], result[1], result[2]};
367 }
368 case 4:
369 {
370 Eigen::Matrix<double, 2, 2> m;
371 m = Eigen::Map<Eigen::Matrix<double, 2, 2> const>(result.data(), 2,
372 2);
373 return m;
374 }
375 case 9:
376 {
377 Eigen::Matrix<double, 3, 3> m;
378 m = Eigen::Map<Eigen::Matrix<double, 3, 3> const>(result.data(), 3,
379 3);
380 return m;
381 }
382 }
383 OGS_FATAL("Cannot convert a vector of size {} to a PropertyDataType",
384 result.size());
385}
386
387static std::vector<std::string> collectVariables(
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)
391{
392 std::vector<std::string> variables;
393
394 auto collect_variables = [&](auto string_expressions)
395 {
396 for (auto const& string_expression : string_expressions)
397 {
398 if (!exprtk::collect_variables(string_expression, variables))
399 {
400 OGS_FATAL(
401 "Collecting variables from expression '{}' didn't work.",
402 string_expression);
403 }
404 }
405 };
406
407 collect_variables(value_string_expressions);
408 for (auto const& var_name_string_expression : dvalue_string_expressions)
409 {
410 collect_variables(var_name_string_expression.second);
411 }
412
413 BaseLib::makeVectorUnique(variables);
414 return variables;
415}
416
418 std::string name,
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&
424 curves)
425{
426 name_ = std::move(name);
427
428 auto const variables =
429 collectVariables(value_string_expressions, dvalue_string_expressions);
430
431 // filter the strings unrelated to time or spatial position
432 static std::unordered_set<std::string> filter_not_variables = {"t", "x",
433 "y", "z"};
434 for (auto const& curve : curves)
435 {
436 filter_not_variables.insert(curve.first);
437 }
438
439 variables_ =
440 variables |
441 std::views::filter([](const std::string& s)
442 { return !filter_not_variables.contains(s); }) |
443 std::views::transform([](std::string const& s)
444 { return convertStringToVariable(s); }) |
445 ranges::to<std::vector>;
446
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);
451}
452
453std::variant<Function::Implementation<2>*, Function::Implementation<3>*>
455 VariableArray const& variable_array) const
456{
457 if (variable_array.is2D())
458 {
459 return impl2_.get();
460 }
461 if (variable_array.is3D())
462 {
463 return impl3_.get();
464 }
465
466 OGS_FATAL(
467 "Variable array has vectors for 2 and 3 dimensions simultaneously. "
468 "Mixed dimensions cannot be dealt within Function evaluation.");
469}
470
473 double const t, double const /*dt*/) const
474{
475 return std::visit(
476 [&](auto&& impl_ptr)
477 {
478 return evaluateExpressions(variables_, variable_array, pos, t,
479 impl_ptr->value_expressions,
480 impl_ptr->variable_array, mutex_,
481 impl_ptr->spatial_position_is_required);
482 },
484}
485
487 Variable const variable,
489 double const t, double const /*dt*/) const
490{
491 return std::visit(
492 [&](auto&& impl_ptr)
493 {
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; });
498
499 if (it == end(impl_ptr->dvalue_expressions))
500 {
501 OGS_FATAL(
502 "Requested derivative with respect to the variable {:s} "
503 "not "
504 "provided for Function-type property {:s}.",
505 variable_enum_to_string[static_cast<int>(variable)], name_);
506 }
507
508 return evaluateExpressions(variables_, variable_array, pos, t,
509 it->second, impl_ptr->variable_array,
510 mutex_,
511 impl_ptr->spatial_position_is_required);
512 },
514}
515
516Function::~Function() = default;
517
518} // namespace MaterialPropertyLib
#define OGS_FATAL(...)
Definition Error.h:26
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)
Definition Function.cpp:108
std::vector< Expression > value_expressions
Definition Function.cpp:91
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)
Definition Function.cpp:190
std::vector< std::pair< Variable, std::vector< Expression > > > dvalue_expressions
Definition Function.cpp:96
exprtk::expression< double > Expression
Definition Function.cpp:67
std::map< std::string, CurveWrapper > _curve_wrappers
Definition Function.cpp:102
std::unique_ptr< Implementation< 3 > > impl3_
Definition Function.h:60
std::vector< Variable > variables_
Variables used in the exprtk expressions.
Definition Function.h:67
std::variant< Function::Implementation< 2 > *, Function::Implementation< 3 > * > getImplementationForDimensionOfVariableArray(VariableArray const &variable_array) const
Definition Function.cpp:454
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)
Definition Function.cpp:417
PropertyDataType dValue(VariableArray const &variable_array, Variable const variable, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const override
Definition Function.cpp:486
std::unique_ptr< Implementation< 2 > > impl2_
Definition Function.h:59
virtual PropertyDataType value() const
Definition Property.cpp:76
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)
Definition Algorithm.h:180
static std::vector< exprtk::expression< T > > compileExpressions(exprtk::symbol_table< T > &symbol_table, std::vector< std::string > const &string_expressions)
Definition Function.cpp:43
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)
Definition Function.cpp:387
static void updateVariableArrayValues(std::vector< Variable > const &variables, VariableArray const &new_variable_array, VariableArray &variable_array)
Definition Function.cpp:215
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)
Definition Function.cpp:313
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
Definition Property.h:31
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
Definition Function.cpp:35
CurveWrapper(const MathLib::PiecewiseLinearInterpolation &curve)
Definition Function.cpp:29
const MathLib::PiecewiseLinearInterpolation & _curve
Definition Function.cpp:38