PortHamiltonianModelReduction.TimeDomainData
PortHamiltonianModelReduction.bt
PortHamiltonianModelReduction.chirp
PortHamiltonianModelReduction.ephminreal
PortHamiltonianModelReduction.estimate_ham
PortHamiltonianModelReduction.hdss
PortHamiltonianModelReduction.irka
PortHamiltonianModelReduction.klap
PortHamiltonianModelReduction.klap_inital_guess
PortHamiltonianModelReduction.matchnrg
PortHamiltonianModelReduction.matchnrg_are
PortHamiltonianModelReduction.matchnrg_barrier
PortHamiltonianModelReduction.matchnrg_sdp
PortHamiltonianModelReduction.opinf
PortHamiltonianModelReduction.passivate
PortHamiltonianModelReduction.passivate_lmi
PortHamiltonianModelReduction.passivate_lmi_tp
PortHamiltonianModelReduction.pgprojection
PortHamiltonianModelReduction.phdmd
PortHamiltonianModelReduction.phdmd_initial_guess
PortHamiltonianModelReduction.phdmd_sdp
PortHamiltonianModelReduction.phirka
PortHamiltonianModelReduction.pod
PortHamiltonianModelReduction.podbasis
PortHamiltonianModelReduction.prbs
PortHamiltonianModelReduction.prbt
PortHamiltonianModelReduction.sawtooth
PortHamiltonianModelReduction.skew_proc
PortHamiltonianModelReduction.step
PortHamiltonianModelReduction.tddata
PortHamiltonianModelReduction.TimeDomainData
— TypeTimeDomainData{T}
Data structure to hold time-domain simulation data. See tddata
for a convenient constructor.
Fields
Ẋ
: State derivative matrix.X
: State matrix.U
: Input matrix.Y
: Output matrix.t
: Time vector.
PortHamiltonianModelReduction.bt
— MethodΣr = bt(Σ::StateSpace, r; Lx=grampd(Σ, :o), Ly=grampd(Σ, :c)')
Reduces the state dimension of the system Σ
to r
using standard square root balanced truncation (see for instance [1]). The cholesky factors of the Gramians can passed as optional arguments Lx
and Ly
. The default values are the cholesky factors of the observability and controllability Gramians.
PortHamiltonianModelReduction.chirp
— Functionchirp(A=1.0, f0=1e-2, f1=1e1, T=2.0)
Generates a chirp signal with amplitude A
, starting frequency f0
, ending frequency f1
, and duration T
.
PortHamiltonianModelReduction.ephminreal
— MethodΣph = ephminreal(Σph::PortHamiltonianStateSpace; kwargs...)
Computes a minimal realization of the (extended) pH system (see [2, Cor. 4.6])
PortHamiltonianModelReduction.estimate_ham
— MethodQ, err = estimate_ham(data::TimeDomainData, H::AbstractVector; kwargs...)
Estimates the Hamiltonian Hessian Q
from time-domain data
and measurements H
of the Hamiltonian.
PortHamiltonianModelReduction.hdss
— MethodΣqo = hdss(J, R, Q, G, P, S, N)
Σqo = hdss(Σph)
Converts a PortHamiltonianStateSpace
to a QuadraticOutputStateSpace
(Hamiltonian dynamic).
PortHamiltonianModelReduction.irka
— MethodΣr = irka(Σ::StateSpace, r; tol=1e-3, max_iter=200)
Reduces the state dimension of the system Σ
to r
using the iterative rational Krylov algorithm (IRKA) [3].
PortHamiltonianModelReduction.klap
— Functionklap(Σ::StateSpace; L0=L0(Σ), M=M(Σ), P=gram(Σ, :c); recycl=:schur, restart=false, α=1e-8, ε=1e-4, verbose=true, kwargs...) -> Σp, res
Passivates a system Σ
using KLAP [4]. The optimization problem is solved using LBFGS.
PortHamiltonianModelReduction.klap_inital_guess
— Functionklap_inital_guess(Σ, ΔD=0.0; ε=1e-8) -> L0, ΔD
Computes an initial guess for KLAP [4]. The initial guess is computed by perturbing the feedthrough matrix to achieve a passive realization. Then the perturbed system is used to compute the initial guess. The perturbation ΔD
can be specified, otherwise it is computed using ΔD(Σ)
.
PortHamiltonianModelReduction.matchnrg
— MethodΣrem = matchnrg(Σ::PortHamiltonianStateSpace, Σr::PortHamiltonianStateSpace; solver=:BFGS, kwargs...)
Applies energy matching [2] to the ROM Σr
to match the Hamiltonian dynamics of the original system Σ
. If solver is :Hypatia
or :COSMO
, it uses semidefinite programming with the specified optimizer. If solver is :BFGS
, it uses the barrier method with BFGS optimization. If solver is :ARE
, the best solution of the ARE is returned.
PortHamiltonianModelReduction.matchnrg_are
— MethodΣphr = matchnrg_are(Σ::QuadraticOutputStateSpace, Σr::StateSpace; kwargs...)
Solves the energy matching problem where the solution set is replaced with solutions of the ARE.
PortHamiltonianModelReduction.matchnrg_barrier
— MethodΣph = matchnrg_barrier(Σ, Σr; Σr0=nothing, kwargs...)
Solves the energy matching problem using the barrier method.
PortHamiltonianModelReduction.matchnrg_sdp
— FunctionΣphr = matchnrg_sdp(Σ::QuadraticOutputStateSpace, Σr::StateSpace; optimizer=COSMO.Optimizer, ε=1e-8, kwargs...)
Solves the energy matching problem using semidefinite programming.
PortHamiltonianModelReduction.opinf
— MethodΣr, res = opinf(data::TimeDomainData, Wr::Matrix)
Computes a ROM Σr
from time-domain data via operator inference (OpInf) [5]. The error of the OpInf problem is returned as res
. The data is projected using with the matrix Wr
, which is typically a POD basis.
PortHamiltonianModelReduction.passivate
— Functionpassivate(Σ::StateSpace, method=:klap, args...; kwargs...) -> Σp, res
Passivate the system Σ
using the method method
. The available methods are:
:klap
: KLAP optimization [4]:lmi
: LMI optimization [6]:lmi_tp
: LMI optimization with trace parametrization [7, 8]
The remaining arguments args
and keyword arguments kwargs
are passed to the corresponding passivation function.
PortHamiltonianModelReduction.passivate_lmi
— Methodpassivate_lmi(Σ::StateSpace; kwargs...) -> Σp, model
Solves the passivation problem for a given state-space system Σ
using the positive real LMI constraints [6].
PortHamiltonianModelReduction.passivate_lmi_tp
— Methodpassivate_lmi_tp(Σ::StateSpace; kwargs...) -> Σp, model
Solves the passivation problem for a given state-space system Σ
using the positive real LMI constraints with trace parametrization [7, 8].
PortHamiltonianModelReduction.pgprojection
— Functionpgprojection(Σph::PortHamiltonianStateSpace, V::Matrix, W::Matrix=Σph.Q * V / (V' * Σph.Q * V))
pgprojection(Σ::StateSpace, V::Matrix, W::Matrix=V)
Applies a (Petrov-)Galerkin projection to a pH or LTI. In the pH case, the W
matrix has to be chosen such that the structure is preserved. The default choice is W = Σph.Q * V / (V' * Σph.Q * V)
, which is proposed in [9].
PortHamiltonianModelReduction.phdmd
— MethodΣph, err = phdmd(data::TimeDomainData, H::AbstractVector; kwargs...)
Σph, err = phdmd(data::TimeDomainData, Q::AbstractMatrix; kwargs...)
Computes a pH system Σph
that approximates given time-domain data
and a candidate Q
for the Hessian of the Hamiltonian or measurements H
of the Hamiltonian using the pHDMD method [10].
PortHamiltonianModelReduction.phdmd_initial_guess
— MethodΓ, W, err = phdmd_initial_guess(data::TimeDomainData, Q::AbstractMatrix; kwargs...)
Computes an initial guess Γ, W
for the pHDMD problem from time-domain data data
and a candidate Q
for the Hessian of the Hamiltonian (see [10, Thm. 3.7]).
PortHamiltonianModelReduction.phdmd_sdp
— MethodΣph = phdmd_sdp(data::TimeDomainData; kwargs...)
Σph = phdmd_sdp(data::TimeDomainData, Q::AbstractMatrix; kwargs...)
Computes a pH system Σph
that approximates given time-domain data
using semidefinite programming. Either the Hessian of the Hamiltonian Q
is provided or also learned from the data.
PortHamiltonianModelReduction.phirka
— MethodΣr = phirka(Σph::PortHamiltonianStateSpace, r; tol=1e-3, max_iter=50)
Reduces the state dimension of the port-Hamiltonian system Σph
to r
using pHIRKA [9].
PortHamiltonianModelReduction.pod
— MethodΣr = pod(Σ, X, r)
Σr = pod(Σ, V, W)
Returns a ROM Σr
of size r
for the system Σ
using the POD (Petrov-)Galerkin method.
PortHamiltonianModelReduction.podbasis
— Methodpodbasis(X::AbstractMatrix, r::Int; kwargs...)
podbasis(F::SVD, r::Int; kwargs...)
Returns the POD basis via the truncated singular value decomposition of a data matrix X
. Alternatively, it can be cast on a F::SVD
object, such that the SVD is not recomputed.
PortHamiltonianModelReduction.prbs
— Functionprbs(t, A=1.0, order=7)
Generates a pseudo-random binary sequence (PRBS) signal with amplitude A
and specified order
.
PortHamiltonianModelReduction.prbt
— MethodΣr = prbt(Σ, r; Lx=prgrampd(Σ, :o), Ly=prgrampd(Σ, :c))
Reduces the state dimension of the system Σ
to r
using positive real balanced truncation (PRBT) [11]. The cholesky factors of the positive real Gramians can passed as optional arguments Lx
and Ly
, in order to avoid recomputation. Σ
can be a StateSpace
or a PortHamiltonianStateSpace
, which will then also be the return type. In the case of a PortHamiltonianStateSpace
, for the conversion to the ROM the minimal solution of the KYP inequality is used as the Hessian of the Hamiltonian.
PortHamiltonianModelReduction.sawtooth
— Functionsawtooth(A=1.0, T=1.0)
Generates a sawtooth wave signal with amplitude A
and period T
.
PortHamiltonianModelReduction.skew_proc
— MethodΓ, err = skew_proc(T::AbstractMatrix, Z::AbstractMatrix; kwargs...)
Computes the analytic solution of the skew-symmetric Procrustes problem (see for instance [10, Thm. 3.4]).
PortHamiltonianModelReduction.step
— Functionstep(A=1.0, ts=1.0)
Generates a step signal with amplitude A
and step time ts
.
PortHamiltonianModelReduction.tddata
— Methoddata = tddata(Ẋ::Matrix{T},X::Matrix{T},U::Matrix{T},Y::Matrix{T},t::Vector{T})
data = tddata(res::SimResult)
Creates a TimeDomainData
object with element type T
.