ModelingToolkitStandardLibrary: Blocks
Index
ModelingToolkitStandardLibrary.Blocks.AddModelingToolkitStandardLibrary.Blocks.Add3ModelingToolkitStandardLibrary.Blocks.Atan2ModelingToolkitStandardLibrary.Blocks.CeilModelingToolkitStandardLibrary.Blocks.ConstantModelingToolkitStandardLibrary.Blocks.DeadZoneModelingToolkitStandardLibrary.Blocks.DerivativeModelingToolkitStandardLibrary.Blocks.DivisionModelingToolkitStandardLibrary.Blocks.FeedbackModelingToolkitStandardLibrary.Blocks.FirstOrderModelingToolkitStandardLibrary.Blocks.FloorModelingToolkitStandardLibrary.Blocks.GainModelingToolkitStandardLibrary.Blocks.IntegratorModelingToolkitStandardLibrary.Blocks.MatrixGainModelingToolkitStandardLibrary.Blocks.ModuloModelingToolkitStandardLibrary.Blocks.PIModelingToolkitStandardLibrary.Blocks.PowerModelingToolkitStandardLibrary.Blocks.ProductModelingToolkitStandardLibrary.Blocks.SISOModelingToolkitStandardLibrary.Blocks.SecondOrderModelingToolkitStandardLibrary.Blocks.StaticNonLinearityModelingToolkitStandardLibrary.Blocks.SumModelingToolkitStandardLibrary.Blocks.UnaryMinusModelingToolkitStandardLibrary.Blocks.AbsModelingToolkitStandardLibrary.Blocks.AcosModelingToolkitStandardLibrary.Blocks.AsinModelingToolkitStandardLibrary.Blocks.AtanModelingToolkitStandardLibrary.Blocks.ContinuousClockModelingToolkitStandardLibrary.Blocks.CosModelingToolkitStandardLibrary.Blocks.CoshModelingToolkitStandardLibrary.Blocks.CosineModelingToolkitStandardLibrary.Blocks.ExpModelingToolkitStandardLibrary.Blocks.ExpSineModelingToolkitStandardLibrary.Blocks.LimPIModelingToolkitStandardLibrary.Blocks.LimPIDModelingToolkitStandardLibrary.Blocks.LimiterModelingToolkitStandardLibrary.Blocks.LogModelingToolkitStandardLibrary.Blocks.Log10ModelingToolkitStandardLibrary.Blocks.MIMOModelingToolkitStandardLibrary.Blocks.PIDModelingToolkitStandardLibrary.Blocks.RampModelingToolkitStandardLibrary.Blocks.RealInputModelingToolkitStandardLibrary.Blocks.RealInputArrayModelingToolkitStandardLibrary.Blocks.RealOutputModelingToolkitStandardLibrary.Blocks.RealOutputArrayModelingToolkitStandardLibrary.Blocks.SignModelingToolkitStandardLibrary.Blocks.SinModelingToolkitStandardLibrary.Blocks.SineModelingToolkitStandardLibrary.Blocks.SinhModelingToolkitStandardLibrary.Blocks.SlewRateLimiterModelingToolkitStandardLibrary.Blocks.SqrtModelingToolkitStandardLibrary.Blocks.StateSpaceModelingToolkitStandardLibrary.Blocks.StepModelingToolkitStandardLibrary.Blocks.TanModelingToolkitStandardLibrary.Blocks.TanhModelingToolkitStandardLibrary.Blocks.TransferFunction
Utility Blocks
ModelingToolkitStandardLibrary.Blocks.RealInput — FunctionRealInput(;name, guess)Connector with one input signal of type Real.
Parameters:
guess=0: Guess value foru.
States:
u: Value of the connector which is a scalar.
ModelingToolkitStandardLibrary.Blocks.RealOutput — FunctionRealOutput(;name, guess)Connector with one output signal of type Real.
Parameters:
guess=0: Guess value foru.
States:
u: Value of the connector which is a scalar.
ModelingToolkitStandardLibrary.Blocks.RealInputArray — FunctionRealInputArray(;name, nin, guess)Connector with an array of input signals of type Real.
Parameters:
nin: Number of inputs.guess=zeros(nin): Guess value foru.
States:
u: Value of the connector which is an array.
ModelingToolkitStandardLibrary.Blocks.RealOutputArray — FunctionRealOutputArray(;name, nout, guess)Connector with an array of output signals of type Real.
Parameters:
nout: Number of outputs.guess=zeros(nout): Guess value foru.
States:
u: Value of the connector which is an array.
ModelingToolkitStandardLibrary.Blocks.SISO — ConstantSISO(;name, u_start = 0.0, y_start = 0.0)Single input single output (SISO) continuous system block.
Parameters:
u_start: Initial value for the inputy_start: Initial value for the output
ModelingToolkitStandardLibrary.Blocks.MIMO — FunctionMIMO(; name, nin = 1, nout = 1, u_start = zeros(nin), y_start = zeros(nout))Base class for a multiple input multiple output (MIMO) continuous system block.
Parameters:
nin: Input dimensionnout: Output dimensionu_start: Initial value for the inputy_start: Initial value for the output
Math Blocks
ModelingToolkitStandardLibrary.Blocks.Gain — ConstantGain(; name, k)Output the product of a gain value with the input signal.
Parameters:
k: Scalar gain
Connectors:
inputoutput
ModelingToolkitStandardLibrary.Blocks.MatrixGain — ConstantMatrixGain(; K::AbstractArray, name)Output the product of a gain matrix with the input signal vector.
Structural parameters:
K: Matrix gain
Connectors:
inputoutput
ModelingToolkitStandardLibrary.Blocks.Sum — ConstantSum(; input__nin::Int, name)Output the sum of the elements of the input port vector. Input port dimension can be set with input__nin
Connectors:
inputoutput
ModelingToolkitStandardLibrary.Blocks.Feedback — ConstantFeedback(; name)Output difference between reference input (input1) and feedback input (input2).
Connectors:
input1input2output
ModelingToolkitStandardLibrary.Blocks.Add — ConstantAdd(; name, k1 = 1.0, k2 = 1.0)Output the sum of the two scalar inputs.
Parameters:
k1: Gain for first inputk2: Gain for second input
Connectors:
input1input2output
ModelingToolkitStandardLibrary.Blocks.Add3 — ConstantAdd(; name, k1 = 1.0, k2 = 1.0, k3 = 1.0)Output the sum of the three scalar inputs.
Parameters:
k1: Gain for first inputk2: Gain for second inputk3: Gain for third input
Connectors:
input1input2input3output
ModelingToolkitStandardLibrary.Blocks.Product — ConstantProduct(; name)Output product of the two inputs.
Connectors:
input1input2output
ModelingToolkitStandardLibrary.Blocks.Division — ConstantDivision(; name)Output first input divided by second input.
Connectors:
input1input2output
ModelingToolkitStandardLibrary.Blocks.UnaryMinus — ConstantUnaryMinus(; name)Output the product of -1 and the input.
Connectors:
inputoutput
ModelingToolkitStandardLibrary.Blocks.Power — ConstantPower(; name)Output the exponential with base as the first input and exponent as second input i.e u1^u2
Connectors:
baseexponentoutput
ModelingToolkitStandardLibrary.Blocks.Modulo — ConstantModulo(; name)Output the remainder when the first input is divided by second input.
Connectors:
dividenddivisorremainder
ModelingToolkitStandardLibrary.Blocks.Floor — ConstantFloor(; name)Output the floor rounding of the input.
Connectors:
inputoutput
ModelingToolkitStandardLibrary.Blocks.Ceil — ConstantCeil(; name)Output the ceiling rounding of the input.
Connectors:
inputoutput
ModelingToolkitStandardLibrary.Blocks.StaticNonLinearity — ConstantStaticNonLinearity(func; name)Applies the given function to the input.
If the given function is not composed of simple core methods (e.g. sin, abs, ...), it has to be registered via @register_symbolic func(u)
Connectors:
inputoutput
ModelingToolkitStandardLibrary.Blocks.Abs — FunctionModelingToolkitStandardLibrary.Blocks.Sign — FunctionModelingToolkitStandardLibrary.Blocks.Sqrt — FunctionSqrt(; name)Output the square root of the input (input >= 0 required).
Connectors:
ModelingToolkitStandardLibrary.Blocks.Sin — FunctionModelingToolkitStandardLibrary.Blocks.Cos — FunctionModelingToolkitStandardLibrary.Blocks.Tan — FunctionModelingToolkitStandardLibrary.Blocks.Asin — FunctionModelingToolkitStandardLibrary.Blocks.Acos — FunctionModelingToolkitStandardLibrary.Blocks.Atan — FunctionModelingToolkitStandardLibrary.Blocks.Atan2 — ConstantAtan2(; name)Output the arc tangent of the input.
Connectors:
input1input2output
ModelingToolkitStandardLibrary.Blocks.Sinh — FunctionModelingToolkitStandardLibrary.Blocks.Cosh — FunctionModelingToolkitStandardLibrary.Blocks.Tanh — FunctionModelingToolkitStandardLibrary.Blocks.Exp — FunctionModelingToolkitStandardLibrary.Blocks.Log — FunctionModelingToolkitStandardLibrary.Blocks.Log10 — FunctionSource Blocks
ModelingToolkitStandardLibrary.Blocks.Constant — ConstantConstant(; name, k = 0.0)Generate constant signal.
Parameters:
k: Constant output value
Connectors:
output
ModelingToolkitStandardLibrary.Blocks.Sine — FunctionSine(; name, frequency, amplitude = 1, phase = 0, offset = 0, start_time = 0,
smooth = false)Generate sine signal.
Parameters:
frequency: [Hz] Frequency of sine waveamplitude: Amplitude of sine wavephase: [rad] Phase of sine waveoffset: Offset of output signalstart_time: [s] Outputy = offsetfort < start_timesmooth: Iftrue, returns a smooth wave. Defaults tofalseIt uses a default smoothing factor ofδ=1e-5, but this can be changed by supplyingsmooth=δ.
Connectors:
output
ModelingToolkitStandardLibrary.Blocks.Cosine — FunctionCosine(; name, frequency, amplitude = 1, phase = 0, offset = 0, start_time = 0,
smooth = false)Generate cosine signal.
Parameters:
frequency: [Hz] Frequency of cosine waveamplitude: Amplitude of cosine wavephase: [rad] Phase of cosine waveoffset: Offset of output signalstart_time: [s] Outputy = offsetfort < start_timesmooth: Iftrue, returns a smooth wave. Defaults tofalseIt uses a default smoothing factor ofδ=1e-5, but this can be changed by supplyingsmooth=δ.
Connectors:
output
ModelingToolkitStandardLibrary.Blocks.ContinuousClock — FunctionContinuousClock(; name, offset = 0, start_time = 0)Generate current time signal.
Parameters:
offset: Offset of output signalstart_time: [s] Outputy = offsetfort < start_time
Connectors:
output
ModelingToolkitStandardLibrary.Blocks.Ramp — FunctionRamp(; name, height = 1, duration = 1, offset = 0, start_time = 0, smooth = false)
Generate ramp signal.
Parameters:
height: Height of rampduration: [s] Duration of ramp (= 0.0 gives a Step)offset: Offset of output signalstart_time: [s] Outputy = offsetfort < start_timesmooth: Iftrue, returns a smooth wave. Defaults tofalseIt uses a default smoothing factor ofδ=1e-5, but this can be changed by supplyingsmooth=δ.
Connectors:
output
ModelingToolkitStandardLibrary.Blocks.Step — FunctionStep(;name, height=1, offset=0, start_time=0, duration=Inf, smooth=true)Generate step signal.
Parameters:
height: Height of stepoffset: Offset of output signalstart_time: [s] Outputy = offsetfort < start_timeand thereafteroffset+height.duration: [s] Ifduration < Infis supplied, the output will revert tooffsetafterdurationseconds.smooth: Iftrue, returns a smooth wave. Defaults totrueIt uses a default smoothing factor ofδ=1e-5, but this can be changed by supplyingsmooth=δ.
Connectors:
output
ModelingToolkitStandardLibrary.Blocks.ExpSine — FunctionExpSine(; name, frequency, amplitude = 1, damping = 0.1, phase = 0, offset = 0, start_time = 0, smooth = false)Exponentially damped sine signal.
Parameters:
frequency: [Hz] Frequency of sine waveamplitude: Amplitude of sine wavedamping: [1/s] Damping coefficient of sine wavephase: [rad] Phase of sine waveoffset: Offset of output signalstart_time: [s] Outputy = offsetfort < start_timesmooth: Iftrue, returns a smooth wave. Defaults tofalseIt uses a default smoothing factor ofδ=1e-5, but this can be changed by supplyingsmooth=δ.
Connectors:
output
Nonlinear Blocks
ModelingToolkitStandardLibrary.Blocks.Limiter — FunctionLimiter(;name, y_max, y_min = y_max > 0 ? -y_max : -Inf)Limit the range of a signal.
Parameters:
y_max: Maximum of output signaly_min: Minimum of output signal
Connectors:
inputoutput
ModelingToolkitStandardLibrary.Blocks.DeadZone — ConstantDeadZone(; name, u_max, u_min = -u_max)The DeadZone block defines a region of zero output. If the input is within u_min ... u_max, the output is zero. Outside of this zone, the output is a linear function of the input with a slope of 1.
y▲
│ /
│ /
u_min │ /
─────|──┼──|───────► u
/ │ u_max
/ │
/ │Parameters:
u_max: Upper limit of dead zoneu_min: Lower limit of dead zone
Connectors:
inputoutput
ModelingToolkitStandardLibrary.Blocks.SlewRateLimiter — FunctionSlewRateLimiter(; name, y_start, rising = 1.0, falling = -rising, Td = 0.001)Limits the slew rate of a signal. Initial value of state Y can be set with int.y
Parameters:
rising: Maximum rising slew ratefalling: Maximum falling slew rateTd: [s] Derivative time constanty_start: Initial value ofystate of SISO
Connectors:
inputoutput
Continuous Blocks
ModelingToolkitStandardLibrary.Blocks.Integrator — ConstantIntegrator(;name, k = 1, x = 0.0)Outputs y = ∫k*u dt, corresponding to the transfer function $1/s$. Initial value of integrator state $x$ can be set with x
Connectors:
inputoutput
Parameters:
k: Gain of integrator
Unknowns:
x: State of Integrator. Defaults to 0.0.
ModelingToolkitStandardLibrary.Blocks.Derivative — ConstantDerivative(; name, k = 1, T, x = 0.0)Outputs an approximate derivative of the input. The transfer function of this block is
k k ks
─ - ─────── = ──────
T sT² + T sT + 1and a state-space realization is given by ss(-1/T, 1/T, -k/T, k/T) where T is the time constant of the filter. A smaller T leads to a more ideal approximation of the derivative.
Initial value of the state $x$ can be set with x.
Parameters:
k: GainT: [s] Time constant (T>0 required; T=0 is ideal derivative block)
Unknowns:
x: Unknown of Derivative. Defaults to 0.0.
Connectors:
inputoutput
ModelingToolkitStandardLibrary.Blocks.FirstOrder — ConstantFirstOrder(; name, k = 1.0, T, x = 0.0, lowpass = true)A first-order filter with a single real pole at s = -1/T and gain k. If lowpass=true (default), the transfer function is given by $Y(s)/U(s) =$
k
───────
sT + 1and if lowpass=false, by
sT + 1 - k
──────────
sT + 1Initial value of the state x can be set with x
Parameters:
k: GainT: [s] Time constant (T>0 required)
Connectors:
inputoutput
See also SecondOrder
ModelingToolkitStandardLibrary.Blocks.SecondOrder — ConstantSecondOrder(; name, k = 1.0, w = 1.0, d = 1.0, x = 0.0, xd = 0.0)A second-order filter with gain k, a bandwidth of w rad/s and relative damping d. The transfer function is given by Y(s)/U(s) =
k*w^2
─────────────────
s² + 2d*w*s + w^2Critical damping corresponds to d=1, which yields the fastest step response without overshoot, d < 1 results in an underdamped filter while d > 1 results in an overdamped filter. d = 1/√2 corresponds to a Butterworth filter of order 2 (maximally flat frequency response). Initial value of the state x can be set with x, and of derivative state xd with xd.
Parameters:
k: Gainw: [rad/s] Angular frequencyd: Damping
Connectors:
inputoutput
ModelingToolkitStandardLibrary.Blocks.StateSpace — FunctionStateSpace(A, B, C, D = 0; x = zeros(size(A,1)), u0 = zeros(size(B,2)), y0 = zeros(size(C,1)), name)A linear, time-invariant state-space system on the form.
\[\begin{aligned} ẋ &= Ax + Bu \\ y &= Cx + Du \end{aligned}\]
Transfer functions can also be simulated by converting them to a StateSpace form.
y0 and u0 can be used to set an operating point, providing them changes the dynamics from an LTI system to the affine system
\[\begin{aligned} ẋ &= Ax + B(u - u0) \\ y &= Cx + D(u - u0) + y0 \end{aligned}\]
For a nonlinear system
\[\begin{aligned} ẋ &= f(x, u) \\ y &= h(x, u) \end{aligned}\]
linearized around the operating point x₀, u₀, we have y0, u0 = h(x₀, u₀), u₀.
ModelingToolkitStandardLibrary.Blocks.TransferFunction — FunctionTransferFunction(; b, a, name)A single input, single output, linear time-invariant system provided as a transfer-function.
Y(s) = b(s) / a(s) U(s)where b and a are vectors of coefficients of the numerator and denominator polynomials, respectively, ordered such that the coefficient of the highest power of s is first.
The internal state realization is on controller canonical form, with state variable x, output variable y and input variable u. For numerical robustness, the realization used by the integrator is scaled by the last entry of the a parameter. The internally scaled state variable is available as x_scaled.
To set the initial state, it's recommended to set the initial condition for x, and let that of x_scaled be computed automatically.
Parameters:
b: Numerator polynomial coefficients, e.g.,2s + 3is specified as[2, 3]a: Denominator polynomial coefficients, e.g.,s² + 2ωs + ω^2is specified as[1, 2ω, ω^2]
Connectors:
inputoutput
See also StateSpace which handles MIMO systems, as well as ControlSystemsMTK.jl for an interface between ControlSystems.jl and ModelingToolkit.jl for advanced manipulation of transfer functions and linear statespace systems. For linearization, see linearize and Linear Analysis.
ModelingToolkitStandardLibrary.Blocks.PI — ConstantPI(;name, k = 1.0, T = 1.0, int.x = 0.0)Textbook version of a PI-controller without actuator saturation and anti-windup measure. The proportional gain can be set with k Initial value of integrator state x can be set with int.x
The PI controller is implemented on standard form:
\[U(s) = k (1 + \dfrac{1}{sT}) E(S)\]
Parameters:
k: Proportional gainT: [s] Integrator time constant (T>0 required)
Connectors:
err_inputctr_output
See also LimPI
ModelingToolkitStandardLibrary.Blocks.LimPI — FunctionLimPI(; name, k = 1.0, T, Ta, int__x = 0.0, u_max = 1.0, u_min = -u_max)Text-book version of a PI-controller with actuator saturation and anti-windup measure.
The PI controller is implemented on standard form
\[u(t) = sat(k (e(t) + ∫\dfrac{1}{T}e(t) dt) )\]
The simplified expression above is given without the anti-windup protection.
Parameters:
k: Proportional gainT: [s] Integrator time constant (T>0 required)Ta: [s] Tracking time constant (Ta>0 required)
Connectors:
err_inputctr_output
ModelingToolkitStandardLibrary.Blocks.PID — FunctionPID(;name, k=1, Ti=false, Td=false, Nd=10, int__x=0, der__x=0)Text-book version of a PID-controller without actuator saturation and anti-windup measure.
Parameters:
k: GainTi: [s] Integrator time constant (Ti>0 required). If set to false, no integral action is used.Td: [s] Derivative time constant (Td>0 required). If set to false, no derivative action is used.Nd: [s] Time constant for the derivative approximation (Nd>0 required; Nd=0 is ideal derivative).int__x: Initial value for the integrator.der__x: Initial value for the derivative state.
Connectors:
err_inputctr_output
See also LimPID
ModelingToolkitStandardLibrary.Blocks.LimPID — FunctionLimPID(; k, Ti=false, Td=false, wp=1, wd=1, Ni, Nd=12, u_max=Inf, u_min=-u_max, gains = false, name)Proportional-Integral-Derivative (PID) controller with output saturation, set-point weighting and integrator anti-windup.
The equation for the control signal is roughly
k(ep + 1/Ti * ∫e + Td * d/dt(ed))
e = u_r - u_y
ep = wp*u_r - u_y
ed = wd*u_r - u_ywhere the transfer function for the derivative includes additional filtering, see ? Derivative for more details.
Parameters:
k: Proportional gainTi: [s] Integrator time constant. Set tofalseto turn off integral action.Td: [s] Derivative time constant. Set tofalseto turn off derivative action.wp: [0,1] Set-point weighting in the proportional part.wd: [0,1] Set-point weighting in the derivative part.Nd: [1/s] Derivative limit, limits the derivative gain to Nd/Td. Reasonable values are ∈ [8, 20]. A higher value gives a better approximation of an ideal derivative at the expense of higher noise amplification.Ni:Ni*Ticontrols the time constantTaof anti-windup tracking. A common (default) choice isTa = √(Ti*Td)which is realized byNi = √(Td / Ti). Anti-windup can be effectively turned off by settingNi = Inf.gains: Ifgains = true,TiandTdwill be interpreted as gains with a fundamental PID transfer function on parallel formki=Ti, kd=Td, k + ki/s + kd*s.
Connectors:
referencemeasurementctr_output