#
Time Varying Behavior
So far, all the components we've created have the same behavior at all points in time. But what if we want to create a component whose behavior changes over time? Let's walk through a few examples of models that have time varying behavior.
#
Step
The first case we'll consider is a step change in the input voltage to our circuit. The follow model acts as a step voltage source:
@mtkmodel StepVoltage begin
@extend (v,) = twopin = TwoPin()
@parameters begin
V0, [unit = u"V"]
Vf, [unit = u"V"]
t_step, [unit = u"s"]
end
@equations begin
v ~ ifelse(t > t_step, Vf, V0)
end
end
export StepVoltage
In this model, the V0 parameter is the initial voltage provided by this
source, t_step is the time when the voltage makes a step change and Vf is
the voltage after the step occurs.
We can build a system model
using ModelingToolkit
using ExampleLibrary
using DifferentialEquations
using Plots
@mtkmodel RLCModel begin
@components begin
resistor = Resistor(R = 100.0)
capacitor = Capacitor(C = 1.0e-3)
inductor = Inductor(L = 1.0)
source = StepVoltage(V0 = 0.0, Vf = 24.0, t_step = 0.5)
ground = Ground()
end
@equations begin
connect(source.p, inductor.p)
connect(inductor.n, resistor.p)
connect(inductor.n, capacitor.p)
connect(resistor.n, ground.g)
connect(capacitor.n, ground.g)
connect(source.n, ground.g)
end
end
@mtkbuild rlc_model = RLCModel()
u0 = [
rlc_model.capacitor.v => 0.0,
rlc_model.inductor.i => 0.0
]
prob = ODEProblem(rlc_model, u0, (0, 1.0))
sol = solve(prob)
display(plot(sol, idxs=[rlc_model.capacitor.v, rlc_model.inductor.i]))
Simulating this system leads to the following response:
#
Sine Wave
The following model implements a sine wave voltage source:
@mtkmodel SineVoltage begin
@extend (v,) = twopin = TwoPin()
@parameters begin
t_start, [unit = u"s"]
Voffset, [unit = u"V"]
amplitude, [unit = u"V"]
frequency, [unit = u"Hz"]
end
@equations begin
v ~ ifelse(t > t_start, amplitude*sin((t-t_start)*2*pi*frequency)+Voffset, Voffset)
end
end
export SineVoltage
In this case, t_start is the time when the sine wave starts, Voffset is the offset voltage of the sine wave, amplitude is the amplitude of the sine wave and frequency is the frequency of the sine wave.
Creating a model with this component is nearly identical to the previous case:
using ModelingToolkit
using ExampleLibrary
using DifferentialEquations
using Plots
@mtkmodel RLCModel begin
@components begin
resistor = Resistor(R = 100.0)
capacitor = Capacitor(C = 1.0e-3)
inductor = Inductor(L = 1.0)
source = SineVoltage(amplitude=24, Voffset=0, frequency=5, t_start = 0.5)
ground = Ground()
end
@equations begin
connect(source.p, inductor.p)
connect(inductor.n, resistor.p)
connect(inductor.n, capacitor.p)
connect(resistor.n, ground.g)
connect(capacitor.n, ground.g)
connect(source.n, ground.g)
end
end
@mtkbuild rlc_model = RLCModel()
u0 = [
rlc_model.capacitor.v => 0.0,
rlc_model.inductor.i => 0.0
]
prob = ODEProblem(rlc_model, u0, (0, 2))
sol = solve(prob)
display(plot(sol, idxs=[rlc_model.capacitor.v, rlc_model.source.v, rlc_model.inductor.i]))
And the results from this simulation look just as we would expect: