Maximimum likelihood estimate plant model parameters from noisy input-output data. The following cases can be handled: 1. 1 MIMO experiment with random excitations; or the concatenation of several data records (see the ArbLocalPolyAnal function). 2. 1 MIMO experiment with periodic excitations (see the FastLocalPolyAnal function). 3. nexp >=1 MIMO experiments with periodic excitations (see the RobustLocalPolyAnal function for the case nu = nexp). [Theta, Cost, smax, smin, wscale] = MIMO_ML(data, Sel, Theta0, ModelVar, IterVar) Stochastic framework: errors-in-variables Y = B/A U0 + Ig/A + NY U = U0 + NU System with nu inputs and ny outputs Y: ny x 1 observed output U: nu x 1 observed input NY: ny x 1 output noise NU: nu x 1 input noise Model class: common denominator model G = B/A: ny x nu plant transfer function Tg = Ig/A: ny x 1 plant transient term A: polynomial of order OrderA B: ny x nu matrix polynomial of order OrderB Ig: ny x 1 vector polynomial of order OrderIg 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., 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 Theta = estimated value 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 wscale = angular frequency scaling 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 nexp MIMO experiments: ny x nexp x F with nexp >= 1 data.U = input DFT spectra of 1 or nu independent MIMO experiments 1 MIMO experiment: nu x F nexp MIMO experiments: nu x nexp x F with nexp >= 1 data.freq = vector of frequency values (Hz), size: F x 1 or 1 x F data.Ts = sampling time (s) data.CY = (sample) noise covariance matrix of Y 1 MIMO experiment: ny x ny x F nexp MIMO experiments: ny x ny x nexp x F with nexp >= 1 data.CU = (sample) noise covariance matrix of U 1 MIMO experiment: nu x nu x F nexp MIMO experiments: nu x nu x nexp x F with nexp >= 1 data.CYU = (sample) noise covariance matrix of U 1 MIMO experiment: ny x nu x F nexp MIMO experiments: ny x nu x nexp x F with nexp >= 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 Theta0 = starting value plant, and initial conditions parameters structure with fields 'A', 'B', 'Ig' Theta0 = struct('A',[],'B',[], 'Ig', []) Theta0.A = 1 x (OrderA+1) Theta0.A(r) = coefficient a(r-1) of Omega^(r-1) Theta0.B = ny x nu x (OrderB+1) Theta0.B(i,j,r) = coefficient b(i,j,r-1) of Omega^(r-1) Theta0.Ig = ny x (OrderIg+1) Theta0.Ig(i,r) = coefficient ig(i,r-1) of Omega^(r-1) ModelVar = contains the information about the model to be identified structure with fields 'Transient', 'ThePlane', 'TheModel', 'RecipPlant' ModelVar = struct('Transient', [], 'PlantPlane', [], 'NoisePlane', [], '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) IterVar = contains the information about the minimization procedure structure with fields 'LM', 'MaxIter', 'TolParam', 'TolCost', 'TraceOn' IterVar = struct('LM', [], 'MaxIter', [], 'TolParam', [], 'TolCost', [], 'TraceOn', []) IterVar.LM = 1 then the Levenberg-Marquardt minimization scheme is used IterVar.MaxIter = maximum number of itterations of the minimization procedure IterVar.TolParam = relative precision on parameter estimates IterVar.TolCost = relative precision on cost function IterVar.TraceOn = 1 then output the iterations 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
0001 function [Theta, Cost, smax, smin, wscale] = MIMO_ML(data, Sel, Theta0, ModelVar, IterVar) 0002 % 0003 % Maximimum likelihood estimate plant model parameters from noisy 0004 % input-output data. The following cases can be handled: 0005 % 0006 % 1. 1 MIMO experiment with random excitations; or the concatenation of several data records 0007 % (see the ArbLocalPolyAnal function). 0008 % 2. 1 MIMO experiment with periodic excitations 0009 % (see the FastLocalPolyAnal function). 0010 % 3. nexp >=1 MIMO experiments with periodic excitations 0011 % (see the RobustLocalPolyAnal function for the case nu = nexp). 0012 % 0013 % 0014 % [Theta, Cost, smax, smin, wscale] = MIMO_ML(data, Sel, Theta0, ModelVar, IterVar) 0015 % 0016 % 0017 % Stochastic framework: errors-in-variables 0018 % Y = B/A U0 + Ig/A + NY 0019 % U = U0 + NU 0020 % 0021 % System with nu inputs and ny outputs 0022 % Y: ny x 1 observed output 0023 % U: nu x 1 observed input 0024 % NY: ny x 1 output noise 0025 % NU: nu x 1 input noise 0026 % 0027 % Model class: common denominator model 0028 % G = B/A: ny x nu plant transfer function 0029 % Tg = Ig/A: ny x 1 plant transient term 0030 % A: polynomial of order OrderA 0031 % B: ny x nu matrix polynomial of order OrderB 0032 % Ig: ny x 1 vector polynomial of order OrderIg 0033 % 0034 % Coefficients polynomials in raising powers of Omega, where 0035 % s-domain Omega = j*2*pi*freq 0036 % sqrt(s)-domain Omega = sqrt(j*2*pi*freq) 0037 % z-domain Omega = exp(-j*2*pi*freq*Ts) 0038 % 0039 % References: 0040 % 0041 % Pintelon, R., J. Schoukens, G. Vandersteen, and K. Barb�(2010). Estimation of nonparametric noise and 0042 % FRF models for multivariable systems - Part II: extensions, applications, Mechanical Systems and 0043 % Signal Processing, vol. 24, no. 3, pp. 596-616. 0044 % 0045 % Pintelon, R., G. Vandersteen, J. Schoukens, and Y. Rolain (2011). Improved (non-)parametric identification of dynamic 0046 % systems excited by periodic signals - The multivariate case, Mechanical Systems and Signal Processing, vol. 25, no. 8, 0047 % pp. 2892-2922. 0048 % 0049 % Pintelon, R., and J. Schoukens (2012). System Identification: A Frequency Domain Approach, second edition, 0050 % IEEE Press-Wiley, Piscataway (USA). 0051 % 0052 % 0053 % Output parameters 0054 % 0055 % Theta = estimated value plant, noise, and initial conditions parameters 0056 % structure with fields 'A', 'B', 'Ig' 0057 % Theta = struct('A',[],'B',[], 'Ig', []) 0058 % Theta.A = 1 x (OrderA+1) 0059 % Theta.A(r) = coefficient a(r-1) of Omega^(r-1) 0060 % Theta.B = ny x nu x (OrderB+1) 0061 % Theta.B(i,j,r) = coefficient b(i,j,r-1) of Omega^(r-1) 0062 % Theta.Ig = ny x (OrderIg+1) 0063 % Theta.Ig(i,r) = coefficient ig(i,r-1) of Omega^(r-1) 0064 % Note: all coefficients (except those for which Sel = 0) are free 0065 % during the minimization + in each iteration step the following 0066 % constraints are imposed: 0067 % norm([a, vec(b), vec(ig)] = 1 0068 % 0069 % Cost = value of the cost function in the last iteration step 0070 % 0071 % smax = largest singular value of the Jacobian matrix 0072 % 0073 % smin = smallest singular value of the Jacobian matrix 0074 % 0075 % wscale = angular frequency scaling 0076 % 0077 % 0078 % Input parameters 0079 % 0080 % data = structure containing the non-parametric data required for the identification 0081 % data.Y = output DFT spectra of 1 or nu independent MIMO experiments 0082 % 1 MIMO experiment: ny x F 0083 % nexp MIMO experiments: ny x nexp x F with nexp >= 1 0084 % data.U = input DFT spectra of 1 or nu independent MIMO experiments 0085 % 1 MIMO experiment: nu x F 0086 % nexp MIMO experiments: nu x nexp x F with nexp >= 1 0087 % data.freq = vector of frequency values (Hz), size: F x 1 or 1 x F 0088 % data.Ts = sampling time (s) 0089 % data.CY = (sample) noise covariance matrix of Y 0090 % 1 MIMO experiment: ny x ny x F 0091 % nexp MIMO experiments: ny x ny x nexp x F with nexp >= 1 0092 % data.CU = (sample) noise covariance matrix of U 0093 % 1 MIMO experiment: nu x nu x F 0094 % nexp MIMO experiments: nu x nu x nexp x F with nexp >= 1 0095 % data.CYU = (sample) noise covariance matrix of U 0096 % 1 MIMO experiment: ny x nu x F 0097 % nexp MIMO experiments: ny x nu x nexp x F with nexp >= 1 0098 % 0099 % Sel = structure with fields 'A', 'B', 'Ig' 0100 % Sel = struct('A',[],'B',[], 'Ig', []) 0101 % Sel.A = 1 x (OrderA+1) 0102 % Sel.A(r) = 1 if coeff. a(r-1) is unknown 0103 % Sel.A(r) = 0 if coeff. a(r-1) = 0 0104 % Sel.B = ny x nu x (OrderB+1) 0105 % Sel.B(i,j,r) = 1 if coeff. b(i,j,r-1) is unknown 0106 % Sel.B(i,j,r) = 0 if coeff. b(i,j,r-1) = 0 0107 % Sel.Ig = ny x (OrderIg+1) 0108 % Sel.Ig(i,r) = 1 if coeff. ig(i,r-1) is unknown 0109 % Sel.Ig(i,r) = 0 if coeff. ig(i,r-1) = 0 0110 % 0111 % Theta0 = starting value plant, and initial conditions parameters 0112 % structure with fields 'A', 'B', 'Ig' 0113 % Theta0 = struct('A',[],'B',[], 'Ig', []) 0114 % Theta0.A = 1 x (OrderA+1) 0115 % Theta0.A(r) = coefficient a(r-1) of Omega^(r-1) 0116 % Theta0.B = ny x nu x (OrderB+1) 0117 % Theta0.B(i,j,r) = coefficient b(i,j,r-1) of Omega^(r-1) 0118 % Theta0.Ig = ny x (OrderIg+1) 0119 % Theta0.Ig(i,r) = coefficient ig(i,r-1) of Omega^(r-1) 0120 % 0121 % ModelVar = contains the information about the model to be identified 0122 % structure with fields 'Transient', 'ThePlane', 'TheModel', 'RecipPlant' 0123 % ModelVar = struct('Transient', [], 'PlantPlane', [], 'NoisePlane', [], 'Struct', [], 'RecipPlant', []) 0124 % ModelVar.Transient = 1 then the initial conditions of the plant and/or noise are estimated 0125 % ModelVar.PlantPlane = plane of the plant model 0126 % 's': continuous-time; 0127 % 'w': sqrt(s)-domain 0128 % 'z': discrete-time; 0129 % '': plane not defined 0130 % ModelVar.Struct = model structure 0131 % 'EIV': errors-in-variables (noisy input-output data) 0132 % 'OE': generalised output error (known input, noisy output) 0133 % ModelVar.RecipPlant = 1 if plant model is reciprocal: G(i,j) = G(j,i) 0134 % 0135 % IterVar = contains the information about the minimization procedure 0136 % structure with fields 'LM', 'MaxIter', 'TolParam', 'TolCost', 'TraceOn' 0137 % IterVar = struct('LM', [], 'MaxIter', [], 'TolParam', [], 'TolCost', [], 'TraceOn', []) 0138 % IterVar.LM = 1 then the Levenberg-Marquardt minimization scheme is used 0139 % IterVar.MaxIter = maximum number of itterations of the minimization procedure 0140 % IterVar.TolParam = relative precision on parameter estimates 0141 % IterVar.TolCost = relative precision on cost function 0142 % IterVar.TraceOn = 1 then output the iterations 0143 % 0144 % 0145 % Copyright (c) Rik Pintelon, Vrije Universiteit Brussel - dept. ELEC, November 2009 0146 % All rights reserved. 0147 % Software can be used freely for non-commercial applications only. 0148 % Version 24 October 2011 0149 % 0150 0151 0152 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% 0153 % initialisation of the variables, and compatibility check of the input arguments % 0154 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% 0155 0156 % model structure variables 0157 ModelVar.PlantPlane = lower(ModelVar.PlantPlane); 0158 ModelVar.Struct = upper(ModelVar.Struct); 0159 0160 % add number of inputs and outputs to ModelVar 0161 ModelVar.ny = size(data.Y,1); 0162 ModelVar.nu = size(data.U,1); 0163 0164 % 1. imposes the compatibility of the parameter vector Theta0 and 0165 % the free model parameters with the model structure 0166 % 2. puts the order of the polynomials in ModelVar 0167 [Theta0, Sel, ModelVar] = MIMO_ML_ModelCompatibility(Theta0, Sel, ModelVar); 0168 0169 % check if DC and Nyquist belong to the frequency set 0170 if data.freq(1) == 0, data.DC = 1; else data.DC = 0; end 0171 if data.freq(end) == 1/(2*data.Ts), data.Nyquist = 1; else data.Nyquist = 0; end 0172 data.freq = data.freq(:); 0173 0174 % determine the number of MIMO experiments 0175 NumberDim = length(size(data.U)); % number of matrix dimensions 0176 if NumberDim == 2 0177 data.NumberExp = 1; % number of MIMO experiments 0178 elseif NumberDim == 3 0179 data.NumberExp = size(data.U, 2); % number of MIMO experiments 0180 end % if 0181 0182 % hermitian symmetric square root inverse output covariance matrix 0183 F = length(data.freq); 0184 try 0185 Weighting = ~isempty(data.CY); 0186 catch 0187 Weighting = 0; 0188 end 0189 if Weighting 0190 data.sqrtCYinv = zeros(size(data.CY)); 0191 if data.NumberExp == 1 0192 for kk = 1:F 0193 [uu, ss, vv] = svd(squeeze(data.CY(:,:,kk)), 0); 0194 data.sqrtCYinv(:,:,kk) = vv*diag(diag(ss).^(-0.5))*vv'; 0195 end % kk 0196 else % more than 1 MIMO experiment 0197 for ee = 1:data.NumberExp 0198 for kk = 1:F 0199 [uu, ss, vv] = svd(squeeze(data.CY(:,:,ee,kk)), 0); 0200 data.sqrtCYinv(:,:,ee,kk) = vv*diag(diag(ss).^(-0.5))*vv'; 0201 end % kk 0202 end % ee, MIMO experiments 0203 end % if 1 MIMO experiment 0204 else 0205 ny = ModelVar.ny; 0206 nu = ModelVar.nu; 0207 if data.NumberExp == 1 0208 % if empty put eye(ny, ny) for every frequency in data.CY 0209 data.CY = repmat(eye(ny),[1,1,F]); 0210 data.sqrtCYinv = data.CY; 0211 else % more than 1 MIMO experiment 0212 % if empty put eye(ny, ny) for every experiment and frequency in data.CY 0213 data.CY = repmat(eye(ny),[1,1,nu,F]); 0214 data.sqrtCYinv = data.CY; 0215 end % if 1 MIMO experiment 0216 end 0217 0218 % jw, sqrt(jw), or exp(-jwTs) values for plant and noise model 0219 x = struct('Plant', []); 0220 % matrix of powers of x 0221 xMat = struct('Plant', []); 0222 0223 0224 %%%%%%%%%%%%%%%%%%%%%%%%%%%%% 0225 % domain of the plant model % 0226 %%%%%%%%%%%%%%%%%%%%%%%%%%%%% 0227 0228 % the vector s represents z^-1, s, or sqrt(s) of the plant model 0229 switch ModelVar.PlantPlane 0230 case {'s','w'} 0231 if ModelVar.PlantPlane == 's' 0232 x.Plant = sqrt(-1)*2*pi*data.freq; 0233 elseif ModelVar.PlantPlane == 'w' 0234 x.Plant = (sqrt(-1)*2*pi*data.freq).^(0.5); 0235 end; 0236 wscale = median(abs(x.Plant)); 0237 x.Plant = x.Plant/wscale; 0238 case 'z' 0239 x.Plant = exp(-sqrt(-1)*2*pi*data.freq*data.Ts); 0240 wscale = 1; 0241 case '' 0242 x.Plant = ones(size(data.freq)); 0243 wscale = 1; 0244 otherwise, disp('Invalid plant plane ...'), return 0245 end 0246 0247 nmax = max([ModelVar.na, ModelVar.nb, ModelVar.nig]); 0248 xMat.Plant = MIMO_ML_CalcOmegaMat(x.Plant, nmax); 0249 0250 % normalize the starting values for s-, and sqrt(s)-domains 0251 Theta0 = MIMO_ML_Normalise(Theta0, wscale, ModelVar); 0252 0253 0254 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% 0255 % maximum likelihood estimate of the model parameters % 0256 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% 0257 0258 % cost function starting values 0259 Cost = MIMO_ML_Cost(data, Theta0, x, ModelVar); 0260 0261 % initialisation itertation parameters 0262 MaxRelErr = NaN; 0263 RelVarCost = NaN; 0264 DeltaCost = NaN; 0265 OldTheta = Theta0; 0266 OldCost = Cost; 0267 OldMaxRelErr = MaxRelErr; 0268 OldRelVarCost = RelVarCost; 0269 OldDeltaCost = DeltaCost; 0270 0271 if IterVar.TraceOn == 1 0272 home 0273 disp(['ML iteration # ',int2str(0)]); 0274 fprintf('condition number: NaN \n'); 0275 fprintf('maximum relative variation parameters: %e \n',MaxRelErr); 0276 fprintf('cost function (delta cost): %e (%e ) \n',Cost,DeltaCost); 0277 end 0278 0279 % starting values only 0280 if IterVar.MaxIter < 1 0281 Theta = OldTheta; 0282 end 0283 0284 % Newton-Gauss or Levenberg-Marquardt iteration loop 0285 Lambda = 0; 0286 StartLM = 1; % initialization for Levenberg-Marquardt minimization scheme 0287 0288 iii = 0; 0289 while iii <= IterVar.MaxIter 0290 0291 iii = iii + 1; 0292 0293 [Theta, Cost, smax, smin] = MIMO_ML_NewtonGaussStep(data, x, xMat, OldTheta, Sel, ModelVar, Lambda); 0294 0295 if ~((IterVar.LM == 1) & (Cost > OldCost)) 0296 0297 % estimated parameters polynomial A 0298 A = Theta.A(Sel.A == 1); 0299 OldA = OldTheta.A(Sel.A == 1); 0300 MaxRelErr = max([abs((A - OldA)./A)]); 0301 0302 % estimated parameters ny x nu matrix polynomial B 0303 for ii = 1:ModelVar.ny 0304 for jj = 1:ModelVar.nu 0305 B = squeeze(Theta.B(ii, jj, Sel.B(ii,jj,:) == 1)).'; 0306 OldB = squeeze(OldTheta.B(ii, jj, Sel.B(ii,jj,:) == 1)).'; 0307 MaxRelErr = max([MaxRelErr, abs((B - OldB)./B)]); 0308 end % jj 0309 end % ii 0310 0311 % estimated parameters ny x 1 vector polynomial Ig 0312 for ii = 1:ModelVar.ny 0313 Ig = Theta.Ig(ii, Sel.Ig(ii,:) == 1); 0314 OldIg = OldTheta.Ig(ii, Sel.Ig(ii,:) == 1); 0315 MaxRelErr = max([MaxRelErr, abs((Ig - OldIg)./Ig)]); 0316 end % ii 0317 0318 DeltaCost = Cost - OldCost; 0319 RelVarCost = abs(DeltaCost)/Cost; 0320 OldTheta = Theta; 0321 OldCost = Cost; 0322 OldMaxRelErr = MaxRelErr; 0323 OldRelVarCost = RelVarCost; 0324 OldDeltaCost = DeltaCost; 0325 0326 end % if not Levenberg-Marquardt or if Cost <= OldCost 0327 0328 if IterVar.LM 0329 0330 if Cost > OldCost 0331 0332 if StartLM == 1 0333 Lambda = smax^2; 0334 StartLM = 0; % initialization Lambda parameter is done 0335 else 0336 Lambda = 10*Lambda; 0337 end % if start Levenberg-Marquardt 0338 0339 Theta = OldTheta; 0340 Cost = OldCost; 0341 MaxRelErr = OldMaxRelErr; 0342 RelVarCost = OldRelVarCost; 0343 DeltaCost = OldDeltaCost; 0344 0345 else 0346 0347 Lambda = 0.4*Lambda; 0348 0349 end % if Cost > OldCost 0350 0351 end % if Levenberg-Marquardt 0352 0353 if IterVar.TraceOn 0354 0355 home 0356 disp(['ML iteration # ',int2str(iii)]); 0357 disp(['condition number:',num2str(smax/smin)]); 0358 disp(['maximum relative variation parameters: ',num2str(MaxRelErr)]); 0359 disp(['cost function (delta cost): ',num2str(Cost),' (',num2str(DeltaCost),')']); 0360 if IterVar.LM == 1 0361 disp(['Lambda factor Levenberg-Marquardt iteration scheme: ',num2str(Lambda)]); 0362 end 0363 0364 end % if trace on 0365 0366 if (MaxRelErr < IterVar.TolParam) | (RelVarCost < IterVar.TolCost) 0367 iii = IterVar.MaxIter + 1; 0368 end % if relative variation parameters or cost function is small enough 0369 0370 end % while main iteration loop 0371 0372 0373 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% 0374 % denormalize the estimated parameters for s-, and sqrt(s)-domains % 0375 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% 0376 0377 Theta = MIMO_ML_DeNormalise(Theta, wscale, ModelVar);