OGS
Function.cpp
Go to the documentation of this file.
1
11
12#ifdef _OPENMP
13#include <omp.h>
14#endif
15
16#include <exprtk.hpp>
17#include <numeric>
18#include <range/v3/range/conversion.hpp>
19#include <range/v3/view/transform.hpp>
20#include <ranges>
21#include <typeinfo>
22#include <unordered_set>
23
24#include "BaseLib/Algorithm.h"
29
30namespace MaterialPropertyLib
31{
32struct CurveWrapper : public exprtk::ifunction<double>
33{
35 : exprtk::ifunction<double>(1), _curve(curve)
36 {
37 exprtk::disable_has_side_effects(*this);
38 }
39
40 double operator()(const double& t) override { return _curve.getValue(t); }
41
42private:
44};
45// Passing symbol table by reference as required by register_symbol_table()
46// call.
47template <typename T>
48static std::vector<exprtk::expression<T>> compileExpressions(
49 exprtk::symbol_table<T>& symbol_table,
50 std::vector<std::string> const& string_expressions)
51{
52 exprtk::parser<T> parser;
53
54 std::vector<exprtk::expression<T>> expressions(string_expressions.size());
55 for (unsigned i = 0; i < string_expressions.size(); ++i)
56 {
57 expressions[i].register_symbol_table(symbol_table);
58 if (!parser.compile(string_expressions[i], expressions[i]))
59 {
60 OGS_FATAL("Error: {:s}\tExpression: {:s}\n",
61 parser.error(),
62 string_expressions[i]);
63 }
64 }
65 return expressions;
66}
67
68template <int D>
69class Function::Implementation
70{
71public:
72 using Expression = exprtk::expression<double>;
73
74public:
76 int num_threads,
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,
81 std::map<std::string,
82 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
83 curves);
84
85private:
88 exprtk::symbol_table<double> createSymbolTable(
89 std::vector<std::string> const& variables,
90 std::map<std::string,
91 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
92 curves,
93 VariableArray& variable_array);
94
96 std::vector<exprtk::symbol_table<double>> createSymbolTables(
97 int num_threads,
98 std::vector<std::string> const& variables,
99 std::map<std::string,
100 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
101 curves);
102
103public:
106 std::vector<std::vector<Expression>> value_expressions;
107
110 std::vector<std::vector<std::pair<Variable, std::vector<Expression>>>>
112
115 std::vector<VariableArray> variable_arrays;
116
117 std::map<std::string, CurveWrapper> _curve_wrappers;
118
120};
121
122template <int D>
124 std::vector<std::string> const& variables,
125 std::map<std::string,
126 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
127 curves,
128 VariableArray& variable_array)
129{
130 exprtk::symbol_table<double> symbol_table;
131 symbol_table.add_constants();
132
133 std::unordered_set<std::string> curve_names;
134 for (auto const& curve : curves)
135 {
136 curve_names.insert(curve.first);
137 }
138
139 std::unordered_set<std::string> used_curves;
140
141 symbol_table.create_variable("t");
142
143 for (auto const& v : variables)
144 {
145 if (v == "t")
146 {
147 continue;
148 }
149 else if (v == "x" || v == "y" || v == "z")
150 {
151 symbol_table.create_variable("x");
152 symbol_table.create_variable("y");
153 symbol_table.create_variable("z");
154 spatial_position_is_required = true;
155 }
156 else if (curve_names.contains(v))
157 {
158 used_curves.insert(v);
159 }
160 else
161 {
162 auto add_scalar = [&v, &symbol_table](double& value)
163 { symbol_table.add_variable(v, value); };
164
165 auto add_vector =
166 [&v, &symbol_table](double* ptr, std::size_t const size)
167 { symbol_table.add_vector(v, ptr, size); };
168
169 auto add_any_variable = BaseLib::Overloaded{
170 [&add_scalar](VariableArray::Scalar* address)
171 { add_scalar(*address); },
172 [&add_vector](VariableArray::KelvinVector* address)
173 {
174 auto constexpr size =
176 auto& result = address->template emplace<
177 Eigen::Matrix<double, size, 1>>();
178 add_vector(result.data(), size);
179 },
180 [&add_vector](VariableArray::DeformationGradient* address)
181 {
182 auto constexpr size = MathLib::VectorizedTensor::size(D);
183 auto& result = address->template emplace<
184 Eigen::Matrix<double, size, 1>>();
185 add_vector(result.data(), size);
186 }};
187
188 Variable const variable = convertStringToVariable(v);
189 variable_array.visitVariable(add_any_variable, variable);
190 }
191 }
192
193 for (const auto& name : used_curves)
194 {
195 const auto& curve_ptr = curves.at(name);
196 _curve_wrappers.emplace(name, CurveWrapper(*curve_ptr));
197 }
198 for (auto& [name, wrapper] : _curve_wrappers)
199 {
200 symbol_table.add_function(name, wrapper);
201 }
202 return symbol_table;
203}
204template <int D>
205std::vector<exprtk::symbol_table<double>>
207 int num_threads,
208 std::vector<std::string> const& variables,
209 std::map<std::string,
210 std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
211 curves)
212{
213 std::vector<exprtk::symbol_table<double>> symbol_tables(num_threads);
214
215 for (int thread_id = 0; thread_id < num_threads; ++thread_id)
216 {
217 symbol_tables[thread_id] =
218 createSymbolTable(variables, curves, variable_arrays[thread_id]);
219 }
220
221 return symbol_tables;
222}
223
224template <int D>
226 int num_threads,
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&
233 curves)
234{
235 variable_arrays.resize(num_threads);
236 value_expressions.resize(num_threads);
237 dvalue_expressions.resize(num_threads);
238
239 auto symbol_tables = createSymbolTables(num_threads, variables, curves);
240
241 // value expressions per thread.
242 for (int thread_id = 0; thread_id < num_threads; ++thread_id)
243 {
244 value_expressions[thread_id] = compileExpressions(
245 symbol_tables[thread_id], value_string_expressions);
246 }
247
248 // dValue expressions per thread.
249 for (int thread_id = 0; thread_id < num_threads; ++thread_id)
250 {
251 for (auto const& [variable_name, string_expressions] :
252 dvalue_string_expressions)
253 {
254 dvalue_expressions[thread_id].emplace_back(
255 convertStringToVariable(variable_name),
256 compileExpressions(symbol_tables[thread_id],
257 string_expressions));
258 }
259 }
260}
261
262static void updateVariableArrayValues(std::vector<Variable> const& variables,
263 VariableArray const& new_variable_array,
264 VariableArray& variable_array)
265{
266 for (auto const& variable : variables)
267 {
268 auto assign_scalar =
269 [&variable, &new_variable_array](VariableArray::Scalar* address)
270 {
271 double const value = *std::get<VariableArray::Scalar const*>(
272 new_variable_array.address_of(variable));
273
274 if (std::isnan(value))
275 {
276 OGS_FATAL(
277 "Function property: Scalar variable '{:s}' is not "
278 "initialized.",
279 variable_enum_to_string[static_cast<int>(variable)]);
280 }
281
282 *address = value;
283 };
284 auto assign_kelvin_vector = [&variable, &new_variable_array](
286 {
287 auto assign_value = [&destination = *address,
288 &variable]<typename S>(S const& source)
289 {
290 if constexpr (std::is_same_v<S, std::monostate>)
291 {
292 OGS_FATAL(
293 "Function property: Kelvin vector variable '{:s}' is "
294 "not initialized.",
295 variable_enum_to_string[static_cast<int>(variable)]);
296 }
297 else
298 {
299 if (std::holds_alternative<S>(destination))
300 {
301 std::get<S>(destination) = MathLib::KelvinVector::
303 }
304 else
305 {
306 OGS_FATAL(
307 "Function property: Mismatch of Kelvin vector "
308 "sizes for variable {:s}.",
309 variable_enum_to_string[static_cast<int>(
310 variable)]);
311 }
312 }
313 };
314 std::visit(assign_value,
315 *std::get<VariableArray::KelvinVector const*>(
316 new_variable_array.address_of(variable)));
317 };
318 auto assign_deformation_gradient =
319 [&variable,
320 &new_variable_array](VariableArray::DeformationGradient* address)
321 {
322 auto assign_value = [&destination = *address,
323 &variable]<typename S>(S const& source)
324 {
325 if constexpr (std::is_same_v<S, std::monostate>)
326 {
327 OGS_FATAL(
328 "Function property: Vectorized tensor variable '{:s}' "
329 "is not initialized.",
330 variable_enum_to_string[static_cast<int>(variable)]);
331 }
332 else
333 {
334 if (std::holds_alternative<S>(destination))
335 {
336 std::get<S>(destination) = source;
337 }
338 else
339 {
340 OGS_FATAL(
341 "Function property: Mismatch of vectorized tensor "
342 "sizes for variable {:s}.",
343 variable_enum_to_string[static_cast<int>(
344 variable)]);
345 }
346 }
347 };
348 std::visit(assign_value,
349 *std::get<VariableArray::DeformationGradient const*>(
350 new_variable_array.address_of(variable)));
351 };
352
353 variable_array.visitVariable(
354 BaseLib::Overloaded{assign_scalar, assign_kelvin_vector,
355 assign_deformation_gradient},
356 variable);
357 }
358}
359
361 std::vector<Variable> const& variables,
362 VariableArray const& new_variable_array,
363 ParameterLib::SpatialPosition const& pos, double const t,
364 std::vector<exprtk::expression<double>> const& expressions,
365 VariableArray& variable_array, bool const spatial_position_is_required)
366{
367 std::vector<double> result(expressions.size());
368
369 updateVariableArrayValues(variables, new_variable_array, variable_array);
370
371 std::transform(begin(expressions), end(expressions), begin(result),
372 [t, &pos, spatial_position_is_required](auto const& e)
373 {
374 auto& symbol_table = e.get_symbol_table();
375 symbol_table.get_variable("t")->ref() = t;
376
377 if (spatial_position_is_required)
378 {
379 if (!pos.getCoordinates())
380 {
381 OGS_FATAL(
382 "FunctionParameter: The spatial position "
383 "has to be set by "
384 "coordinates.");
385 }
386 auto const coords = pos.getCoordinates().value();
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];
390 }
391
392 return e.value();
393 });
394
395 switch (result.size())
396 {
397 case 1:
398 {
399 return result[0];
400 }
401 case 2:
402 {
403 return Eigen::Vector2d{result[0], result[1]};
404 }
405 case 3:
406 {
407 return Eigen::Vector3d{result[0], result[1], result[2]};
408 }
409 case 4:
410 {
411 Eigen::Matrix<double, 2, 2> m;
412 m = Eigen::Map<Eigen::Matrix<double, 2, 2> const>(result.data(), 2,
413 2);
414 return m;
415 }
416 case 9:
417 {
418 Eigen::Matrix<double, 3, 3> m;
419 m = Eigen::Map<Eigen::Matrix<double, 3, 3> const>(result.data(), 3,
420 3);
421 return m;
422 }
423 }
424 OGS_FATAL("Cannot convert a vector of size {} to a PropertyDataType",
425 result.size());
426}
427
428static std::vector<std::string> collectVariables(
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)
432{
433 std::vector<std::string> variables;
434
435 auto collect_variables = [&](auto string_expressions)
436 {
437 for (auto const& string_expression : string_expressions)
438 {
439 if (!exprtk::collect_variables(string_expression, variables))
440 {
441 OGS_FATAL(
442 "Collecting variables from expression '{}' didn't work.",
443 string_expression);
444 }
445 }
446 };
447
448 collect_variables(value_string_expressions);
449 for (auto const& var_name_string_expression : dvalue_string_expressions)
450 {
451 collect_variables(var_name_string_expression.second);
452 }
453
454 BaseLib::makeVectorUnique(variables);
455 return variables;
456}
457
459 std::string name,
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&
465 curves)
466{
467 name_ = std::move(name);
468
469 auto const variables =
470 collectVariables(value_string_expressions, dvalue_string_expressions);
471
472 // filter the strings unrelated to time or spatial position
473 static std::unordered_set<std::string> filter_not_variables = {"t", "x",
474 "y", "z"};
475 for (auto const& curve : curves)
476 {
477 filter_not_variables.insert(curve.first);
478 }
479
480 variables_ =
481 variables |
482 std::views::filter([](const std::string& s)
483 { return !filter_not_variables.contains(s); }) |
484 std::views::transform([](std::string const& s)
485 { return convertStringToVariable(s); }) |
486 ranges::to<std::vector>;
487
488 auto const get_number_omp_threads = []()
489 {
490#ifdef _OPENMP
491 return omp_get_max_threads();
492#else
493 return 1;
494#endif
495 };
496
497 int const num_threads = std::max(BaseLib::getNumberOfAssemblyThreads(),
498 get_number_omp_threads());
499
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);
506}
507
508std::variant<Function::Implementation<2>*, Function::Implementation<3>*>
510 VariableArray const& variable_array) const
511{
512 if (variable_array.is2D())
513 {
514 return impl2_.get();
515 }
516 if (variable_array.is3D())
517 {
518 return impl3_.get();
519 }
520
521 OGS_FATAL(
522 "Variable array has vectors for 2 and 3 dimensions simultaneously. "
523 "Mixed dimensions cannot be dealt within Function evaluation.");
524}
525
528 double const t, double const /*dt*/) const
529{
530#ifdef _OPENMP
531 int const thread_id = omp_get_thread_num();
532#else
533 int const thread_id = 0;
534#endif
535 return std::visit(
536 [&](auto&& impl_ptr)
537 {
538 if (thread_id >=
539 static_cast<int>(impl_ptr->value_expressions.size()))
540 {
541 OGS_FATAL(
542 "In Function-type property '{:s}' evaluation the "
543 "OMP-thread with id {:d} exceeds the number of allocated "
544 "threads {:d}.",
545 name_, thread_id, impl_ptr->value_expressions.size());
546 }
547 return evaluateExpressions(variables_, variable_array, pos, t,
548 impl_ptr->value_expressions[thread_id],
549 impl_ptr->variable_arrays[thread_id],
550 impl_ptr->spatial_position_is_required);
551 },
553}
554
556 Variable const variable,
558 double const t, double const /*dt*/) const
559{
560#ifdef _OPENMP
561 int const thread_id = omp_get_thread_num();
562#else
563 int const thread_id = 0;
564#endif
565 return std::visit(
566 [&](auto&& impl_ptr)
567 {
568 if (thread_id >=
569 static_cast<int>(impl_ptr->dvalue_expressions.size()))
570 {
571 OGS_FATAL(
572 "In Function-type property '{:s}' evaluation the "
573 "OMP-thread with id {:d} exceeds the number of allocated "
574 "threads {:d}.",
575 name_, thread_id, impl_ptr->value_expressions.size());
576 }
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; });
581
582 if (it == end(impl_ptr->dvalue_expressions[thread_id]))
583 {
584 OGS_FATAL(
585 "Requested derivative with respect to the variable {:s} "
586 "not "
587 "provided for Function-type property {:s}.",
588 variable_enum_to_string[static_cast<int>(variable)], name_);
589 }
590
591 return evaluateExpressions(variables_, variable_array, pos, t,
592 it->second,
593 impl_ptr->variable_arrays[thread_id],
594 impl_ptr->spatial_position_is_required);
595 },
597}
598
599Function::~Function() = default;
600
601} // namespace MaterialPropertyLib
#define OGS_FATAL(...)
Definition Error.h:26
Definition of the PiecewiseLinearInterpolation class.
std::vector< VariableArray > variable_arrays
Definition Function.cpp:115
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.
Definition Function.cpp:206
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)
Definition Function.cpp:225
exprtk::expression< double > Expression
Definition Function.cpp:72
std::map< std::string, CurveWrapper > _curve_wrappers
Definition Function.cpp:117
std::vector< std::vector< std::pair< Variable, std::vector< Expression > > > > dvalue_expressions
Definition Function.cpp:111
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)
Definition Function.cpp:123
std::vector< std::vector< Expression > > value_expressions
Definition Function.cpp:106
std::unique_ptr< Implementation< 3 > > impl3_
Definition Function.h:59
std::vector< Variable > variables_
Variables used in the exprtk expressions.
Definition Function.h:66
std::variant< Function::Implementation< 2 > *, Function::Implementation< 3 > * > getImplementationForDimensionOfVariableArray(VariableArray const &variable_array) const
Definition Function.cpp:509
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:458
PropertyDataType dValue(VariableArray const &variable_array, Variable const variable, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const override
Definition Function.cpp:555
std::unique_ptr< Implementation< 2 > > impl2_
Definition Function.h:58
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
int getNumberOfAssemblyThreads()
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:48
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:428
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)
Definition Function.cpp:360
static void updateVariableArrayValues(std::vector< Variable > const &variables, VariableArray const &new_variable_array, VariableArray &variable_array)
Definition Function.cpp:262
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:40
CurveWrapper(const MathLib::PiecewiseLinearInterpolation &curve)
Definition Function.cpp:34
const MathLib::PiecewiseLinearInterpolation & _curve
Definition Function.cpp:43