-
-
Notifications
You must be signed in to change notification settings - Fork 215
/
Copy pathlinearization.jl
680 lines (588 loc) · 23.6 KB
/
linearization.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
"""
lin_fun, simplified_sys = linearization_function(sys::AbstractSystem, inputs, outputs; simplify = false, initialize = true, initialization_solver_alg = TrustRegion(), kwargs...)
Return a function that linearizes the system `sys`. The function [`linearize`](@ref) provides a higher-level and easier to use interface.
`lin_fun` is a function `(variables, p, t) -> (; f_x, f_z, g_x, g_z, f_u, g_u, h_x, h_z, h_u)`, i.e., it returns a NamedTuple with the Jacobians of `f,g,h` for the nonlinear `sys` (technically for `simplified_sys`) on the form
```math
\\begin{aligned}
ẋ &= f(x, z, u) \\\\
0 &= g(x, z, u) \\\\
y &= h(x, z, u)
\\end{aligned}
```
where `x` are differential unknown variables, `z` algebraic variables, `u` inputs and `y` outputs. To obtain a linear statespace representation, see [`linearize`](@ref). The input argument `variables` is a vector defining the operating point, corresponding to `unknowns(simplified_sys)` and `p` is a vector corresponding to the parameters of `simplified_sys`. Note: all variables in `inputs` have been converted to parameters in `simplified_sys`.
The `simplified_sys` has undergone [`structural_simplify`](@ref) and had any occurring input or output variables replaced with the variables provided in arguments `inputs` and `outputs`. The unknowns of this system also indicate the order of the unknowns that holds for the linearized matrices.
# Arguments:
- `sys`: An [`ODESystem`](@ref). This function will automatically apply simplification passes on `sys` and return the resulting `simplified_sys`.
- `inputs`: A vector of variables that indicate the inputs of the linearized input-output model.
- `outputs`: A vector of variables that indicate the outputs of the linearized input-output model.
- `simplify`: Apply simplification in tearing.
- `initialize`: If true, a check is performed to ensure that the operating point is consistent (satisfies algebraic equations). If the op is not consistent, initialization is performed.
- `initialization_solver_alg`: A NonlinearSolve algorithm to use for solving for a feasible set of state and algebraic variables that satisfies the specified operating point.
- `kwargs`: Are passed on to `find_solvables!`
See also [`linearize`](@ref) which provides a higher-level interface.
"""
function linearization_function(sys::AbstractSystem, inputs,
outputs; simplify = false,
initialize = true,
initializealg = nothing,
initialization_abstol = 1e-5,
initialization_reltol = 1e-3,
op = Dict(),
p = DiffEqBase.NullParameters(),
zero_dummy_der = false,
initialization_solver_alg = TrustRegion(),
eval_expression = false, eval_module = @__MODULE__,
warn_initialize_determined = true,
guesses = Dict(),
kwargs...)
op = Dict(op)
inputs isa AbstractVector || (inputs = [inputs])
outputs isa AbstractVector || (outputs = [outputs])
inputs = mapreduce(vcat, inputs; init = []) do var
symbolic_type(var) == ArraySymbolic() ? collect(var) : [var]
end
outputs = mapreduce(vcat, outputs; init = []) do var
symbolic_type(var) == ArraySymbolic() ? collect(var) : [var]
end
ssys, diff_idxs, alge_idxs, input_idxs = io_preprocessing(sys, inputs, outputs;
simplify,
kwargs...)
if zero_dummy_der
dummyder = setdiff(unknowns(ssys), unknowns(sys))
defs = Dict(x => 0.0 for x in dummyder)
@set! ssys.defaults = merge(defs, defaults(ssys))
op = merge(defs, op)
end
sys = ssys
if initializealg === nothing
initializealg = initialize ? OverrideInit() : NoInit()
end
prob = ODEProblem{true, SciMLBase.FullSpecialize}(
sys, op, (nothing, nothing), p; allow_incomplete = true,
algebraic_only = true, guesses)
u0 = state_values(prob)
ps = parameters(sys)
h = build_explicit_observed_function(sys, outputs; eval_expression, eval_module)
initialization_kwargs = (;
abstol = initialization_abstol, reltol = initialization_reltol,
nlsolve_alg = initialization_solver_alg)
lin_fun = LinearizationFunction(
diff_idxs, alge_idxs, input_idxs, length(unknowns(sys)),
prob, h, u0 === nothing ? nothing : similar(u0),
ForwardDiff.Chunk(input_idxs), initializealg, initialization_kwargs)
return lin_fun, sys
end
"""
$(TYPEDEF)
A callable struct which linearizes a system.
# Fields
$(TYPEDFIELDS)
"""
struct LinearizationFunction{
DI <: AbstractVector{Int}, AI <: AbstractVector{Int}, II, P <: ODEProblem,
H, C, Ch, IA <: SciMLBase.DAEInitializationAlgorithm, IK}
"""
The indexes of differential equations in the linearized system.
"""
diff_idxs::DI
"""
The indexes of algebraic equations in the linearized system.
"""
alge_idxs::AI
"""
The indexes of parameters in the linearized system which represent
input variables.
"""
input_idxs::II
"""
The number of unknowns in the linearized system.
"""
num_states::Int
"""
The `ODEProblem` of the linearized system.
"""
prob::P
"""
A function which takes `(u, p, t)` and returns the outputs of the linearized system.
"""
h::H
"""
Any required cache buffers.
"""
caches::C
# TODO: Use DI?
"""
A `ForwardDiff.Chunk` for taking the jacobian with respect to the inputs.
"""
chunk::Ch
"""
The initialization algorithm to use.
"""
initializealg::IA
"""
Keyword arguments to be passed to `SciMLBase.get_initial_values`.
"""
initialize_kwargs::IK
end
SymbolicIndexingInterface.symbolic_container(f::LinearizationFunction) = f.prob
SymbolicIndexingInterface.state_values(f::LinearizationFunction) = state_values(f.prob)
function SymbolicIndexingInterface.parameter_values(f::LinearizationFunction)
parameter_values(f.prob)
end
SymbolicIndexingInterface.current_time(f::LinearizationFunction) = current_time(f.prob)
function Base.show(io::IO, mime::MIME"text/plain", lf::LinearizationFunction)
printstyled(io, "LinearizationFunction"; bold = true, color = :blue)
println(io, " which wraps:")
show(io, mime, lf.prob)
end
"""
$(TYPEDSIGNATURES)
Linearize the wrapped system at the point given by `(unknowns, p, t)`.
"""
function (linfun::LinearizationFunction)(u, p, t)
if eltype(p) <: Pair
p = todict(p)
newps = copy(parameter_values(linfun.prob))
for (k, v) in p
if is_parameter(linfun, k)
v = fixpoint_sub(v, p)
setp(linfun, k)(newps, v)
end
end
p = newps
end
fun = linfun.prob.f
if u !== nothing # Handle systems without unknowns
linfun.num_states == length(u) ||
error("Number of unknown variables ($(linfun.num_states)) does not match the number of input unknowns ($(length(u)))")
integ_cache = (linfun.caches,)
integ = MockIntegrator{true}(u, p, t, integ_cache, nothing)
u, p, success = SciMLBase.get_initial_values(
linfun.prob, integ, fun, linfun.initializealg, Val(true);
linfun.initialize_kwargs...)
if !success
error("Initialization algorithm $(linfun.initializealg) failed with `unknowns = $u` and `p = $p`.")
end
uf = SciMLBase.UJacobianWrapper(fun, t, p)
fg_xz = ForwardDiff.jacobian(uf, u)
h_xz = ForwardDiff.jacobian(
let p = p, t = t, h = linfun.h
xz -> h(xz, p, t)
end, u)
pf = SciMLBase.ParamJacobianWrapper(fun, t, u)
fg_u = jacobian_wrt_vars(pf, p, linfun.input_idxs, linfun.chunk)
else
linfun.num_states == 0 ||
error("Number of unknown variables (0) does not match the number of input unknowns ($(length(u)))")
fg_xz = zeros(0, 0)
h_xz = fg_u = zeros(0, length(linfun.input_idxs))
end
hp = let u = u, t = t, h = linfun.h
_hp(p) = h(u, p, t)
_hp
end
h_u = jacobian_wrt_vars(hp, p, linfun.input_idxs, linfun.chunk)
(f_x = fg_xz[linfun.diff_idxs, linfun.diff_idxs],
f_z = fg_xz[linfun.diff_idxs, linfun.alge_idxs],
g_x = fg_xz[linfun.alge_idxs, linfun.diff_idxs],
g_z = fg_xz[linfun.alge_idxs, linfun.alge_idxs],
f_u = fg_u[linfun.diff_idxs, :],
g_u = fg_u[linfun.alge_idxs, :],
h_x = h_xz[:, linfun.diff_idxs],
h_z = h_xz[:, linfun.alge_idxs],
h_u = h_u,
x = u,
p,
t)
end
"""
$(TYPEDEF)
Mock `DEIntegrator` to allow using `CheckInit` without having to create a new integrator
(and consequently depend on `OrdinaryDiffEq`).
# Fields
$(TYPEDFIELDS)
"""
struct MockIntegrator{iip, U, P, T, C, O} <: SciMLBase.DEIntegrator{Nothing, iip, U, T}
"""
The state vector.
"""
u::U
"""
The parameter object.
"""
p::P
"""
The current time.
"""
t::T
"""
The integrator cache.
"""
cache::C
"""
Integrator "options" for `CheckInit`.
"""
opts::O
end
function MockIntegrator{iip}(u::U, p::P, t::T, cache::C, opts::O) where {iip, U, P, T, C, O}
return MockIntegrator{iip, U, P, T, C, O}(u, p, t, cache, opts)
end
SymbolicIndexingInterface.state_values(integ::MockIntegrator) = integ.u
SymbolicIndexingInterface.parameter_values(integ::MockIntegrator) = integ.p
SymbolicIndexingInterface.current_time(integ::MockIntegrator) = integ.t
SciMLBase.get_tmp_cache(integ::MockIntegrator) = integ.cache
"""
$(TYPEDEF)
A struct representing a linearization operation to be performed. Can be symbolically
indexed to efficiently update the operating point for multiple linearizations in a loop.
The value of the independent variable can be set by mutating the `.t` field of this struct.
"""
mutable struct LinearizationProblem{F <: LinearizationFunction, T}
"""
The wrapped `LinearizationFunction`
"""
const f::F
t::T
end
function Base.show(io::IO, mime::MIME"text/plain", prob::LinearizationProblem)
printstyled(io, "LinearizationProblem"; bold = true, color = :blue)
println(io, " at time ", prob.t, " which wraps:")
show(io, mime, prob.f.prob)
end
"""
$(TYPEDSIGNATURES)
Construct a `LinearizationProblem` for linearizing the system `sys` with the given
`inputs` and `outputs`.
# Keyword arguments
- `t`: The value of the independent variable
All other keyword arguments are forwarded to `linearization_function`.
"""
function LinearizationProblem(sys::AbstractSystem, inputs, outputs; t = 0.0, kwargs...)
linfun, _ = linearization_function(sys, inputs, outputs; kwargs...)
return LinearizationProblem(linfun, t)
end
SymbolicIndexingInterface.symbolic_container(p::LinearizationProblem) = p.f
SymbolicIndexingInterface.state_values(p::LinearizationProblem) = state_values(p.f)
SymbolicIndexingInterface.parameter_values(p::LinearizationProblem) = parameter_values(p.f)
SymbolicIndexingInterface.current_time(p::LinearizationProblem) = p.t
function Base.getindex(prob::LinearizationProblem, idx)
getu(prob, idx)(prob)
end
function Base.setindex!(prob::LinearizationProblem, val, idx)
setu(prob, idx)(prob, val)
end
function Base.getproperty(prob::LinearizationProblem, x::Symbol)
if x == :ps
return ParameterIndexingProxy(prob)
end
return getfield(prob, x)
end
function CommonSolve.solve(prob::LinearizationProblem; allow_input_derivatives = false)
u0 = state_values(prob)
p = parameter_values(prob)
t = current_time(prob)
linres = prob.f(u0, p, t)
f_x, f_z, g_x, g_z, f_u, g_u, h_x, h_z, h_u, x, p, t = linres
nx, nu = size(f_u)
nz = size(f_z, 2)
ny = size(h_x, 1)
D = h_u
if isempty(g_z)
A = f_x
B = f_u
C = h_x
@assert iszero(g_x)
@assert iszero(g_z)
@assert iszero(g_u)
else
gz = lu(g_z; check = false)
issuccess(gz) ||
error("g_z not invertible, this indicates that the DAE is of index > 1.")
gzgx = -(gz \ g_x)
A = [f_x f_z
gzgx*f_x gzgx*f_z]
B = [f_u
gzgx * f_u] # The cited paper has zeros in the bottom block, see derivation in https://github.com/SciML/ModelingToolkit.jl/pull/1691 for the correct formula
C = [h_x h_z]
Bs = -(gz \ g_u) # This equation differ from the cited paper, the paper is likely wrong since their equaiton leads to a dimension mismatch.
if !iszero(Bs)
if !allow_input_derivatives
der_inds = findall(vec(any(!=(0), Bs, dims = 1)))
error("Input derivatives appeared in expressions (-g_z\\g_u != 0), the following inputs appeared differentiated: $(inputs(sys)[der_inds]). Call `linearize` with keyword argument `allow_input_derivatives = true` to allow this and have the returned `B` matrix be of double width ($(2nu)), where the last $nu inputs are the derivatives of the first $nu inputs.")
end
B = [B [zeros(nx, nu); Bs]]
D = [D zeros(ny, nu)]
end
end
(; A, B, C, D), (; x, p, t)
end
"""
(; A, B, C, D), simplified_sys = linearize_symbolic(sys::AbstractSystem, inputs, outputs; simplify = false, allow_input_derivatives = false, kwargs...)
Similar to [`linearize`](@ref), but returns symbolic matrices `A,B,C,D` rather than numeric. While `linearize` uses ForwardDiff to perform the linearization, this function uses `Symbolics.jacobian`.
See [`linearize`](@ref) for a description of the arguments.
# Extended help
The named tuple returned as the first argument additionally contains the jacobians `f_x, f_z, g_x, g_z, f_u, g_u, h_x, h_z, h_u` of
```math
\\begin{aligned}
ẋ &= f(x, z, u) \\\\
0 &= g(x, z, u) \\\\
y &= h(x, z, u)
\\end{aligned}
```
where `x` are differential unknown variables, `z` algebraic variables, `u` inputs and `y` outputs.
"""
function linearize_symbolic(sys::AbstractSystem, inputs,
outputs; simplify = false, allow_input_derivatives = false,
eval_expression = false, eval_module = @__MODULE__,
kwargs...)
sys, diff_idxs, alge_idxs, input_idxs = io_preprocessing(
sys, inputs, outputs; simplify,
kwargs...)
sts = unknowns(sys)
t = get_iv(sys)
ps = parameters(sys; initial_parameters = true)
p = reorder_parameters(sys, ps)
fun_expr = generate_function(sys, sts, ps; expression = Val{true})[1]
fun = eval_or_rgf(fun_expr; eval_expression, eval_module)
dx = fun(sts, p, t)
h = build_explicit_observed_function(sys, outputs; eval_expression, eval_module)
y = h(sts, p, t)
fg_xz = Symbolics.jacobian(dx, sts)
fg_u = Symbolics.jacobian(dx, inputs)
h_xz = Symbolics.jacobian(y, sts)
h_u = Symbolics.jacobian(y, inputs)
f_x = fg_xz[diff_idxs, diff_idxs]
f_z = fg_xz[diff_idxs, alge_idxs]
g_x = fg_xz[alge_idxs, diff_idxs]
g_z = fg_xz[alge_idxs, alge_idxs]
f_u = fg_u[diff_idxs, :]
g_u = fg_u[alge_idxs, :]
h_x = h_xz[:, diff_idxs]
h_z = h_xz[:, alge_idxs]
nx, nu = size(f_u)
nz = size(f_z, 2)
ny = size(h_x, 1)
D = h_u
if isempty(g_z) # ODE
A = f_x
B = f_u
C = h_x
else
gz = lu(g_z; check = false)
issuccess(gz) ||
error("g_z not invertible, this indicates that the DAE is of index > 1.")
gzgx = -(gz \ g_x)
A = [f_x f_z
gzgx*f_x gzgx*f_z]
B = [f_u
gzgx * f_u] # The cited paper has zeros in the bottom block, see derivation in https://github.com/SciML/ModelingToolkit.jl/pull/1691 for the correct formula
C = [h_x h_z]
Bs = -(gz \ g_u) # This equation differ from the cited paper, the paper is likely wrong since their equaiton leads to a dimension mismatch.
if !iszero(Bs)
if !allow_input_derivatives
der_inds = findall(vec(any(!iszero, Bs, dims = 1)))
@show typeof(der_inds)
error("Input derivatives appeared in expressions (-g_z\\g_u != 0), the following inputs appeared differentiated: $(ModelingToolkit.inputs(sys)[der_inds]). Call `linearize_symbolic` with keyword argument `allow_input_derivatives = true` to allow this and have the returned `B` matrix be of double width ($(2nu)), where the last $nu inputs are the derivatives of the first $nu inputs.")
end
B = [B [zeros(nx, nu); Bs]]
D = [D zeros(ny, nu)]
end
end
(; A, B, C, D, f_x, f_z, g_x, g_z, f_u, g_u, h_x, h_z, h_u), sys
end
function markio!(state, orig_inputs, inputs, outputs; check = true)
fullvars = get_fullvars(state)
inputset = Dict{Any, Bool}(i => false for i in inputs)
outputset = Dict{Any, Bool}(o => false for o in outputs)
for (i, v) in enumerate(fullvars)
if v in keys(inputset)
if v in keys(outputset)
v = setio(v, true, true)
outputset[v] = true
else
v = setio(v, true, false)
end
inputset[v] = true
fullvars[i] = v
elseif v in keys(outputset)
v = setio(v, false, true)
outputset[v] = true
fullvars[i] = v
else
if isinput(v)
push!(orig_inputs, v)
end
v = setio(v, false, false)
fullvars[i] = v
end
end
if check
ikeys = keys(filter(!last, inputset))
if !isempty(ikeys)
error(
"Some specified inputs were not found in system. The following variables were not found ",
ikeys)
end
end
check && (all(values(outputset)) ||
error(
"Some specified outputs were not found in system. The following Dict indicates the found variables ",
outputset))
state, orig_inputs
end
"""
(; A, B, C, D), simplified_sys, extras = linearize(sys, inputs, outputs; t=0.0, op = Dict(), allow_input_derivatives = false, zero_dummy_der=false, kwargs...)
(; A, B, C, D), extras = linearize(simplified_sys, lin_fun; t=0.0, op = Dict(), allow_input_derivatives = false, zero_dummy_der=false)
Linearize `sys` between `inputs` and `outputs`, both vectors of variables. Return a NamedTuple with the matrices of a linear statespace representation
on the form
```math
\\begin{aligned}
ẋ &= Ax + Bu\\\\
y &= Cx + Du
\\end{aligned}
```
The first signature automatically calls [`linearization_function`](@ref) internally,
while the second signature expects the outputs of [`linearization_function`](@ref) as input.
`op` denotes the operating point around which to linearize. If none is provided,
the default values of `sys` are used.
If `allow_input_derivatives = false`, an error will be thrown if input derivatives (``u̇``) appear as inputs in the linearized equations. If input derivatives are allowed, the returned `B` matrix will be of double width, corresponding to the input `[u; u̇]`.
`zero_dummy_der` can be set to automatically set the operating point to zero for all dummy derivatives.
The return value `extras` is a NamedTuple `(; x, p, t)` containing the result of the initialization problem that was solved to determine the operating point.
See also [`linearization_function`](@ref) which provides a lower-level interface, [`linearize_symbolic`](@ref) and [`ModelingToolkit.reorder_unknowns`](@ref).
See extended help for an example.
The implementation and notation follows that of
["Linear Analysis Approach for Modelica Models", Allain et al. 2009](https://ep.liu.se/ecp/043/075/ecp09430097.pdf)
# Extended help
This example builds the following feedback interconnection and linearizes it from the input of `F` to the output of `P`.
```
r ┌─────┐ ┌─────┐ ┌─────┐
───►│ ├──────►│ │ u │ │
│ F │ │ C ├────►│ P │ y
└─────┘ ┌►│ │ │ ├─┬─►
│ └─────┘ └─────┘ │
│ │
└─────────────────────┘
```
```julia
using ModelingToolkit
using ModelingToolkit: t_nounits as t, D_nounits as D
function plant(; name)
@variables x(t) = 1
@variables u(t)=0 y(t)=0
eqs = [D(x) ~ -x + u
y ~ x]
ODESystem(eqs, t; name = name)
end
function ref_filt(; name)
@variables x(t)=0 y(t)=0
@variables u(t)=0 [input = true]
eqs = [D(x) ~ -2 * x + u
y ~ x]
ODESystem(eqs, t, name = name)
end
function controller(kp; name)
@variables y(t)=0 r(t)=0 u(t)=0
@parameters kp = kp
eqs = [
u ~ kp * (r - y),
]
ODESystem(eqs, t; name = name)
end
@named f = ref_filt()
@named c = controller(1)
@named p = plant()
connections = [f.y ~ c.r # filtered reference to controller reference
c.u ~ p.u # controller output to plant input
p.y ~ c.y]
@named cl = ODESystem(connections, t, systems = [f, c, p])
lsys0, ssys = linearize(cl, [f.u], [p.x])
desired_order = [f.x, p.x]
lsys = ModelingToolkit.reorder_unknowns(lsys0, unknowns(ssys), desired_order)
@assert lsys.A == [-2 0; 1 -2]
@assert lsys.B == [1; 0;;]
@assert lsys.C == [0 1]
@assert lsys.D[] == 0
## Symbolic linearization
lsys_sym, _ = ModelingToolkit.linearize_symbolic(cl, [f.u], [p.x])
@assert substitute(lsys_sym.A, ModelingToolkit.defaults(cl)) == lsys.A
```
"""
function linearize(sys, lin_fun::LinearizationFunction; t = 0.0,
op = Dict(), allow_input_derivatives = false,
p = DiffEqBase.NullParameters())
prob = LinearizationProblem(lin_fun, t)
op = anydict(op)
evaluate_varmap!(op, unknowns(sys))
for (k, v) in op
v === nothing && continue
if is_parameter(prob, Initial(k))
setu(prob, Initial(k))(prob, v)
else
setu(prob, k)(prob, v)
end
end
p = anydict(p)
for (k, v) in p
setu(prob, k)(prob, v)
end
return solve(prob; allow_input_derivatives)
end
function linearize(sys, inputs, outputs; op = Dict(), t = 0.0,
allow_input_derivatives = false,
zero_dummy_der = false,
kwargs...)
lin_fun, ssys = linearization_function(sys,
inputs,
outputs;
zero_dummy_der,
op,
kwargs...)
mats, extras = linearize(ssys, lin_fun; op, t, allow_input_derivatives)
mats, ssys, extras
end
"""
(; Ã, B̃, C̃, D̃) = similarity_transform(sys, T; unitary=false)
Perform a similarity transform `T : Tx̃ = x` on linear system represented by matrices in NamedTuple `sys` such that
```
à = T⁻¹AT
B̃ = T⁻¹ B
C̃ = CT
D̃ = D
```
If `unitary=true`, `T` is assumed unitary and the matrix adjoint is used instead of the inverse.
"""
function similarity_transform(sys::NamedTuple, T; unitary = false)
if unitary
A = T'sys.A * T
B = T'sys.B
else
Tf = lu(T)
A = Tf \ sys.A * T
B = Tf \ sys.B
end
C = sys.C * T
D = sys.D
(; A, B, C, D)
end
"""
reorder_unknowns(sys::NamedTuple, old, new)
Permute the state representation of `sys` obtained from [`linearize`](@ref) so that the state unknown is changed from `old` to `new`
Example:
```
lsys, ssys = linearize(pid, [reference.u, measurement.u], [ctr_output.u])
desired_order = [int.x, der.x] # Unknowns that are present in unknowns(ssys)
lsys = ModelingToolkit.reorder_unknowns(lsys, unknowns(ssys), desired_order)
```
See also [`ModelingToolkit.similarity_transform`](@ref)
"""
function reorder_unknowns(sys::NamedTuple, old, new)
nx = length(old)
length(new) == nx || error("old and new must have the same length")
perm = [findfirst(isequal(n), old) for n in new]
issorted(perm) && return sys # shortcut return, no reordering
P = zeros(Int, nx, nx)
for i in 1:nx # Build permutation matrix
P[i, perm[i]] = 1
end
similarity_transform(sys, P; unitary = true)
end