Skip to content
Snippets Groups Projects
Commit 51c1252e authored by Dmitri Naumov's avatar Dmitri Naumov
Browse files

Merge branch 'AllowCurveUsageInFunctionParameters' into 'master'

Allow curve usage in function parameters

See merge request ogs/ogs!3058
parents cad22edf a3bda72c
No related branches found
No related tags found
No related merge requests found
Mathematical expression of the function (currently x,y,z, and t are supported variables). For non-scalar values, an expression should be given for each component.
Mathematical expression of a function.
For non-scalar values, an expression should be given for each component.
Currently x,y,z, and t are supported variables, and curves from the `<curves>`
section can be called using their names. A curve is a single argument function
and can be used in an expression like `curveA(sin(t))`.
......@@ -15,7 +15,10 @@
namespace ParameterLib
{
std::unique_ptr<ParameterBase> createFunctionParameter(
std::string const& name, BaseLib::ConfigTree const& config)
std::string const& name, BaseLib::ConfigTree const& config,
std::map<std::string,
std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
curves)
{
//! \ogs_file_param{prj__parameters__parameter__type}
config.checkConfigParameter("type", "Function");
......@@ -29,7 +32,8 @@ std::unique_ptr<ParameterBase> createFunctionParameter(
vec_expressions.emplace_back(expression_str);
}
return std::make_unique<FunctionParameter<double>>(name, vec_expressions);
return std::make_unique<FunctionParameter<double>>(name, vec_expressions,
curves);
}
} // namespace ParameterLib
......@@ -28,6 +28,23 @@ namespace ParameterLib
template <typename T>
struct FunctionParameter final : public Parameter<T>
{
class CurveWrapper : public exprtk::ifunction<T>
{
public:
CurveWrapper(MathLib::PiecewiseLinearInterpolation const& curve)
: exprtk::ifunction<T>(1), _curve(curve)
{
exprtk::disable_has_side_effects(*this);
}
double operator()(double const& t) override
{
return _curve.getValue(t);
}
private:
MathLib::PiecewiseLinearInterpolation const& _curve;
};
using symbol_table_t = exprtk::symbol_table<T>;
using expression_t = exprtk::expression<T>;
using parser_t = exprtk::parser<T>;
......@@ -40,16 +57,34 @@ struct FunctionParameter final : public Parameter<T>
* @param vec_expression_str a vector of mathematical expressions
* The vector size specifies the number of components of the parameter.
*/
FunctionParameter(std::string const& name,
std::vector<std::string> const& vec_expression_str)
FunctionParameter(
std::string const& name,
std::vector<std::string> const& vec_expression_str,
std::map<std::string,
std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
curves)
: Parameter<T>(name, nullptr), _vec_expression_str(vec_expression_str)
{
// Convert curves to function objects callable by the exprtk.
_curves.reserve(curves.size());
std::transform(
begin(curves), end(curves), std::back_inserter(_curves),
[](auto const& curve) -> std::pair<std::string, CurveWrapper> {
return {curve.first, CurveWrapper(*curve.second)};
});
// Create symbol table for variables and functions.
_symbol_table.add_constants();
_symbol_table.create_variable("x");
_symbol_table.create_variable("y");
_symbol_table.create_variable("z");
_symbol_table.create_variable("t");
for (auto& curve : _curves)
{
_symbol_table.add_function(curve.first, curve.second);
}
// Compile expressions.
_vec_expression.resize(_vec_expression_str.size());
for (unsigned i = 0; i < _vec_expression_str.size(); i++)
{
......@@ -108,6 +143,7 @@ private:
std::vector<std::string> const _vec_expression_str;
symbol_table_t _symbol_table;
std::vector<expression_t> _vec_expression;
std::vector<std::pair<std::string, CurveWrapper>> _curves;
};
std::unique_ptr<ParameterBase> createFunctionParameter(
......
......@@ -59,7 +59,7 @@ std::unique_ptr<ParameterBase> createParameter(
if (type == "Function")
{
INFO("FunctionParameter: {:s}", name);
return createFunctionParameter(name, config);
return createFunctionParameter(name, config, curves);
}
if (type == "Group")
{
......
......@@ -161,14 +161,8 @@
</parameter>
<parameter>
<name>h1</name>
<type>CurveScaled</type>
<curve>h1_curve</curve>
<parameter>h1_spatial</parameter>
</parameter>
<parameter>
<name>h1_spatial</name>
<type>Constant</type>
<value>1</value>
<type>Function</type>
<expression>h1_curve(t)</expression>
</parameter>
</parameters>
<curves>
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment