Skip to content

Commit 3a62cdf

Browse files
committed
add generic transfer function
1 parent 1c7eb69 commit 3a62cdf

File tree

3 files changed

+142
-1
lines changed

3 files changed

+142
-1
lines changed

src/Blocks/Blocks.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -25,7 +25,7 @@ include("sources.jl")
2525
export Limiter, DeadZone, SlewRateLimiter
2626
include("nonlinear.jl")
2727

28-
export Integrator, Derivative, FirstOrder, SecondOrder, StateSpace
28+
export Integrator, Derivative, FirstOrder, SecondOrder, StateSpace, TransferFunction
2929
export PI, LimPI, PID, LimPID
3030
include("continuous.jl")
3131

src/Blocks/continuous.jl

Lines changed: 83 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -549,3 +549,86 @@ linearized around the operating point `x₀, u₀`, we have `y0, u0 = h(x₀, u
549549
end
550550

551551
StateSpace(A, B, C, D = nothing; kwargs...) = StateSpace(; A, B, C, D, kwargs...)
552+
553+
symbolic_eps(t) = eps(t)
554+
@register_symbolic symbolic_eps(t)
555+
556+
557+
"""
558+
TransferFunction(; b, a, name)
559+
560+
A linear time-invariant system provided as a transfer-function.
561+
```
562+
Y(s) = b(s) / a(s) U(s)
563+
```
564+
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.
565+
566+
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`.
567+
568+
To set the initial state, it's recommended to set the initial condition for `x`, and let that of `x_scaled` be computed automatically.
569+
570+
# Parameters:
571+
- `b`: Numerator polynomial coefficients, e.g., `2s + 3` is specified as `[2, 3]`
572+
- `a`: Denomenator polynomial coefficients, e.g., `s + 2ws + w^2` is specified as `[1, 2w, w^2]`
573+
574+
# Connectors:
575+
- `input`
576+
- `output`
577+
578+
See also [`StateSpace`](@ref) as well as [ControlSystemsMTK.jl](https://juliacontrol.github.io/ControlSystemsMTK.jl/stable/) for an interface between [ControlSystems.jl](https://juliacontrol.github.io/ControlSystems.jl/stable/) and ModelingToolkit.jl for advanced manipulation of transfer functions and linear statespace systems. For linearization, see [`linearize`](@ref) and [Linear Analysis](https://docs.sciml.ai/ModelingToolkitStandardLibrary/stable/API/linear_analysis/).
579+
"""
580+
@component function TransferFunction(; b = [1], a = [1, 1], name)
581+
nb = length(b)
582+
na = length(a)
583+
nb <= na ||
584+
error("Transfer function is not proper, the numerator must not be longer than the denominator")
585+
nx = na - 1
586+
nbb = max(0, na - nb)
587+
588+
@named begin
589+
input = RealInput()
590+
output = RealOutput()
591+
end
592+
593+
@parameters begin
594+
b[1:nb] = b,
595+
[
596+
description = "Numerator coefficients of transfer function (e.g., 2s + 3 is specified as [2,3])",
597+
]
598+
a[1:na] = a,
599+
[
600+
description = "Denominator coefficients of transfer function (e.g., s + 2ws + w^2 is specified as [1, 2w, w^2])",
601+
]
602+
bb[1:(nbb + nb)] = [zeros(nbb); b]
603+
d = bb[1] / a[1]
604+
end
605+
606+
a = collect(a)
607+
@parameters a_end = ifelse(a[end] > 100 * symbolic_eps(sqrt(a' * a)), a[end], 1.0)
608+
609+
pars = [collect(b); a; collect(bb); d; a_end]
610+
@variables begin
611+
x(t)[1:nx] = zeros(nx),
612+
[description = "State of transfer function on controller canonical form"]
613+
x_scaled(t)[1:nx] = collect(x) * a_end, [description = "Scaled vector x"]
614+
u(t), [description = "Input of transfer function"]
615+
y(t), [description = "Output of transfer function"]
616+
end
617+
618+
x = collect(x)
619+
x_scaled = collect(x_scaled)
620+
621+
sts = [x; x_scaled; y; u]
622+
623+
if nx == 0
624+
eqs = [y ~ d * u]
625+
else
626+
eqs = [D(x_scaled[1]) ~ (-a[2:na]'x_scaled + a_end * u) / a[1]
627+
D.(x_scaled[2:nx]) .~ x_scaled[1:(nx - 1)]
628+
y ~ ((bb[2:na] - d * a[2:na])'x_scaled) / a_end + d * u
629+
x .~ x_scaled ./ a_end]
630+
end
631+
push!(eqs, input.u ~ u)
632+
push!(eqs, output.u ~ y)
633+
compose(ODESystem(eqs, t, sts, pars; name = name), input, output)
634+
end

test/Blocks/continuous.jl

Lines changed: 58 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -408,4 +408,62 @@ end
408408
@test sol[plant.output.u][end]re_val atol=1e-3 # zero control error after 100s
409409
@test all(-1.5 .<= sol[pid_controller.ctr_output.u] .<= 1.5) # test limit
410410
end
411+
412+
@testset "TransferFunction" begin
413+
pt1_func(t, k, T) = k * (1 - exp(-t / T)) # Known solution to first-order system
414+
415+
@named c = Constant(; k = 1)
416+
@named pt1 = TransferFunction(b = [1.2], a = [3.14, 1])
417+
@named iosys = ODESystem(connect(c.output, pt1.input), t, systems = [pt1, c])
418+
sys = structural_simplify(iosys)
419+
prob = ODEProblem(sys, Pair[pt1.a_end => 1], (0.0, 100.0))
420+
sol = solve(prob, Rodas4())
421+
@test sol.retcode == Success
422+
@test sol[pt1.output.u]pt1_func.(sol.t, 1.2, 3.14) atol=1e-3
423+
424+
# Test logic for a_end by constructing an integrator
425+
@named c = Constant(; k = 1)
426+
@named pt1 = TransferFunction(b = [1.2], a = [3.14, 0])
427+
@named iosys = ODESystem(connect(c.output, pt1.input), t, systems = [pt1, c])
428+
sys = structural_simplify(iosys)
429+
prob = ODEProblem(sys, Pair[pt1.a_end => 1], (0.0, 100.0))
430+
sol = solve(prob, Rodas4())
431+
@test sol.retcode == Success
432+
@test sol[pt1.output.u] sol.t .* (1.2 / 3.14)
433+
@test sol[pt1.x[1]] sol.t .* (1 / 3.14) # Test that scaling of state works properly
434+
435+
# Test higher order
436+
437+
function pt2_func(t, k, w, d)
438+
y = if d == 0
439+
-k * (-1 + cos(t * w))
440+
else
441+
d = complex(d)
442+
real(k * (1 +
443+
(-cosh(sqrt(-1 + d^2) * t * w) -
444+
(d * sinh(sqrt(-1 + d^2) * t * w)) / sqrt(-1 + d^2)) /
445+
exp(d * t * w)))
446+
end
447+
end
448+
449+
k, w, d = 1.0, 1.0, 0.5
450+
@named pt1 = TransferFunction(b = [w^2], a = [1, 2d * w, w^2])
451+
@named iosys = ODESystem(connect(c.output, pt1.input), t, systems = [pt1, c])
452+
sys = structural_simplify(iosys)
453+
prob = ODEProblem(sys, Pair[pt1.a_end => 1], (0.0, 100.0))
454+
sol = solve(prob, Rodas4())
455+
@test sol.retcode == Success
456+
@test sol[pt1.output.u]pt2_func.(sol.t, k, w, d) atol=1e-3
457+
458+
# test zeros (high-pass version of first test)
459+
@named c = Constant(; k = 1)
460+
@named pt1 = TransferFunction(b = [1, 0], a = [1, 1])
461+
@named iosys = ODESystem(connect(c.output, pt1.input), t, systems = [pt1, c])
462+
sys = structural_simplify(iosys)
463+
prob = ODEProblem(sys, Pair[pt1.a_end => 1], (0.0, 100.0))
464+
sol = solve(prob, Rodas4())
465+
@test sol.retcode == Success
466+
@test sol[pt1.output.u]1 .- pt1_func.(sol.t, 1, 1) atol=1e-3
467+
@test sol[pt1.x[1]]pt1_func.(sol.t, 1, 1) atol=1e-3 # Test that scaling of state works properly
468+
end
411469
end

0 commit comments

Comments
 (0)