comparison m-toolbox/classes/@ao/sDomainFit.m @ 0:f0afece42f48

Import.
author Daniele Nicolodi <nicolodi@science.unitn.it>
date Wed, 23 Nov 2011 19:22:13 +0100
parents
children
comparison
equal deleted inserted replaced
-1:000000000000 0:f0afece42f48
1 % sDomainFit performs a fitting loop to identify model order and
2 % parameters.
3 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
4 %
5 % DESCRIPTION: sDomainFit fit a partial fraction model to frequency
6 % response data using the function utils.math.vcfit.
7 %
8 % The function performs a fitting loop to automatically identify model
9 % order and parameters in s-domain. Output is a s-domain model expanded
10 % in partial fractions:
11 %
12 % r1 rN
13 % f(s) = ------- + ... + ------- + d
14 % s - p1 s - pN
15 %
16 % The function attempt to perform first the identification of a model
17 % with d = 0, then if the operation do not succeed, it try the
18 % identification with d different from zero.
19 % % Identification loop stop when the stop condition is reached.
20 % Stop criterion is based on three different approachs:
21 %
22 % 1) Mean Squared Error and variation
23 % Check if the normalized mean squared error is lower than the value specified in
24 % FITTOL and if the relative variation of the mean squared error is lower
25 % than the value specified in MSEVARTOL.
26 % E.g. FITTOL = 1e-3, MSEVARTOL = 1e-2 search for a fit with
27 % normalized magnitude error lower than 1e-3 and and MSE relative
28 % variation lower than 1e-2.
29 %
30 % 1) Log residuals difference and root mean squared error
31 % Log Residuals difference
32 % Check if the minimum of the logarithmic difference between data and
33 % residuals is larger than a specified value. ie. if the conditioning
34 % value is 2, the function ensures that the difference between data and
35 % residuals is at lest 2 order of magnitude lower than data itsleves.
36 % Root Mean Squared Error
37 % Check that the variation of the root mean squared error is lower than
38 % 10^(-1*value).
39 %
40 % 2) Residuals spectral flatness and root mean squared error
41 % Residuals Spectral Flatness
42 % In case of a fit on noisy data, the residuals from a good fit are
43 % expected to be as much as possible similar to a white noise. This
44 % property can be used to test the accuracy of a fit procedure. In
45 % particular it can be tested that the spectral flatness coefficient of
46 % the residuals is larger than a certain qiantity sf such that 0<sf<1.
47 % Root Mean Squared Error
48 % Check that the variation of the root mean squared error is lower than
49 % 10^(-1*value).
50 %
51 % Both in the first and second approach the fitting loop stops when the
52 % two stopping conditions are satisfied.
53 % The output are AOs containing the frequency response of the fitted
54 % model, while the Model parameters are output as a parfrac model
55 % in the output AOs procinfo filed.
56 %
57 % The function can also perform a single loop without taking care of
58 % the stop conditions. This happens when 'AutoSearch' parameter is
59 % setted to 'off'.
60 %
61 % If you provide more than one AO as input, they will be fitted
62 % together with a common set of poles.
63 %
64 % CALL: mod = sDomainFit(a, pl)
65 %
66 % INPUTS: a - input AOs to fit to. If you provide more than one AO as
67 % input, they will be fitted together with a common set
68 % of poles. Only frequency domain (fsdata) data can be
69 % fitted. Each non fsdata object will be ignored. Input
70 % objects must have the same number of elements.
71 % pl - parameter list (see below)
72 %
73 % OUTPUTS:
74 % mod - matrix of one parfrac object for each input AO.
75 % Usseful fit information are stored in the procinfoi
76 % field:
77 % FIT_RESP - model frequency response.
78 % FIT_RESIDUALS - analysis object containing the fit
79 % residuals.
80 % FIT_MSE - analysis object containing the mean squared
81 % error progression during the fitting loop.
82 %
83 %
84 % Note: all the input objects are assumed to caontain the same X
85 % (frequencies) values
86 %
87 %
88 % EXAMPLES:
89 %
90 % 1) Fit to a frequency-series using Mean Squared Error and variation stop
91 % criterion
92 %
93 % % Create a frequency-series AO
94 % pl_data = plist('fsfcn', '0.01./(0.0001+f)', 'f1', 1e-5, 'f2', 5, 'nf', 1000);
95 % a = ao(pl_data);
96 %
97 % % Fitting parameter list
98 % pl_fit = plist('AutoSearch','on',...
99 % 'StartPoles',[],...
100 % 'StartPolesOpt','clog',...
101 % 'maxiter',5,...
102 % 'minorder',2,...
103 % 'maxorder',20,...
104 % 'weights',[],...
105 % 'CONDTYPE','MSE',...
106 % 'FITTOL',1e-3,...
107 % 'MSEVARTOL',1e-2,...
108 % 'Plot','off',...
109 % 'ForceStability','off',...
110 % 'direct term','off',...
111 % 'CheckProgress','off');
112 %
113 % % Do fit
114 % b = sDomainFit(a, pl_fit);
115 %
116 % 2) Fit to a frequency-series using Log residuals difference and mean
117 % squared error variation stop criterion
118 %
119 % % Create a frequency-series AO
120 % pl_data = plist('fsfcn', '0.01./(0.0001+f)', 'f1', 1e-5, 'f2', 5, 'nf', 1000);
121 % a = ao(pl_data);
122 %
123 % % Fitting parameter list
124 % pl_fit = plist('FS',[],...
125 % 'AutoSearch','on',...
126 % 'StartPoles',[],...
127 % 'StartPolesOpt','clog',...
128 % 'maxiter',5,...
129 % 'minorder',2,...
130 % 'maxorder',20,...
131 % 'weights',[],...
132 % 'weightparam','abs',...
133 % 'CONDTYPE','RLD',...
134 % 'FITTOL',1e-3,...
135 % 'MSEVARTOL',1e-2,...
136 % 'Plot','off',...
137 % 'ForceStability','off',...
138 % 'CheckProgress','off');
139 %
140 % % Do fit
141 % b = sDomainFit(a, pl_fit);
142 %
143 % 3) Fit to a frequency-series using Residuals spectral flatness and mean
144 % squared error variation stop criterion
145 %
146 % % Create a frequency-series AO
147 % pl_data = plist('fsfcn', '0.01./(0.0001+f)', 'f1', 1e-5, 'f2', 5, 'nf', 1000);
148 % a = ao(pl_data);
149 %
150 % % Fitting parameter list
151 % pl_fit = plist('FS',[],...
152 % 'AutoSearch','on',...
153 % 'StartPoles',[],...
154 % 'StartPolesOpt','clog',...
155 % 'maxiter',5,...
156 % 'minorder',2,...
157 % 'maxorder',20,...
158 % 'weights',[],...
159 % 'weightparam','abs',...
160 % 'CONDTYPE','RSF',...
161 % 'FITTOL',0.5,...
162 % 'MSEVARTOL',1e-2,...
163 % 'Plot','off',...
164 % 'ForceStability','off',...
165 % 'CheckProgress','off');
166 %
167 % % Do fit
168 % b = sDomainFit(a, pl_fit);
169 %
170 %
171 % <a href="matlab:utils.helper.displayMethodInfo('ao', 'sDomainFit')">Parameters Description</a>
172 %
173 % VERSION: $Id: sDomainFit.m,v 1.32 2011/08/15 09:46:44 hewitson Exp $
174 %
175 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
176
177 function varargout = sDomainFit(varargin)
178
179 % Check if this is a call for parameters
180 if utils.helper.isinfocall(varargin{:})
181 varargout{1} = getInfo(varargin{3});
182 return
183 end
184
185 import utils.const.*
186 utils.helper.msg(msg.PROC3, 'running %s/%s', mfilename('class'), mfilename);
187
188 % Collect input variable names
189 in_names = cell(size(varargin));
190 for ii = 1:nargin,in_names{ii} = inputname(ii);end
191
192 % Collect all AOs and plists
193 [as, ao_invars] = utils.helper.collect_objects(varargin(:), 'ao', in_names);
194 pl = utils.helper.collect_objects(varargin(:), 'plist', in_names);
195
196 if nargout == 0
197 error('### sDomainFit cannot be used as a modifier. Please give an output variable.');
198 end
199
200 %%% Decide on a deep copy or a modify
201 bs = copy(as, nargout);
202 inhists = [as.hist];
203
204 % combine plists
205 pl = parse(pl, getDefaultPlist());
206
207 %%%%% Extract necessary parameters %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
208
209 extpoles = find(pl, 'StartPoles'); % Check if external poles are providied
210 spolesopt = 0;
211 if isempty(extpoles) % if no external poles set them internally
212 splopt = find(pl, 'StartPolesOpt');
213 switch lower(splopt)
214 case 'real'
215 spolesopt = 1;
216 case 'clog'
217 spolesopt = 2;
218 case 'clin'
219 spolesopt = 3;
220 end
221 end
222
223 maxiter = find(pl, 'maxiter'); % set the maximum number of iterations
224 minorder = find(pl, 'minorder'); % set the minimum function order
225 maxoredr = find(pl, 'maxorder');% set the maximum function order
226
227 extweights = find(pl, 'weights'); % check if external weights are provided
228 weightparam = 0;
229 if isempty(extweights) % set internally the weights on the basis of the input options
230 wtparam = find(pl, 'weightparam');
231 switch lower(wtparam)
232 case 'ones'
233 weightparam = 1;
234 case 'abs'
235 weightparam = 2;
236 case 'sqrt'
237 weightparam = 3;
238 end
239 end
240
241 % decide to plot or not
242 plt = find(pl, 'plot');
243 switch lower(plt)
244 case 'on'
245 showplot = 1;
246 case 'off'
247 showplot = 0;
248 end
249
250 % Make a decision between Fit conditioning type
251 condtype = find(pl, 'CONDTYPE');
252 condtype = upper(condtype);
253 switch condtype
254 case 'MSE'
255 ctp = 'chivar'; % use normalized mean squared error value and relative variation
256 lrscond = find(pl, 'FITTOL');
257 % give an error for strange values of lrscond
258 if lrscond<0
259 error('!!! Negative values for FITTOL are not allowed !!!')
260 end
261 % handling data
262 lrscond = -1*log10(lrscond);
263 % give a warning for strange values of lrscond
264 if lrscond<0
265 warning('You are searching for a MSE lower than %s', num2str(10^(-1*lrscond)))
266 end
267 case 'RLD'
268 ctp = 'lrsmse'; % use residuals log difference and MSE relative variation
269 lrscond = find(pl, 'FITTOL');
270 % give a warning for strange values of lrscond
271 if lrscond<0
272 error('!!! Negative values for FITTOL are not allowed !!!')
273 end
274 if lrscond<1
275 warning('You are searching for a frequency by frequency residuals log difference of %s', num2str(lrscond))
276 end
277 case 'RSF'
278 ctp = 'rftmse'; % use residuals spectral flatness and MSE relative variation
279 lrscond = find(pl, 'FITTOL');
280 % give a warning for strange values of lrscond
281 if lrscond<0 || lrscond>1
282 error('!!! Values <0 or >1 for FITTOL are not allowed when CONDTYPE is RSF !!!')
283 end
284 end
285
286 % Tolerance for the MSE relative variation
287 msevar = find(pl, 'MSEVARTOL');
288 % handling data
289 msevar = -1*log10(msevar);
290 % give a warning for strange values of msevar
291 if msevar<0
292 warning('You are searching for MSE relative variation lower than %s', num2str(10^(-1*msevar)))
293 end
294
295 % decide to stabilize or not the model
296 stab = find(pl, 'ForceStability');
297 switch lower(stab)
298 case 'on'
299 stabfit = 1;
300 case 'off'
301 stabfit = 0;
302 end
303
304 % decide to fit with or whitout direct term
305 dtm = find(pl, 'direct term');
306 switch lower(dtm)
307 case 'on'
308 dterm = 1;
309 case 'off'
310 dterm = 0;
311 end
312
313 % decide to disp or not the fitting progress in matlab command window
314 prg = find(pl, 'CheckProgress');
315 switch lower(prg)
316 case 'on'
317 spy = 1;
318 case 'off'
319 spy = 0;
320 end
321
322 % decide to perform or not a full automatic model search
323 autos = find(pl, 'AutoSearch');
324 switch lower(autos)
325 case 'on'
326 fullauto = 1;
327 case 'off'
328 fullauto = 0;
329 end
330
331 % extract delay
332 delay = find(pl, 'delay');
333
334 %%%%% End Extract necessary parameters %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
335
336
337
338 %%%%% Fitting %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
339
340 % Fit parameters
341 params = struct('spolesopt',spolesopt,...
342 'extpoles', extpoles,...
343 'Nmaxiter',maxiter,...
344 'minorder',minorder,...
345 'maxorder',maxoredr,...
346 'weightparam',weightparam,...
347 'extweights', extweights,...
348 'plot',showplot,...
349 'ctp',ctp,...
350 'lrscond',lrscond,...
351 'msevar',msevar,...
352 'stabfit',stabfit,...
353 'dterm',dterm,...
354 'spy',spy,...
355 'fullauto',fullauto);
356
357 %%% extracting elements from AOs
358
359 % Finding the index of the first fsdata
360 for gg = 1:numel(bs)
361 if isa(bs(gg).data, 'fsdata')
362 prm = gg;
363 break
364 end
365 end
366
367 y = zeros(length(bs(prm).data.getY),numel(bs)); % initialize input vector
368 k = numel(bs(prm).data.getY); % setting a comparison constant
369 idx = true(numel(bs),1); % initialize the control index
370 for jj=1:numel(bs)
371 % checking that AOs are fsdata and skipping non fsdata objects
372 if ~isa(bs(jj).data, 'fsdata')
373 % skipping data if non fsdata
374 warning('!!! %s expects ao/fsdata objects. Skipping AO %s', mfilename, ao_invars{jj});
375 idx(jj) = false; % set the corresponding value of the control index to false
376 else
377 % preparing data for fit
378 yt = bs(jj).data.getY;
379 if numel(yt)~=k
380 error('Input AOs must have the same number of elements')
381 end
382 if size(yt,2)>1 % wish to work with columns
383 y(:,jj) = yt.';
384 else
385 y(:,jj) = yt;
386 end
387 end
388 end
389 %%% extracting frequencies
390 % Note: all the objects are assumed to caontain the same X (frequencies) values
391 f = bs(prm).data.getX;
392
393 % reshaping y to contain only Y from fsdata, subtract delay if given by
394 % user
395 if ~isempty(delay)
396 y = y(:,idx)./exp(-2*pi*1i*f*delay);
397 else
398 y = y(:,idx);
399 end
400
401 % Fitting loop
402 [res,poles,dterm,mresp,rdl,mse] = utils.math.autocfit(y,f,params);
403
404 %%%%% Building output AOs with model responses, model parameters are %%%%
405
406 for kk = 1:numel(bs)
407 if idx(kk) % set the corresponding Y values of fitted data
408
409 % if delay is input we return a pzmodel with the corresponding delay
410 if isempty(delay)
411 mdl(kk) = parfrac(plist('res', res(:,kk),'poles', poles, 'dir',...
412 dterm(:,kk), 'name', sprintf('fit(%s)', ao_invars{kk})));
413 else
414 mdl_aux = parfrac(plist('res', res(:,kk),'poles', poles, 'dir',...
415 dterm(:,kk), 'name', sprintf('fit(%s)', ao_invars{kk})));
416 mdl(kk) = pzmodel(mdl_aux);
417 mdl(kk).setDelay(delay);
418 end
419
420 % Output also response, residuals and mse progression in the procinfo
421
422 rsp = mresp(:,kk);
423 bs(kk).data.setY(rsp);
424
425 % Set output AO name
426 bs(kk).name = sprintf('fit(%s)', ao_invars{kk});
427
428 res_ao = copy(bs(kk),1);
429 trdl = rdl(:,kk);
430 res_ao.data.setY(trdl);
431
432 % Set output AO name
433 res_ao.name = sprintf('fit_residuals(%s)', ao_invars{kk});
434
435 d = cdata();
436 tmse = mse(:,kk);
437 d.setY(tmse);
438 mse_ao = ao(d);
439
440 % Set output AO name
441 mse_ao.name = sprintf('fit_mse(%s)', ao_invars{kk});
442
443 procpl = plist('fit_resp',bs(kk),...
444 'fit_residuals',res_ao,...
445 'fit_mse',mse_ao);
446
447 mdl(kk).setProcinfo(procpl);
448
449 else
450 mdl(kk) = parfrac();
451 end
452
453 end
454
455 % set output as matrix if multiple inputs
456 if numel(mdl) ~= 1
457 mmdl = matrix(mdl);
458 else
459 mmdl = mdl;
460 end
461
462 mmdl.setName(sprintf('fit(%s)', ao_invars{:}));
463
464 mmdl.addHistory(getInfo('None'), pl, [ao_invars(:)], [inhists(:)]);
465
466 % ----- Set outputs -----
467 if nargout == 1
468 varargout{1} = mmdl;
469 else
470 % multiple output is not supported
471 error('### Multiple output is not supported ###')
472 end
473
474 end
475
476 %--------------------------------------------------------------------------
477 % Get Info Object
478 %--------------------------------------------------------------------------
479 function ii = getInfo(varargin)
480 if nargin == 1 && strcmpi(varargin{1}, 'None')
481 sets = {};
482 pl = [];
483 else
484 sets = {'Default'};
485 pl = getDefaultPlist;
486 end
487 % Build info object
488 ii = minfo(mfilename, 'ao', 'ltpda', utils.const.categories.sigproc, '$Id: sDomainFit.m,v 1.32 2011/08/15 09:46:44 hewitson Exp $', sets, pl);
489 ii.setModifier(false);
490 end
491
492 %--------------------------------------------------------------------------
493 % Get Default Plist
494 %--------------------------------------------------------------------------
495 function plout = getDefaultPlist()
496 persistent pl;
497 if exist('pl', 'var')==0 || isempty(pl)
498 pl = buildplist();
499 end
500 plout = pl;
501 end
502
503 function pl = buildplist()
504
505 pl = plist();
506
507 % AutoSearch
508 p = param({'AutoSearch', ['''on'': Parform a full automatic search for the<br>'...
509 'transfer function order. The fitting<br>'...
510 'procedure will stop when stop conditions<br>'...
511 'defined are satisfied.<br>'...
512 '''off'': Perform a fitting loop as long as the<br>'...
513 'number of iteration reach ''maxiter''. The order<br>'...
514 'of the fitting function will be that<br>'...
515 'specified in ''minorder''.']}, ...
516 {1, {'on', 'off'}, paramValue.SINGLE});
517 pl.append(p);
518
519 % StartPoles
520 p = param({'StartPoles', ['A vector of starting poles. Providing a fixed<br>'...
521 'set of starting poles fixes the function<br>'...
522 'order. If it is left empty starting poles are<br>'...
523 'internally assigned.']}, paramValue.EMPTY_DOUBLE);
524 pl.append(p);
525
526 % StartPolesOpt
527 p = param({'StartPolesOpt', ['Define the characteristics of internally<br>'...
528 'assigned starting poles. Admitted values<br>'...
529 'are:<ul>'...
530 '<li>''real'' linear-spaced real poles</li>'...
531 '<li>''clog'' log-spaced complex poles</li>'...
532 '<li>''clin'' linear-spaced complex poles<li></ul>']}, ...
533 {2, {'real', 'clog', 'clin'}, paramValue.SINGLE});
534 pl.append(p);
535
536 % MaxIter
537 p = param({'MaxIter', 'Maximum number of iterations in fit routine.'}, paramValue.DOUBLE_VALUE(50));
538 pl.append(p);
539
540 % MinOrder
541 p = param({'MinOrder', 'Minimum order to fit with.'}, paramValue.DOUBLE_VALUE(2));
542 pl.append(p);
543
544 % MaxOrder
545 p = param({'MaxOrder', 'Maximum order to fit with.'}, paramValue.DOUBLE_VALUE(20));
546 pl.append(p);
547
548 % Weights
549 p = param({'Weights', ['A vector with the desired weights. If a single<br>'...
550 'Ao is input weights must be a Nx1 vector where<br>'...
551 'N is the number of elements in the input Ao. If<br>'...
552 'M Aos are passed as input, then weights must<br>'...
553 'be a NxM matrix. If it is leaved empty weights<br>'...
554 'are internally assigned basing on the input<br>'...
555 'parameters']}, paramValue.EMPTY_DOUBLE);
556 pl.append(p);
557
558 % Weightparam
559 p = param({'weightparam', ['Specify the characteristics of the internally<br>'...
560 'assigned weights. Admitted values are:<ul>'...
561 '<li>''ones'' assigns weights equal to 1 to all data.<li>'...
562 '<li>''abs'' weights data with <tt>1./abs(y)</tt></li>'...
563 '<li>''sqrt'' weights data with <tt>1./sqrt(abs(y))</tt></li>']}, ...
564 {2, {'ones', 'abs', 'sqrt'}, paramValue.SINGLE});
565 pl.append(p);
566
567 % CONDTYPE
568 p = param({'CONDTYPE', ['Fit conditioning type. Admitted values are:<ul>'...
569 '<li>''MSE'' Mean Squared Error and variation</li>'...
570 '<li>''RLD'' Log residuals difference and mean squared error variation<li>'...
571 '<li>''RSF'' Residuals spectral flatness and mean squared error variation<li></ul>']}, ...
572 {1, {'MSE', 'RLD', 'RSF'}, paramValue.SINGLE});
573 pl.append(p);
574
575 % FITTOL
576 p = param({'FITTOL', 'Fit tolerance.'}, paramValue.DOUBLE_VALUE(1e-3));
577 pl.append(p);
578
579 % MSEVARTOL
580 p = param({'MSEVARTOL', ['Mean Squared Error Variation - Check if the<br>'...
581 'realtive variation of the mean squared error is<br>'...
582 'smaller than the value specified. This<br>'...
583 'option is useful for finding the minimum of Chi-squared.']}, ...
584 paramValue.DOUBLE_VALUE(1e-2));
585 pl.append(p);
586
587 % Plot
588 p = param({'Plot', 'Plot results of each fitting step.'}, {2, {'on', 'off'}, paramValue.SINGLE});
589 p.val.setValIndex(2);
590 pl.append(p);
591
592 % ForceStability
593 p = param({'ForceStability', 'Force poles to be stable'}, ...
594 {2, {'on', 'off'}, paramValue.SINGLE});
595 pl.append(p);
596
597 % direct term
598 p = param({'direct term', 'Fit with direct term.'}, {2, {'on', 'off'}, paramValue.SINGLE});
599 pl.append(p);
600
601 % CheckProgress
602 p = param({'CheckProgress', 'Display the status of the fit iteration.'}, ...
603 {2, {'on', 'off'}, paramValue.SINGLE});
604 pl.append(p);
605
606 % Delay
607 p = param({'delay', 'Innput a delay that will be subtracted from the fit.<br>'...
608 'The output is a pzmodel which includes the inputted delay.'},paramValue.EMPTY_DOUBLE);
609 pl.append(p);
610 end
611 % END
612
613
614 % PARAMETERS:
615 % 'AutoSearch' - 'on': Parform a full automatic search for the
616 % transfer function order. The fitting
617 % procedure will stop when stop conditions
618 % defined are satisfied. [Default]
619 % 'off': Perform a fitting loop as long as the
620 % number of iteration reach 'maxiter'. The order
621 % of the fitting function will be that
622 % specified in 'minorder'.
623 % 'StartPoles' - A vector of starting poles. Providing a fixed
624 % set of starting poles fixes the function
625 % order. If it is left empty starting poles are
626 % internally assigned. [Default []]
627 % 'StartPolesOpt' - Define the characteristics of internally
628 % assigned starting poles. Admitted values
629 % are:
630 % 'real' linspaced real poles
631 % 'clog' logspaced complex poles [Default]
632 % 'clin' linspaced complex poles
633 % 'maxiter' - Maximum number of allowed iteration. [Deafult
634 % 50].
635 % 'minorder' - Minimum model function order. [Default 2]
636 % 'maxorder' - Maximum model function order. [Default 20]
637 % 'weights' - A vector with the desired weights. If a single
638 % Ao is input weights must be a Nx1 vector where
639 % N is the number of elements in the input Ao. If
640 % M Aos are passed as input, then weights must
641 % be a NxM matrix. If it is leaved empty weights
642 % are internally assigned basing on the input
643 % parameters. [Default []]
644 % 'weightparam' - Specify the characteristics of the internally
645 % assigned weights. Admitted values are:
646 % 'ones' assigns weights equal to 1 to all data.
647 % 'abs' weights data with 1./abs(y) [Default]
648 % 'sqrt' weights data with 1./sqrt(abs(y))
649 % 'CONDTYPE' - Fit conditioning type. Admitted values are:
650 % - 'MSE' Mean Squared Error and variation
651 % [Default]
652 % - 'RLD' Log residuals difference and mean
653 % squared error variation
654 % - 'RSF' Residuals spectral flatness and mean
655 % squared error variation
656 % 'FITTOL' - Fit tolerance [Default, 1e-3]
657 % 'MSEVARTOL' - This allow to check if the relative variation
658 % of mean squared error is lower than the value
659 % sepcified. [Default 1e-2]
660 % 'Plot' - Plot fit result: 'on' or 'off' [default]
661 % 'ForceStability' - Force poles to be stable, values are
662 % 'on' or 'off'. [Default 'off']
663 % 'direct term' - Fit with direct term if 'on', without if
664 % 'off'. [Default 'off']
665 % 'CheckProgress' - Disply the status of the fit iteration.
666 % Values are 'on and 'off'. [Default 'off']
667 %
668 %
669