|  | 
|  | 1 | +# Linear Analysis | 
|  | 2 | + | 
|  | 3 | +Linear analysis refers to the process of linearizing a nonlinear model and analysing the resulting linear dynamical system. To facilitate linear analysis, ModelingToolkit provides the concept of an [`AnalysisPoint`](@ref), which can be inserted in-between two causal blocks (such as those from `ModelingToolkitStandardLibrary.Blocks` sub module). Once a model containing analysis points is built, several operations are available: | 
|  | 4 | + | 
|  | 5 | +  - [`get_sensitivity`](@ref) get the [sensitivity function (wiki)](https://en.wikipedia.org/wiki/Sensitivity_(control_systems)), $S(s)$, as defined in the field of control theory. | 
|  | 6 | +  - [`get_comp_sensitivity`](@ref) get the complementary sensitivity function $T(s) : S(s)+T(s)=1$. | 
|  | 7 | +  - [`get_looptransfer`](@ref) get the (open) loop-transfer function where the loop starts and ends in the analysis point. For a typical simple feedback connection with a plant $P(s)$ and a controller $C(s)$, the loop-transfer function at the plant output is $P(s)C(s)$. | 
|  | 8 | +  - [`linearize`](@ref) can be called with two analysis points denoting the input and output of the linearized system. | 
|  | 9 | +  - [`open_loop`](@ref) return a new (nonlinear) system where the loop has been broken in the analysis point, i.e., the connection the analysis point usually implies has been removed. | 
|  | 10 | + | 
|  | 11 | +An analysis point can be created explicitly using the constructor [`AnalysisPoint`](@ref), or automatically when connecting two causal components using `connect`: | 
|  | 12 | + | 
|  | 13 | +```julia | 
|  | 14 | +connect(comp1.output, :analysis_point_name, comp2.input) | 
|  | 15 | +``` | 
|  | 16 | + | 
|  | 17 | +A single output can also be connected to multiple inputs: | 
|  | 18 | + | 
|  | 19 | +```julia | 
|  | 20 | +connect(comp1.output, :analysis_point_name, comp2.input, comp3.input, comp4.input) | 
|  | 21 | +``` | 
|  | 22 | + | 
|  | 23 | +!!! warning "Causality" | 
|  | 24 | +     | 
|  | 25 | +    Analysis points are *causal*, i.e., they imply a directionality for the flow of information. The order of the connections in the connect statement is thus important, i.e., `connect(out, :name, in)` is different from `connect(in, :name, out)`. | 
|  | 26 | + | 
|  | 27 | +The directionality of an analysis point can be thought of as an arrow in a block diagram, where the name of the analysis point applies to the arrow itself. | 
|  | 28 | + | 
|  | 29 | +``` | 
|  | 30 | +┌─────┐         ┌─────┐ | 
|  | 31 | +│     │  name   │     │ | 
|  | 32 | +│  out├────────►│in   │ | 
|  | 33 | +│     │         │     │ | 
|  | 34 | +└─────┘         └─────┘ | 
|  | 35 | +``` | 
|  | 36 | + | 
|  | 37 | +This is signified by the name being the middle argument to `connect`. | 
|  | 38 | + | 
|  | 39 | +Of the above mentioned functions, all except for [`open_loop`](@ref) return the output of [`ModelingToolkit.linearize`](@ref), which is | 
|  | 40 | + | 
|  | 41 | +```julia | 
|  | 42 | +matrices, simplified_sys = linearize(...) | 
|  | 43 | +# matrices = (; A, B, C, D) | 
|  | 44 | +``` | 
|  | 45 | + | 
|  | 46 | +i.e., `matrices` is a named tuple containing the matrices of a linear state-space system on the form | 
|  | 47 | + | 
|  | 48 | +```math | 
|  | 49 | +\begin{aligned} | 
|  | 50 | +\dot x &= Ax + Bu\\ | 
|  | 51 | +y &= Cx + Du | 
|  | 52 | +\end{aligned} | 
|  | 53 | +``` | 
|  | 54 | + | 
|  | 55 | +## Example | 
|  | 56 | + | 
|  | 57 | +The following example builds a simple closed-loop system with a plant $P$ and a controller $C$. Two analysis points are inserted, one before and one after $P$. We then derive a number of sensitivity functions and show the corresponding code using the package ControlSystemBase.jl | 
|  | 58 | + | 
|  | 59 | +```@example LINEAR_ANALYSIS | 
|  | 60 | +using ModelingToolkitStandardLibrary.Blocks, ModelingToolkit | 
|  | 61 | +@named P = FirstOrder(k = 1, T = 1) # A first-order system with pole in -1 | 
|  | 62 | +@named C = Gain(-1)             # A P controller | 
|  | 63 | +t = ModelingToolkit.get_iv(P) | 
|  | 64 | +eqs = [connect(P.output, :plant_output, C.input)  # Connect with an automatically created analysis point called :plant_output | 
|  | 65 | +       connect(C.output, :plant_input, P.input)] | 
|  | 66 | +sys = ODESystem(eqs, t, systems = [P, C], name = :feedback_system) | 
|  | 67 | +
 | 
|  | 68 | +matrices_S = get_sensitivity(sys, :plant_input)[1] # Compute the matrices of a state-space representation of the (input)sensitivity function. | 
|  | 69 | +matrices_T = get_comp_sensitivity(sys, :plant_input)[1] | 
|  | 70 | +``` | 
|  | 71 | + | 
|  | 72 | +Continued linear analysis and design can be performed using ControlSystemsBase.jl. | 
|  | 73 | +We create `ControlSystemsBase.StateSpace` objects using | 
|  | 74 | + | 
|  | 75 | +```@example LINEAR_ANALYSIS | 
|  | 76 | +using ControlSystemsBase, Plots | 
|  | 77 | +S = ss(matrices_S...) | 
|  | 78 | +T = ss(matrices_T...) | 
|  | 79 | +bodeplot([S, T], lab = ["S" "" "T" ""], plot_title = "Bode plot of sensitivity functions", | 
|  | 80 | +    margin = 5Plots.mm) | 
|  | 81 | +``` | 
|  | 82 | + | 
|  | 83 | +The sensitivity functions obtained this way should be equivalent to the ones obtained with the code below | 
|  | 84 | + | 
|  | 85 | +```@example LINEAR_ANALYSIS_CS | 
|  | 86 | +using ControlSystemsBase | 
|  | 87 | +P = tf(1.0, [1, 1]) |> ss | 
|  | 88 | +C = 1                      # Negative feedback assumed in ControlSystems | 
|  | 89 | +S = sensitivity(P, C)      # or feedback(1, P*C) | 
|  | 90 | +T = comp_sensitivity(P, C) # or feedback(P*C) | 
|  | 91 | +``` | 
|  | 92 | + | 
|  | 93 | +We may also derive the loop-transfer function $L(s) = P(s)C(s)$ using | 
|  | 94 | + | 
|  | 95 | +```@example LINEAR_ANALYSIS | 
|  | 96 | +matrices_L = get_looptransfer(sys, :plant_output)[1] | 
|  | 97 | +L = ss(matrices_L...) | 
|  | 98 | +``` | 
|  | 99 | + | 
|  | 100 | +which is equivalent to the following with ControlSystems | 
|  | 101 | + | 
|  | 102 | +```@example LINEAR_ANALYSIS_CS | 
|  | 103 | +L = P * (-C) # Add the minus sign to build the negative feedback into the controller | 
|  | 104 | +``` | 
|  | 105 | + | 
|  | 106 | +To obtain the transfer function between two analysis points, we call `linearize` | 
|  | 107 | + | 
|  | 108 | +```@example LINEAR_ANALYSIS | 
|  | 109 | +using ModelingToolkit # hide | 
|  | 110 | +matrices_PS = linearize(sys, :plant_input, :plant_output)[1] | 
|  | 111 | +``` | 
|  | 112 | + | 
|  | 113 | +this particular transfer function should be equivalent to the linear system `P(s)S(s)`, i.e., equivalent to | 
|  | 114 | + | 
|  | 115 | +```@example LINEAR_ANALYSIS_CS | 
|  | 116 | +feedback(P, C) | 
|  | 117 | +``` | 
|  | 118 | + | 
|  | 119 | +### Obtaining transfer functions | 
|  | 120 | + | 
|  | 121 | +A statespace system from [ControlSystemsBase](https://juliacontrol.github.io/ControlSystems.jl/stable/man/creating_systems/) can be converted to a transfer function using the function `tf`: | 
|  | 122 | + | 
|  | 123 | +```@example LINEAR_ANALYSIS_CS | 
|  | 124 | +tf(S) | 
|  | 125 | +``` | 
|  | 126 | + | 
|  | 127 | +## Gain and phase margins | 
|  | 128 | + | 
|  | 129 | +Further linear analysis can be performed using the [analysis methods from ControlSystemsBase](https://juliacontrol.github.io/ControlSystems.jl/stable/lib/analysis/). For example, calculating the gain and phase margins of a system can be done using | 
|  | 130 | + | 
|  | 131 | +```@example LINEAR_ANALYSIS_CS | 
|  | 132 | +margin(P) | 
|  | 133 | +``` | 
|  | 134 | + | 
|  | 135 | +(they are infinite for this system). A Nyquist plot can be produced using | 
|  | 136 | + | 
|  | 137 | +```@example LINEAR_ANALYSIS_CS | 
|  | 138 | +nyquistplot(P) | 
|  | 139 | +``` | 
|  | 140 | + | 
|  | 141 | +## Index | 
|  | 142 | + | 
|  | 143 | +```@index | 
|  | 144 | +Pages = ["linear_analysis.md"] | 
|  | 145 | +``` | 
|  | 146 | + | 
|  | 147 | +```@autodocs | 
|  | 148 | +Modules = [ModelingToolkit] | 
|  | 149 | +Pages   = ["systems/analysis_points.jl"] | 
|  | 150 | +Order   = [:function, :type] | 
|  | 151 | +Private = false | 
|  | 152 | +``` | 
0 commit comments