ml_dt_dss_bt

Balanced truncation for discrete-time descriptor systems.

Contents

Syntax

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

Description

This function computes the generalized balanced truncation for a discrete-time descriptor system of the form

   E*x(t+1) = 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 discrete-time Lyapunov equations

   Ap*Pp*Ap' - Ep*Pp*Ep' + Bp*Bp' = 0,
   Ap'*Qp*Ap - Ep'*Qp*Ep + Cp'*Cp = 0

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+1) = Ar*x(t) + Br*u(t),                                  (3)
        y(t) = Cr*x(t) + Dr*u(t)                                   (4)

is computed, such that for the original transfer function G and the reduced-order transfer function Gr with an r-th order strictly proper part it holds

   ||G - Gr||_{\infty} <= 2 * (Hsvp(r+1) + ... + Hsvp(n)),

with Hsvp, a vector containing the proper 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

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

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
dlyapdlopts
structure, containing the optional parameters for the computation of the discrete-time Lyapunov equations, see ml_dlyapdl_smith_fac
default: struct()
gdlyapdlopts
structure, containing the optional parameters for the computation of the generalized discrete-time Lyapunov equations, see ml_gdlyapdl_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()
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'
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 absolute error bound
default: 'tolerance'
stabdecopts
structure, containing the optional parameters for the decomposition of the stable and unstable parts of the system using the sign function and subspace extraction method, see ml_signm and ml_getqz
default: struct()
StoreProjection
{0, 1}, used to disable/enable storing of the computed projection matrices W and T
default: 0
Tolerance
{!}
nonnegative scalar, tolerance used for the computation of the size of the reduced-order model by an absolute error bound if 'tolerance' is set for OrderComputation
default: 1.0e-02

Note: Parameters marked with {!} may also be a cell array containing multiple arguments. In this case an cell array of the same size is returned with one entry computed for each input argument and the marked fields of the info struct are cells as well. When multiple arguments are given as cells, they are expected to have the same length.

Output

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

Entry
Meaning
AbsErrBound
{!}
computed error bound for the absolute error of the reduced-order model in H-infinity norm
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
infoDLYAPDL
structure, containing information about the discrete-time dual Lyapunov equations solver, see ml_dlyapdl_smith_fac
infoADTF
structure, containing information about the additive decomposition of the system into its infinite, finite stable and finite anti-stable parts, see ml_dt_dss_adtf
infoGDLYAPDL
structure, containing information about the generalized discrete-time Lyapunov equation solver for the improper Gramians, see ml_gdlyapdl_smith_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
T
{!}
projection matrix used as right state-space transformation to obtain the resulting block system, if opts.StoreProjection == 1
W
{!}
projection matrix used as left state-space transformation to obtain the resulting block system, if opts.StoreProjection == 1

Reference

V. Mehrmann, T. Stykel, Balanced truncation model reduction for large-scale systems in descriptor form, in: P. Benner, V. Mehrmann, D. C. Sorensen (Eds.), Dimension Reduction of Large-Scale Systems, Vol. 45 of Lect. Notes Comput. Sci. Eng., Springer-Verlag, Berlin/Heidelberg, Germany, 2005, pp. 83--115. https://doi.org/10.1007/3-540-27909-1_3

See Also

ml_dt_ss_bt | ml_morlabopts