Model estimation

The estimation of the measurand in the analysis of dynamic measurements typically corresponds to a deconvolution problem. Therefore, a digital filter can be designed whose input is the measured system output signal and whose output is an estimate of the measurand. The package Model estimation implements methods for the design of such filters given an array of frequency response values or the reciprocal of frequency response values with associated uncertainties for the measurement system.

The package Model estimation also contains a function for the identification of transfer function models.

The package consists of the following modules:

Fitting filters to frequency response or reciprocal

This module assists in carrying out least-squares IIR and FIR filter fits

It is possible to carry out a least-squares fit of digital, time-discrete IIR and FIR filters to a given complex frequency response and the design of digital deconvolution filters by least-squares fitting to the reciprocal of a given frequency response each with propagation of associated uncertainties.

This module contains the following functions:

  • LSIIR(): Least-squares (time-discrete) IIR filter fit to a given frequency response or its reciprocal optionally propagating uncertainties.

  • LSFIR(): Least-squares (time-discrete) FIR filter fit to a given frequency response or its reciprocal optionally propagating uncertainties either via Monte Carlo or via a singular-value decomposition and linear matrix propagation.

PyDynamic.model_estimation.fit_filter.LSFIR(H: ndarray, N: int, f: ndarray, Fs: float, tau: int, weights: Optional[ndarray] = None, verbose: Optional[bool] = True, inv: Optional[bool] = False, UH: Optional[ndarray] = None, mc_runs: Optional[int] = None, trunc_svd_tol: Optional[float] = None) Tuple[ndarray, ndarray][source]

Design of FIR filter as fit to freq. resp. or its reciprocal with uncertainties

Least-squares fit of a (time-discrete) digital FIR filter to the reciprocal of the frequency response values or actual frequency response values for which associated uncertainties are given for its real and imaginary part. Uncertainties are propagated either using a Monte Carlo method if mc_runs is provided as integer greater than one or otherwise using a truncated singular-value decomposition and linear matrix propagation. The Monte Carlo approach may help in cases where the weighting matrix or the Jacobian are ill-conditioned, resulting in false uncertainties associated with the filter coefficients.

Note

Uncertainty propagation via singular-value decomposition is not yet implemented, when fitting to the actual frequency response and not its reciprocal. Alternatively specify the number mc_runs of runs to propagate the uncertainties via the Monte Carlo method.

Parameters:
  • H (array_like of shape (M,) or (2M,)) – (Complex) frequency response values in dtype complex or as a vector containing the real parts in the first half followed by the imaginary parts

  • N (int) – FIR filter order

  • f (array_like of shape (M,)) – frequencies at which H is given

  • Fs (float) – sampling frequency of digital FIR filter

  • tau (int) – time delay in samples for improved fitting

  • weights (array_like of shape (2M,), optional) – vector of weights for a weighted least-squares method (default results in no weighting)

  • verbose (bool, optional) – If True (default) verbose output is printed to the command line

  • inv (bool, optional) – If False (default) apply the fit to the frequency response values directly, otherwise fit to the reciprocal of the frequency response values

  • UH (array_like of shape (2M,2M), optional) – uncertainties associated with the real and imaginary part of H

  • mc_runs (int, optional) – Number of Monte Carlo runs greater than one. Only used, if uncertainties associated with the real and imaginary part of H are provided. Only one of mc_runs and trunc_svd_tol can be provided.

  • trunc_svd_tol (float, optional) – Lower bound for singular values to be considered for pseudo-inverse. Values smaller than this threshold are considered zero. Defaults to zero. Only one of mc_runs and trunc_svd_tol can be provided.

Returns:

  • b (array_like of shape (N+1,)) – The FIR filter coefficient vector in a 1-D sequence

  • Ub (array_like of shape (N+1,N+1)) – Uncertainties associated with b. Will be None if UH is not provided or is None.

Raises:

NotImplementedError – The least-squares fitting of a digital FIR filter to a frequency response H with propagation of associated uncertainties using a truncated singular-value decomposition and linear matrix propagation is not yet implemented. Alternatively specify the number mc_runs of runs to propagate the uncertainties via the Monte Carlo method.

References

PyDynamic.model_estimation.fit_filter.LSIIR(H: ndarray, Nb: int, Na: int, f: ndarray, Fs: float, tau: Optional[int] = 0, verbose: Optional[bool] = True, return_rms: Optional[bool] = False, max_stab_iter: Optional[int] = 50, inv: Optional[bool] = False, UH: Optional[ndarray] = None, mc_runs: Optional[int] = 1000) Union[Tuple[ndarray, ndarray, int, Optional[ndarray], float], Tuple[ndarray, ndarray, int, Optional[ndarray]]][source]

Least-squares (time-discrete) IIR filter fit to frequency response or reciprocal

For fitting an IIR filter model to the reciprocal of the frequency response values or directly to the frequency response values provided by the user, this method uses a least-squares fit to determine an estimate of the filter coefficients. The filter then optionally is stabilized by pole mapping and introduction of a time delay. Associated uncertainties are optionally propagated when provided using the GUM S2 Monte Carlo method.

Parameters:
  • H (array_like of shape (M,)) – (Complex) frequency response values.

  • Nb (int) – Order of IIR numerator polynomial.

  • Na (int) – Order of IIR denominator polynomial.

  • f (array_like of shape (M,)) – Frequencies at which H is given.

  • Fs (float) – Sampling frequency for digital IIR filter.

  • tau (int, optional) – Initial estimate of time delay for obtaining a stable filter (default = 0).

  • verbose (bool, optional) – If True (default) be more talkative on stdout. Otherwise no output is written anywhere.

  • return_rms (bool, optional) – If True (default is False), the root-mean-square (rms) value of the fit is returned as an additional output value.

  • max_stab_iter (int, optional) – Maximum count of iterations for stabilizing the resulting filter. If no stabilization should be carried out, this parameter can be set to 0 (default = 50). This parameter replaced the previous justFit which was dropped in PyDynamic 2.0.0.

  • inv (bool, optional) – If False (default) apply the fit to the frequency response values directly, otherwise fit to the reciprocal of the frequency response values.

  • UH (array_like of shape (2M, 2M), optional) – Uncertainties associated with real and imaginary part of H.

  • mc_runs (int, optional) – Number of Monte Carlo runs (default = 1000). Only used if uncertainties UH are provided.

Returns:

  • b (np.ndarray) – The IIR filter numerator coefficient vector in a 1-D sequence.

  • a (np.ndarray) – The IIR filter denominator coefficient vector in a 1-D sequence.

  • tau (int) – Filter time delay (in samples).

  • Uab (np.ndarray of shape (Nb+Na+1, Nb+Na+1)) – Uncertainties associated with [a[1:],b]. Will be None if UH is not provided or is None.

  • rms (float) – The root-mean-square error of the fit. Only returned, if return_rms == True.

References

Identification of transfer function models

This module contains a function for the identification of transfer function models:

  • fit_som(): Fit second-order model to complex-valued frequency response

PyDynamic.model_estimation.fit_transfer.fit_som(f: ndarray, H: ndarray, UH: Optional[ndarray] = None, weighting: Optional[ndarray] = None, MCruns: Optional[int] = 10000, scaling: Optional[float] = 0.001, verbose: Optional[bool] = False)[source]

Fit second-order model to complex-valued frequency response

Fit second-order model (spring-damper model) with parameters \(S_0, delta\) and \(f_0\) to complex-valued frequency response with uncertainty associated with real and imaginary parts.

For a transformation of an uncertainty associated with amplitude and phase to one associated with real and imaginary parts, see PyDynamic.uncertainty.propagate_DFT.AmpPhase2DFT.

Parameters:
  • f ((M,) np.ndarray) – vector of frequencies

  • H ((2M,) np.ndarray) – real and imaginary parts of measured frequency response values at frequencies f

  • UH ((2M,) or (2M,2M) np.ndarray, optional) – uncertainties associated with real and imaginary parts When UH is one-dimensional, it is assumed to contain standard uncertainties; otherwise it is taken as covariance matrix. When UH is not specified no uncertainties assoc. with the fit are calculated, which is the default behaviour.

  • weighting (str or (2M,) np.ndarray, optional) – Type of weighting associated with frequency responses, can be (‘diag’, ‘cov’) if UH is given, or Numpy array of weights, defaults to None, which means all values are considered equally important

  • MCruns (int, optional) – Number of Monte Carlo trials for propagation of uncertainties, defaults to 10000. When MCruns is set to ‘None’, matrix multiplication is used for the propagation of uncertainties. However, in some cases this can cause trouble.

  • scaling (float, optional) – scaling of least-squares design matrix for improved fit quality, defaults to 1e-3

  • verbose (bool, optional) – if True a progressbar will be printed to console during the Monte Carlo simulations, if False nothing will be printed out, defaults to False

Returns:

  • p (np.ndarray) – vector of estimated model parameters [S0, delta, f0]

  • Up (np.ndarray) – covariance associated with parameter estimate