diff --git a/Documentation/ProjectFile/prj/parameters/parameter/Function/t_expression.md b/Documentation/ProjectFile/prj/parameters/parameter/Function/t_expression.md
index f9eaacbf26b6efa4cedc4acaa0d3cc050bd96ea7..823536ab6268685d048feda0048b66cd905384cc 100644
--- a/Documentation/ProjectFile/prj/parameters/parameter/Function/t_expression.md
+++ b/Documentation/ProjectFile/prj/parameters/parameter/Function/t_expression.md
@@ -1 +1,6 @@
-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))`.
diff --git a/ParameterLib/FunctionParameter.cpp b/ParameterLib/FunctionParameter.cpp
index 11da80a75809b894742d1f35b191db2e1d6877ac..2c3511be564dbde6c146694dcde55ac38e704c93 100644
--- a/ParameterLib/FunctionParameter.cpp
+++ b/ParameterLib/FunctionParameter.cpp
@@ -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
diff --git a/ParameterLib/FunctionParameter.h b/ParameterLib/FunctionParameter.h
index 812fdeac365dfe348e6f4079bf2c813b9547ac73..8a316c5f1fe741e08ce885cbbf44b0be288137d0 100644
--- a/ParameterLib/FunctionParameter.h
+++ b/ParameterLib/FunctionParameter.h
@@ -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(
diff --git a/ParameterLib/Parameter.cpp b/ParameterLib/Parameter.cpp
index cd39a7309d7a42994d227920429a07dba5f0d4b1..3bafd8cde755651bb807f0143b62d0c5291beb0d 100644
--- a/ParameterLib/Parameter.cpp
+++ b/ParameterLib/Parameter.cpp
@@ -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")
     {
diff --git a/Tests/Data/Parabolic/T/t2_1D2bt/t2_1D2bt.prj b/Tests/Data/Parabolic/T/t2_1D2bt/t2_1D2bt.prj
index 3a3617bda01725f06646ee2b87c48eb50856e462..51b340f9e194cdcef36128605003654ccaa7dc40 100644
--- a/Tests/Data/Parabolic/T/t2_1D2bt/t2_1D2bt.prj
+++ b/Tests/Data/Parabolic/T/t2_1D2bt/t2_1D2bt.prj
@@ -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>