MIMO_ML_CRtheta

PURPOSE ^

SYNOPSIS ^

function [CRbound, sqrtCRbound, TheCond] = MIMO_ML_CRtheta(PolyTrans, Deriv, Sel, ModelVar, data);

DESCRIPTION ^

 function [CRbound, sqrtCRbound, TheCond] = MIMO_ML_CRtheta(PolyTrans, Deriv, Sel, ModelVar, data);


   Output parameter

        CRbound        =    Cramer-Rao bound of the estimated model parameters, the estimated plant model, and the estimated noise model
                            CRbound = struct('A', [], 'vecB',[], 'AvecB', [], 'all', [])
                            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 free parameters vecB(i) and vecB(j)
                            CRbound.all = dim(Theta) x dim(Theta), where dim(Theta) = FreeParam.A + FreeParam.B + FreeParam.C + FreeParam.D
                                CRbound.all(i,j) = covariance between free parameters Theta(i) and Theta(j)
                            Notes:
                               - in s-, sqrt(s) domains the CR-bound of the normalised parameters is calculated
                               - the (normalised) model parameters satisfy the following constraints
                                   z-domain:               a(0)      = 1
                                   s-, sqrt(s)-domains:    a(OrderA) = 1

       sqrtCRbound =   square root of Cramer-Rao lower bound on all the model parameters to guarantee a numerical stable calculation of
                           dim(Theta) x dim(Theta); and CRbound.all = sqrtCRbound * sqrtCRbound.'

       TheCond     =   condition number Jacobian matrix used to calculate the Cramer-Rao bound (cond(CRbound) = TheCond^2)


   Input parameters

        PolyTrans    =    structure containing the polynomials and transfer functions evaluated in x
                            PolyTrans.A             =    denominator polynomial plant transfer function evaluated in x.Plant 
                                                       size 1 x F 
                            PolyTrans.G             =    plant transfer matrix evaluated in x.Plant
                                                       size ny x nu x F 
                            PolyTrans.Tg            =    plant transient term evaluated in x.Plant
                                                       size ny x F 
                           PolyTrans.sqrtCEinv     =   hermitian symmetric square root of the inverse of the covariance of the 
                                                       output error (Cov(NY-G*NU)) 

        Deriv        =    structure containing the derivative of vec(G), vec(H), vec(G'), vec(H') w.r.t.
                        all the plant model parameters a, b, c, d
                            Deriv.vecGa             =    derivative vec(G) w.r.t. a;    size ny*nu x (na+1) x F 
                            Deriv.vecGb             =    derivative vec(G) w.r.t. b;    size ny*nu x ny*nu*(nb+1) x F
                            Deriv.vecGHa            =    derivative vec(G') w.r.t. a; size ny*nu x (na+1) x F
                            Deriv.vecGHb            =    derivative vec(G') w.r.t. b; size ny*nu x ny*nu*(nb+1) x F

        Sel            =    structure with fields 'A', 'B', 'Ig'
                            Sel = struct('A',[],'B',[], 'Ig', [])
                            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
                            Sel.Ig = ny x (OrderIg+1)
                                Sel.Ig(i,r) = 1 if coeff. ig(i,r-1) is unknown
                                Sel.Ig(i,r) = 0 if coeff. ig(i,r-1) = 0

        ModelVar    =    contains the information about the model to be identified structure with the following fields
                            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)
                            ModelVar.nu                =    number of inputs
                            ModelVar.ny                =     number of outputs
                            ModelVar.na                =    order polynomial A
                            ModelVar.nb                =     order matrix polynomial B

        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 nu MIMO experiments:    nu x nu x F 
                                    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 
                                                           - nu MIMO experiments:              ny x ny x nu 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 
                                                           - nu MIMO experiments:              nu x nu x nu 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 
                                                           - nu MIMO experiments:              ny x nu x nu x F 
                                   data.dof        =   degrees of freedom of the sample covariance estimates


 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, sqrtCRbound, TheCond] = MIMO_ML_CRtheta(PolyTrans, Deriv, Sel, ModelVar, data);
0002 %
0003 % function [CRbound, sqrtCRbound, TheCond] = MIMO_ML_CRtheta(PolyTrans, Deriv, Sel, ModelVar, data);
0004 %
0005 %
0006 %   Output parameter
0007 %
0008 %        CRbound        =    Cramer-Rao bound of the estimated model parameters, the estimated plant model, and the estimated noise model
0009 %                            CRbound = struct('A', [], 'vecB',[], 'AvecB', [], 'all', [])
0010 %                            CRbound.A = FreeParam.A x FreeParam.A
0011 %                                CRbound.A(i,j) = covariance between free coefficients a(i-1) and a(j-1)
0012 %                            CRbound.AvecB = FreeParam.A x FreeParam.B
0013 %                                CRbound.AvecB(i,j) = covariance between free coefficients a(i-1) and vecB(j)
0014 %                            CRbound.vecB = FreeParam.B x FreeParam.B
0015 %                                CRbound.vecB(i,j) = covariance between free parameters vecB(i) and vecB(j)
0016 %                            CRbound.all = dim(Theta) x dim(Theta), where dim(Theta) = FreeParam.A + FreeParam.B + FreeParam.C + FreeParam.D
0017 %                                CRbound.all(i,j) = covariance between free parameters Theta(i) and Theta(j)
0018 %                            Notes:
0019 %                               - in s-, sqrt(s) domains the CR-bound of the normalised parameters is calculated
0020 %                               - the (normalised) model parameters satisfy the following constraints
0021 %                                   z-domain:               a(0)      = 1
0022 %                                   s-, sqrt(s)-domains:    a(OrderA) = 1
0023 %
0024 %       sqrtCRbound =   square root of Cramer-Rao lower bound on all the model parameters to guarantee a numerical stable calculation of
0025 %                           dim(Theta) x dim(Theta); and CRbound.all = sqrtCRbound * sqrtCRbound.'
0026 %
0027 %       TheCond     =   condition number Jacobian matrix used to calculate the Cramer-Rao bound (cond(CRbound) = TheCond^2)
0028 %
0029 %
0030 %   Input parameters
0031 %
0032 %        PolyTrans    =    structure containing the polynomials and transfer functions evaluated in x
0033 %                            PolyTrans.A             =    denominator polynomial plant transfer function evaluated in x.Plant
0034 %                                                       size 1 x F
0035 %                            PolyTrans.G             =    plant transfer matrix evaluated in x.Plant
0036 %                                                       size ny x nu x F
0037 %                            PolyTrans.Tg            =    plant transient term evaluated in x.Plant
0038 %                                                       size ny x F
0039 %                           PolyTrans.sqrtCEinv     =   hermitian symmetric square root of the inverse of the covariance of the
0040 %                                                       output error (Cov(NY-G*NU))
0041 %
0042 %        Deriv        =    structure containing the derivative of vec(G), vec(H), vec(G'), vec(H') w.r.t.
0043 %                        all the plant model parameters a, b, c, d
0044 %                            Deriv.vecGa             =    derivative vec(G) w.r.t. a;    size ny*nu x (na+1) x F
0045 %                            Deriv.vecGb             =    derivative vec(G) w.r.t. b;    size ny*nu x ny*nu*(nb+1) x F
0046 %                            Deriv.vecGHa            =    derivative vec(G') w.r.t. a; size ny*nu x (na+1) x F
0047 %                            Deriv.vecGHb            =    derivative vec(G') w.r.t. b; size ny*nu x ny*nu*(nb+1) x F
0048 %
0049 %        Sel            =    structure with fields 'A', 'B', 'Ig'
0050 %                            Sel = struct('A',[],'B',[], 'Ig', [])
0051 %                            Sel.A = 1 x (OrderA+1)
0052 %                                Sel.A(r) = 1 if coeff. a(r-1) is unknown
0053 %                                Sel.A(r) = 0 if coeff. a(r-1) = 0
0054 %                            Sel.B = ny x nu x (OrderB+1)
0055 %                                Sel.B(i,j,r) = 1 if coeff. b(i,j,r-1) is unknown
0056 %                                Sel.B(i,j,r) = 0 if coeff. b(i,j,r-1) = 0
0057 %                            Sel.Ig = ny x (OrderIg+1)
0058 %                                Sel.Ig(i,r) = 1 if coeff. ig(i,r-1) is unknown
0059 %                                Sel.Ig(i,r) = 0 if coeff. ig(i,r-1) = 0
0060 %
0061 %        ModelVar    =    contains the information about the model to be identified structure with the following fields
0062 %                            ModelVar.Transient        =    1 then the initial conditions of the plant and/or noise are estimated
0063 %                            ModelVar.PlantPlane        =    plane of the plant model
0064 %                                                            's':    continuous-time
0065 %                                                            'w':    sqrt(s)-domain
0066 %                                                            'z':    discrete-time
0067 %                                                            '':     plane not defined
0068 %                            ModelVar.Struct            =    model structure
0069 %                                                           'EIV':    errors-in-variables (noisy input-output data)
0070 %                                                           'OE':    generalised output error (known input, noisy output)
0071 %                            ModelVar.RecipPlant        =    1 if plant model is reciprocal: G(i,j) = G(j,i)
0072 %                            ModelVar.nu                =    number of inputs
0073 %                            ModelVar.ny                =     number of outputs
0074 %                            ModelVar.na                =    order polynomial A
0075 %                            ModelVar.nb                =     order matrix polynomial B
0076 %
0077 %        data        =    structure containing the non-parametric data required for the identification
0078 %                                    data.U          =    For random signals one of the two following cases:
0079 %                                                           - input DFT spectrum 1 MIMO experiment:     nu x F
0080 %                                                           - power spectrum 1 MIMO experiment:         nu x nu x F
0081 %                                                       For periodic signals one of the two following cases:
0082 %                                                           - input DFT spectrum of 1 MIMO experiment:    nu x F
0083 %                                                           - input DFT specta of nu MIMO experiments:    nu x nu x F
0084 %                                    data.freq       =    vector of frequency values (Hz), size: F x 1
0085 %                                    data.Ts         =    sampling time (s)
0086 %                                    data.CY         =    (sample) noise covariance matrix of Y.
0087 %                                                       For random signals:
0088 %                                                           - covariance of 1 MIMO experiment:  ny x ny x F
0089 %                                                       For periodic signals one of the two following cases:
0090 %                                                           - 1 MIMO experiment:                ny x ny x F
0091 %                                                           - nu MIMO experiments:              ny x ny x nu x F
0092 %                                   data.CU         =   (sample) noise covariance matrix of U
0093 %                                                       For random signals:
0094 %                                                           - covariance of 1 MIMO experiment:  nu x nu x F
0095 %                                                       For periodic signals one of the two following cases:
0096 %                                                           - 1 MIMO experiment:                nu x nu x F
0097 %                                                           - nu MIMO experiments:              nu x nu x nu x F
0098 %                                   data.CYU        =   (sample) noise covariance matrix of U
0099 %                                                       For random signals:
0100 %                                                           - covariance of 1 MIMO experiment:  ny x nu x F
0101 %                                                       For periodic signals one of the two following cases:
0102 %                                                           - 1 MIMO experiment:                ny x nu x F
0103 %                                                           - nu MIMO experiments:              ny x nu x nu x F
0104 %                                   data.dof        =   degrees of freedom of the sample covariance estimates
0105 %
0106 %
0107 % Copyright (c) Rik Pintelon, Vrije Universiteit Brussel - dept. ELEC, November 2009
0108 % All rights reserved.
0109 % Software can be used freely for non-commercial applications only.
0110 % Version 24 October 2011
0111 %
0112 
0113 
0114 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0115 % initialisation of the variables %
0116 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0117 
0118 na = ModelVar.na;                               % order polynomial A
0119 nb = ModelVar.nb;                               % order matrix polynomial B
0120 nig = ModelVar.nig;                             % order vector polynomial Ig
0121 nu = ModelVar.nu;                               % number of inputs
0122 ny = ModelVar.ny;                               % number of outputs
0123 ntheta = (na+1) + (nb+1)*nu*ny + (nig+1)*ny;    % total number of model parameters
0124 F = length(data.freq);                          % number of frequencies
0125 
0126 % check the number of MIMO experiments via the size
0127 % of the noise covariance matrix
0128 sizeCov = size(data.CY);
0129 if length(sizeCov) == 4
0130     NumberExp = size(data.CY, 3);               % nexp MIMO experiments
0131 else
0132     NumberExp = 1;                              % 1 MIMO experiment
0133 end % if
0134 
0135 % intermediate data structure for calculating the CR-bound
0136 dataee = struct('CY', [], 'CU', [], 'CYU', []);
0137 dataee.CY = zeros(ny, ny, F);
0138 switch ModelVar.Struct
0139     case 'EIV'
0140         dataee.CU = zeros(nu, nu, F);           % input covariance of 1 MIMO experiment
0141         dataee.CYU = zeros(ny, nu, F);          % output-input covariance of 1 MIMO experiment
0142     case 'OE'
0143         dataee.sqrtCYinv = zeros(ny, ny, F);
0144 end % switch
0145 
0146 CRbound = struct('A', [], 'AvecB', [], 'vecB', [], 'all', []);
0147 % free parameters in each (matrix or vector) polynomial
0148 FreeParam.A = sum(Sel.A);
0149 FreeParam.B = sum(sum(sum(Sel.B)));
0150 FreeParam.Theta = FreeParam.A + FreeParam.B;
0151 CRbound.A = zeros(FreeParam.A, FreeParam.A);
0152 CRbound.AvecB = zeros(FreeParam.A, FreeParam.B);
0153 CRbound.vecB = zeros(FreeParam.B, FreeParam.B);
0154 CRbound.all = zeros(FreeParam.Theta, FreeParam.Theta);
0155 
0156 % input (power) spectrum
0157 DimU = length(size(data.U));    % dimension matrix data.U
0158 switch DimU
0159     case 2,
0160         LL = 1;                 % DFT spectrum nu x 1 reference signal is given, size nu x 1 x F
0161     case 3,
0162         if (NumberExp == 1) && (length(sizeCov) == 3)
0163             LL = nu;            % power spectrum nu x 1 reference signal is given, size nu x nu x F
0164         else % more than 1 MIMO experiment
0165             LL = 1;             % 1 MIMO experiment at the time should be handled
0166         end % if 1 MIMO experiment
0167 end
0168 
0169 % calculates the square root input power spectrum for random inputs only
0170 % data.U = nu x nu x F matrix and data.CY = ny x ny x F matrix (1 experiment with a random signal)
0171 if (DimU == 3) && (length(sizeCov) == 3)
0172     % calculate a hermitian symmetric square root of data.U
0173     UR = zeros(nu, nu, F);
0174     for kk = 1:F
0175         [uR, sR, vR] = svd(data.U(:,:,kk),0);           % data.U is a hermitian symmetric positive definite matrix
0176         SqrtSR = uR * diag(diag(sR).^(0.5)) * uR';      % hermitian symmetric square root
0177         UR(:, :, kk) = SqrtSR;
0178     end % kk, frequencies
0179 else % periodic signal
0180     UR = zeros(nu, F);
0181 end % if input power spectrum is provided
0182 
0183 % Jacobian matrix for all MIMO experiments
0184 % If LL = nu then NumberExp = 1
0185 TheJacob = zeros(LL*NumberExp*ny*F, ntheta);
0186 
0187 
0188 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0189 % Jacobian matrix w.r.t. all model parameters and all MIMO experiments %
0190 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0191 
0192 for ee = 1:NumberExp
0193     
0194     % data of MIMO experiment no. ee
0195     if length(sizeCov) == 4 % periodic excitations only
0196         UR(:,:) = data.U(:,ee,:);
0197         dataee.CY(:,:,:) = data.CY(:,:,ee,:);
0198         if strcmp(ModelVar.Struct, 'EIV')
0199             dataee.CU(:,:,:) = data.CU(:,:,ee,:);
0200             dataee.CYU(:,:,:) = data.CYU(:,:,ee,:);              
0201         end % if errors-in-variables
0202     else % 1 MIMO experiment; random or periodic excitations
0203         if LL == 1
0204             UR(:,:) = data.U;
0205         end % if not the input power spectrum is provided as input data
0206         dataee.CY(:,:,:) = data.CY;
0207          if strcmp(ModelVar.Struct, 'EIV')
0208             dataee.CU(:,:,:) = data.CU;
0209             dataee.CYU(:,:,:) = data.CYU;  
0210          end % if errors-in-variables
0211     end % if more than 1 MIMO experiment
0212 
0213 
0214     %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0215     % calculate a hermitian square root of the inverse %
0216     % of the covariance of the output error (NY-G*NU)  %
0217     %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0218 
0219     switch ModelVar.Struct
0220 
0221         case 'EIV'
0222             PolyTrans = MIMO_ML_InvCovOutputError(dataee, PolyTrans);
0223 
0224         case 'OE'
0225             % relative upper limit under which the singular values are set to zero
0226             % in the pseudo-inverse of the square root of the weighting matrix
0227             delta = 1e-6;
0228             % pseudo-inverse square root covariance matrix
0229             PolyTrans.sqrtCEinv = Sqrt_Inv(dataee.CY, delta);
0230 
0231     end % switch
0232 
0233     % inverse of a hermitian symmetric square root
0234     % of the covariance of the output error (NY-G*NU)
0235     sqrtCYinv = PolyTrans.sqrtCEinv;
0236 
0237 
0238     %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0239     % Jacobian vec(W) w.r.t. all plant model parameters where          %
0240     % W = CY^(-1/2) * G  * U                                           %
0241     % size W:                                                          %
0242     %              ny x F for LL = 1                                   %
0243     %              ny x nu x F for LL = nu                             %
0244     % JacW = Jacobian vec(W) w.r.t. Theta, dimension ny x ntheta x F   %
0245     %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0246 
0247     JacW = zeros(LL*ny, ntheta, F);     
0248 
0249     for kk = 1:F
0250 
0251         % derivatives w.r.t. a
0252         Low = 1;
0253         Upp = (na+1);
0254         switch DimU
0255             case 2,
0256                 MatW = kron(UR(:, kk).', sqrtCYinv(:,:,kk));
0257             case 3,
0258                 if (NumberExp == 1) && (length(sizeCov) == 3) % random input, input power spectrum provided
0259                     MatW = kron(UR(:,:,kk).', sqrtCYinv(:,:,kk)); 
0260                 else % periodic input, nu MIMO experiments
0261                     MatW = kron(UR(:,kk).', sqrtCYinv(:,:,kk));                    
0262                 end % if 1 MIMO experiment
0263         end
0264         JacW(:, Low:Upp, kk) = MatW * Deriv.vecGa(:,:,kk);
0265 
0266         % derivatives w.r.t. b
0267         Low = Upp + 1;
0268         Upp = Low + (nb+1)*ny*nu - 1;
0269         JacW(:, Low:Upp, kk) = MatW * Deriv.vecGb(:,:,kk);
0270 
0271     end % frequencies kk
0272 
0273     % DC and Nyquist count for 1/2 in the sums => as 1/sqrt(2) in the Jacobian matrices
0274     if data.DC == 1
0275         JacW(:,:,1) = JacW(:,:,1)/sqrt(2);
0276     end % if DC
0277     if data.Nyquist == 1
0278         JacW(:,:,end) = JacW(:,:,end)/sqrt(2);
0279     end % if Nyquist
0280 
0281 
0282     %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0283     % 1. Put the Jacobian in the following size:                           %
0284     %       JacW: LL*ny*F x ntheta                                         %
0285     % 2. impose common parameter structure and eliminate excess parameters %
0286     % 3. put real and imaginary parts on top of each other                 %
0287     %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0288 
0289     % put all frequency contributions on top of each other
0290     JacW = reshape(permute(JacW, [1,3,2]), LL*ny*F, ntheta);
0291     
0292     % put experiment no. ee in the Jacobian matrix and prediction error vector
0293     SelectRows = [(ee-1)*LL*ny*F+1:ee*LL*ny*F];
0294     TheJacob(SelectRows, :) = JacW;
0295     
0296 end % ee, MIMO experiments
0297 
0298 % impose common parameter structure and eliminate excess parameters
0299 TheJacob = MIMO_ML_AddSelectColumns(TheJacob, Sel, ModelVar);
0300 
0301 % put real and imaginary parts on top of each other
0302 % factor sqrt(2) needed because the noise cov. of the complex values is used
0303 TheJacob = sqrt(2)*[real(TheJacob); imag(TheJacob)]; 
0304 
0305 
0306 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0307 % CR-bound free parameters Theta %
0308 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0309 
0310 [uj, sj, vj] = svd(TheJacob, 0);
0311 
0312 % final CR-bound estimated model parameters
0313 vj = vj * diag(1./diag(sj));
0314 sqrtCRbound = vj;
0315 CRbound.all = vj * vj.';
0316 TheCond = sj(1,1)/sj(end,end);
0317 
0318 % scaling CRbound for the variability of the estimated noise covariances
0319 if isinf(data.dof)
0320     ScaleCRML = 1;                              % exactly known covariances
0321 else
0322     dof = data.dof;
0323     ScaleCRML = (dof)^2/(dof+1-ny)/(dof-ny-1);  % estimated covariances
0324 end % if
0325 CRbound.all = ScaleCRML * CRbound.all;
0326 sqrtCRbound = sqrt(ScaleCRML) * sqrtCRbound;
0327 
0328 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0329 % CR-bound free parameters (matrix or vector) polynomials %
0330 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0331 
0332 % free parameters A-polynomial
0333 Low = 1;
0334 Upp = Low + FreeParam.A - 1;
0335 CRbound.A = CRbound.all(Low:Upp, Low:Upp);       % CR-bound
0336 LowA = Low; UppA = Upp;
0337 
0338 % free parameters B-polynomial matrix
0339 Low = Upp + 1;
0340 Upp = Low + FreeParam.B - 1;
0341 CRbound.vecB = CRbound.all(Low:Upp, Low:Upp);         % CR-bound
0342 
0343 % covariance between free A- and B-coefficients
0344 CRbound.AvecB = CRbound.all(LowA:UppA, Low:Upp);      % CR-bound
0345

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