Contents
MORLAB Demo: Model Reduction Default Usage
This demo script contains the general application of model reduction methods for dynamical systems in MORLAB.
See also ml_ct_brbt, ml_ct_bst, ml_ct_bt, ml_ct_flbt, ml_ct_hinfbt, ml_ct_tlbt, ml_ct_hna, ml_ct_krylov, ml_ct_lqgbt, ml_ct_mt, ml_ct_prbt ml_ct_tlbt, ml_ct_twostep_mor, ml_dt_bt, ml_dt_krylov, ml_dt_lqgbt, ml_dt_mt, ml_dt_twostep_mor.
% % This file is part of the MORLAB toolbox % (https://www.mpi-magdeburg.mpg.de/projects/morlab). % Copyright (C) 2006-2023 Peter Benner, Jens Saak, and Steffen W. R. Werner % All rights reserved. % License: BSD 2-Clause License (see COPYING) %
General Handling
The MORLAB toolbox implements model reduction methods for continuous-time and discrete-time systems in standard form
xr'(t) = Ar*xr(t) + Br*u(t), yr(t) = Cr*xr(t) + Dr*u(t),
in descriptor form
E*x'(t) = A*x(t) + B*u(t), y(t) = C*x(t) + D*u(t),
or in second-order form
M*x''(t) = -K*x(t) - E*x'(t) + Bu*u(t),
y(t) = Cp*x(t) + Cv*x'(t) + D*u(t).In practice, those systems are described by a large number of differential(-algebraic) equations, which makes the usage in optimization and controller design difficult. The goal in model order reduction is the construction of a suitable approximation with less differential(-algebraic) equations than in the original system. The new system has to approximate the input-output behavior of the original system, i.e., the same inputs u(t) the outputs are close in a certain sense |y(t) - yr(t)| < tol.
To get started with the model reduction in MORLAB, we load a prepared data file containing a stable standard system.
if exist('OCTAVE_VERSION', 'builtin') orig_warn = warning('off', 'Octave:data-file-in-path'); load morlab_data_std_br.mat; warning(orig_warn); else load morlab_data_std_br.mat; end
For the application of a model reduction methods, we just need to give the system to the model reduction function we want to use. This can be done in different ways. First, MORLAB supports the single matrices given to the function. As example we want to apply the classical modal truncation to the loaded matrices.
[Ar, Br, Cr, Dr, info] = ml_ct_mt(A, B, C, D);
The output matrices Ar, Br, Cr, Dr are now the system matrices of the approximation. The info struct contains general information about the underlying used methods.
disp(info);
infoSIGNM: [1×1 struct]
infoSYLV: [1×1 struct]
N: 1
SystemType: 'continuous-time dense standard state-space system'
V: []
W: []
A more advanced way of using the model reduction function is to construct first the system as a struct. This is recommended, since most of the evaluation and analysis tools in MORLAB assume the system to be in a struct.
sys = struct( ... 'A', A, ... 'B', B, ... 'C', C, ... 'D', D);
Now we want to apply the balanced truncation method to this system.
rom = ml_ct_bt(sys);
The reduced-order model is now given in the rom struct with the same naming of the matrices as in sys.
disp(rom);
A: -0.4978
B: [0.3681 0.3710 0.4207]
C: [2×1 double]
D: [2×3 double]
To take a closer look on the result, we can use the analysis tools of MORLAB, e.g., let's compute the transfer functions of the original system and the reduction in the frequency domain by the sigma plot function.
ml_sigmaplot(sys, rom);
Or we can compare the two systems in a time simulation.
ml_ct_ss_simulate_ss11(sys, rom);
All model reduction methods allow customization by optional parameters. See also the morlab_demo_morlabopts for a get started guide with the option contructor. Here we will just call the function for simplicity and adjust the parameters as we want to.
opts = ml_morlabopts('ml_ct_bt');
The top-level model reduction routines in MORLAB automatically determine the structure of the system and use the appropriate specialized method. For demonstration, we load another system data but this time a descriptor system with E matrix.
if exist('OCTAVE_VERSION', 'builtin') orig_warn = warning('off', 'Octave:data-file-in-path'); sys = load('morlab_data_desc_br.mat'); warning(orig_warn); else sys = load('morlab_data_desc_br.mat'); end disp(sys);
A: [100×100 double]
E: [100×100 double]
B: [100×3 double]
C: [2×100 double]
D: [2×3 double]
datainfo: [1×1 struct]
If we apply now a model reduction method to the new system, for example, the Hankel-norm approximation, the reduced-order system will also be a descriptor system.
[rom, info] = ml_ct_hna(sys); disp(rom);
A: [4×4 double]
B: [4×3 double]
C: [2×4 double]
D: [2×3 double]
E: [4×4 double]
Remarks
- All model reduction routines are more or less equivalent to each other, which means that all usage examples here also apply to the other routines.
- The shown functionality also applies to the case of sparse dynamical systems. Note that the M-M.E.S.S. toolbox needs to be installed for sparse model reduction.
- This demo script only contains continuous-time examples for demonstration but MORLAB also contains discrete-time model reduction methods, e.g., ml_dt_bt is the discrete-time balanced truncation.