Discrete block library
Index
ModelingToolkitSampledData.ClockChanger
ModelingToolkitSampledData.Delay
ModelingToolkitSampledData.Difference
ModelingToolkitSampledData.DiscreteDerivative
ModelingToolkitSampledData.DiscreteIntegrator
ModelingToolkitSampledData.DiscreteOnOffController
ModelingToolkitSampledData.DiscretePIDParallel
ModelingToolkitSampledData.DiscretePIDStandard
ModelingToolkitSampledData.DiscreteSlewRateLimiter
ModelingToolkitSampledData.ExponentialFilter
ModelingToolkitSampledData.MovingAverageFilter
ModelingToolkitSampledData.NormalNoise
ModelingToolkitSampledData.Quantization
ModelingToolkitSampledData.SampleWithADEffects
ModelingToolkitSampledData.Sampler
ModelingToolkitSampledData.UniformNoise
ModelingToolkitSampledData.ZeroOrderHold
ModelingToolkitSampledData.seeded_rand
ModelingToolkitSampledData.seeded_randn
Basic blocks
Noise and measurement corruption
Controllers
Discrete-time filters
Docstrings
ModelingToolkitSampledData.ClockChanger
— ConstantClockChanger()
ClockChanger(; to::AbstractClock, from::AbstractClock)
Connectors:
input
(continuous-time signal)output
(discrete-time signal)
The ClockChanger component is experimental and has known correctness issues. Please use with caution.
ModelingToolkitSampledData.Delay
— ConstantDelay(; n = 1)
A discrete-time delay of n
samples, corresponding to the transfer function $z^{-n}$.
Connectors:
input
output
Structural Parameters:
n
: Number of delay samples
ModelingToolkitSampledData.Difference
— ConstantDifference()
A discrete-time difference operator, corresponding to the transfer function $1 - z^{-1}$.
Connectors:
input
output
For a discrete-time finite-difference derivative approximation, see DiscreteDerivative
.
ModelingToolkitSampledData.DiscreteDerivative
— ConstantDiscreteDerivative(; k = 1)
A discrete-time derivative filter, corresponding to the transfer function $k (z-1) / (T_s z)$
where $T_s$ is the sample time of the derivative filter.
Connectors:
input
output
Parameters:
k
: Gain of derivative filter
ModelingToolkitSampledData.DiscreteIntegrator
— ConstantDiscreteIntegrator(;name, k = 1, x = 0.0, method = :backward)
Outputs y = ∫k*u dt
, discretized to correspond to the one of the discrete-time transfer functions
method = :forward
: $T_s / (z - 1)$method = :backward
: $T_s z / (z - 1)$method = :trapezoidal
: $(T_s / 2) (z + 1) / (z - 1)$
where $T_s$ is the sample time of the integrator.
Initial value of integrator state $x$ can be set with x
Connectors:
input
output
Parameters:
k
: Gain of integrator
ModelingToolkitSampledData.DiscreteOnOffController
— ConstantDiscreteOnOffController(b = 0.1, bool = true)
Discrete-time On-Off controller with hysteresis. The controller switches between two states based on the error signal reference-input
. The controller is in the on-state if the error signal is within the bandwidth b
around the reference signal, and in the off-state otherwise.
Connectors:
reference
: The reference signal to the controllerinput
: The measurement feedbackoutput
: The control signal output
Parameters:
b
: Bandwidth around reference signal within which the controller does not reactbool
: (structural) If true (default), the controller switches between 0 and 1. If false, the controller switches between -1 and 1.k
: Controller gain. The output of the controller is scaled by this gain, i.e.,k = 2, bool = false
will result in an output of -2 or 2.
ModelingToolkitSampledData.DiscretePIDParallel
— ConstantDiscretePIDParallel(;name, kp = 1, ki = 1, kd = 1, Ni = √(max(kd * ki, 1e-6)), Nd = 10kp, u_max = Inf, u_min = -u_max, wp = 1, wd = 1, Ts = 1, with_I = true, with_D = true, Imethod = :forward, Dmethod = :backward)
Discrete-time PID controller on parallel form with anti-windup and set-point weighting.
The controller is implemented on parallel form:
Simplified:
\[u = k_p e + \int k_i e dt + k_d \dfrac{de}{dt} \]
Detailed:
\[u = k_p(w_p r - y) + \int \big( k_i (r - y) + N_i e_s \big ) dt + k_d \dfrac{d}{dt}(w_d r - y)\]
where e_s = u - v
is the saturated error signal, v
is the unsaturated control signal and u
is the saturated control signal.
The derivative is filtered to allow a maximum gain of $N_d$.
The integrator is discretized using the method specified by Imethod
, options include
Imethod = :forward
(default): Corresponding to the transfer function $T_s / (z - 1)$Imethod = :backward
: Corresponding to the transfer function $T_s z / (z - 1)$Imethod = :trapezoidal
: Corresponding to the transfer function $(T_s / 2) (z + 1) / (z - 1)$
The derivative is discretized using the method specified by Dmethod
, options include
Dmethod = :forward
: Corresponding to the transfer function $\dfrac{N (z-1)}{z - \dfrac{k_d-N T_s}{k_d}}$.Dmethod = :backward
(default): Corresponding to the transfer function $\dfrac{\dfrac{Nk_d}{k_d + N T_s}(z-1)}{z - \dfrac{k_d}{k_d + N T_s}}$
Anti windup is realized by tracking using the gain $N_i$ on the error signal $e_s$ when the output is saturated.
To use the controller in 1DOF mode, i.e., with only the control error as input, connect the error signal to the reference
connector, connect a Constant(; k = 0)
to the measurement
connector and set wp = wd = 1
.
Connectors:
reference
: The reference signal to the controller (or the error signal if used in 1DOF mode)measurement
: The measurement feedbackctr_output
: The control signal output
Parameters:
kp
: Proportional gainki
: Integral gain (only active ifwith_I = true
)kd
: Derivative gain (only active ifwith_D = true
)Ni
: Anti-windup gain (only active ifwith_I = true
)Nd
: Maximum derivative gain (only active ifwith_D = true
). Typically set to 10-100 times the proportional gain.u_max
: Maximum output above which the output is saturatedu_min
: Minimum output below which the output is saturated. This defaults to-u_max
ifu_max > 0
and-Inf
otherwise.wp
:[0, 1]
Set-point weighting in the proportional part. Set to0
to prevent step changes in the output due to step changes in the reference.wd
:[0, 1]
Set-point weighting in the derivative part. Set to0
to prevent very large impulsive changes in the output due to step changes in the reference.with_I
: Whether or not to include the integral partwith_D
: Whether or not to include the derivative partImethod
: Discretization method for the integrator (see details above)Dmethod
: Discretization method for the derivative (see details above)
Extended help:
Internal variables:
I
: State of integratorD
: State of filtered derivativer
: Reference signal internal variabley
: Measurement signal internal variablewde
: Setpoint-weighted error for derivativev
: Un-saturated output of the controlleru
: Saturated output of the controllereI
: Error signal input to integrator including anit-windup tracking signale
: Error signal
ModelingToolkitSampledData.DiscretePIDStandard
— ConstantDiscretePIDStandard(;name, K = 1, Ti = 1, Td = 1, Ni = √(max(kd * ki, 1e-6)), Nd = 10, u_max = Inf, u_min = -u_max, wp = 1, wd = 1, Ts = 1, with_I = true, with_D = true, Imethod = :forward, Dmethod = :backward)
Discrete-time PID controller on standard (ideal) form with anti-windup and set-point weighting.
The controller is implemented on standard form
Simplified:
\[u = K \left( e + \int \frac{1}{T_i} e dt + T_d \dfrac{de}{dt} \right)\]
Detailed:
\[u = K \left( (w_p r - y) + \int \big( \frac{1}{T_i} (r - y) + N_i e_s \big ) dt + T_d \dfrac{d}{dt}(w_d r - y) \right)\]
where e_s = u - v
is the saturated error signal, v
is the unsaturated control signal and u
is the saturated control signal.
The derivative is filtered to allow a maximum gain of $N_d$.
The integrator is discretized using the method specified by Imethod
, options include
Imethod = :forward
(default): Corresponding to the transfer function $T_s / (z - 1)$Imethod = :backward
: Corresponding to the transfer function $T_s z / (z - 1)$Imethod = :trapezoidal
: Corresponding to the transfer function $(T_s / 2) (z + 1) / (z - 1)$
The derivative is discretized using the method specified by Dmethod
, options include
Dmethod = :forward
: Corresponding to the transfer function $\dfrac{N (z-1)}{z - \dfrac{T_d-N T_s}{T_d}}$.Dmethod = :backward
(default): Corresponding to the transfer function $\dfrac{\dfrac{NT_d}{T_d + N T_s}(z-1)}{z - \dfrac{T_d}{T_d + N T_s}}$
Anti windup is realized by tracking using the gain $N_i$ on the error signal $e_s$ when the output is saturated.
To use the controller in 1DOF mode, i.e., with only the control error as input, connect the error signal to the reference
connector, connect a Constant(; k = 0)
to the measurement
connector and set wp = wd = 1
.
Connectors:
reference
: The reference signal to the controller (or the error signal if used in 1DOF mode)measurement
: The measurement feedbackctr_output
: The control signal output
Parameters:
K
: Proportional gainTi
: Integral time constant (only active ifwith_I = true
)Td
: Derivative time (only active ifwith_D = true
)Ni
: Anti-windup gain (only active ifwith_I = true
)Nd
: Maximum derivative gain (only active ifwith_D = true
). Typically set to 10-100.u_max
: Maximum output above which the output is saturatedu_min
: Minimum output below which the output is saturated. This defaults to-u_max
ifu_max > 0
and-Inf
otherwise.wp
:[0, 1]
Set-point weighting in the proportional part. Set to0
to prevent step changes in the output due to step changes in the reference.wd
:[0, 1]
Set-point weighting in the derivative part. Set to0
to prevent very large impulsive changes in the output due to step changes in the reference.with_I
: Whether or not to include the integral partwith_D
: Whether or not to include the derivative partImethod
: Discretization method for the integrator (see details above)Dmethod
: Discretization method for the derivative (see details above)
Extended help:
Internal variables:
I
: State of integratorD
: State of filtered derivativer
: Reference signal internal variabley
: Measurement signal internal variablewde
: Setpoint-weighted error for derivativev
: Un-saturated output of the controlleru
: Saturated output of the controllereI
: Error signal input to integrator including anit-windup tracking signale
: Error signal
ModelingToolkitSampledData.DiscreteSlewRateLimiter
— ConstantDiscreteSlewRateLimiter(rate = 1.0, rate_negative = rate)
A discrete-time slew rate limiter that limits the rate of change of the input signal.
Note, the sample interval is not taken into account when computing the rate of change, the difference between two consequetive samples is saturated.
Parameters:
rate
: Slew rate limit (in positive/increasing direction). Must be a positive number.rate_negative
: Negative slew rate limit, defaults torate
. Must be a positive number.
Variables
d
: Unsaturated rate of change of the input signalu
: Input signaly
: Output signal (saturated slew rate)
Connectors:
input
output
Example
cl = Clock(0.1)
z = ShiftIndex(cl)
@mtkmodel SlweRateLimiterModel begin
@components begin
input = Sine(amplitude=1, frequency=0.8)
limiter = DiscreteSlewRateLimiter(; z, rate=0.4, rate_negative = 0.3)
end
@variables begin
x(t) = 0 # Dummy variable to workaround JSCompiler bug
end
@equations begin
connect(input.output, limiter.input)
D(x) ~ 0.1x + Hold(limiter.y)
end
end
@named m = SlweRateLimiterModel()
m = complete(m)
ssys = structural_simplify(IRSystem(m))
prob = ODEProblem(ssys, [m.limiter.y(z-1) => 0], (0.0, 2.0))
sol = solve(prob, Tsit5(), dtmax=0.01)
plot(sol, idxs=[m.input.output.u, m.limiter.y], title="Slew rate limited sine wave")
ModelingToolkitSampledData.ExponentialFilter
— ConstantExponentialFilter(a = 0.1)
Exponential filtering (first-order filter) with input-output relation $y(z) = (1 - a) y(z-1) + a u(z-1)$, transfer function
\[Y(z) = \dfrac{a}{1 - (1 - a) z^{-1}} U(z)\]
Parameters:
a
: Filter parameter[0, 1]
, a small value implies stronger filtering.
Variables:
u
: Input signaly
: Output signal
Connectors:
input::RealInput
: Input signaloutput::RealOutput
: Output signal
ModelingToolkitSampledData.MovingAverageFilter
— ConstantMovingAverageFilter(N = 3)
Exponential filtering with input-output relation $y(z) = \dfrac{1}{N} \sum_{i=0}^{N-1} u(z-i)$.
Please note: this implementation of a moving average filter is not optimized for very large number of filter taps N
.
Parameters:
N
: (structural) Number of samples to average over
Variables:
u
: Input signaly
: Output signal
Connectors:
input::RealInput
: Input signaloutput::RealOutput
: Output signal
ModelingToolkitSampledData.NormalNoise
— ConstantNormalNoise()
A discrete-time noise source that returns a normally-distributed value at each clock tick.
Parameters
mu = 0
: Meansigma = 1
: Standard deviationseed = 1
: Seed for the random number generator
Structural parameters
z
: TheShiftIndex
used to indicate clock partition.
Connectors:
output
ModelingToolkitSampledData.Quantization
— ConstantQuantization
A quantization block that quantizes the input signal to a specified number of bits.
Parameters:
y_max
: Upper limit of outputy_min
: Lower limit of outputbits
: Number of bits of quantizationquantized
: If quantization effects shall be computed. If false, the output is equal to the input, which may be useful for, e.g., linearization.midrise
: (structural) If true (default), the quantizer is a midrise quantizer, otherwise it is a midtread quantizer. See Docs: Quantization for more details.
Connectors:
input
output
Variables
y
: Output signal, equal tooutput.u
u
: Input signal, equal toinput.u
Automatic differentiation
This block is not differentiable, the derivative is zero everywhere exect for at the level transition where it is ill defined. To use in a differentiable context, set quantized = false
which turns this block into the identity function.
ModelingToolkitSampledData.SampleWithADEffects
— ConstantSampleWithADEffects(quantized = false, noisy = false)
A sampler with additional effects that appear in practical systems, such as measurement noise and quantization.
The operations occur in the order
- Sampling
- Noise addition
- Quantization
Structural parameters:
quantized
: If true, the output is quantized. When this option is used, the output is quantized to the number of bits specified by thebits
parameter. The quantization is midrise ifmidrise = true
, otherwise it is midtread. The output is also limited to the range[y_min, y_max]
.noisy
: If true, the output is corrupted by additive white Gaussian noise with standard deviationsigma
(defaults to 0.1). Ifnoisy = false
, the noise block is a unit gain.dt
: Sample interval of the sampler. If not specified, the sample interval is inferred from the clock of the system.clock
: Clock signal of the system. If not specified, the sample interval is inferred from the clock of the system. Ifclock
is specified, the parameterdt
has no effect.
Parameters:
y_min
: Lower limit of output, defaults to -1. Only used ifquantized = true
.y_max
: Upper limit of output, defaults to 1. Only used ifquantized = true
.bits
: Number of bits of quantization, defaults to 8 (256 output levels betweeny_min
andy_max
). Only used ifquantized = true
.sigma
: Standard deviation of the additive Gaussian noise, defaults to 0.1. Only used ifnoisy = true
.
ModelingToolkitSampledData.Sampler
— ConstantSampler()
Sampler(; dt::Real)
Sampler(; clock::AbstractClock)
Sampler
transforms a continuous-time signal into discrete time by sampling the input signal every time the associated clock ticks. The clock can be specified explicitly using the clock
keyword argument, or implicitly by providing a sample time dt
, in which case a standard periodic Clock
is used.
Connectors:
input
(continuous-time signal)output
(discrete-time signal)
ModelingToolkitSampledData.UniformNoise
— ConstantUniformNoise()
A discrete-time noise source that returns a uniformly distributed value at each clock tick.
Parameters
l = 0
: Lower boundu = 1
: Upper boundseed = 1
: Seed for the random number generator
Structural parameters
rng
: A random number generator, defaults toRandom.Xoshiro()
.z
: TheShiftIndex
used to indicate clock partition.
Connectors:
output
ModelingToolkitSampledData.ZeroOrderHold
— ConstantZeroOrderHold()
Zero-order-Hold translates a discrete time (clocked) signal into continuous time by holding the value of the discrete signal constant until the next sample.
Connectors:
input
(discrete-time signal)output
(continuous-time signal)
ModelingToolkitSampledData.seeded_rand
— Methodseeded_rand(seed, t)
Internal function. This function seeds the seed parameter as well as the current simulation time.
ModelingToolkitSampledData.seeded_randn
— Methodseeded_randn(seed, t)
Internal function. This function seeds the seed parameter as well as the current simulation time.