MIMO_ML_Jacob

PURPOSE ^

SYNOPSIS ^

function TheJacob = MIMO_ML_Jacob(data, xMat, Error, PolyTrans, ModelVar);

DESCRIPTION ^

      TheJacob = MIMO_ML_Jacob(data, xMat, Error, PolyTrans, ModelVar)

   Calculates the jacobian matrix w.r.t. ALL the plant and transient model
   parameters.
    The selection of the free model parameters, and the imposition of the model constraints
    (norm one, reciprocity) is done in the routine MIMO_ML_Add_SelectColumns


 Output parameters

     TheJacob         =    jacobian matrix, size: ny x ntheta x F 

 Input parameters

        data        =    structure containing the non-parametric data
                            data.Y                  =    DFT spectrum ny x 1 output signal, size: ny x F 
                            data.U                  =    DFT spectrum nu x 1 input signal, size: 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, size: ny x ny x F 
                           data.CU                 =   (sample) noise covariance matrix of U, size: nu x nu x F 
                           data.CYU                =   (sample) noise covariance matrix of U, size: ny x nu x F 
                            data.sqrtCYinv          =    CY^(-0.5), size: ny x ny x F 
                            data.DC                 =    1 if DC present otherwise 0
                            data.Nyquist            =    1 if Nyquist frequency present otherwise 0

        xMat        =    structure with tables of powers of (jwk)^r or (zk^-r)
                            xMat.Plant              =    plant model, size: F x max order

        Error       =    prediction error model equations, size ny x F 

        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)) 

        ModelVar    =    contains the information about the model to be identified
                        structure with fields 'Transient', 'ThePlane', 'TheModel', 'Reciprocal', ...
                            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
                            ModelVar.nig            =    order vector polynomial Ig


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

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function TheJacob = MIMO_ML_Jacob(data, xMat, Error, PolyTrans, ModelVar);
0002 %
0003 %      TheJacob = MIMO_ML_Jacob(data, xMat, Error, PolyTrans, ModelVar)
0004 %
0005 %   Calculates the jacobian matrix w.r.t. ALL the plant and transient model
0006 %   parameters.
0007 %    The selection of the free model parameters, and the imposition of the model constraints
0008 %    (norm one, reciprocity) is done in the routine MIMO_ML_Add_SelectColumns
0009 %
0010 %
0011 % Output parameters
0012 %
0013 %     TheJacob         =    jacobian matrix, size: ny x ntheta x F
0014 %
0015 % Input parameters
0016 %
0017 %        data        =    structure containing the non-parametric data
0018 %                            data.Y                  =    DFT spectrum ny x 1 output signal, size: ny x F
0019 %                            data.U                  =    DFT spectrum nu x 1 input signal, size: nu x F
0020 %                            data.freq               =    vector of frequency values (Hz), size: F x 1
0021 %                            data.Ts                 =    sampling time (s)
0022 %                            data.CY                 =    (sample) noise covariance matrix of Y, size: ny x ny x F
0023 %                           data.CU                 =   (sample) noise covariance matrix of U, size: nu x nu x F
0024 %                           data.CYU                =   (sample) noise covariance matrix of U, size: ny x nu x F
0025 %                            data.sqrtCYinv          =    CY^(-0.5), size: ny x ny x F
0026 %                            data.DC                 =    1 if DC present otherwise 0
0027 %                            data.Nyquist            =    1 if Nyquist frequency present otherwise 0
0028 %
0029 %        xMat        =    structure with tables of powers of (jwk)^r or (zk^-r)
0030 %                            xMat.Plant              =    plant model, size: F x max order
0031 %
0032 %        Error       =    prediction error model equations, size ny x F
0033 %
0034 %        PolyTrans    =    structure containing the polynomials and transfer functions evaluated in x
0035 %                            PolyTrans.A             =    denominator polynomial plant transfer function evaluated in x.Plant
0036 %                                                       size 1 x F
0037 %                            PolyTrans.G             =    plant transfer matrix evaluated in x.Plant
0038 %                                                       size ny x nu x F
0039 %                            PolyTrans.Tg            =    plant transient term evaluated in x.Plant
0040 %                                                       size ny x F
0041 %                           PolyTrans.sqrtCEinv     =   hermitian symmetric square root of the inverse of the covariance of the
0042 %                                                       output error (Cov(NY-G*NU))
0043 %
0044 %        ModelVar    =    contains the information about the model to be identified
0045 %                        structure with fields 'Transient', 'ThePlane', 'TheModel', 'Reciprocal', ...
0046 %                            ModelVar.Transient        =    1 then the initial conditions of the plant and/or noise are estimated
0047 %                            ModelVar.PlantPlane        =    plane of the plant model
0048 %                                                            's':    continuous-time;
0049 %                                                            'w':    sqrt(s)-domain
0050 %                                                            'z':    discrete-time;
0051 %                                                            '':        plane not defined
0052 %                            ModelVar.Struct            =    model structure
0053 %                                                            'EIV':  errors-in-variables (noisy input-output data)
0054 %                                                            'OE':    generalised output error (known input, noisy output)
0055 %                            ModelVar.RecipPlant        =    1 if plant model is reciprocal: G(i,j) = G(j,i)
0056 %                            ModelVar.nu                =    number of inputs
0057 %                            ModelVar.ny                =     number of outputs
0058 %                            ModelVar.na                =    order polynomial A
0059 %                            ModelVar.nb                =     order matrix polynomial B
0060 %                            ModelVar.nig            =    order vector polynomial Ig
0061 %
0062 %
0063 % Copyright (c) Rik Pintelon, Vrije Universiteit Brussel - dept. ELEC, 27 November 2009
0064 % All rights reserved.
0065 % Software can be used freely for non-commercial applications only.
0066 %
0067 
0068 % note that DC and Nyquist have a contribution 1/2 to the cost function
0069 % therefore the appropriate variables must be scaled by 1/sqrt(2) at DC
0070 % and Nyquist; this is already done for the variable Error
0071 
0072 
0073 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0074 % calculate the derivative of the transfer functions w.r.t. the model parameters %
0075 %   a = Theta.A(:)                                                               %
0076 %   b = reshape(permute(Theta.B, [3,1,2]), [ny*nu*(nb+1),1])                     %
0077 %   ig = reshape(permute(Theta.Ig, [2,1]), [ny*(nig+1),1])                       %
0078 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0079 
0080 Deriv = MIMO_ML_CalcDeriv(xMat, PolyTrans, ModelVar);
0081 % derivatives of the hermitian transpose of the plant transfer function w.r.t. the model parameters
0082 if strcmp(ModelVar.Struct, 'EIV')
0083     Deriv.vecGHa = MIMO_ML_DvecZ2DvecZH(Deriv.vecGa, ModelVar.ny, ModelVar.nu);
0084     Deriv.vecGHb = MIMO_ML_DvecZ2DvecZH(Deriv.vecGb, ModelVar.ny, ModelVar.nu);    
0085 end % if errors-in-variables
0086 
0087 
0088 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0089 % derivative prediction error w.r.t. ALL model parameters %
0090 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0091 
0092 % number of frequencies
0093 F = size(xMat.Plant,1);
0094 
0095 % order polynomials
0096 na = ModelVar.na;
0097 nb = ModelVar.nb;
0098 nig = ModelVar.nig;
0099 
0100 % number of outputs and inputs
0101 ny = ModelVar.ny;
0102 nu = ModelVar.nu;
0103 
0104 % total number of model parameters
0105 ntheta = (na+1) + (nb+1)*nu*ny + (nig+1)*ny;
0106 
0107 TheJacob = zeros(ny, ntheta, F);
0108 if strcmp(ModelVar.Struct, 'EIV')
0109     
0110     % EIV part Jacobian matrix
0111     Jac_eiv = zeros(ny, ntheta, F);
0112            
0113     %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0114     % Fast calculation of matrices needed for calculating Jac_eiv. The lines below are equivalent with %
0115     % M_vecG = zeros(ny, ny*nu, F);                                                                    %
0116     % M_vecGH = zeros(ny, ny*nu, F);                                                                   %
0117     % for kk = 1:F                                                                                     %
0118     %    M1(:,:,kk) = PolyTrans.G(:,:,kk) * data.CU(:,:,kk) - data.CYU(:,:,kk);                        %
0119     %    V1(:,kk) = PolyTrans.sqrtCEinv(:,:,kk) * Error(:,kk);                                         %
0120     %    M_vecG(:,:,kk) = kron(V1(:,kk).'*conj(M1(:,:,kk)), PolyTrans.sqrtCEinv(:,:,kk));              %
0121     %    M_vecGH(:,:,kk) = kron(V1(:,kk).', PolyTrans.sqrtCEinv(:,:,kk)*M1(:,:,kk));                   %
0122     % end % kk, frequencies                                                                            %
0123     %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0124     
0125     M1 = Mat_Mult(PolyTrans.G, data.CU) - data.CYU;                     % Mat_Mult = matrix multiplication
0126     ErrorT = zeros(ny, 1, F);
0127     ErrorT(:,1,:) = Error;
0128     V1 = Mat_Mult(PolyTrans.sqrtCEinv, ErrorT);
0129     V1 = permute(V1, [2,1,3]);                                          % transpose of V1
0130     M_vecG = Kron_Row_Mat(Mat_Mult(V1, conj(M1)), PolyTrans.sqrtCEinv); % Kron_Row_Mat = Kronecker product of a row and a matrix
0131     M_vecGH = Kron_Row_Mat(V1, Mat_Mult(PolyTrans.sqrtCEinv, M1));
0132     
0133 end % if errors-in-variables
0134 
0135 % hermitian square root of the inverse of
0136 % the covariance matrix Cov(NY-G*NU)
0137 sqrtCEinv = PolyTrans.sqrtCEinv;
0138 
0139 
0140 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0141 % Fast calculation of the derivative w.r.t. a. The lines below are equivalent with               %
0142 % Low = 1;                                                                                       %
0143 % Upp = (na+1);                                                                                  %
0144 % for kk = 1:F                                                                                   %
0145 %    Mat = kron(data.U(:,kk).', sqrtCEinv(:,:,kk));                                               %
0146 %    TheJacob(:,Low:Upp,kk) = - Mat * Deriv.vecGa(:,:,kk) - sqrtCEinv(:,:,kk) * Deriv.Tga(:,:,kk) %
0147 % end                                                                                            %
0148 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0149 
0150 % fast calculation Kronecker product
0151 dataUT = zeros(1,nu,F);
0152 dataUT(1,:,:) = data.U;
0153 Mat = Kron_Row_Mat(dataUT, sqrtCEinv);
0154 
0155 % fast calculation Jacobian matrix w.r.t. a
0156 Low = 1;
0157 Upp = (na+1);
0158 TheJacob(:,Low:Upp,:) = - Mat_Mult(Mat, Deriv.vecGa) - Mat_Mult(sqrtCEinv, Deriv.Tga);
0159 
0160 
0161 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0162 % Fast calculation of the derivative w.r.t. a. The lines below are equivalent with                         %
0163 % Low = 1;                                                                                                 %
0164 % Upp = (na+1);                                                                                            %
0165 % for kk = 1:F                                                                                             %
0166 %    Jac_eiv(:,Low:Upp,kk) = M_vecG(:,:,kk) * Deriv.vecGa(:,:,kk) + M_vecGH(:,:,kk) * Deriv.vecGHa(:,:,kk); %
0167 % end                                                                                                      %
0168 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0169 
0170 if strcmp(ModelVar.Struct, 'EIV')
0171     
0172     % fast calculation EIV part Jacobian matrix w.r.t a
0173     Jac_eiv(:,Low:Upp,:) = Mat_Mult(M_vecG, Deriv.vecGa) + Mat_Mult(M_vecGH, Deriv.vecGHa);
0174     
0175 end % if errors-in-variables
0176 
0177 
0178 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0179 % Fast calculation of the derivative w.r.t. b. The lines below are equivalent with              %
0180 % Low = Upp + 1;                                                                                %
0181 % Upp = Upp = Low + (nb+1)*ny*nu - 1;                                                           %
0182 % for kk = 1:F                                                                                  %
0183 %    TheJacob(:,Low:Upp,kk) = - Mat * Deriv.vecGb(:,:,kk)                                        %
0184 % end                                                                                           %
0185 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0186 
0187 % fast calculation derivative w.r.t. b
0188 Low = Upp + 1;
0189 Upp = Low + (nb+1)*ny*nu - 1;
0190 TheJacob(:,Low:Upp,:) = - Mat_Mult(Mat, Deriv.vecGb);
0191 
0192 
0193 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0194 % Fast calculation of the derivative w.r.t. b. The lines below are equivalent with                         %
0195 % Low = Upp+1;                                                                                             %
0196 % Upp = Upp = Low + (nb+1)*ny*nu - 1;                                                                      %
0197 % for kk = 1:F                                                                                             %
0198 %    Jac_eiv(:,Low:Upp,kk) = M_vecG(:,:,kk) * Deriv.vecGb(:,:,kk) + M_vecGH(:,:,kk) * Deriv.vecGHb(:,:,kk); %
0199 % end                                                                                                      %
0200 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0201 
0202 if strcmp(ModelVar.Struct, 'EIV')
0203 
0204     % fast calculation EIV part Jacobian matrix w.r.t b
0205     Jac_eiv(:,Low:Upp,:) = Mat_Mult(M_vecG, Deriv.vecGb) + Mat_Mult(M_vecGH, Deriv.vecGHb);
0206     
0207 end % if errors-in-variables
0208 
0209 
0210 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0211 % Fast calculation of the derivative w.r.t. ig. The lines below are equivalent with             %
0212 % for kk = 1:F                                                                                  %
0213 %    Low = Upp+1;                                                                                %
0214 %    Upp = Low + (nig+1)*ny - 1;                                                                 %
0215 %    TheJacob(:,Low:Upp,kk) = - sqrtCEinv(:,:,kk) * Deriv.Tgig(:,:,kk)                           %
0216 % end                                                                                           %
0217 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0218 
0219 % fast calculation derivative w.r.t. ig
0220 Low = Upp+1;
0221 Upp = Low + (nig+1)*ny - 1;
0222 TheJacob(:,Low:Upp,:) = - Mat_Mult(sqrtCEinv, Deriv.Tgig);
0223 
0224 
0225 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0226 % Set the scaling for DC and Nyquist and add the EIV part %
0227 % DC and Nyquist count for 1/sqrt(2)                      %
0228 % This is already done in Error and hence also in Jac_eiv %
0229 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0230 
0231 if data.DC == 1
0232     TheJacob(:,:,1) = TheJacob(:,:,1)/sqrt(2);
0233 end
0234 if data.Nyquist == 1
0235     TheJacob(:,:,end) = TheJacob(:,:,end)/sqrt(2);
0236 end
0237 
0238 if strcmp(ModelVar.Struct, 'EIV')
0239     TheJacob = TheJacob - 0.5*Jac_eiv;
0240 end % if errors-in-variables

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