MIMO_ML_CRbound

PURPOSE ^

SYNOPSIS ^

function [CRbound, G, Theta, CovThetan, Thetan, Seln, wscale, TheCond] = MIMO_ML_CRbound(data, Sel, Theta, ModelVar);

DESCRIPTION ^

 [CRbound, G, Theta, CovThetan, Thetan, Seln, wscale, TheCond] = MIMO_ML_CRbound(data, Sel, Theta, ModelVar);

   Calculates the CR-bound of the physical plant model parameters
   (asymptotically the influence of the plant and noise transients is zero)
                    
                Stochastic framework: errors-in-variables
                    Y = B/A U0 + NY
                    U = U0 + NU    

                System with nu inputs and ny outputs
                    Y:                    ny x 1 observed output
                    U:                    nu x 1 observed input

                Model class: common denominator model
                    G = B/A:            ny x nu plant transfer function
                    A:                    polynomial of order OrderA
                    B:                    ny x nu matrix polynomial of order OrderB

                Coefficients polynomials in raising powers of Omega, where
                    s-domain            Omega = j*2*pi*freq
                    sqrt(s)-domain        Omega = sqrt(j*2*pi*freq)
                    z-domain            Omega = exp(-j*2*pi*freq*Ts)

       References:

                   Pintelon, R., P. Guillaume, and J. Schoukens (2007). Uncertainty calculation in (operational) modal analysis, 
                   Mechanical Systems and Signal Processing, vol. 21, no. 6, pp. 2359-2373.

                   Pintelon, R., J. Schoukens, G. Vandersteen, and K. Barb�(2010). Estimation of nonparametric noise and 
                   FRF models for multivariable systems - Part II: extensions, applications, Mechanical Systems and 
                   Signal Processing, vol. 24, no. 3, pp. 596-616.

                   Pintelon, R., G. Vandersteen, J. Schoukens, and Y. Rolain (2011). Improved (non-)parametric identification of dynamic 
                   systems excited by periodic signals - The multivariate case, Mechanical Systems and Signal Processing, vol. 25, no. 8, 
                   pp. 2892-2922.

                   Pintelon, R., and J. Schoukens (2012). System Identification: A Frequency Domain Approach, second edition, 
                   IEEE Press-Wiley, Piscataway (USA). 

    Output parameters

        CRbound                =    Cramer-Rao bound of the estimated plant model parameters 
                                CRbound = struct('A', [], 'AvecB', [], 'vecB', [], 'vecG', [], 'res', [], 'poles')
                                    CRbound.A               =   FreeParam.A x FreeParam.A
                                                                   CRbound.A(i,j)      =    covariance between free coefficients  
                                                                                           a(i-1) and a(j-1) 
                                    CRbound.AvecB           =   FreeParam.A x FreeParam.B
                                                                   CRbound.AvecB(i,j)    =   covariance between free coefficients 
                                                                                           a(i-1) and vecB(j) 
                                    CRbound.vecB            =   FreeParam.B x FreeParam.B
                                                                   CRbound.vecB(i,j)   =   covariance between vecB(i) and vecB(j) 
                                                                                           where vecB = permute(B, [3, 1, 2]); 
                                                                                           vecB = vecB(:)
                                    CRbound.vecG            =   (ny*nu) x (ny*nu) x F
                                                                   CRbound.vecG(i,j,r) =   covariance between vecG(i,r) and 
                                                                                           vecG(j,r)
                                    CRbound.G               =   ny x nu x F
                                                                   CRbound.G(i,j,r)    =   variance G(i,j,r)

                                   CRbound.res             =   struct{'all', 'sv', 'lsv', 'rsv'}
                                                                   CRbound.res.all     =   covariance matrix of (vec(Res))re; where ()re puts the real
                                                                                           and imaginary parts of the matrix on top of each other; vec()
                                                                                           stacks the columns of the matrix on top of each other; and Res
                                                                                           is a residue matrix of the ny x nu transfer function matrix G = B/A;
                                                                                           size: (2*ny*nu) x (2*ny*nu) x na
                                                                   CRbound.res.sv      =   variance singular values of the residues;
                                                                                           size: min(ny, nu) x na
                                                                   CRbound.res.lsv     =   covariance left singular vectors of the residues [real(ur); imag(ur)];
                                                                                           size: 2*ny x 2*ny x min(ny, nu) x na
                                                                   CRbound.res.rsv     =   covariance right singular vectors of the residues [real(vr); imag(vr)];
                                                                                           size: 2*nu x 2*nu x min(ny, nu) x na
                                   CRbound.poles           =   struct{'root', 'all', 'damp', 'freq', 'time'}
                                                                   CRbound.poles.root    =   cov((root)re); where ()re stacks the real and imaginary part of the root
                                                                                           on top of each other; ; size 2 x 2 x number of roots
                                                                   CRbound.poles.all    =   Cov(roots.all) = cov((roots)re); where the real and imaginary parts of
                                                                                           the vector roots are stacked on top of each other;
                                                                                           size 2*(number of roots) x 2*(number of roots)
                                                                   CRbound.poles.damp    =   variance damping complex roots; entry is NaN for real roots
                                                                   CRbound.poles.freq    =   variance frequency complex roots; entry is NaN for real roots
                                                                   CRbound.poles.time    =   variance time constant real roots; entry is NaN for complex roots
                                   Notes:
                                       - Theta = [A; vec(B)]; size: Number of free parameters x 1
                                       - vec(B) is calculated as follows: vecB = permute(B, [3, 1, 2]); vecB = vecB(:);
                                       - for discrete-time systems the damping (damp), resonance frequency (freq), and time constants 
                                         (time) are calculated via the inverse impulse invariant transformation s = log(z)/Ts 
                                       - for w = sqrt(s) systems the damping (damp), resonance frequency (freq), and time constants 
                                         (time) are calculated via the transformation s = w^2 

       G                   =    plant transfer matrix evaluated at data.freq, size ny x nu x F

       Theta               =   constrained physical model parameters, where the free input parameters are defined in Seln,
                               and where the poles, residue matrices, and singular value decomposition of the residue matrices
                               have been added
                               struct('A', [], 'B', [], 'res', [], 'poles', [], 'sv') 
                                    Theta.A         =   1 x (OrderA+1)
                                                           Theta.A(r) = coefficient a(r-1) of Omega^(r-1) 
                                    Theta.B         =   ny x nu x (OrderB+1)
                                                           Theta.B(i,j,r) = coefficient b(i,j,r-1) of Omega^(r-1)
                                   Theta.res       =   struct{'all', 'sv', 'lsv', 'rsv'}
                                                           Theta.res.all       =   residue matrices of the transfer function matrix G = B/A;
                                                                                   size: ny x nu x na 
                                                           Theta.res.sv        =   singular values of the residue matrices;
                                                                                   size: min(ny, nu) x na
                                                           Theta.res.lsv       =   left singular vectors of the residue matrices;
                                                                                   size: ny x min(ny, nu) x na
                                                           Theta.res.rsv       =   right singular vectors of the residue matrices;
                                                                                   size: nu x min(ny, nu) x na
                                                           Notes:
                                                                 - the residues are sorted in the same order as the poles 
                                                                 - the largest element of the right singular vector is made positive real
                                                                   by multiplying the left and right singular vectors by exp(-j*fi) where
                                                                   fi is the phase of the largest element of the right singular vector
                                   Theta.poles     =   struct{'root', 'damp', 'freq', 'time}
                                                           Theta.poles.root    =   poles of the transfer function matrix G = B/A;
                                                                                   size na x 1
                                                           Theta.poles.damp    =   damping ratio of the poles;
                                                                                   size na x 1
                                                           Theta.poles.freq    =   resonance frequency in Hz of the poles;
                                                                                   size na x 1
                                                           Theta.poles.time    =   time constant in seconds of the real poles;
                                                                                   size na x 1
                                   Notes:
                                       - for discrete-time systems the damping (damp), resonance frequency (freq), and time constants 
                                         (time) are calculated via the inverse impulse invariant transformation s = log(z)/Ts 
                                       - for w = sqrt(s) systems the damping (damp), resonance frequency (freq), and time constants 
                                         (time) are calculated via the transformation s = w^2 

        CovThetan           =    Cramer-Rao bound of the normalised estimated model parameters, the estimated plant model, and the estimated noise model
                                CovThetan = struct('A', [], 'AvecB', [], 'vecB', [], 'all', [], 'poles', [], 'res', [])
                               See CRbound for the defintion of the structure. 
                               Should be used for numerical stable calculation of uncertainty bounds; to guarantee the
                               positive semi-definiteness of X.' * CovThetan * X, one should calculate
                                                 (sqrtCovThetan * X).' * (sqrtCovThetan * X) 
                               where sqrtCovThetan is a square root of CovThetan calculated via an SVD. 

       Seln                =   structure containing the information about the estimated parameters of Thetan; 
                               see Sel in the input parameters

        Thetan                =    (normalised) estimated plant and noise parameters (see input parameters) were the
                                following constraints have been imposed:
                                       - in s-, sqrt(s) domains the parameters are normalised with wscale,
                                         for example, ar * s^r => (ar*wscale^r) * (s/wscale)^r
                                         where wscale.Plant and wscale.Noise are used for the plant and noise model
                                         parameters respectively
                                       - the following constraints on the (normalised) model parameters are imposed
                                           z-domain:               a(0)  =  1,
                                           s-, sqrt(s)-domains:    a(na) =  1

        wscale                =    angular frequency scaling

        TheCond                =    condition number Jacobian matrix used to calculate the CR-bound; cond(CR-bound) = TheCond^2


    Input parameters

        data                =    structure containing the non-parametric data required for the identification
                                    data.U          =    For random signals one of the two following cases: 
                                                           - input DFT spectrum 1 MIMO experiment:         nu x F 
                                                           - power spectrum 1 MIMO experiment:             nu x nu x F 
                                                       For periodic signals one of the two following cases: 
                                                           - input DFT spectrum of 1 MIMO experiment:      nu x F 
                                                           - input DFT specta of nexp MIMO experiments:    nu x nexp x F 
                                                             with nexp >= 1 
                                    data.freq       =    vector of frequency values (Hz), size: F x 1
                                    data.Ts         =    sampling time (s)
                                    data.CY         =    (sample) noise covariance matrix of Y.
                                                       For random signals:
                                                           - covariance of 1 MIMO experiment:              ny x ny x F 
                                                       For periodic signals one of the two following cases:  
                                                           - 1 MIMO experiment:                            ny x ny x F 
                                                           - nexp MIMO experiments (nexp >= 1):            ny x ny x nexp x F 
                                   data.CU         =   (sample) noise covariance matrix of U  
                                                       For random signals:
                                                           - covariance of 1 MIMO experiment:  nu x nu x F 
                                                       For periodic signals one of the two following cases:  
                                                           - 1 MIMO experiment:                            nu x nu x F 
                                                           - nexp MIMO experiments (nexp >= 1):            nu x nu x nexp x F 
                                   data.CYU        =   (sample) noise covariance matrix of U 
                                                       For random signals:
                                                           - covariance of 1 MIMO experiment:  ny x nu x F 
                                                       For periodic signals one of the two following cases:  
                                                           - 1 MIMO experiment:                            ny x nu x F 
                                                           - nexp MIMO experiments (nexp >= 1):            ny x nu x nexp x F 
                                   data.dof        =   degrees of freedom of the sample covariance estimates
                                                       (optional, default value dof = inf) 

        Sel                    =    struct('A', [], 'B', [])
                                    Sel.A           =   1 x (OrderA+1)
                                                           Sel.A(r) = 1 if coeff. a(r-1) is unknown
                                                           Sel.A(r) = 0 if coeff. a(r-1) = 0
                                    Sel.B           =   ny x nu x (OrderB+1)
                                                           Sel.B(i,j,r) = 1 if coeff. b(i,j,r-1) is unknown
                                                           Sel.B(i,j,r) = 0 if coeff. b(i,j,r-1) = 0

        Theta                =    estimated plant model parameters
                                Theta = struct('A', [], 'B', [])
                                    Theta.A         =   1 x (OrderA+1)
                                                           Theta.A(r) = coefficient a(r-1) of Omega^(r-1) 
                                    Theta.B         =   ny x nu x (OrderB+1)
                                                           Theta.B(i,j,r) = coefficient b(i,j,r-1) of Omega^(r-1)

        ModelVar            =    contains the information about the model to be identified
                                struct('Transient', [], 'PlantPlane', [], 'Struct', [], 'RecipPlant', [])
                                    ModelVar.Transient        =    1 then the initial conditions of the plant and/or noise are estimated
                                    ModelVar.PlantPlane        =    plane of the plant model
                                                                    's':    continuous-time;
                                                                    'w':    sqrt(s)-domain
                                                                    'z':    discrete-time;
                                                                    '':        plane not defined
                                    ModelVar.Struct            =    model structure
                                                                   'EIV':  errors-in-variables (noisy input-output data)
                                                                   'OE':    generalised output error (known input, noisy output)
                                    ModelVar.RecipPlant        =    1 if plant model is reciprocal: G(i,j) = G(j,i)


 Copyright (c) Rik Pintelon, Vrije Universiteit Brussel - dept. ELEC, November 2009
 All rights reserved.
 Software can be used freely for non-commercial applications only.
 Version 24 October 2011

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function [CRbound, G, Theta, CovThetan, Thetan, Seln, wscale, TheCond] = MIMO_ML_CRbound(data, Sel, Theta, ModelVar);
0002 %
0003 % [CRbound, G, Theta, CovThetan, Thetan, Seln, wscale, TheCond] = MIMO_ML_CRbound(data, Sel, Theta, ModelVar);
0004 %
0005 %   Calculates the CR-bound of the physical plant model parameters
0006 %   (asymptotically the influence of the plant and noise transients is zero)
0007 %
0008 %                Stochastic framework: errors-in-variables
0009 %                    Y = B/A U0 + NY
0010 %                    U = U0 + NU
0011 %
0012 %                System with nu inputs and ny outputs
0013 %                    Y:                    ny x 1 observed output
0014 %                    U:                    nu x 1 observed input
0015 %
0016 %                Model class: common denominator model
0017 %                    G = B/A:            ny x nu plant transfer function
0018 %                    A:                    polynomial of order OrderA
0019 %                    B:                    ny x nu matrix polynomial of order OrderB
0020 %
0021 %                Coefficients polynomials in raising powers of Omega, where
0022 %                    s-domain            Omega = j*2*pi*freq
0023 %                    sqrt(s)-domain        Omega = sqrt(j*2*pi*freq)
0024 %                    z-domain            Omega = exp(-j*2*pi*freq*Ts)
0025 %
0026 %       References:
0027 %
0028 %                   Pintelon, R., P. Guillaume, and J. Schoukens (2007). Uncertainty calculation in (operational) modal analysis,
0029 %                   Mechanical Systems and Signal Processing, vol. 21, no. 6, pp. 2359-2373.
0030 %
0031 %                   Pintelon, R., J. Schoukens, G. Vandersteen, and K. Barb�(2010). Estimation of nonparametric noise and
0032 %                   FRF models for multivariable systems - Part II: extensions, applications, Mechanical Systems and
0033 %                   Signal Processing, vol. 24, no. 3, pp. 596-616.
0034 %
0035 %                   Pintelon, R., G. Vandersteen, J. Schoukens, and Y. Rolain (2011). Improved (non-)parametric identification of dynamic
0036 %                   systems excited by periodic signals - The multivariate case, Mechanical Systems and Signal Processing, vol. 25, no. 8,
0037 %                   pp. 2892-2922.
0038 %
0039 %                   Pintelon, R., and J. Schoukens (2012). System Identification: A Frequency Domain Approach, second edition,
0040 %                   IEEE Press-Wiley, Piscataway (USA).
0041 %
0042 %    Output parameters
0043 %
0044 %        CRbound                =    Cramer-Rao bound of the estimated plant model parameters
0045 %                                CRbound = struct('A', [], 'AvecB', [], 'vecB', [], 'vecG', [], 'res', [], 'poles')
0046 %                                    CRbound.A               =   FreeParam.A x FreeParam.A
0047 %                                                                   CRbound.A(i,j)      =    covariance between free coefficients
0048 %                                                                                           a(i-1) and a(j-1)
0049 %                                    CRbound.AvecB           =   FreeParam.A x FreeParam.B
0050 %                                                                   CRbound.AvecB(i,j)    =   covariance between free coefficients
0051 %                                                                                           a(i-1) and vecB(j)
0052 %                                    CRbound.vecB            =   FreeParam.B x FreeParam.B
0053 %                                                                   CRbound.vecB(i,j)   =   covariance between vecB(i) and vecB(j)
0054 %                                                                                           where vecB = permute(B, [3, 1, 2]);
0055 %                                                                                           vecB = vecB(:)
0056 %                                    CRbound.vecG            =   (ny*nu) x (ny*nu) x F
0057 %                                                                   CRbound.vecG(i,j,r) =   covariance between vecG(i,r) and
0058 %                                                                                           vecG(j,r)
0059 %                                    CRbound.G               =   ny x nu x F
0060 %                                                                   CRbound.G(i,j,r)    =   variance G(i,j,r)
0061 %
0062 %                                   CRbound.res             =   struct{'all', 'sv', 'lsv', 'rsv'}
0063 %                                                                   CRbound.res.all     =   covariance matrix of (vec(Res))re; where ()re puts the real
0064 %                                                                                           and imaginary parts of the matrix on top of each other; vec()
0065 %                                                                                           stacks the columns of the matrix on top of each other; and Res
0066 %                                                                                           is a residue matrix of the ny x nu transfer function matrix G = B/A;
0067 %                                                                                           size: (2*ny*nu) x (2*ny*nu) x na
0068 %                                                                   CRbound.res.sv      =   variance singular values of the residues;
0069 %                                                                                           size: min(ny, nu) x na
0070 %                                                                   CRbound.res.lsv     =   covariance left singular vectors of the residues [real(ur); imag(ur)];
0071 %                                                                                           size: 2*ny x 2*ny x min(ny, nu) x na
0072 %                                                                   CRbound.res.rsv     =   covariance right singular vectors of the residues [real(vr); imag(vr)];
0073 %                                                                                           size: 2*nu x 2*nu x min(ny, nu) x na
0074 %                                   CRbound.poles           =   struct{'root', 'all', 'damp', 'freq', 'time'}
0075 %                                                                   CRbound.poles.root    =   cov((root)re); where ()re stacks the real and imaginary part of the root
0076 %                                                                                           on top of each other; ; size 2 x 2 x number of roots
0077 %                                                                   CRbound.poles.all    =   Cov(roots.all) = cov((roots)re); where the real and imaginary parts of
0078 %                                                                                           the vector roots are stacked on top of each other;
0079 %                                                                                           size 2*(number of roots) x 2*(number of roots)
0080 %                                                                   CRbound.poles.damp    =   variance damping complex roots; entry is NaN for real roots
0081 %                                                                   CRbound.poles.freq    =   variance frequency complex roots; entry is NaN for real roots
0082 %                                                                   CRbound.poles.time    =   variance time constant real roots; entry is NaN for complex roots
0083 %                                   Notes:
0084 %                                       - Theta = [A; vec(B)]; size: Number of free parameters x 1
0085 %                                       - vec(B) is calculated as follows: vecB = permute(B, [3, 1, 2]); vecB = vecB(:);
0086 %                                       - for discrete-time systems the damping (damp), resonance frequency (freq), and time constants
0087 %                                         (time) are calculated via the inverse impulse invariant transformation s = log(z)/Ts
0088 %                                       - for w = sqrt(s) systems the damping (damp), resonance frequency (freq), and time constants
0089 %                                         (time) are calculated via the transformation s = w^2
0090 %
0091 %       G                   =    plant transfer matrix evaluated at data.freq, size ny x nu x F
0092 %
0093 %       Theta               =   constrained physical model parameters, where the free input parameters are defined in Seln,
0094 %                               and where the poles, residue matrices, and singular value decomposition of the residue matrices
0095 %                               have been added
0096 %                               struct('A', [], 'B', [], 'res', [], 'poles', [], 'sv')
0097 %                                    Theta.A         =   1 x (OrderA+1)
0098 %                                                           Theta.A(r) = coefficient a(r-1) of Omega^(r-1)
0099 %                                    Theta.B         =   ny x nu x (OrderB+1)
0100 %                                                           Theta.B(i,j,r) = coefficient b(i,j,r-1) of Omega^(r-1)
0101 %                                   Theta.res       =   struct{'all', 'sv', 'lsv', 'rsv'}
0102 %                                                           Theta.res.all       =   residue matrices of the transfer function matrix G = B/A;
0103 %                                                                                   size: ny x nu x na
0104 %                                                           Theta.res.sv        =   singular values of the residue matrices;
0105 %                                                                                   size: min(ny, nu) x na
0106 %                                                           Theta.res.lsv       =   left singular vectors of the residue matrices;
0107 %                                                                                   size: ny x min(ny, nu) x na
0108 %                                                           Theta.res.rsv       =   right singular vectors of the residue matrices;
0109 %                                                                                   size: nu x min(ny, nu) x na
0110 %                                                           Notes:
0111 %                                                                 - the residues are sorted in the same order as the poles
0112 %                                                                 - the largest element of the right singular vector is made positive real
0113 %                                                                   by multiplying the left and right singular vectors by exp(-j*fi) where
0114 %                                                                   fi is the phase of the largest element of the right singular vector
0115 %                                   Theta.poles     =   struct{'root', 'damp', 'freq', 'time}
0116 %                                                           Theta.poles.root    =   poles of the transfer function matrix G = B/A;
0117 %                                                                                   size na x 1
0118 %                                                           Theta.poles.damp    =   damping ratio of the poles;
0119 %                                                                                   size na x 1
0120 %                                                           Theta.poles.freq    =   resonance frequency in Hz of the poles;
0121 %                                                                                   size na x 1
0122 %                                                           Theta.poles.time    =   time constant in seconds of the real poles;
0123 %                                                                                   size na x 1
0124 %                                   Notes:
0125 %                                       - for discrete-time systems the damping (damp), resonance frequency (freq), and time constants
0126 %                                         (time) are calculated via the inverse impulse invariant transformation s = log(z)/Ts
0127 %                                       - for w = sqrt(s) systems the damping (damp), resonance frequency (freq), and time constants
0128 %                                         (time) are calculated via the transformation s = w^2
0129 %
0130 %        CovThetan           =    Cramer-Rao bound of the normalised estimated model parameters, the estimated plant model, and the estimated noise model
0131 %                                CovThetan = struct('A', [], 'AvecB', [], 'vecB', [], 'all', [], 'poles', [], 'res', [])
0132 %                               See CRbound for the defintion of the structure.
0133 %                               Should be used for numerical stable calculation of uncertainty bounds; to guarantee the
0134 %                               positive semi-definiteness of X.' * CovThetan * X, one should calculate
0135 %                                                 (sqrtCovThetan * X).' * (sqrtCovThetan * X)
0136 %                               where sqrtCovThetan is a square root of CovThetan calculated via an SVD.
0137 %
0138 %       Seln                =   structure containing the information about the estimated parameters of Thetan;
0139 %                               see Sel in the input parameters
0140 %
0141 %        Thetan                =    (normalised) estimated plant and noise parameters (see input parameters) were the
0142 %                                following constraints have been imposed:
0143 %                                       - in s-, sqrt(s) domains the parameters are normalised with wscale,
0144 %                                         for example, ar * s^r => (ar*wscale^r) * (s/wscale)^r
0145 %                                         where wscale.Plant and wscale.Noise are used for the plant and noise model
0146 %                                         parameters respectively
0147 %                                       - the following constraints on the (normalised) model parameters are imposed
0148 %                                           z-domain:               a(0)  =  1,
0149 %                                           s-, sqrt(s)-domains:    a(na) =  1
0150 %
0151 %        wscale                =    angular frequency scaling
0152 %
0153 %        TheCond                =    condition number Jacobian matrix used to calculate the CR-bound; cond(CR-bound) = TheCond^2
0154 %
0155 %
0156 %    Input parameters
0157 %
0158 %        data                =    structure containing the non-parametric data required for the identification
0159 %                                    data.U          =    For random signals one of the two following cases:
0160 %                                                           - input DFT spectrum 1 MIMO experiment:         nu x F
0161 %                                                           - power spectrum 1 MIMO experiment:             nu x nu x F
0162 %                                                       For periodic signals one of the two following cases:
0163 %                                                           - input DFT spectrum of 1 MIMO experiment:      nu x F
0164 %                                                           - input DFT specta of nexp MIMO experiments:    nu x nexp x F
0165 %                                                             with nexp >= 1
0166 %                                    data.freq       =    vector of frequency values (Hz), size: F x 1
0167 %                                    data.Ts         =    sampling time (s)
0168 %                                    data.CY         =    (sample) noise covariance matrix of Y.
0169 %                                                       For random signals:
0170 %                                                           - covariance of 1 MIMO experiment:              ny x ny x F
0171 %                                                       For periodic signals one of the two following cases:
0172 %                                                           - 1 MIMO experiment:                            ny x ny x F
0173 %                                                           - nexp MIMO experiments (nexp >= 1):            ny x ny x nexp x F
0174 %                                   data.CU         =   (sample) noise covariance matrix of U
0175 %                                                       For random signals:
0176 %                                                           - covariance of 1 MIMO experiment:  nu x nu x F
0177 %                                                       For periodic signals one of the two following cases:
0178 %                                                           - 1 MIMO experiment:                            nu x nu x F
0179 %                                                           - nexp MIMO experiments (nexp >= 1):            nu x nu x nexp x F
0180 %                                   data.CYU        =   (sample) noise covariance matrix of U
0181 %                                                       For random signals:
0182 %                                                           - covariance of 1 MIMO experiment:  ny x nu x F
0183 %                                                       For periodic signals one of the two following cases:
0184 %                                                           - 1 MIMO experiment:                            ny x nu x F
0185 %                                                           - nexp MIMO experiments (nexp >= 1):            ny x nu x nexp x F
0186 %                                   data.dof        =   degrees of freedom of the sample covariance estimates
0187 %                                                       (optional, default value dof = inf)
0188 %
0189 %        Sel                    =    struct('A', [], 'B', [])
0190 %                                    Sel.A           =   1 x (OrderA+1)
0191 %                                                           Sel.A(r) = 1 if coeff. a(r-1) is unknown
0192 %                                                           Sel.A(r) = 0 if coeff. a(r-1) = 0
0193 %                                    Sel.B           =   ny x nu x (OrderB+1)
0194 %                                                           Sel.B(i,j,r) = 1 if coeff. b(i,j,r-1) is unknown
0195 %                                                           Sel.B(i,j,r) = 0 if coeff. b(i,j,r-1) = 0
0196 %
0197 %        Theta                =    estimated plant model parameters
0198 %                                Theta = struct('A', [], 'B', [])
0199 %                                    Theta.A         =   1 x (OrderA+1)
0200 %                                                           Theta.A(r) = coefficient a(r-1) of Omega^(r-1)
0201 %                                    Theta.B         =   ny x nu x (OrderB+1)
0202 %                                                           Theta.B(i,j,r) = coefficient b(i,j,r-1) of Omega^(r-1)
0203 %
0204 %        ModelVar            =    contains the information about the model to be identified
0205 %                                struct('Transient', [], 'PlantPlane', [], 'Struct', [], 'RecipPlant', [])
0206 %                                    ModelVar.Transient        =    1 then the initial conditions of the plant and/or noise are estimated
0207 %                                    ModelVar.PlantPlane        =    plane of the plant model
0208 %                                                                    's':    continuous-time;
0209 %                                                                    'w':    sqrt(s)-domain
0210 %                                                                    'z':    discrete-time;
0211 %                                                                    '':        plane not defined
0212 %                                    ModelVar.Struct            =    model structure
0213 %                                                                   'EIV':  errors-in-variables (noisy input-output data)
0214 %                                                                   'OE':    generalised output error (known input, noisy output)
0215 %                                    ModelVar.RecipPlant        =    1 if plant model is reciprocal: G(i,j) = G(j,i)
0216 %
0217 %
0218 % Copyright (c) Rik Pintelon, Vrije Universiteit Brussel - dept. ELEC, November 2009
0219 % All rights reserved.
0220 % Software can be used freely for non-commercial applications only.
0221 % Version 24 October 2011
0222 %
0223 
0224 
0225 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0226 % initialisation of the variables, and compatibility check of the input arguments %
0227 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0228 
0229 % degrees of freedom of the noise covariance estimates
0230 try
0231     if isempty(data.dof)
0232         data.dof = inf;         % known noise covariances
0233     end
0234 catch
0235     data.dof = inf;             % known noise covariances
0236 end % try
0237 
0238 % model structure variables
0239 ModelVar.PlantPlane = lower(ModelVar.PlantPlane);
0240 ModelVar.Struct = upper(ModelVar.Struct);
0241 
0242 ModelVar.ny = size(data.Y,1);
0243 ModelVar.nu = size(data.U,1);
0244 
0245 % 1. imposes the compatibility of the parameter vector Theta and
0246 %    the free model parameters with the model structure
0247 % 2. puts the order of the polynomials in ModelVar
0248 [Theta, Sel, ModelVar] = MIMO_ML_ModelCompatibility(Theta, Sel, ModelVar);
0249 
0250 % Check if DC and Nyquist belong to the frequency set
0251 if data.freq(1) == 0; data.DC = 1; else data.DC = 0; end
0252 if data.freq(end) == 1/(2*data.Ts); data.Nyquist = 1; else data.Nyquist = 0; end
0253 data.freq = data.freq(:);
0254 
0255 % jw, sqrt(jw), or exp(-jwTs) values for plant and noise model
0256 x = struct('Plant', []);
0257 % matrix of powers of x
0258 xMat = struct('Plant', []);
0259 
0260 CovThetan = struct('A', [], 'vecB', [], 'AvecB', [], 'all', [], 'poles', [], 'res', []);
0261 
0262 
0263 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0264 % Domain of the plant model %
0265 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0266 
0267 % the vector s represents z^-1, s, or sqrt(s) of the plant model
0268 switch ModelVar.PlantPlane
0269     case {'s','w'}
0270         if ModelVar.PlantPlane == 's'
0271             x.Plant = sqrt(-1)*2*pi*data.freq;
0272         elseif ModelVar.PlantPlane == 'w'
0273             x.Plant = (sqrt(-1)*2*pi*data.freq).^(0.5);
0274         end;
0275         wscale = median(abs(x.Plant));
0276         x.Plant = x.Plant/wscale;
0277     case 'z'
0278         x.Plant = exp(-sqrt(-1)*2*pi*data.freq*data.Ts);
0279         wscale = 1;
0280     case ''
0281         x.Plant = ones(size(data.freq));
0282         wscale = 1;
0283     otherwise, disp('Invalid plant plane ...'), return
0284 end
0285 
0286 nmax = max([ModelVar.na, ModelVar.nb, ModelVar.nig]);
0287 xMat.Plant = MIMO_ML_CalcOmegaMat(x.Plant, nmax);
0288 
0289 
0290 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0291 % Normalise the model parameters for s-, and sqrt(s)-domains %
0292 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0293 
0294 Thetan = MIMO_ML_Normalise(Theta, wscale, ModelVar);
0295 
0296 
0297 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0298 % Impose constraints on the (normalised) model parameters      %
0299 %   z-domain:           a0 = 1                                 %
0300 %   s-, sqrt(s)-domain: ana = 1                                %
0301 % and modify the free model parameter selection accordingly    %
0302 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0303 
0304 % impose the constraints on the model parameters
0305 [Thetan, Seln] = MIMO_ML_MonicModel(Thetan, Sel, ModelVar);
0306 % store the normalised plant transient term because later on it is deleted
0307 SaveIg = Thetan.Ig;   
0308 
0309 
0310 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0311 % Cramer-Rao lower bound of the free model parameters %
0312 %   CRBound.A                                         %
0313 %   CRbound.vecB                                      %
0314 %   CRbound.AvecB                                     %
0315 %   CRbound.all                                       %
0316 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0317 
0318 % calculation polynomials and transfer functions
0319 PolyTrans = MIMO_ML_CalcPolyTrans(Thetan, x);
0320 G = PolyTrans.G;
0321 
0322 % Since they are asymptotically zero, the Cramer-Rao lower bound is not
0323 % calculated for the plant and noise transient terms
0324 ModelVar.Transient = 0;
0325 ModelVar.nig = 0;
0326 Seln.Ig = zeros(ModelVar.ny, 1);
0327 Thetan.Ig = zeros(ModelVar.ny, 1);
0328 PolyTrans.Tg = zeros(size(PolyTrans.Tg));
0329 
0330 % calculate the derivative of the transfer functions w.r.t. the model parameters
0331 Deriv = MIMO_ML_CalcDeriv(xMat, PolyTrans, ModelVar);
0332 
0333 % derivatives of the hermitian transpose of the plant transfer function w.r.t. the model parameters
0334 Deriv.vecGHa = MIMO_ML_DvecZ2DvecZH(Deriv.vecGa, ModelVar.ny, ModelVar.nu);
0335 Deriv.vecGHb = MIMO_ML_DvecZ2DvecZH(Deriv.vecGb, ModelVar.ny, ModelVar.nu);
0336 
0337 [CRbound, sqrtCRbound, TheCond] = MIMO_ML_CRtheta(PolyTrans, Deriv, Seln, ModelVar, data);
0338 
0339 
0340 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0341 % Cramer-Rao lower bound of the poles and residue matrices %
0342 %   CRbound.poles                                          %
0343 %   CRbound.res                                            %
0344 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0345 
0346 % calculate poles and residue matrices and their covariances
0347 [CovParam, Param] = CovResidues(Thetan, CRbound, Seln, ModelVar.PlantPlane, data.Ts);
0348 
0349 % add poles and residue matrices to Thetan
0350 Thetan.poles = Param.poles;
0351 Thetan.res = [];
0352 Thetan.res.all = Param.res;
0353 Thetan.res.sv = Param.sv;
0354 Thetan.res.lsv = Param.lsv;
0355 Thetan.res.rsv = Param.rsv;
0356 
0357 % add poles and residue matrices to CRbound
0358 CRbound.poles = CovParam.poles;
0359 CRbound.res = CovParam.res;
0360 
0361 
0362 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0363 % Cramer-Rao lower bound of the plant transfer function matrix %
0364 %   CRBound.vecG                                               %
0365 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0366 
0367 CRbound = MIMO_ML_CRtf(CRbound, sqrtCRbound, PolyTrans, Deriv, Seln, ModelVar);
0368 
0369 
0370 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0371 % CR-bound of the estimated normalised parameters %
0372 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0373 
0374 CovThetan.A = CRbound.A;
0375 CovThetan.vecB = CRbound.vecB;
0376 CovThetan.AvecB = CRbound.AvecB;
0377 CovThetan.all = CRbound.all;
0378 CovThetan.poles = CRbound.poles;
0379 CovThetan.res = CRbound.res;
0380 
0381 
0382 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0383 % CR-bound of the estimated physical model parameters %
0384 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0385 
0386 % frequency denormalisation CR-bounds plant model parameters
0387 CRbound = MIMO_ML_CRboundPhysicalParam(CRbound, Seln, wscale, ModelVar);
0388 
0389 
0390 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0391 % Restore the estimated plant transient parameters %
0392 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0393 
0394 Thetan.Ig = SaveIg;
0395 Seln.Ig = Sel.Ig;
0396 
0397 
0398 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0399 % Calculation of the constrained physical model parameters %
0400 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0401 
0402 % Thetan: constrained frequency normalised parameters A, B, Ig, res, poles
0403 % Theta: constrained physical parameters A, B, Ig
0404 Theta = MIMO_ML_DeNormalise(Thetan, wscale, ModelVar);
0405 
0406 if ~strcmp(ModelVar.PlantPlane, 'z')
0407     
0408     switch ModelVar.PlantPlane 
0409         case 's'
0410             wscale2 = wscale;
0411         case 'w'
0412             % poles are squared in sqrt(s)-domain before calculating
0413             % the resonance frequencies and damping ratios
0414             wscale2 = wscale^2;   
0415     end % if
0416     
0417     % denormalisation of the poles, frequency and time constants,
0418     % but not the damping because it has no units
0419     Theta.poles = Thetan.poles;
0420     Theta.poles.root = Theta.poles.root * wscale;
0421     Theta.poles.freq = Theta.poles.freq * wscale2;
0422     Theta.poles.time = Theta.poles.time / wscale2;
0423 
0424     % denormalisation of the residue matrices and the singular values,
0425     % but not the left and right singular vectors because they have no units
0426     Theta.res = Thetan.res;
0427     Theta.res.all = Theta.res.all * wscale;
0428     Theta.res.sv = Theta.res.sv * wscale;
0429 
0430 end % if not z-domain
0431 
0432

Generated on Thu 07-Jun-2012 11:58:58 by m2html © 2005