MIMO_ML_NewtonGaussStep

PURPOSE ^

SYNOPSIS ^

function [Theta, Cost, smax, smin] = MIMO_ML_NewtonGaussStep(data, x, xMat, OldTheta, Sel, ModelVar, Lambda);

DESCRIPTION ^


    [Theta, Cost, smax, smin] = MIMO_ML_NewtonGaussStep(data, x, xMat, OldTheta, Sel, ModelVar, Lambda)


    Output parameters

        Theta        =    new estimate plant, noise, and initial conditions parameters
                        structure with fields 'A', 'B', 'Ig'
                            Theta = struct('A',[],'B',[], 'Ig', [])
                            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.Ig = ny x (OrderIg+1)
                                Theta.Ig(i,r) = coefficient ig(i,r-1) of Omega^(r-1) 
                        Note:    all coefficients (except those for which Sel = 0) are free
                                during the minimization + in each iteration step the following
                                constraints are imposed:
                                    norm([a, vec(b), vec(ig)] = 1

        Cost        =    value of the cost function in the last iteration step

        smax        =    largest singular value of the Jacobian matrix

        smin        =    smallest singular value of the Jacobian matrix


    Input parameters

        data        =    structure containing the non-parametric data required for the identification
                            data.Y          =    output DFT spectra of 1 or nu independent MIMO experiments 
                                                   1 MIMO experiment:   ny x F 
                                                    nu MIMO experiments: ny x nu x F 
                            data.U          =    input DFT spectra of 1 or nu independent MIMO experiments 
                                                    1 MIMO experiment:   nu x F 
                                                   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 
                                                   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  
                                                   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 
                                                    1 MIMO experiment:   ny x nu x F 
                                                   nu MIMO experiments: ny x nu x nu x F 
                            data.sqrtCYinv    =    CY^(-0.5); if OE 
                                                    1 MIMO experiment:   ny x ny x F 
                                                   nu MIMO experiments: ny x ny x nu x F 
                            data.DC            =    1 if DC present otherwise 0
                            data.Nyquist    =    1 if Nyquist frequency present otherwise 0
                           data.NumberExp  =   number of independent MIMO experiments (1 or nu) 

        x            =    structure containing (jwk) or (zk^-1) values
                            x.Plant    =    plant model, size: F x 1

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

        OldTheta    =    previous estimate plant, noise, and initial conditions parameters
                        structure with fields 'A', 'B', 'Ig'
                        OldTheta = struct('A',[],'B',[], 'Ig', [])
                            OldTheta.A = 1 x (OrderA+1)
                                OldTheta.A(r) = coefficient a(r-1) of Omega^(r-1) 
                            OldTheta.B = ny x nu x (OrderB+1)
                                OldTheta.B(i,j,r) = coefficient b(i,j,r-1) of Omega^(r-1)
                            OldTheta.Ig = ny x (OrderIg+1)
                                OldTheta.Ig(i,r) = coefficient ig(i,r-1) of Omega^(r-1) 
                        Note:    all coefficients (except those for which Sel = 0) are free
                                during the minimization + in each iteration step the following
                                constraints are imposed:
                                    norm([a, vec(b), vec(ig)] = 1

        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 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
        Lambda            =    Levenberg-Marquardt parameter


 Copyright (c) Rik Pintelon, Vrije Universiteit Brussel - dept. ELEC, 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 [Theta, Cost, smax, smin] = MIMO_ML_NewtonGaussStep(data, x, xMat, OldTheta, Sel, ModelVar, Lambda);
0002 %
0003 %
0004 %    [Theta, Cost, smax, smin] = MIMO_ML_NewtonGaussStep(data, x, xMat, OldTheta, Sel, ModelVar, Lambda)
0005 %
0006 %
0007 %    Output parameters
0008 %
0009 %        Theta        =    new estimate plant, noise, and initial conditions parameters
0010 %                        structure with fields 'A', 'B', 'Ig'
0011 %                            Theta = struct('A',[],'B',[], 'Ig', [])
0012 %                            Theta.A = 1 x (OrderA+1)
0013 %                                Theta.A(r) = coefficient a(r-1) of Omega^(r-1)
0014 %                            Theta.B = ny x nu x (OrderB+1)
0015 %                                Theta.B(i,j,r) = coefficient b(i,j,r-1) of Omega^(r-1)
0016 %                            Theta.Ig = ny x (OrderIg+1)
0017 %                                Theta.Ig(i,r) = coefficient ig(i,r-1) of Omega^(r-1)
0018 %                        Note:    all coefficients (except those for which Sel = 0) are free
0019 %                                during the minimization + in each iteration step the following
0020 %                                constraints are imposed:
0021 %                                    norm([a, vec(b), vec(ig)] = 1
0022 %
0023 %        Cost        =    value of the cost function in the last iteration step
0024 %
0025 %        smax        =    largest singular value of the Jacobian matrix
0026 %
0027 %        smin        =    smallest singular value of the Jacobian matrix
0028 %
0029 %
0030 %    Input parameters
0031 %
0032 %        data        =    structure containing the non-parametric data required for the identification
0033 %                            data.Y          =    output DFT spectra of 1 or nu independent MIMO experiments
0034 %                                                   1 MIMO experiment:   ny x F
0035 %                                                    nu MIMO experiments: ny x nu x F
0036 %                            data.U          =    input DFT spectra of 1 or nu independent MIMO experiments
0037 %                                                    1 MIMO experiment:   nu x F
0038 %                                                   nu MIMO experiments: nu x nu x F
0039 %                            data.freq       =    vector of frequency values (Hz), size: F x 1
0040 %                            data.Ts         =    sampling time (s)
0041 %                            data.CY         =    (sample) noise covariance matrix of Y
0042 %                                                   1 MIMO experiment:   ny x ny x F
0043 %                                                    nu MIMO experiments: ny x ny x nu x F
0044 %                           data.CU         =   (sample) noise covariance matrix of U
0045 %                                                   1 MIMO experiment:   nu x nu x F
0046 %                                                    nu MIMO experiments: nu x nu x nu x F
0047 %                           data.CYU        =   (sample) noise covariance matrix of U
0048 %                                                    1 MIMO experiment:   ny x nu x F
0049 %                                                   nu MIMO experiments: ny x nu x nu x F
0050 %                            data.sqrtCYinv    =    CY^(-0.5); if OE
0051 %                                                    1 MIMO experiment:   ny x ny x F
0052 %                                                   nu MIMO experiments: ny x ny x nu x F
0053 %                            data.DC            =    1 if DC present otherwise 0
0054 %                            data.Nyquist    =    1 if Nyquist frequency present otherwise 0
0055 %                           data.NumberExp  =   number of independent MIMO experiments (1 or nu)
0056 %
0057 %        x            =    structure containing (jwk) or (zk^-1) values
0058 %                            x.Plant    =    plant model, size: F x 1
0059 %
0060 %        xMat        =    structure with tables of powers of (jwk)^r or (zk^-r)
0061 %                            xMat.Plant        =    plant model, size: F x max order
0062 %
0063 %        OldTheta    =    previous estimate plant, noise, and initial conditions parameters
0064 %                        structure with fields 'A', 'B', 'Ig'
0065 %                        OldTheta = struct('A',[],'B',[], 'Ig', [])
0066 %                            OldTheta.A = 1 x (OrderA+1)
0067 %                                OldTheta.A(r) = coefficient a(r-1) of Omega^(r-1)
0068 %                            OldTheta.B = ny x nu x (OrderB+1)
0069 %                                OldTheta.B(i,j,r) = coefficient b(i,j,r-1) of Omega^(r-1)
0070 %                            OldTheta.Ig = ny x (OrderIg+1)
0071 %                                OldTheta.Ig(i,r) = coefficient ig(i,r-1) of Omega^(r-1)
0072 %                        Note:    all coefficients (except those for which Sel = 0) are free
0073 %                                during the minimization + in each iteration step the following
0074 %                                constraints are imposed:
0075 %                                    norm([a, vec(b), vec(ig)] = 1
0076 %
0077 %        Sel            =    structure with fields 'A', 'B', 'Ig'
0078 %                            Sel = struct('A',[],'B',[], 'Ig', [])
0079 %                            Sel.A = 1 x (OrderA+1)
0080 %                                Sel.A(r) = 1 if coeff. a(r-1) is unknown
0081 %                                Sel.A(r) = 0 if coeff. a(r-1) = 0
0082 %                            Sel.B = ny x nu x (OrderB+1)
0083 %                                Sel.B(i,j,r) = 1 if coeff. b(i,j,r-1) is unknown
0084 %                                Sel.B(i,j,r) = 0 if coeff. b(i,j,r-1) = 0
0085 %                            Sel.Ig = ny x (OrderIg+1)
0086 %                                Sel.Ig(i,r) = 1 if coeff. ig(i,r-1) is unknown
0087 %                                Sel.Ig(i,r) = 0 if coeff. ig(i,r-1) = 0
0088 %
0089 %        ModelVar        =    contains the information about the model to be identified
0090 %                            structure with fields 'Transient', 'ThePlane', 'TheModel', 'Reciprocal', ...
0091 %                                ModelVar.Transient        =    1 then the initial conditions of the plant and/or noise are estimated
0092 %                                ModelVar.PlantPlane        =    plane of the plant model
0093 %                                                                's':    continuous-time;
0094 %                                                                'w':    sqrt(s)-domain
0095 %                                                                'z':    discrete-time;
0096 %                                                                '':        plane not defined
0097 %                                ModelVar.Struct            =    model structure
0098 %                                                               'EIV':  errors-in-variables (noisy input-output data)
0099 %                                                               'OE':    generalised output error (known input, noisy output)
0100 %                                ModelVar.RecipPlant        =    1 if plant model is reciprocal: G(i,j) = G(j,i)
0101 %                                ModelVar.nu                =    number of inputs
0102 %                                ModelVar.ny                =     number of outputs
0103 %                                ModelVar.na                =    order polynomial A
0104 %                                ModelVar.nb                =     order matrix polynomial B
0105 %                                ModelVar.nig            =    order vector polynomial Ig
0106 %        Lambda            =    Levenberg-Marquardt parameter
0107 %
0108 %
0109 % Copyright (c) Rik Pintelon, Vrije Universiteit Brussel - dept. ELEC, November 2009
0110 % All rights reserved.
0111 % Software can be used freely for non-commercial applications only.
0112 %
0113 
0114 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0115 % initialisation of the variables %
0116 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0117 
0118 F = length(data.freq);                          % number of frequencies
0119 NumberExp = data.NumberExp;                     % number of independent MIMO experiments
0120 ny = ModelVar.ny;                               % number of outputs
0121 nu = ModelVar.nu;                               % number of inputs
0122 na = ModelVar.na;                               % order A polynomial
0123 nb = ModelVar.nb;                               % order B matrix polynomial
0124 nig = ModelVar.nig;                             % order Ig vector polynomial
0125 ntheta = (na+1) + (nb+1)*nu*ny + (nig+1)*ny;    % total number of parameters
0126 TheJacob = zeros(NumberExp*ny*F, ntheta);       % the jacobian matrix of all MIMO experiments
0127 TheError = zeros(NumberExp*ny*F, 1);            % the prediction error vector of all MIMO experiments
0128 dataee = data;                                  % information about 1 MIMO experiment (see below)
0129 dataee.Y = zeros(ny, F);                        % output of 1 MIMO experiment
0130 dataee.U = zeros(nu, F);                        % input of 1 MIMO experiment
0131 dataee.CY = zeros(ny, ny, F);                   % output covariance of 1 MIMO experiment
0132 switch ModelVar.Struct
0133     case 'EIV'
0134         dataee.CU = zeros(nu, nu, F);           % input covariance of 1 MIMO experiment
0135         dataee.CYU = zeros(ny, nu, F);          % output-input covariance of 1 MIMO experiment
0136     case 'OE'
0137         dataee.sqrtCYinv = zeros(ny, ny, F);
0138 end % switch
0139 
0140 
0141 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0142 % calculate the variation of the free model parameters %
0143 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0144 
0145 % transfer functions and polynomials
0146 PolyTrans = MIMO_ML_CalcPolyTrans(OldTheta, x);
0147 
0148 for ee = 1:NumberExp
0149     
0150     % data of MIMO experiment no. ee
0151     if NumberExp > 1
0152         dataee.Y(:,:) = data.Y(:,ee,:);
0153         dataee.U(:,:) = data.U(:,ee,:);
0154         dataee.CY(:,:,:) = data.CY(:,:,ee,:);
0155          switch ModelVar.Struct
0156             case 'EIV'
0157                 dataee.CU(:,:,:) = data.CU(:,:,ee,:);
0158                 dataee.CYU(:,:,:) = data.CYU(:,:,ee,:);              
0159             case 'OE'
0160                 dataee.sqrtCYinv(:,:,:) = data.sqrtCYinv(:,:,ee,:);
0161         end % switch
0162     else % 1 MIMO experiment
0163         dataee.Y(:,:) = data.Y;
0164         dataee.U(:,:) = data.U;
0165         dataee.CY(:,:,:) = data.CY;
0166         switch ModelVar.Struct
0167             case 'EIV'
0168                 dataee.CU(:,:,:) = data.CU;
0169                 dataee.CYU(:,:,:) = data.CYU;              
0170             case 'OE'
0171                 dataee.sqrtCYinv(:,:,:) = data.sqrtCYinv;
0172         end % switch
0173     end % if more than 1 MIMO experiment
0174     
0175     % calculate a hermitian square root of the inverse
0176     % of the covariance of the output error (NY-G*NU).
0177     switch ModelVar.Struct
0178         case 'EIV'
0179             PolyTrans = MIMO_ML_InvCovOutputError(dataee, PolyTrans);
0180         case 'OE'
0181             PolyTrans.sqrtCEinv = dataee.sqrtCYinv;
0182     end % switch
0183 
0184     % prediction Error, and sample covariance matrix CovE of the prediction error
0185     Error = MIMO_ML_PredError(dataee, PolyTrans);
0186 
0187     % Jacobian matrix of prediction error w.r.t. ALL model parameters
0188     Jacob = MIMO_ML_Jacob(dataee, xMat, Error, PolyTrans, ModelVar);
0189 
0190     % put the different frequency contributions on top of each other
0191     % frequency contributions on top of each other
0192     % size Jacob in: ny x ntheta x F
0193     % size Jacob out: ny*F x ntheta
0194     Jacob = reshape(permute(Jacob, [1,3,2]), [ny*F, ntheta]);
0195     
0196     % put experiment no. ee in the Jacobian matrix and prediction error vector
0197     SelectRows = [(ee-1)*ny*F+1:ee*ny*F];
0198     TheError(SelectRows) = Error(:);
0199     TheJacob(SelectRows, :) = Jacob;
0200    
0201 end % ee, number of MIMO experiments
0202 
0203 % impose the common parameter structure and eliminate the excess parameters
0204 % size Jacob in: ny*F x ntheta
0205 % size Jacob out: ny*F x number of free model parameters
0206 TheJacob = MIMO_ML_AddSelectColumns(TheJacob, Sel, ModelVar);
0207 
0208 % scaling of the columns of the Jacobian matrix to improve the condition number
0209 TheScale = sum(abs(TheJacob.^2), 1).^0.5;
0210 IndexZeroes = find(TheScale == 0);
0211 TheScale(IndexZeroes) = 1;
0212 TheJacob = TheJacob ./ repmat(TheScale, [size(TheJacob, 1), 1]);
0213 
0214 % calculate the parameter variation DeltaParam
0215 [u, ss, v] = svd([real(TheJacob); imag(TheJacob)], 0);
0216 
0217 nn =1;
0218 smax = ss(1,1);
0219 smin = ss(end-nn, end-nn);
0220 diags = diag(ss);
0221 diags = diags(1:end-nn);
0222 DeltaParam = -v*diag([diags./(diags.^2+Lambda); zeros(nn,1)])*(u.'*[real(TheError);imag(TheError)]);
0223 
0224 % denormalisation parameter variation
0225 DeltaParam = DeltaParam ./ TheScale.';
0226 
0227 
0228 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0229 % extract the plant model parameters from DeltaParam %
0230 % order in DeltaParam: a, b, ig                      %
0231 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0232 
0233 Theta = MIMO_ML_ExtractParam(DeltaParam, OldTheta, Sel);    
0234 
0235 
0236 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0237 % Impose constraints on updated parameters %
0238 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0239 
0240 Theta = MIMO_ML_Constrain(Theta, ModelVar);
0241 
0242 
0243 %%%%%%%%%%%%%%%%%
0244 % cost function %
0245 %%%%%%%%%%%%%%%%%
0246 
0247 Cost = MIMO_ML_Cost(data, Theta, x, ModelVar);
0248

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