| Name | Description |
|---|---|
| Bicubic function | |
| Biquadratic function | |
| Polynomial function | |
| Power function that is linearized below a user-defined threshold | |
| Function that is quadratic in first argument and linear in second argument | |
| Power function, regularized near zero, but nonzero value for x=0 | |
| Once continuously differentiable approximation to exp(-|x|) in interval |x| < delta | |
| Once continuously differentiable approximation to the Heaviside function | |
| Once continuously differentiable approximation to the maximum function | |
| Once continuously differentiable approximation to the minimum function | |
| Once continuously differentiable approximation to the limit function | |
| Integration using the trapezoidal rule | |
| Collection of models that illustrate model use and test models | |
| Package with base classes for Buildings.Utilities.Math.Functions | |
| Function that approximates 1/x by a twice continuously differentiable function |
y = a1 + a2 x1 + a3 x12 + a4 x2 + a5 x22 + a6 x1 x2 + a7 x1^3 + a8 x2^3 + a9 x12 x2 + a10 x1 x22
| Type | Name | Default | Description |
|---|---|---|---|
| Real | a[10] | Coefficients | |
| Real | x1 | Independent variable | |
| Real | x2 | Independent variable |
| Type | Name | Description |
|---|---|---|
| Real | y | Result |
function bicubic "Bicubic function"
input Real a[10] "Coefficients";
input Real x1 "Independent variable";
input Real x2 "Independent variable";
output Real y "Result";
protected
Real x1Sq "= x1^2";
Real x2Sq "= x2^2";
algorithm
x1Sq :=x1*x1;
x2Sq :=x2*x2;
y := a[1] + a[2] * x1 + a[3] * x1^2
+ a[4] * x2 + a[5] * x2^2
+ a[6] * x1 * x2
+ a[7] * x1Sq * x1
+ a[8] * x2Sq * x2
+ a[9] * x1Sq * x2
+ a[10] * x1 * x2Sq;
end bicubic;
y = a1 + a2 x1 + a3 x12 + a4 x2 + a5 x22 + a6 x1 x2
| Type | Name | Default | Description |
|---|---|---|---|
| Real | a[6] | Coefficients | |
| Real | x1 | Independent variable | |
| Real | x2 | Independent variable |
| Type | Name | Description |
|---|---|---|
| Real | y | Result |
function biquadratic "Biquadratic function" input Real a[6] "Coefficients"; input Real x1 "Independent variable"; input Real x2 "Independent variable"; output Real y "Result"; algorithm y :=a[1] + x1*(a[2] + a[3]*x1) + x2*(a[4]+ a[5]*x2) + a[6]*x1*x2;end biquadratic;
y = a1 + a2 x + a3 x2 + ...
| Type | Name | Default | Description |
|---|---|---|---|
| Real | x | Independent variable | |
| Real | a[:] | Coefficients |
| Type | Name | Description |
|---|---|---|
| Real | y | Result |
function polynomial "Polynomial function"
annotation(derivative=Buildings.Utilities.Math.Functions.BaseClasses.der_polynomial);
input Real x "Independent variable";
input Real a[:] "Coefficients";
output Real y "Result";
protected
parameter Integer n = size(a, 1)-1;
Real xp[n+1] "Powers of x";
algorithm
xp[1] :=1;
for i in 1:n loop
xp[i+1] :=xp[i]*x;
end for;
y :=a*xp;
end polynomial;
Function that approximates y=xn where 0 < n so that
For x < x0, this function replaces y=xn by a linear function that is continuously differentiable everywhere.
A typical use of this function is to replace T = T4(1/4) in a radiation balance to ensure that the function is defined everywhere. This can help solving the initialization problem when a solver may be far from a solution and hence T4 < 0.
See the package Examples for the graph.
| Type | Name | Default | Description |
|---|---|---|---|
| Real | x | Abscissa value | |
| Real | n | Exponent | |
| Real | x0 | Abscissa value below which linearization occurs |
| Type | Name | Description |
|---|---|---|
| Real | y | Function value |
function powerLinearized "Power function that is linearized below a user-defined threshold" input Real x "Abscissa value"; input Real n "Exponent"; input Real x0 "Abscissa value below which linearization occurs"; output Real y "Function value"; algorithm if x > x0 then y := x^n; else y := x0^n * (1-n) + n * x0^(n-1) * x; end if;end powerLinearized;
y = a1 + a2 x1 + a3 x12 + (a4 + a5 x1 + a6 x12) x2
| Type | Name | Default | Description |
|---|---|---|---|
| Real | a[6] | Coefficients | |
| Real | x1 | Independent variable for quadratic part | |
| Real | x2 | Independent variable for linear part |
| Type | Name | Description |
|---|---|---|
| Real | y | Result |
function quadraticLinear "Function that is quadratic in first argument and linear in second argument" input Real a[6] "Coefficients"; input Real x1 "Independent variable for quadratic part"; input Real x2 "Independent variable for linear part"; output Real y "Result"; protected Real x1Sq; algorithm x1Sq :=x1*x1; y :=a[1] + a[2]*x1 + a[3]*x1Sq + (a[4] + a[5]*x1 + a[6]*x1Sq)*x2;end quadraticLinear;
Function that approximates y=|x|n where n > 0 so that
This function replaces y=|x|n in the interval -δ...+δ by a 4-th order polynomial that has the same function value and the first and second derivative at x=± δ.
A typical use of this function is to replace the function for the convective heat transfer coefficient for forced or free convection that is of the form h=c |dT|n for some constant c and exponent 0 ≤ n ≤ 1. By using this function, the original function that has an infinite derivative near zero and that takes on zero at the origin is replaced by a function with a bounded derivative and a non-zero value at the origin. Physically, the region -δ...+δ may be interpreted as the region where heat conduction dominates convection in the boundary layer.
See the package Examples for the graph.
| Type | Name | Default | Description |
|---|---|---|---|
| Real | x | Abscissa value | |
| Real | n | Exponent | |
| Real | delta | 0.01 | Abscissa value where transition occurs |
| Type | Name | Description |
|---|---|---|
| Real | y | Function value |
function regNonZeroPower "Power function, regularized near zero, but nonzero value for x=0" annotation(derivative=BaseClasses.der_regNonZeroPower); input Real x "Abscissa value"; input Real n "Exponent"; input Real delta = 0.01 "Abscissa value where transition occurs"; output Real y "Function value"; protected Real a1; Real a3; Real a5; Real delta2; Real x2; Real y_d "=y(delta)"; Real yP_d "=dy(delta)/dx"; Real yPP_d "=d^2y(delta)/dx^2"; algorithm if abs(x) > delta then y := abs(x)^n; else delta2 :=delta*delta; x2 :=x*x; y_d :=delta^n; yP_d :=n*delta^(n - 1); yPP_d :=n*(n - 1)*delta^(n - 2); a1 := -(yP_d/delta - yPP_d)/delta2/8; a3 := (yPP_d - 12 * a1 * delta2)/2; a5 := (y_d - delta2 * (a3 + delta2 * a1)); y := a5 + x2 * (a3 + x2 * a1); assert(a5>0, "Delta is too small for this exponent."); end if;end regNonZeroPower;
Function to provide a once continuously differentiable approximation to exp(- |x| ) in the interval |x| < δ for some positive δ
| Type | Name | Default | Description |
|---|---|---|---|
| Real | x | Input argument | |
| Real | delta | Transition point where approximation occurs |
| Type | Name | Description |
|---|---|---|
| Real | y | Output argument |
function smoothExponential
"Once continuously differentiable approximation to exp(-|x|) in interval |x| < delta"
input Real x "Input argument";
input Real delta "Transition point where approximation occurs";
output Real y "Output argument";
protected
Real absX "Absolute value of x";
Real a2 "Coefficient for approximating function";
Real a3 "Coefficient for approximating function";
Real e;
Real d2;
Real x2;
algorithm
absX :=abs(x);
if absX > delta then
y := exp(-absX);
else
d2 := delta*delta;
e := exp(-delta);
a2 := (delta*e-4*(1-e))/2/d2;
a3 := (e-1-a2*d2)/d2/d2;
x2 := x*x;
y := 1+x2*(a2+x2*a3);
end if;
end smoothExponential;
Once Lipschitz continuously differentiable approximation to the
Heaviside(.,.) function.
| Type | Name | Default | Description |
|---|---|---|---|
| Real | x | Argument | |
| Real | delta | Parameter used for scaling |
| Type | Name | Description |
|---|---|---|
| Real | y | Result |
function smoothHeaviside "Once continuously differentiable approximation to the Heaviside function" input Real x "Argument"; input Real delta "Parameter used for scaling"; output Real y "Result"; algorithm y := spliceFunction(1, 0, x, delta);end smoothHeaviside;
Once continuously differentiable approximation to the max(.,.) function.
| Type | Name | Default | Description |
|---|---|---|---|
| Real | x1 | First argument | |
| Real | x2 | Second argument | |
| Real | deltaX | Width of transition interval |
| Type | Name | Description |
|---|---|---|
| Real | y | Result |
function smoothMax
"Once continuously differentiable approximation to the maximum function"
input Real x1 "First argument";
input Real x2 "Second argument";
input Real deltaX "Width of transition interval";
output Real y "Result";
algorithm
y := Buildings.Utilities.Math.Functions.spliceFunction(
pos=x1, neg=x2, x=x1-x2, deltax=deltaX);
end smoothMax;
Once continuously differentiable approximation to the max(.,.) function.
| Type | Name | Default | Description |
|---|---|---|---|
| Real | x1 | First argument | |
| Real | x2 | Second argument | |
| Real | deltaX | Width of transition interval |
| Type | Name | Description |
|---|---|---|
| Real | y | Result |
function smoothMin
"Once continuously differentiable approximation to the minimum function"
input Real x1 "First argument";
input Real x2 "Second argument";
input Real deltaX "Width of transition interval";
output Real y "Result";
algorithm
y := Buildings.Utilities.Math.Functions.spliceFunction(
pos=x1, neg=x2, x=x2-x1, deltax=deltaX);
end smoothMin;
Once continuously differentiable approximation to the limit(.,.) function.
The output is bounded to be in [0, 1].
| Type | Name | Default | Description |
|---|---|---|---|
| Real | x | Variable | |
| Real | l | Low limit | |
| Real | u | Upper limit | |
| Real | deltaX | Width of transition interval |
| Type | Name | Description |
|---|---|---|
| Real | y | Result |
function smoothLimit "Once continuously differentiable approximation to the limit function" input Real x "Variable"; input Real l "Low limit"; input Real u "Upper limit"; input Real deltaX "Width of transition interval"; output Real y "Result"; protected Real cor; algorithm cor :=deltaX/10; y := Buildings.Utilities.Math.Functions.smoothMax(x,l+deltaX,cor); y := Buildings.Utilities.Math.Functions.smoothMin(y,u-deltaX,cor);end smoothLimit;
Function to provide a once continuously differentialbe transition between to arguments.
The function is adapted from Modelica.Media.Air.MoistAir.Utilities.spliceFunction and provided here for easier accessability to model developers.
| Type | Name | Default | Description |
|---|---|---|---|
| Real | pos | Argument of x > 0 | |
| Real | neg | Argument of x < 0 | |
| Real | x | Independent value | |
| Real | deltax | Half width of transition interval |
| Type | Name | Description |
|---|---|---|
| Real | out | Smoothed value |
function spliceFunction
annotation(derivative=BaseClasses.der_spliceFunction);
input Real pos "Argument of x > 0";
input Real neg "Argument of x < 0";
input Real x "Independent value";
input Real deltax "Half width of transition interval";
output Real out "Smoothed value";
protected
Real scaledX;
Real scaledX1;
Real y;
algorithm
scaledX1 := x/deltax;
scaledX := scaledX1*Modelica.Math.asin(1);
if scaledX1 <= -0.999999999 then
y := 0;
elseif scaledX1 >= 0.999999999 then
y := 1;
else
y := (Modelica.Math.tanh(Modelica.Math.tan(scaledX)) + 1)/2;
end if;
out := pos*y + (1 - y)*neg;
end spliceFunction;
This function computes a definite integral using the trapezoidal rule.
| Type | Name | Default | Description |
|---|---|---|---|
| Integer | N | Number of integrand points | |
| Real | f[:] | Integrands | |
| Real | deltaX | Width of interval for Trapezoidal integration |
| Type | Name | Description |
|---|---|---|
| Real | result | Result |
function trapezoidalIntegration
"Integration using the trapezoidal rule"
input Integer N "Number of integrand points";
input Real[:] f "Integrands";
input Real deltaX "Width of interval for Trapezoidal integration";
output Real result "Result";
algorithm
assert(N >= 2, "N must be no less than 2.");
result := 0;
for i in 1:N loop
result := result + f[i];
end for;
result := 2*result;
result := result - f[1] - f[N];
result := result*deltaX/2;
end trapezoidalIntegration;
Function that approximates y=1 ⁄ x inside the interval -δ ≤ x ≤ δ. The approximation is twice continuously differentiable with a bounded derivative on the whole real line.
See the package Examples for the graph.
| Type | Name | Default | Description |
|---|---|---|---|
| Real | x | Abscissa value | |
| Real | delta | Abscissa value below which approximation occurs |
| Type | Name | Description |
|---|---|---|
| Real | y | Function value |
function inverseXRegularized
"Function that approximates 1/x by a twice continuously differentiable function"
input Real x "Abscissa value";
input Real delta(min=0) "Abscissa value below which approximation occurs";
output Real y "Function value";
protected
Real delta2 "Delta^2";
Real x2_d2 "=x^2/delta^2";
algorithm
delta2 :=delta*delta;
x2_d2 := x*x/delta2;
y :=smooth(2, if (abs(x) > delta) then 1/x else
x/delta2 + x*abs(x/delta2/delta*(2 - x2_d2*(3 - x2_d2))));
end inverseXRegularized;