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