Buildings.Utilities.Math.Functions

Package with mathematical functions

Information

This package contains functions for commonly used mathematical operations. The functions are used in the blocks Buildings.Utilities.Math.

Extends from Modelica.Icons.VariantsPackage (Icon for package containing variants).

Package Content

Name Description
Buildings.Utilities.Math.Functions.average average Average of a vector
Buildings.Utilities.Math.Functions.bicubic bicubic Bicubic function
Buildings.Utilities.Math.Functions.biquadratic biquadratic Biquadratic function
Buildings.Utilities.Math.Functions.booleanReplicator booleanReplicator Replicates Boolean signals
Buildings.Utilities.Math.Functions.cubicHermiteLinearExtrapolation cubicHermiteLinearExtrapolation Interpolate using a cubic Hermite spline with linear extrapolation
Buildings.Utilities.Math.Functions.integerReplicator integerReplicator Replicates integer signals
Buildings.Utilities.Math.Functions.inverseXRegularized inverseXRegularized Function that approximates 1/x by a twice continuously differentiable function
Buildings.Utilities.Math.Functions.isMonotonic isMonotonic Returns true if the argument is a monotonic sequence
Buildings.Utilities.Math.Functions.polynomial polynomial Polynomial function
Buildings.Utilities.Math.Functions.powerLinearized powerLinearized Power function that is linearized below a user-defined threshold
Buildings.Utilities.Math.Functions.quadraticLinear quadraticLinear Function that is quadratic in first argument and linear in second argument
Buildings.Utilities.Math.Functions.quinticHermite quinticHermite Quintic Hermite spline function for interpolation between two functions with a continuous second derivative
Buildings.Utilities.Math.Functions.regNonZeroPower regNonZeroPower Power function, regularized near zero, but nonzero value for x=0
Buildings.Utilities.Math.Functions.regStep regStep Approximation of a general step, such that the approximation is continuous and differentiable
Buildings.Utilities.Math.Functions.smoothExponential smoothExponential Once continuously differentiable approximation to exp(-|x|) in interval |x| < delta
Buildings.Utilities.Math.Functions.smoothHeaviside smoothHeaviside Once continuously differentiable approximation to the Heaviside function
Buildings.Utilities.Math.Functions.smoothInterpolation smoothInterpolation Interpolate using a cubic Hermite spline with linear extrapolation for a vector xSup[], ySup[] and independent variable x
Buildings.Utilities.Math.Functions.smoothLimit smoothLimit Once continuously differentiable approximation to the limit function
Buildings.Utilities.Math.Functions.smoothMax smoothMax Once continuously differentiable approximation to the maximum function
Buildings.Utilities.Math.Functions.smoothMin smoothMin Once continuously differentiable approximation to the minimum function
Buildings.Utilities.Math.Functions.spliceFunction spliceFunction  
Buildings.Utilities.Math.Functions.splineDerivatives splineDerivatives Function to compute the derivatives for cubic hermite spline interpolation
Buildings.Utilities.Math.Functions.trapezoidalIntegration trapezoidalIntegration Integration using the trapezoidal rule
Buildings.Utilities.Math.Functions.Examples Examples Collection of models that illustrate model use and test models
Buildings.Utilities.Math.Functions.BaseClasses BaseClasses Package with base classes for Buildings.Utilities.Math.Functions

Buildings.Utilities.Math.Functions.average Buildings.Utilities.Math.Functions.average

Average of a vector

Information

This function outputs the average of the vector.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Integernin Number of inputs
Realu[nin] Input vector

Outputs

TypeNameDescription
RealyResult

Modelica definition

function average "Average of a vector" extends Modelica.Icons.Function; input Integer nin "Number of inputs"; input Real u[nin] "Input vector"; output Real y "Result"; algorithm y := sum(u)/nin; end average;

Buildings.Utilities.Math.Functions.bicubic Buildings.Utilities.Math.Functions.bicubic

Bicubic function

Information

This function computes

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

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Reala[10] Coefficients
Realx1 Independent variable
Realx2 Independent variable

Outputs

TypeNameDescription
RealyResult

Modelica definition

function bicubic "Bicubic function" extends Modelica.Icons.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;

Buildings.Utilities.Math.Functions.biquadratic Buildings.Utilities.Math.Functions.biquadratic

Biquadratic function

Information

This function computes

y = a1 + a2 x1 + a3 x12 + a4 x2 + a5 x22 + a6 x1 x2

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Reala[6] Coefficients
Realx1 Independent variable
Realx2 Independent variable

Outputs

TypeNameDescription
RealyResult

Modelica definition

function biquadratic "Biquadratic function" extends Modelica.Icons.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;

Buildings.Utilities.Math.Functions.booleanReplicator Buildings.Utilities.Math.Functions.booleanReplicator

Replicates Boolean signals

Information

This function replicates the boolean input signal to an array of nout identical output signals.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Integernout1Number of outouts
Booleanu Boolean input signal

Outputs

TypeNameDescription
Booleany[nout]Boolean output signals

Modelica definition

function booleanReplicator "Replicates Boolean signals" extends Modelica.Icons.Function; input Integer nout=1 "Number of outouts"; input Boolean u "Boolean input signal"; output Boolean y[nout] "Boolean output signals"; algorithm y := fill(u, nout); end booleanReplicator;

Buildings.Utilities.Math.Functions.cubicHermiteLinearExtrapolation Buildings.Utilities.Math.Functions.cubicHermiteLinearExtrapolation

Interpolate using a cubic Hermite spline with linear extrapolation

Information

For x1 < x < x2, this function interpolates using cubic hermite spline. For x outside this interval, the function linearly extrapolates.

For how to use this function, see Buildings.Utilities.Math.Functions.Examples.CubicHermite.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realx Abscissa value
Realx1 Lower abscissa value
Realx2 Upper abscissa value
Realy1 Lower ordinate value
Realy2 Upper ordinate value
Realy1d Lower gradient
Realy2d Upper gradient

Outputs

TypeNameDescription
RealyInterpolated ordinate value

Modelica definition

function cubicHermiteLinearExtrapolation "Interpolate using a cubic Hermite spline with linear extrapolation" extends Modelica.Icons.Function; input Real x "Abscissa value"; input Real x1 "Lower abscissa value"; input Real x2 "Upper abscissa value"; input Real y1 "Lower ordinate value"; input Real y2 "Upper ordinate value"; input Real y1d "Lower gradient"; input Real y2d "Upper gradient"; output Real y "Interpolated ordinate value"; algorithm if (x > x1 and x < x2) then y := Modelica.Fluid.Utilities.cubicHermite( x=x, x1=x1, x2=x2, y1=y1, y2=y2, y1d=y1d, y2d=y2d); elseif x <= x1 then // linear extrapolation y := y1 + (x - x1)*y1d; else y := y2 + (x - x2)*y2d; end if; end cubicHermiteLinearExtrapolation;

Buildings.Utilities.Math.Functions.integerReplicator Buildings.Utilities.Math.Functions.integerReplicator

Replicates integer signals

Information

This function replicates the integer input signal to an array of nout identical output signals.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Integernout Number of outputs
Integeru Integer input signal

Outputs

TypeNameDescription
Integery[nout]Integer output signals

Modelica definition

function integerReplicator "Replicates integer signals" extends Modelica.Icons.Function; input Integer nout "Number of outputs"; input Integer u "Integer input signal"; output Integer y[nout] "Integer output signals"; algorithm y := fill(u, nout); end integerReplicator;

Buildings.Utilities.Math.Functions.inverseXRegularized Buildings.Utilities.Math.Functions.inverseXRegularized

Function that approximates 1/x by a twice continuously differentiable function

Information

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 plot of Buildings.Utilities.Math.Functions.Examples.InverseXRegularized for the graph.

For efficiency, the polynomial coefficients a, b, c, d, e, f and the inverse of the smoothing parameter deltaInv are exposed as arguments to this function. Typically, these coefficients only depend on parameters and hence can be computed once. They must be equal to their default values, otherwise the function is not twice continuously differentiable. By exposing these coefficients as function arguments, models that call this function can compute them as parameters, and assign these parameter values in the function call. This avoids that the coefficients are evaluated for each time step, as they would otherwise be if they were to be computed inside the body of the function. However, assigning the values is optional as otherwise, at the expense of efficiency, the values will be computed each time the function is invoked. See Buildings.Utilities.Math.Functions.Examples.InverseXRegularized for how to efficiently call this function.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realx Abscissa value
Realdelta Abscissa value below which approximation occurs
RealdeltaInv1/deltaInverse value of delta
Reala-15*deltaInvPolynomial coefficient
Realb119*deltaInv^2Polynomial coefficient
Realc-361*deltaInv^3Polynomial coefficient
Reald534*deltaInv^4Polynomial coefficient
Reale-380*deltaInv^5Polynomial coefficient
Realf104*deltaInv^6Polynomial coefficient

Outputs

TypeNameDescription
RealyFunction value

Modelica definition

function inverseXRegularized "Function that approximates 1/x by a twice continuously differentiable function" annotation(derivative=Buildings.Utilities.Math.Functions.BaseClasses.der_inverseXRegularized); extends Modelica.Icons.Function; input Real x "Abscissa value"; input Real delta(min=Modelica.Constants.eps) "Abscissa value below which approximation occurs"; input Real deltaInv=1/delta "Inverse value of delta"; input Real a=-15*deltaInv "Polynomial coefficient"; input Real b=119*deltaInv^2 "Polynomial coefficient"; input Real c=-361*deltaInv^3 "Polynomial coefficient"; input Real d=534*deltaInv^4 "Polynomial coefficient"; input Real e=-380*deltaInv^5 "Polynomial coefficient"; input Real f=104*deltaInv^6 "Polynomial coefficient"; output Real y "Function value"; algorithm y := if (x > delta or x < -delta) then 1/x elseif (x < delta/2 and x > -delta/ 2) then x/(delta*delta) else Buildings.Utilities.Math.Functions.BaseClasses.smoothTransition( x=x, delta=delta, deltaInv=deltaInv, a=a, b=b, c=c, d=d, e=e, f=f); end inverseXRegularized;

Buildings.Utilities.Math.Functions.isMonotonic Buildings.Utilities.Math.Functions.isMonotonic

Returns true if the argument is a monotonic sequence

Information

This function returns true if its argument is monotonic increasing or decreasing, and false otherwise. If strict=true, then strict monotonicity is tested, otherwise weak monotonicity is tested.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realx[:] Sequence to be tested
BooleanstrictfalseSet to true to test for strict monotonicity

Outputs

TypeNameDescription
BooleanmonotonicTrue if x is monotonic increasing or decreasing

Modelica definition

function isMonotonic "Returns true if the argument is a monotonic sequence" extends Modelica.Icons.Function; input Real x[:] "Sequence to be tested"; input Boolean strict=false "Set to true to test for strict monotonicity"; output Boolean monotonic "True if x is monotonic increasing or decreasing"; protected Integer n=size(x, 1) "Number of data points"; algorithm if n == 1 then monotonic := true; else monotonic := true; if strict then if (x[1] >= x[n]) then for i in 1:n - 1 loop if (not x[i] > x[i + 1]) then monotonic := false; end if; end for; else for i in 1:n - 1 loop if (not x[i] < x[i + 1]) then monotonic := false; end if; end for; end if; else // not strict if (x[1] >= x[n]) then for i in 1:n - 1 loop if (not x[i] >= x[i + 1]) then monotonic := false; end if; end for; else for i in 1:n - 1 loop if (not x[i] <= x[i + 1]) then monotonic := false; end if; end for; end if; end if; // strict end if; end isMonotonic;

Buildings.Utilities.Math.Functions.polynomial Buildings.Utilities.Math.Functions.polynomial

Polynomial function

Information

This function computes a polynomial of arbitrary order. The polynomial has the form

y = a1 + a2 x + a3 x2 + ...

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realx Independent variable
Reala[:] Coefficients

Outputs

TypeNameDescription
RealyResult

Modelica definition

function polynomial "Polynomial function" extends Modelica.Icons.Function; 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;

Buildings.Utilities.Math.Functions.powerLinearized Buildings.Utilities.Math.Functions.powerLinearized

Power function that is linearized below a user-defined threshold

Information

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.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realx Abscissa value
Realn Exponent
Realx0 Abscissa value below which linearization occurs

Outputs

TypeNameDescription
RealyFunction value

Modelica definition

function powerLinearized "Power function that is linearized below a user-defined threshold" extends Modelica.Icons.Function; 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;

Buildings.Utilities.Math.Functions.quadraticLinear Buildings.Utilities.Math.Functions.quadraticLinear

Function that is quadratic in first argument and linear in second argument

Information

This function computes

y = a1 + a2 x1 + a3 x12 + (a4 + a5 x1 + a6 x12) x2

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Reala[6] Coefficients
Realx1 Independent variable for quadratic part
Realx2 Independent variable for linear part

Outputs

TypeNameDescription
RealyResult

Modelica definition

function quadraticLinear "Function that is quadratic in first argument and linear in second argument" extends Modelica.Icons.Function; 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;

Buildings.Utilities.Math.Functions.quinticHermite Buildings.Utilities.Math.Functions.quinticHermite

Quintic Hermite spline function for interpolation between two functions with a continuous second derivative

Information

Returns the result y of a quintic Hermite spline, which is a C2 continuous interpolation between two functions. The abscissa value x has to be between x1 and x2. Variables y1, y1d, y1dd are the ordinate, ordinate derivative and ordinate second derivative of the function at x1. Variables y2, y2d, y2dd are respectively the ordinate, ordinate derivative and ordinate second derivative of the function at x2.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realx Abscissa value
Realx1 Lower abscissa value
Realx2 Upper abscissa value
Realy1 Lower ordinate value
Realy2 Upper ordinate value
Realy1d Lower derivative
Realy2d Upper derivative
Realy1dd Lower second derivative
Realy2dd Upper second derivative

Outputs

TypeNameDescription
RealyInterpolated ordinate value

Modelica definition

function quinticHermite "Quintic Hermite spline function for interpolation between two functions with a continuous second derivative" extends Modelica.Icons.Function; input Real x "Abscissa value"; input Real x1 "Lower abscissa value"; input Real x2 "Upper abscissa value"; input Real y1 "Lower ordinate value"; input Real y2 "Upper ordinate value"; input Real y1d "Lower derivative"; input Real y2d "Upper derivative"; input Real y1dd "Lower second derivative"; input Real y2dd "Upper second derivative"; output Real y "Interpolated ordinate value"; protected Real h=x2 - x1; Real hpow2=h*h; Real t=(x - x1)/h; Real tpow2=t*t; Real tpow3=tpow2*t; Real tpow4=tpow3*t; Real tpow5=tpow4*t; Real H0=1 - 10*tpow3 + 15*tpow4 - 6*tpow5; Real H1=t - 6*tpow3 + 8*tpow4 - 3*tpow5; Real H2=0.5*(tpow2 - 3*tpow3 + 3*tpow4 - tpow5); Real H3=0.5*tpow3 - tpow4 + 0.5*tpow5; Real H4=-4*tpow3 + 7*tpow4 - 3*tpow5; Real H5=1 - H0; algorithm y := H0*y1 + H1*y1d*h + H2*y1dd*hpow2 + H3*y2dd*hpow2 + H4*y2d*h + H5*y2; end quinticHermite;

Buildings.Utilities.Math.Functions.regNonZeroPower Buildings.Utilities.Math.Functions.regNonZeroPower

Power function, regularized near zero, but nonzero value for x=0

Information

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.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realx Abscissa value
Realn Exponent
Realdelta0.01Abscissa value where transition occurs

Outputs

TypeNameDescription
RealyFunction value

Modelica definition

function regNonZeroPower "Power function, regularized near zero, but nonzero value for x=0" annotation(derivative=BaseClasses.der_regNonZeroPower); extends Modelica.Icons.Function; 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;

Buildings.Utilities.Math.Functions.regStep Buildings.Utilities.Math.Functions.regStep

Approximation of a general step, such that the approximation is continuous and differentiable

Information

This function is used to approximate the equation

    y = if x > 0 then y1 else y2;

by a smooth characteristic, so that the expression is continuous and differentiable:

   y = smooth(1, if x >  x_small then y1 else
                 if x < -x_small then y2 else f(y1, y2));

In the region -x_small < x < x_small a 2nd order polynomial is used for a smooth transition from y1 to y2.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realx Abscissa value
Realy1 Ordinate value for x > 0
Realy2 Ordinate value for x < 0
Realx_small1e-5Approximation of step for -x_small <= x <= x_small; x_small >= 0 required

Outputs

TypeNameDescription
RealyOrdinate value to approximate y = if x > 0 then y1 else y2

Modelica definition

function regStep "Approximation of a general step, such that the approximation is continuous and differentiable" extends Modelica.Icons.Function; input Real x "Abscissa value"; input Real y1 "Ordinate value for x > 0"; input Real y2 "Ordinate value for x < 0"; input Real x_small(min=0) = 1e-5 "Approximation of step for -x_small <= x <= x_small; x_small >= 0 required"; output Real y "Ordinate value to approximate y = if x > 0 then y1 else y2"; algorithm y := smooth(1, if x > x_small then y1 else if x < -x_small then y2 else if x_small > 0 then (x/x_small)*((x/x_small)^2 - 3)*(y2 - y1)/4 + (y1 + y2)/2 else (y1 + y2)/2); end regStep;

Buildings.Utilities.Math.Functions.smoothExponential Buildings.Utilities.Math.Functions.smoothExponential

Once continuously differentiable approximation to exp(-|x|) in interval |x| < delta

Information

Function to provide a once continuously differentiable approximation to exp(- |x| ) in the interval |x| < δ for some positive δ

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realx Input argument
Realdelta Transition point where approximation occurs

Outputs

TypeNameDescription
RealyOutput argument

Modelica definition

function smoothExponential "Once continuously differentiable approximation to exp(-|x|) in interval |x| < delta" extends Modelica.Icons.Function; 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;

Buildings.Utilities.Math.Functions.smoothHeaviside Buildings.Utilities.Math.Functions.smoothHeaviside

Once continuously differentiable approximation to the Heaviside function

Information

Once Lipschitz continuously differentiable approximation to the Heaviside(.,.) function. See Example Buildings.Utilities.Math.Examples.SmoothHeaviside.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realx Argument
Realdelta Parameter used for scaling

Outputs

TypeNameDescription
RealyResult

Modelica definition

function smoothHeaviside "Once continuously differentiable approximation to the Heaviside function" extends Modelica.Icons.Function; input Real x "Argument"; input Real delta "Parameter used for scaling"; output Real y "Result"; algorithm y := Buildings.Utilities.Math.Functions.regStep( y1=1, y2=0, x=x, x_small=delta); end smoothHeaviside;

Buildings.Utilities.Math.Functions.smoothInterpolation Buildings.Utilities.Math.Functions.smoothInterpolation

Interpolate using a cubic Hermite spline with linear extrapolation for a vector xSup[], ySup[] and independent variable x

Information

For xSup1 ≤ x ≤ xSupn, where n is the size of the support points xSup, which must be strictly monotonically increasing, this function interpolates using cubic hermite spline. For x outside this interval, the function linearly extrapolates.

If n=2, linear interpolation is used an if n=1, the function value y1 is returned.

Note that if xSup and ySup only depend on parameters or constants, then Buildings.Utilities.Math.Functions.cubicHermiteLinearExtrapolation will be more efficient. In contrast to the function Modelica.Math.Vectors.interpolate which provides linear interpolation, this function does not trigger events.

For how to use this function, see Buildings.Utilities.Math.Functions.Examples.SmoothInterpolation.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realx Abscissa value
RealxSup[:] Support points (strictly increasing)
RealySup[size(xSup, 1)] Function values at xSup
BooleanensureMonotonicityisMonotonic(ySup, strict=fal...Set to true to ensure monotonicity of the cubic hermite

Outputs

TypeNameDescription
RealyIntInterpolated ordinate value

Modelica definition

function smoothInterpolation "Interpolate using a cubic Hermite spline with linear extrapolation for a vector xSup[], ySup[] and independent variable x" extends Modelica.Icons.Function; input Real x "Abscissa value"; input Real xSup[:] "Support points (strictly increasing)"; input Real ySup[size(xSup, 1)] "Function values at xSup"; input Boolean ensureMonotonicity=isMonotonic(ySup, strict=false) "Set to true to ensure monotonicity of the cubic hermite"; output Real yInt "Interpolated ordinate value"; protected Integer n=size(xSup, 1) "Number of support points"; Real dy_dx[size(xSup, 1)] "Derivative at xSup"; Integer i "Integer to select data interval"; algorithm if n > 2 then // Most common case with more than 2 data points. // Do cubic spline interpolation. dy_dx := Buildings.Utilities.Math.Functions.splineDerivatives( x=xSup, y=ySup, ensureMonotonicity=ensureMonotonicity); // i is a counter that is used to pick the derivative // that corresponds to the interval that contains x i := 1; for j in 1:n - 1 loop if x > xSup[j] then i := j; end if; end for; assert(xSup[i] < xSup[i + 1], "Support points xSup must be increasing."); yInt := cubicHermiteLinearExtrapolation( x=x, x1=xSup[i], x2=xSup[i + 1], y1=ySup[i], y2=ySup[i + 1], y1d=dy_dx[i], y2d=dy_dx[i + 1]); elseif n == 2 then // Linear interpolation. yInt := ySup[1] + (x - xSup[1])*(ySup[2] - ySup[1])/(xSup[2] - xSup[1]); else // n == 1 yInt := ySup[1]; end if; end smoothInterpolation;

Buildings.Utilities.Math.Functions.smoothLimit Buildings.Utilities.Math.Functions.smoothLimit

Once continuously differentiable approximation to the limit function

Information

Once continuously differentiable approximation to the limit(.,.) function. The output is bounded to be in [l, u].

Note that the limit need not be respected, such as illustrated in Buildings.Utilities.Math.Examples.SmoothMin.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realx Variable
Reall Low limit
Realu Upper limit
RealdeltaX Width of transition interval

Outputs

TypeNameDescription
RealyResult

Modelica definition

function smoothLimit "Once continuously differentiable approximation to the limit function" extends Modelica.Icons.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;

Buildings.Utilities.Math.Functions.smoothMax Buildings.Utilities.Math.Functions.smoothMax

Once continuously differentiable approximation to the maximum function

Information

Once continuously differentiable approximation to the max(.,.) function.

Note that the maximum need not be respected, such as illustrated in Buildings.Utilities.Math.Examples.SmoothMin.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realx1 First argument
Realx2 Second argument
RealdeltaX Width of transition interval

Outputs

TypeNameDescription
RealyResult

Modelica definition

function smoothMax "Once continuously differentiable approximation to the maximum function" extends Modelica.Icons.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.regStep( y1=x1, y2=x2, x=x1 - x2, x_small=deltaX); end smoothMax;

Buildings.Utilities.Math.Functions.smoothMin Buildings.Utilities.Math.Functions.smoothMin

Once continuously differentiable approximation to the minimum function

Information

Once continuously differentiable approximation to the min(.,.) function.

Note that the minimum need not be respected, such as illustrated in Buildings.Utilities.Math.Examples.SmoothMin.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realx1 First argument
Realx2 Second argument
RealdeltaX Width of transition interval

Outputs

TypeNameDescription
RealyResult

Modelica definition

function smoothMin "Once continuously differentiable approximation to the minimum function" extends Modelica.Icons.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.regStep( y1=x1, y2=x2, x=x2 - x1, x_small=deltaX); end smoothMin;

Buildings.Utilities.Math.Functions.spliceFunction Buildings.Utilities.Math.Functions.spliceFunction

Information

Function to provide a once continuously differentiable transition between to arguments.

The function is adapted from Modelica.Media.Air.MoistAir.Utilities.spliceFunction and provided here for easier accessability to model developers.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realpos Argument of x > 0
Realneg Argument of x < 0
Realx Independent value
Realdeltax Half width of transition interval

Outputs

TypeNameDescription
RealoutSmoothed value

Modelica definition

function spliceFunction annotation(derivative=BaseClasses.der_spliceFunction); extends Modelica.Icons.Function; 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 scaledX1; Real y; constant Real asin1=Modelica.Math.asin(1); algorithm scaledX1 := x/deltax; if scaledX1 <= -0.999999999 then out := neg; elseif scaledX1 >= 0.999999999 then out := pos; else y := (Modelica.Math.tanh(Modelica.Math.tan(scaledX1*asin1)) + 1)/2; out := pos*y + (1 - y)*neg; end if; end spliceFunction;

Buildings.Utilities.Math.Functions.splineDerivatives Buildings.Utilities.Math.Functions.splineDerivatives

Function to compute the derivatives for cubic hermite spline interpolation

Information

This function computes the derivatives at the support points xi that can be used as input for evaluating a cubic hermite spline.

If ensureMonotonicity=true, then the support points yi need to be monotone increasing (or increasing), and the computed derivatives di are such that the cubic hermite is monotone increasing (or decreasing). The algorithm to ensure monotonicity is based on the method described in Fritsch and Carlson (1980) for ρ = ρ2.

This function is typically used with Buildings.Utilities.Math.Functions.cubicHermiteLinearExtrapolation which is used to evaluate the cubic spline. Because in many applications, the shape of the spline depends on parameters, this function has been implemented in such a way that all derivatives can be computed at once and then stored for use during the time stepping, in which the above function may be called.

References

F.N. Fritsch and R.E. Carlson, Monotone piecewise cubic interpolation. SIAM J. Numer. Anal., 17 (1980), pp. 238-246.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realx[:] Support point, strict monotone increasing
Realy[size(x, 1)] Function values at x
BooleanensureMonotonicityisMonotonic(y, strict=false)Set to true to ensure monotonicity of the cubic hermite

Outputs

TypeNameDescription
Reald[size(x, 1)]Derivative at the support points

Modelica definition

function splineDerivatives "Function to compute the derivatives for cubic hermite spline interpolation" extends Modelica.Icons.Function; input Real x[:] "Support point, strict monotone increasing"; input Real y[size(x, 1)] "Function values at x"; input Boolean ensureMonotonicity=isMonotonic(y, strict=false) "Set to true to ensure monotonicity of the cubic hermite"; output Real d[size(x, 1)] "Derivative at the support points"; protected Integer n=size(x, 1) "Number of data points"; Real delta[n - 1] "Slope of secant line between data points"; Real alpha "Coefficient to ensure monotonicity"; Real beta "Coefficient to ensure monotonicity"; Real tau "Coefficient to ensure monotonicity"; algorithm if (n > 1) then assert(x[1] < x[n], "x must be strictly increasing. Received x[1] = " + String(x[1]) + " x[" + String(n) + "] = " + String(x[n])); // Check data assert(isMonotonic(x, strict=true), "x-values must be strictly monontone increasing or decreasing."); if ensureMonotonicity then assert(isMonotonic(y, strict=false), "If ensureMonotonicity=true, y-values must be monontone increasing or decreasing."); end if; end if; // Compute derivatives at the support points if n == 1 then // only one data point d[1] := 0; elseif n == 2 then // linear function d[1] := (y[2] - y[1])/(x[2] - x[1]); d[2] := d[1]; else // Slopes of the secant lines between i and i+1 for i in 1:n - 1 loop delta[i] := (y[i + 1] - y[i])/(x[i + 1] - x[i]); end for; // Initial values for tangents at the support points. // End points use one-sided derivatives d[1] := delta[1]; d[n] := delta[n - 1]; for i in 2:n - 1 loop d[i] := (delta[i - 1] + delta[i])/2; end for; end if; // Ensure monotonicity if n > 2 and ensureMonotonicity then for i in 1:n - 1 loop if (abs(delta[i]) < Modelica.Constants.small) then d[i] := 0; d[i + 1] := 0; else alpha := d[i]/delta[i]; beta := d[i + 1]/delta[i]; // Constrain derivative to ensure monotonicity in this interval if (alpha^2 + beta^2) > 9 then tau := 3/(alpha^2 + beta^2)^(1/2); d[i] := delta[i]*alpha*tau; d[i + 1] := delta[i]*beta*tau; end if; end if; end for; end if; end splineDerivatives;

Buildings.Utilities.Math.Functions.trapezoidalIntegration Buildings.Utilities.Math.Functions.trapezoidalIntegration

Integration using the trapezoidal rule

Information

This function computes a definite integral using the trapezoidal rule.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
IntegerN Number of integrand points
Realf[:] Integrands
RealdeltaX Width of interval for Trapezoidal integration

Outputs

TypeNameDescription
RealresultResult

Modelica definition

function trapezoidalIntegration "Integration using the trapezoidal rule" extends Modelica.Icons.Function; 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;