ml_ct_dss_flbt

Frequency-limited balanced truncation for desc. systems.

Contents

Syntax

[Ar, Br, Cr, Dr, Er, info] = ml_ct_dss_flbt(A, B, C, D, E)
[Ar, Br, Cr, Dr, Er, info] = ml_ct_dss_flbt(A, B, C, D, E, opts)
[rom, info] = ml_ct_dss_flbt(sys)
[rom, info] = ml_ct_dss_flbt(sys, opts)

Description

This function computes the generalized frequency-limited balanced
truncation for a descriptor system of the form
    E*x'(t) = A*x(t) + B*u(t),                                      (1)
       y(t) = C*x(t) + D*u(t).                                      (2)
Therefore, first an additive decomposition of the system is performed
using the matrix disk function, such that
         [ Ei  0 ]       [ Ai  0 ]                        [ Ci ]
    E2 = [       ], A2 = [       ], B2 = [ Bi, Bp ], C2 = [    ],
         [ 0  Ep ]       [ 0  Ap ]                        [ Cp ]
with (Ei, Ai, Bi, Ci, D) belonging to the polynomial part and
(Ep, Ap, Bp, Cp, 0) belonging to the strictly proper part.
Now, the two generalized continuous-time Lyapunov equations
    Ap*Pp*Ep' + Ep*Pp*Ap' + BF*Bp' + Bp*BF' = 0,
    Ap'*Qp*Ep + Ep'*Qp*Ap + Cp'*CF + CF'*Cp = 0,
where BF and CF are frequency-dependent matrices, are solved for the
reduction of the strictly proper part, and the two generalized
discrete-time Lyapunov equations
    Ai*Pi*Ai' - Ei*Pi*Ei' - Bi*Bi' = 0,
    Ai'*Qi*Ai - Ei'*Qi*Ei - Ci'*Ci = 0
are solved for the reduction of the polynomial part. As result, a
reduced-order system of the form
    Er*x'(t) = Ar*x(t) + Br*u(t),                                   (3)
        y(t) = Cr*x(t) + Dr*u(t)                                    (4)
is computed. For enforcing stability in the reduced-order model, the
modified Gramian approach can be used, which also gives a global error
bound of the form
    ||G - Gr||_{\infty} <= 2*||JB||*||JC||(Hsvp(r+1) + ... + Hsvp(n)),
with Hsvp, a vector containing the proper frequency-limited Hankel
singular values of the system.
Note: For unstable systems, an additional additive decomposition into
      the stable and anti-stable parts is performed and then only the
      stable part will be reduced.

Input

A    - matrix from (1) with dimensions n x n
B    - matrix from (1) with dimensions n x m
C    - matrix from (2) with dimensions p x n
D    - matrix from (2) with dimensions p x m
E    - matrix from (1) with dimensions n x n
sys  - structure or state-space object, containing the descriptor
       system's matrices:

Entry
Meaning
A
matrix from (1) with dimensions n x n
B
matrix from (1) with dimensions n x m
C
matrix from (2) with dimensions p x n
D
matrix from (2) with dimensions p x m
E
matrix from (1) with dimensions n x n

opts - structure, containing the following optional entries:

Parameter
Meaning
DecompEig
positive scalar, overestimation of the absolute value of the largest finite eigenvalue of s*E - A, if set, replaces the computation with DecompTol
default: []
DecompTol
nonnegative scalar, tolerance multiplied with the largest singular value of E to determine the smallest non-quasi-zero singular value of E
default: log(n)*eps
FreqRange
nonnegative vector, frequency intervals such that [w(1), w(2)] ... [w(2k-1), w(2k)] are approximated
default: [0, 1.0e+03]
gdlyapopts
structure, containing the optional parameters for the computation of the generalized discrete-time Lyapunov equations, see ml_gdlyap_smith_fac
default: struct()
ImproperTrunc
nonnegative scalar, tolerance multiplied with the largest proper Hankel singular value of the system to truncate the improper part
default: log(n)*eps
Index
nonnegative integer, index of the descriptor system used to set an upper bound on the size of the reduced improper part, Inf if unknown
default: Inf
infdecopts
structure, containing the optional parameters for the decomposition of the finite and infinite parts of the system using the disk function and subspace extraction method, see ml_disk and ml_getqz
default: struct()
lyapdlopts
structure, containing the optional parameters for the computation of the generalized continuous-time Lyapunov equations, see ml_lyapdl_sgn_ldl if ModGramian = 0 and ml_lyapdl_sgn_fac if ModGramian = 1
default: struct()
Method
character array, determining algorithm for the computation of the reduced-order model
  • 'sr' - square-root method
  • 'bfsr' - balancing-free square-root method
default: 'sr'
ModGramian
{0, 1}, used to disable/enable the modified Gramian approach
default: 0
Order
positive integer, order of the resulting reduced-order model chosen by the user if 'order' is set for OrderComputation
default: min(10,length(Hsvp)) + Nu + Ni
OrderComputation
character array, determining the method for the computation of the size of the reduced-order model
  • 'order' - take explicit order
  • 'tolerance' - using rel. tolerance for the hsv
  • 'sum' - using rel. tolerance for sum of hsv
default: 'sum'
stabdecopts
structure, containing the optional parameters for the decomposition of the stable and unstable parts of the system using the disk function and subspace extraction method, see ml_disk and ml_getqz
default: struct()
Tolerance
nonnegative scalar, tolerance used for the computation of the size of the reduced-order model if 'tolerance' or 'sum' is set for OrderComputation
default: 1.0e-02

Output

Ar   - matrix of (3) with dimensions r x r
Br   - matrix of (3) with dimensions r x m
Cr   - matrix of (4) with dimensions p x r
Dr   - matrix of (4) with dimensions p x m
Er   - matrix of (3) with dimensions r x r
rom  - structure or state-space object, containing the reduced-order
       descriptor system:

Entry
Meaning
A
matrix from (3) with dimensions r x r
B
matrix from (3) with dimensions r x m
C
matrix from (4) with dimensions p x r
D
matrix from (4) with dimensions p x m
E
matrix from (3) with dimensions r x r

info - structure, containing the following information:

Entry
Meaning
AbsErrBound
computed error bound for the absolute error of the reduced-order model in H-infinity norm, only for the modified Gramian approach
Hsvi
a vector, containing the computed Hankel singular values of the improper part of the system
Hsvp
a vector, containing the computed Hankel singular values of the proper part of the system
infoGADTF
structure, containing information about the additive decomposition of the system into its infinite, finite stable and finite anti-stable parts, see ml_ct_dss_adtf
infoGDLYAP_C
structure, containing information about the generalized discrete-time Lyapunov equation solver for the improper controllability Gramian, see ml_gdlyap_smith_fac
infoGDLYAP_O
structure, containing information about the generalized discrete-time Lyapunov equation solver for the improper observability Gramian, see ml_gdlyap_smith_fac
infoLYAPDL
structure, containing information about the continuous-time dual Lyapunov equations solver, see ml_lyapdl_sgn_ldl or ml_lyapdl_sgn_fac
Ni
Dimension of the improper part in the reduced- order model
Np
Dimension of the proper part in the reduced-order model
Nu
Dimension of the unstable part in the reduced- order model

See Also

ml_ct_dss_bt | ml_ct_ss_flbt | ml_morlabopts