PortHamiltonianSystems.PortHamiltonianStateSpace
ControlSystemsBase.gram
ControlSystemsBase.grampd
ControlSystemsBase.ss
LinearAlgebra.norm
PortHamiltonianSystems.W
PortHamiltonianSystems.compose
PortHamiltonianSystems.decompose
PortHamiltonianSystems.ispassive
PortHamiltonianSystems.ispsd
PortHamiltonianSystems.kyp
PortHamiltonianSystems.kypmat
PortHamiltonianSystems.kypmax
PortHamiltonianSystems.kypmin
PortHamiltonianSystems.lrcholesky
PortHamiltonianSystems.phminreal
PortHamiltonianSystems.phss
PortHamiltonianSystems.phss
PortHamiltonianSystems.popov
PortHamiltonianSystems.prare
PortHamiltonianSystems.prgram
PortHamiltonianSystems.prgrampd
PortHamiltonianSystems.project_psd
PortHamiltonianSystems.sampopov
PortHamiltonianSystems.truncation
PortHamiltonianSystems.tsvd
PortHamiltonianSystems.unvec
PortHamiltonianSystems.unvech
PortHamiltonianSystems.Γ
PortHamiltonianSystems.PortHamiltonianStateSpace
— TypePortHamiltonianStateSpace{T}
An object representing a port-Hamiltonian state-space system.
dx(t)/dt = (J-R)Qx(t) + (G-P)u(t)
y(t) = (G+P)'Qx(t) + (S-N)u(t)
See the function phss
for a user facing constructor.
Fields:
J::SkewHermitian{T}
R::Symmetric{T}
Q::Symmetric{T}
G::Matrix{T}
P::Matrix{T}
S::Symmetric{T}
N::SkewHermitian{T}
ControlSystemsBase.gram
— MethodX = gram(Σph, opt; kwargs...)
Returns the Gramian of the system Σph
. If opt
is :c
or :o
it returns the controllability or observability Gramian, respectively by calling the gram
from ControlSystems package. If opt
is :prc
or :pro
it returns the positive-real controllability or observability Gramian, respectively (see prgram
).
ControlSystemsBase.grampd
— MethodL = grampd(Σph, opt; kwargs...)
Returns a Cholesky factor L
of the Gramian of the system Σph
. If opt
is :c
or :o
it returns the controllability or observability Gramian, respectively by calling the grampd
from ControlSystems.jl package. If opt
is :prc
or :pro
it returns a Cholesky factor L
of the positive-real controllability or observability Gramian, respectively.
ControlSystemsBase.ss
— MethodΣ = ss(J, R, Q, G, P, S, N)
Σ = ss(Σph)
Converts a PortHamiltonianStateSpace
to a standard StateSpace
.
LinearAlgebra.norm
— Functionnorm(Σph, p=2; kwargs...)
Converts a PortHamiltonianStateSpace
to a StateSpace
and calls norm
on it. For more details see ControlSystems.jl package.
PortHamiltonianSystems.W
— MethodW(Σ)
Returns the dissipation matrix of the pH system.
PortHamiltonianSystems.compose
— MethodA, B, C, D = compose(J, R, Q, G, P, S, N)
Composes the port-Hamiltonian matrices to standard state-space matrices according to
A = (J - R) * Q
B = G - P
C = (G + P)' * Q
D = S - N.
PortHamiltonianSystems.decompose
— MethodJ, R, Q, G, P, S, N = decompose(A, B, C, D, X)
Decomposes the standard state-space matrices to port-Hamiltonian matrices for a given Hamiltonian X
according to
Q = X
J = skew(A / X)
R = -sym(A / X)
G = 0.5 * (X \ C' + B)
P = 0.5 * (X \ C' - B)
S = sym(D)
N = skew(D).
PortHamiltonianSystems.ispassive
— Methodispassive(Σ::StateSpace; opt=:lmi, kwargs...)
ispassive(Σ::PortHamiltonianStateSpace; kwargs...)
Checks whether the system Σ
is passive by solving the KYP inequality using kyp
(if opt=:lmi
) or checking the Popov function for passivity violations via sampopov
(if opt=:popov
).
PortHamiltonianSystems.ispsd
— Methodispsd(M)
Returns true
if M
is positive semi-definite, otherwise false
.
PortHamiltonianSystems.kyp
— MethodX = kyp(Σ; kwargs...)
Tries to solve the KYP inequality via semi definite programming.
PortHamiltonianSystems.kypmat
— MethodW = kypmat(Σ, X)
Returns the KYP matrix of the system Σ
for a given matrix X
.
PortHamiltonianSystems.kypmax
— MethodXmin = kypmax(Σ; kwargs...)
Returns the maximal solution to the KYP inequality by solving the ARE equation for the anti-stabilizing solution.
PortHamiltonianSystems.kypmin
— MethodXmin = kypmin(Σ; kwargs...)
Returns the minimal solution to the KYP inequality by solving the ARE for the stabilizing solution.
PortHamiltonianSystems.lrcholesky
— MethodL = lrcholesky(X; trunc_tol=1e-12)
Computes a low-rank approximate Cholesky-like factorization of a symmetric positive semi-definite matrix $X$ s.t. $X = L * L'$ (up to a prescribed tolerance trunc_tol
).
PortHamiltonianSystems.phminreal
— Methodphminreal(Σph::PortHamiltonianStateSpace; trunc_tol=1e-12)
Computes a structure preserving minimal realization of a port-Hamiltonian system [1].
PortHamiltonianSystems.phss
— MethodΣph = phss(J, R, Q, G, P, S, N)
Σph = phss(J, R, Q, G)
Σph = phss(Γ, W, Q)
Creates a port-Hamiltonian state-space model Σph::PortHamiltonianStateSpace{T}
with matrix element type T
. Note that J
, N
must be skew-symmetric and R
, Q
, S
must be symmetric. When using the structure and dissipation matrices Γ
and W
, they must also satisfy these properties, i.e., Γ
must be skew-symmetric and W
must be symmetric.
PortHamiltonianSystems.phss
— MethodΣph = phss(Σ)
Σph = phss(Σ, X)
Converts a passive StateSpace
to a PortHamiltonianStateSpace
by executing decompose
. If X
is not provided, the minimal solution of the KYP inequality is used, which is obtained by solving the ARE.
PortHamiltonianSystems.popov
— Methodpopov(Σ, s)
popov(Σ, ω)
Evaluates the popov function of the system Σ
at the complex variable s
.
Φ(s) = Σ(s) + Σ(-s)'
where Σ(s)
is the frequency response (transfer function) of Σ
at the complex variable s
.
PortHamiltonianSystems.prare
— Methodprare(Σ, opt, X; kwargs...)
Evaluates the positive-real algebraic Riccati equation for system Σ
and candidat solution X
.
If opt
is :o
the positive-real controllability algebraic Riccati equation is evaluated,
A'X + XA + (C' - XB) inv(D + D') (C - B'X).
If opt
is :c
the positive-real observability algebraic Riccati equation is evaluated,
AX + XA' + (B - XC') inv(D + D') (B' - CX).
PortHamiltonianSystems.prgram
— MethodX = prgram(Σ, opt; kwargs...)
Returns the positive-real Gramian of system Σ
. If opt
is :c
or :o
it returns the positive-real controllability or positive-real observability Gramian, respectively, by solving the corresponding positive-real algebraic Riccati equation (see prare
).
PortHamiltonianSystems.prgrampd
— MethodL = prgrampd(Σ, opt; kwargs...)
Returns the Cholesky factor of the positive-real Gramian of system Σ
(see prgram
). If opt
is :c
or :o
it returns the positive-real controllability or positive-real observability Gramian, respectively.
In the case that the solution of the positive-real algebraic Riccati equation is not positive definite (due to numerical errors), it is projected to the set of positive semi-definite matrices calling project_psd
.
PortHamiltonianSystems.project_psd
— MethodMpsd = project_psd(M; eigtol=1e-8)
Returns the nearest positive semi-definite matrix to M
by setting negative eigenvalues to eigtol
.
PortHamiltonianSystems.sampopov
— Methodsampopov(Σ; ω=10 .^ range(-15, stop=5, length=5000))
Samples the Popov function for the system Σ
at ranges of frequencies ω
.
PortHamiltonianSystems.truncation
— Methodtruncation(d, L, trunc_tol) -> (dr, Lr)
Computes a rank revealing factorization for a given LDL-decomposition of $S = L * \mathrm{diag}(d) * L^T$ (up to a prescribed tolerance trunc_tol
) such that $L_r * diag(d_r) * L_r^T \approx S$.
PortHamiltonianSystems.tsvd
— MethodF = tsvd(A; ε=1e-12, kwargs...)
F = tsvd(A, r; kwargs...)
Returns the truncated singular value decomposition of A
by truncating small singular values below ε
or by retaining only the top r
singular values.
PortHamiltonianSystems.unvec
— MethodM = unvec(v)
Returns the matrix M
from the vectorized v
, i.e., the inverse of LinearAlgebra.vec
.
PortHamiltonianSystems.unvech
— MethodM = unvech(v)
Returns the Hermitian matrix M
from the half-vectorized v
, i.e., the inverse of VectorizationTransforms.vech
.
PortHamiltonianSystems.Γ
— MethodΓ(Σ)
Returns the structure matrix of the pH system.