Bootstrapped total least squares 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_BTLS(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) Reference: Pintelon R., P. Guillaume, G. Vandersteen and Y. Rolain (1998). Analyses, development and applications of TLS algorithms in frequency-Domain System Identification, SIAM J. Matrix Anal. Appl., vol. 19, no. 4, pp. 983-1004. 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', [], '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.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 = if 1 then the iteration stops if the corresponding ML cost function increases 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 iterations (optional) Copyright (c) Rik Pintelon, Vrije Universiteit Brussel - dept. ELEC, 30 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_BTLS(data, Sel, Theta0, ModelVar, IterVar); 0002 % 0003 % Bootstrapped total least squares 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_BTLS(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 % Reference: 0040 % 0041 % Pintelon R., P. Guillaume, G. Vandersteen and Y. Rolain (1998). Analyses, development and applications of TLS 0042 % algorithms in frequency-Domain System Identification, SIAM J. Matrix Anal. Appl., vol. 19, no. 4, pp. 983-1004. 0043 % 0044 % 0045 % Output parameters 0046 % 0047 % Theta = estimated value plant, noise, and initial conditions parameters 0048 % structure with fields 'A', 'B', 'Ig' 0049 % Theta = struct('A',[],'B',[], 'Ig', []) 0050 % Theta.A = 1 x (OrderA+1) 0051 % Theta.A(r) = coefficient a(r-1) of Omega^(r-1) 0052 % Theta.B = ny x nu x (OrderB+1) 0053 % Theta.B(i,j,r) = coefficient b(i,j,r-1) of Omega^(r-1) 0054 % Theta.Ig = ny x (OrderIg+1) 0055 % Theta.Ig(i,r) = coefficient ig(i,r-1) of Omega^(r-1) 0056 % Note: all coefficients (except those for which Sel = 0) are free 0057 % during the minimization + in each iteration step the following 0058 % constraints are imposed: 0059 % norm([a, vec(b), vec(ig)] = 1 0060 % 0061 % Cost = value of the cost function in the last iteration step 0062 % 0063 % smax = largest singular value of the Jacobian matrix 0064 % 0065 % smin = smallest singular value of the Jacobian matrix 0066 % 0067 % wscale = angular frequency scaling 0068 % 0069 % 0070 % Input parameters 0071 % 0072 % data = structure containing the non-parametric data required for the identification 0073 % data.Y = output DFT spectra of 1 or nu independent MIMO experiments 0074 % 1 MIMO experiment: ny x F 0075 % nexp MIMO experiments: ny x nexp x F with nexp >= 1 0076 % data.U = input DFT spectra of 1 or nu independent MIMO experiments 0077 % 1 MIMO experiment: nu x F 0078 % nexp MIMO experiments: nu x nexp x F with nexp >= 1 0079 % data.freq = vector of frequency values (Hz), size: F x 1 or 1 x F 0080 % data.Ts = sampling time (s) 0081 % data.CY = (sample) noise covariance matrix of Y 0082 % 1 MIMO experiment: ny x ny x F 0083 % nexp MIMO experiments: ny x ny x nexp x F with nexp >= 1 0084 % data.CU = (sample) noise covariance matrix of U 0085 % 1 MIMO experiment: nu x nu x F 0086 % nexp MIMO experiments: nu x nu x nexp x F with nexp >= 1 0087 % data.CYU = (sample) noise covariance matrix of U 0088 % 1 MIMO experiment: ny x nu x F 0089 % nexp MIMO experiments: ny x nu x nexp x F with nexp >= 1 0090 % 0091 % Sel = structure with fields 'A', 'B', 'Ig' 0092 % Sel = struct('A',[],'B',[], 'Ig', []) 0093 % Sel.A = 1 x (OrderA+1) 0094 % Sel.A(r) = 1 if coeff. a(r-1) is unknown 0095 % Sel.A(r) = 0 if coeff. a(r-1) = 0 0096 % Sel.B = ny x nu x (OrderB+1) 0097 % Sel.B(i,j,r) = 1 if coeff. b(i,j,r-1) is unknown 0098 % Sel.B(i,j,r) = 0 if coeff. b(i,j,r-1) = 0 0099 % Sel.Ig = ny x (OrderIg+1) 0100 % Sel.Ig(i,r) = 1 if coeff. ig(i,r-1) is unknown 0101 % Sel.Ig(i,r) = 0 if coeff. ig(i,r-1) = 0 0102 % 0103 % Theta0 = starting value plant, and initial conditions parameters 0104 % structure with fields 'A', 'B', 'Ig' 0105 % Theta0 = struct('A',[],'B',[], 'Ig', []) 0106 % Theta0.A = 1 x (OrderA+1) 0107 % Theta0.A(r) = coefficient a(r-1) of Omega^(r-1) 0108 % Theta0.B = ny x nu x (OrderB+1) 0109 % Theta0.B(i,j,r) = coefficient b(i,j,r-1) of Omega^(r-1) 0110 % Theta0.Ig = ny x (OrderIg+1) 0111 % Theta0.Ig(i,r) = coefficient ig(i,r-1) of Omega^(r-1) 0112 % 0113 % ModelVar = contains the information about the model to be identified 0114 % structure with fields 'Transient', 'ThePlane', 'TheModel', 'RecipPlant' 0115 % ModelVar = struct('Transient', [], 'PlantPlane', [], 'NoisePlane', [], 'RecipPlant', []) 0116 % ModelVar.Transient = 1 then the initial conditions of the plant and/or noise are estimated 0117 % ModelVar.PlantPlane = plane of the plant model 0118 % 's': continuous-time; 0119 % 'w': sqrt(s)-domain 0120 % 'z': discrete-time; 0121 % '': plane not defined 0122 % ModelVar.RecipPlant = 1 if plant model is reciprocal: G(i,j) = G(j,i) 0123 % 0124 % IterVar = contains the information about the minimization procedure 0125 % structure with fields 'LM', 'MaxIter', 'TolParam', 'TolCost', 'TraceOn' 0126 % IterVar = struct('LM', [], 'MaxIter', [], 'TolParam', [], 'TolCost', [], 'TraceOn', []) 0127 % IterVar.LM = if 1 then the iteration stops if the corresponding ML cost 0128 % function increases 0129 % IterVar.MaxIter = maximum number of itterations of the minimization procedure 0130 % IterVar.TolParam = relative precision on parameter estimates 0131 % IterVar.TolCost = relative precision on cost function 0132 % IterVar.TraceOn = 1 then output iterations (optional) 0133 % 0134 % 0135 % Copyright (c) Rik Pintelon, Vrije Universiteit Brussel - dept. ELEC, 30 November 2009 0136 % All rights reserved. 0137 % Software can be used freely for non-commercial applications only. 0138 % Version 24 October 2011 0139 % 0140 0141 0142 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% 0143 % initialisation of the variables, and compatibility check of the input arguments % 0144 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% 0145 0146 % model structure variables 0147 ModelVar.PlantPlane = lower(ModelVar.PlantPlane); 0148 ModelVar.Struct = 'EIV'; 0149 0150 % add number of inputs and outputs to ModelVar 0151 ModelVar.ny = size(data.Y,1); 0152 ModelVar.nu = size(data.U,1); 0153 0154 % 1. imposes the compatibility of the parameter vector Theta0 and 0155 % the free model parameters with the model structure 0156 % 2. puts the order of the polynomials in ModelVar 0157 [Theta0, Sel, ModelVar] = MIMO_ML_ModelCompatibility(Theta0, Sel, ModelVar); 0158 0159 % check if DC and Nyquist belong to the frequency set 0160 if data.freq(1) == 0, data.DC = 1; else data.DC = 0; end 0161 if data.freq(end) == 1/(2*data.Ts), data.Nyquist = 1; else data.Nyquist = 0; end 0162 data.freq = data.freq(:); 0163 0164 % determine the number of MIMO experiments 0165 NumberDim = length(size(data.U)); % number of matrix dimensions 0166 if NumberDim == 2 0167 data.NumberExp = 1; % number of MIMO experiments 0168 elseif NumberDim == 3 0169 data.NumberExp = size(data.U, 2); % number of MIMO experiments 0170 end % if 0171 0172 % jw, sqrt(jw), or exp(-jwTs) values for plant and noise model 0173 x = struct('Plant', []); 0174 % matrix of powers of x 0175 xMat = struct('Plant', []); 0176 0177 0178 %%%%%%%%%%%%%%%%%%%%%%%%%%%%% 0179 % domain of the plant model % 0180 %%%%%%%%%%%%%%%%%%%%%%%%%%%%% 0181 0182 % the vector s represents z^-1, s, or sqrt(s) of the plant model 0183 switch ModelVar.PlantPlane 0184 case {'s','w'} 0185 if ModelVar.PlantPlane == 's' 0186 x.Plant = sqrt(-1)*2*pi*data.freq; 0187 elseif ModelVar.PlantPlane == 'w' 0188 x.Plant = (sqrt(-1)*2*pi*data.freq).^(0.5); 0189 end; 0190 wscale = median(abs(x.Plant)); 0191 x.Plant = x.Plant/wscale; 0192 case 'z' 0193 x.Plant = exp(-sqrt(-1)*2*pi*data.freq*data.Ts); 0194 wscale = 1; 0195 case '' 0196 x.Plant = ones(size(data.freq)); 0197 wscale = 1; 0198 otherwise, disp('Invalid plant plane ...'), return 0199 end 0200 0201 nmax = max([ModelVar.na, ModelVar.nb, ModelVar.nig]); 0202 xMat.Plant = MIMO_ML_CalcOmegaMat(x.Plant, nmax); 0203 0204 % normalize the starting values for s-, and sqrt(s)-domains 0205 Theta0 = MIMO_ML_Normalise(Theta0, wscale, ModelVar); 0206 0207 0208 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% 0209 % maximum likelihood estimate of the model parameters % 0210 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% 0211 0212 % cost function starting values 0213 Cost = MIMO_ML_Cost(data, Theta0, x, ModelVar); 0214 0215 % initialisation itertation parameters 0216 MaxRelErr = NaN; 0217 RelVarCost = NaN; 0218 DeltaCost = NaN; 0219 OldTheta = Theta0; 0220 OldCost = Cost; 0221 OldMaxRelErr = MaxRelErr; 0222 OldRelVarCost = RelVarCost; 0223 OldDeltaCost = DeltaCost; 0224 0225 if IterVar.TraceOn == 1 0226 home 0227 disp(['BTLS iteration # ',int2str(0)]); 0228 fprintf('condition number: NaN \n'); 0229 fprintf('maximum relative variation parameters: %e \n',MaxRelErr); 0230 fprintf('cost function (delta cost): %e (%e ) \n',Cost,DeltaCost); 0231 end 0232 0233 % starting values only 0234 if IterVar.MaxIter < 1 0235 Theta = OldTheta; 0236 end 0237 0238 0239 iii = 0; 0240 while iii <= IterVar.MaxIter 0241 0242 iii = iii + 1; 0243 0244 % parametric transfer function estimate 0245 PolyTrans = MIMO_ML_CalcPolyTrans(OldTheta, x); 0246 0247 % weighting linear least squares cost function 0248 data.W = MIMO_Calc_IQML_Weight(data, PolyTrans); 0249 0250 % weighted linear least squares estimate 0251 [Theta, smax, smin, wscale] = MIMO_WGTLS(data, Sel, ModelVar); 0252 0253 % normalise the parameters 0254 Theta = MIMO_ML_Normalise(Theta, wscale, ModelVar); 0255 0256 Cost = MIMO_ML_Cost(data, Theta, x, ModelVar); 0257 0258 if ~(IterVar.LM && (Cost > OldCost)) 0259 0260 % estimated parameters polynomial A 0261 A = Theta.A(Sel.A == 1); 0262 OldA = OldTheta.A(Sel.A == 1); 0263 MaxRelErr = max([abs((A - OldA)./A)]); 0264 0265 % estimated parameters ny x nu matrix polynomial B 0266 for ii = 1:ModelVar.ny 0267 for jj = 1:ModelVar.nu 0268 B = squeeze(Theta.B(ii, jj, Sel.B(ii,jj,:) == 1)).'; 0269 OldB = squeeze(OldTheta.B(ii, jj, Sel.B(ii,jj,:) == 1)).'; 0270 MaxRelErr = max([MaxRelErr, abs((B - OldB)./B)]); 0271 end % jj 0272 end % ii 0273 0274 % estimated parameters ny x 1 vector polynomial Ig 0275 for ii = 1:ModelVar.ny 0276 Ig = Theta.Ig(ii, Sel.Ig(ii,:) == 1); 0277 OldIg = OldTheta.Ig(ii, Sel.Ig(ii,:) == 1); 0278 MaxRelErr = max([MaxRelErr, abs((Ig - OldIg)./Ig)]); 0279 end % ii 0280 0281 DeltaCost = Cost - OldCost; 0282 RelVarCost = abs(DeltaCost)/Cost; 0283 OldTheta = Theta; 0284 OldCost = Cost; 0285 OldMaxRelErr = MaxRelErr; 0286 OldRelVarCost = RelVarCost; 0287 OldDeltaCost = DeltaCost; 0288 0289 end 0290 0291 % restore the previous estimates and stop the iteration 0292 if IterVar.LM && (Cost > OldCost) 0293 0294 Theta = OldTheta; 0295 Cost = OldCost; 0296 MaxRelErr = OldMaxRelErr; 0297 RelVarCost = OldRelVarCost; 0298 DeltaCost = OldDeltaCost; 0299 iii = IterVar.MaxIter + 1; 0300 0301 end 0302 0303 if (IterVar.TraceOn && (iii ~= IterVar.MaxIter + 1)) 0304 0305 home 0306 disp(['BTLS iteration # ',int2str(iii)]); 0307 disp(['condition number:',num2str(smax/smin)]); 0308 disp(['maximum relative variation parameters: ',num2str(MaxRelErr)]); 0309 disp(['cost function (delta cost): ',num2str(Cost),' (',num2str(DeltaCost),')']); 0310 0311 end % if trace on 0312 0313 if (MaxRelErr < IterVar.TolParam) | (RelVarCost < IterVar.TolCost) 0314 iii = IterVar.MaxIter + 1; 0315 end % if relative variation parameters or cost function is small enough 0316 0317 end % while main iteration loop 0318 0319 0320 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% 0321 % denormalize the estimated parameters for s-, and sqrt(s)-domains % 0322 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% 0323 0324 Theta = MIMO_ML_DeNormalise(Theta, wscale, ModelVar);