Performance evaluation of FDI filters

  • fditspec Computation of the weak or strong structure matrix.
  • fdisspec Computation of the strong structure matrix.
  • fdiscond Computation of the fault detection sensitivity condition.
  • fdif2ngap Computation of the fault-to-noise gap.
  • fdimmperf Computation of the model-matching performace.
FaultDetectionTools.fditspecFunction
S = fditspec(sysr::FDFilterIF; FDfreq = missing, block = false, poleshift = false, 
             FDtol, FDStol, atol = 0, atol1 = atol, atol2 = atol, rtol, fast = true)

Compute the weak or strong binary structure matrix S of the transfer function matrix Rf(λ) of the transfer channel from the fault inputs to residuals of a fault detection filter internal form object sysr::FDFilterIF. For a filter sysr with q residual outputs and mf fault inputs, Rf(λ) is the q x mf transfer function matrix of the fault inputs channel with the descriptor system representation sysr.sys[:,sysr.faults] := (Af-lambda*Ef,Bf,Cf,Df).

If FDfreq = missing (default), then S contains the weak structure matrix of Rf(λ). For block = false, S is determined as a q x mf binary matrix (BitMatrix), whose (i,j)-th element is S[i,j] = 1, if the (i,j)-th element of Rf(λ) is nonzero, and otherwise, S[i,j] = 0. For block = true, S is determined as a 1 x mf binary matrix, whose (1,j)-th element is S[1,j] = 1, if the j-th column of Rf(λ) is nonzero, and otherwise, S[1,j] = 0.

If FDfreq = freq specifies a vector freq of nf real frequencies which characterize the classes of persistent fault signals, then for a suitable proper and invertible M(λ) (see below), S contains the strong structure matrix of M(λ)*Rf(λ) with respect to a set of nf complex frequencies Ω, defined as follows: if f is a real frequency in freq, then the corresponding complex frequency in Ω is λ := im*f, for a continuous-time system, or λ := exp(im*f*abs(Ts)), for a discrete-time system with sampling-time Ts.

FDtol = tol1 specifies an absolute threshold tol1 for the magnitudes of nonzero elements in the system matrices Bf and Df and is used to determine the weak structure matrix. Its default value is tol1 = 0.0001*max(1, norm(Bf,1), norm(Df,1)).

FDStol = tol2 specifies an absolute threshold tol2 for the magnitudes of nonzero elements in the system matrices Af, Ef, Bf, Cf and Df and is used to determine the strong structure matrix. Its default value is tol2 = epsm*max(1, norm(Ef,1), norm(Af,1), norm(Bf,1), norm(Cf,Inf), norm(Df,1))), where epsm is the working machine precision.

For block = false, then, if poleshift = true, M(λ) is chosen diagonal such that M(λ)*Rf(λ) has no poles in Ω and if poleshift = false (default), M(λ) = I is used and an error is issued if Rf(λ) has poles in Ω. S is determined as a q x mf binary matrix, whose (i,j)-th element is S[i,j] = 1, if the (i,j)-th element of M(λ)*Rf(λ) evaluated for all frequencies in Ω is nonzero, and otherwise, S[i,j] = 0.

For block = true, then, if poleshift = true, M(λ) is chosen such that M(λ)*Rf(λ) as no poles in Ω and if poleshift = false (default), M(λ) = I is used and an error is issued if Rf(λ) has poles in Ω. S is determined as an 1 x mf binary matrix, whose (1,j)-th element is S[1,j] = 1, if the j-th column of M(λ)*Rf(λ) evaluated for all frequencies in Ω is nonzero and otherwise S[1,j] = 0.

The keyword arguments atol1, atol2, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of matrices Af, Bf, Cf, Df, the absolute tolerance for the nonzero elements of Ef, and the relative tolerance for the nonzero elements of Af, Bf, Cf, Df and Ef. The keyword argument atol can be used to simultaneously set atol1 = atol and atol2 = atol.

Method: For the definition of the structure matrix, see [1]. For the determination of the weak structure matrix, minimal realizations are determined for each column of Rf(λ) if block = true or for each element of Rf(λ) if block = false and the nonzero columns or elements in each column are identified (see Corollary 7.1 of [1]). For the determination of the strong structure matrix, minimal realizations are determined for each column of M(λ)*Rf(λ) if block = true or for each element of M(λ)*Rf(λ) if block = false and the full rank of the corresponding system matrix is checked for all frequencies in FDfreq (see Corollary 7.2 in [1]) (i.e., the lack of zeros in all frequencies).

References:

[1] Varga A. Solving Fault Diagnosis Problems - Linear Synthesis Techniques. Springer Verlag, 2017; sec.3.4.

S = fditspec(sysr::FDIFilterIF; FDfreq = missing, poleshift = false, 
             FDtol, FDStol, atol = 0, atol1 = atol, atol2 = atol, rtol, fast = true)

Compute the weak or strong binary structure matrix S of the global transfer function matrix Rf(λ) of the transfer channel from the fault inputs to residuals of a fault detection and isolation filter internal form object sysr::FDIFilterIF. The filter sysr consists of N individual FDI filters sysr.sys[i], for i = 1, ..., N, where the fault to residual channel of the i-th filter sysr.sys[i][:,sysr.faults] has qi residual outputs and mf fault inputs, has the descriptor system representation sysr.sys[i][:,sysr.faults] := (Afi-lambda*Efi,Bfi,Cfi,Dfi) and Rfi(λ) is the corresponding qi x mf transfer function matrix. The global transfer function matrix Rf(λ) is formed by row concatenation of the transfer function matrices of the N individual filters, i.e., Rf(λ) := [ Rf1(λ); Rf2(λ); ...; RfN(λ)]. For the evaluation of the strong structure matrix, the structure matrix of the stable transfer function matrix M(λ)*Rf(λ) is determined, with a M(λ) block-diagonal M(λ) = block-diag(M1(λ), M2(λ), ..., MN(λ)), where Mi(λ) is a suitable square and invertible transfer function matrix (see below).

FDtol = tol1 specifies an absolute threshold tol1 for the magnitudes of nonzero elements in the system matrices Bf and Df and is used to determine the weak structure matrix. Its default value is tol1 = 0.0001*max(1, norm(Bf,1), norm(Df,1)).

FDStol = tol2 specifies an absolute threshold tol2 for the magnitudes of nonzero elements in the system matrices Af, Ef, Bf, Cf and Df and is used to determine the strong structure matrix. Its default value is tol2 = epsm*max(1, norm(Ef,1), norm(Af,1), norm(Bf,1), norm(Cf,Inf), norm(Df,1))), where epsm is the working machine precision.

If FDfreq = missing (default), then S contains the weak structure matrix of Rf(λ). S is determined as a N x mf binary matrix, whose (i,j)-th element is S[i,j] = 1, if the j-th column of Rfi(λ) is nonzero, and otherwise, S[i,j] = 0.

If FDfreq = freq specifies a vector freq of nf real frequencies which characterize the classes of persistent fault signals, then for a suitable proper and invertible M(λ) (see below), S contains the strong structure matrix of M(λ)*Rf(λ) with respect to a set of nf complex frequencies Ω, defined as follows: if f is a real frequency in freq, then the corresponding complex frequency in Ω is λ := im*f, for a continuous-time system, or λ := exp(im*f*abs(Ts)), for a discrete-time system with sampling-time Ts.

S is determined as a N x mf binary matrix, whose (i,j)-th element is S[i,j] = 1, if the j-th column of Mi(λ)*Rfi(λ) is nonzero for all frequencies in Ω, and otherwise, S[i,j] = 0. If poleshift = true, Mi(λ) is chosen such that Mi(λ)*Rfi(λ) has no poles in Ω and if poleshift = false (default), Mi(λ) = I is used and an error is issued if any Rfi(λ) has poles in Ω.

The keyword arguments atol1, atol2, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of matrices Af, Bf, Cf, Df, the absolute tolerance for the nonzero elements of Ef, and the relative tolerance for the nonzero elements of Af, Bf, Cf, Df and Ef. The keyword argument atol can be used to simultaneously set atol1 = atol and atol2 = atol.

Method: For the definition of the structure matrix, see [1]. For the determination of the weak structure matrix, minimal realizations are determined for each column of Rfi(λ) and the nonzero columns are identified (see Corollary 7.1 of [1]). For the determination of the strong structure matrix, minimal realizations are determined for each column of Mi(λ)*Rfi(λ) and the full rank of the corresponding system matrix is checked for all frequencies in Ω (see Corollary 7.2 in [1]) (i.e., the lack of zeros in all frequencies in Ω).

References:

[1] Varga A. Solving Fault Diagnosis Problems - Linear Synthesis Techniques. Springer Verlag, 2017; sec.3.4.

FaultDetectionTools.fdisspecFunction
 S = fdisspec(sysr::FDFilterIF, freq; block = false, stabilize = false, FDGainTol = 0.01, 
                 atol = 0, atol1 = atol, atol2 = atol, rtol, fast = true)

Compute, for a given set of real frequencies freq, the strong binary structure matrix S of the stable transfer function matrix M(λ)*Rf(λ), where Rf(λ) is the transfer function matrix of the transfer channel from the fault inputs to residuals of the fault detection filter internal form object sysr::FDFilterIF and M(λ) is a suitable proper and invertible stabilizing transfer function matrix (see below). For a filter sysr with q residual outputs and mf fault inputs, Rf(λ) is the q x mf transfer function matrix of the fault inputs channel with the descriptor system representation sysr.sys[:,sysr.faults] := (Af-lambda*Ef,Bf,Cf,Df).

freq must contain a real frequency value or a vector of nf real frequencies which characterize the classes of persistent fault signals (default: freq = 0, i.e., characterizing constant faults). S contains the strong structure matrix of M(λ)*Rf(λ) with respect to a set of nf complex frequencies Ω, defined as follows: if f is a real frequency in freq, then the corresponding complex frequency in Ω is λ := im*f, for a continuous-time system, or λ := exp(im*f*abs(Ts)), for a discrete-time system with sampling-time Ts.

FDGainTol = tol specifies an absolute threshold tol for the nonzero magnitudes of the frequency response gains (default: tol = 0.01).

For block = false, then, if stabilize = true, M(λ) is chosen diagonal such that M(λ)*Rf(λ) has only stable poles and if stabilize = false (default), M(λ) = I is used and an error is issued if Rf(λ) has poles in Ω. S is determined as a q x mf binary matrix, whose (i,j)-th element is S[i,j] = 1, if the (i,j)-th element of M(λ)*Rf(λ) evaluated for all frequencies in freq is larger than or equal to tol, and otherwise, S[i,j] = 0.

For block = true, then, if stabilize = true, M(λ) is chosen such that M(λ)*Rf(λ) has only stable poles and if stabilize = false (default), M(λ) = I is used and an error is issued if Rf(λ) has poles in Ω. S is determined as an 1 x mf binary matrix, whose (1,j)-th element is S[1,j] = 1, if the j-th column of M(λ)*Rf(λ) evaluated for all frequencies in Ω is larger than or equal to tol and otherwise, S[1,j] = 0.

The keyword arguments atol1, atol2, atol3, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of matrices Af, Bf, Cf, Df, the absolute tolerance for the nonzero elements of Ef, the absolute tolerance for the nonzero elements of Cf, and the relative tolerance for the nonzero elements of Af, Bf, Cf, Df and Ef. The keyword argument atol can be used to simultaneously set atol1 = atol, atol2 = atol and atol3 = atol.

The computation of minimal realizations of individual input-output channels relies on pencil manipulation algorithms, which employ rank determinations based on either the use of rank revealing QR-decomposition with column pivoting, if fast = true, or the SVD-decomposition. The rank decision based on the SVD-decomposition is generally more reliable, but the involved computational effort is higher.

Method: S is evaluated using the definition of the strong structure matrix in [1].

References:

[1] Varga A. Solving Fault Diagnosis Problems - Linear Synthesis Techniques. Springer Verlag, 2017; sec. 3.4.

 S = fdisspec(sysr::FDIFilterIF, freq; stabilize = false, FDGainTol = 0.01, 
                 atol = 0, atol1 = atol, atol2 = atol, rtol, fast = true)

Compute, for a given set of real frequencies freq, the strong binary structure matrix S of the stable transfer function matrix M(λ)*Rf(λ), where Rf(λ) is the global transfer function matrix of the transfer channel from the fault inputs to residuals of the fault detection and isolation filter internal form object sysr::FDIFilterIF and M(λ) is a suitable block-diagonal proper and invertible stabilizing transfer function matrix (see below). The filter sysr consists of N individual FDI filters sysr.sys[i], for i = 1, ..., N, where the fault to residual channel of the i-th filter sysr.sys[i][:,sysr.faults] has qi residual outputs and mf fault inputs, has the descriptor system representation sysr.sys[i][:,sysr.faults] := (Afi-lambda*Efi,Bfi,Cfi,Dfi) and Rfi(λ) is the corresponding qi x mf transfer function matrix. The global transfer function matrix Rf(λ) is formed by row concatenation of the transfer function matrices of the N individual filters, i.e., Rf(λ) := [ Rf1(λ); Rf2(λ); ...; RfN(λ)]. M(λ) = block-diag(M1(λ), M2(λ), ..., MN(λ)), where Mi(λ) is square and invertible and chosen such that Mi(λ)Rfi(λ) is stable (see below).

freq must contain a real frequency value or a vector of nf real frequencies which characterize the classes of persistent fault signals (default: freq = 0, i.e., characterizing constant faults). S contains the strong structure matrix of M(λ)*Rf(λ) with respect to a set of nf complex frequencies Ω, defined as follows: if f is a real frequency in freq, then the corresponding complex frequency in Ω is λ := im*f, for a continuous-time system, or λ := exp(im*f*abs(Ts)), for a discrete-time system with sampling-time Ts.

FDGainTol = tol specifies an absolute threshold tol for the nonzero magnitudes of the frequency response gains (default: tol = 0.01).

If stabilize = true, Mi(λ) is chosen such that Mi(λ)*Rfi(λ) has only stable poles and if stabilize = false (default), Mi(λ) = I is used and an error is issued if any Rfi(λ) has poles in Ω.

S is determined as a N x mf binary matrix, whose (i,j)-th element is S[i,j] = 1, if the norm of the j-th column of Mi(λ)*Rfi(λ) evaluated for all frequencies in Ω is larger than or equal to tol, and otherwise, S[i,j] = 0.

The keyword arguments atol1, atol2, atol3, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of matrices Afi, Bfi, Cfi, Dfi, the absolute tolerance for the nonzero elements of Efi, the absolute tolerance for the nonzero elements of Cfi, and the relative tolerance for the nonzero elements of Afi, Bfi, Cfi, Dfi and Efi. The keyword argument atol can be used to simultaneously set atol1 = atol, atol2 = atol and atol3 = atol.

The computation of minimal realizations of individual input-output channels relies on pencil manipulation algorithms, which employ rank determinations based on either the use of rank revealing QR-decomposition with column pivoting, if fast = true, or the SVD-decomposition. The rank decision based on the SVD-decomposition is generally more reliable, but the involved computational effort is higher.

Method: S is evaluated using the definition of the strong structure matrix in [1].

References:

[1] Varga A. Solving Fault Diagnosis Problems - Linear Synthesis Techniques. Springer Verlag, 2017; sec. 3.4.

FaultDetectionTools.fdiscondFunction
 fdiscond(sysr::FDFilterIF,freq) -> (scond, β, γ)

Compute for the stable transfer function matrix Rf(λ) of the transfer channel from the fault inputs to residuals of the fault detection filter internal form object sysr::FDFilterIF the quantities: β - the H∞- index of Rf(λ), γ - the maximum of the columns norms of Rf(λ) and the fault detection sensitivity condition scond evaluated as scond := β/γ. If freq is a vector of real frequency values, then β and γ are evaluated over the frequencies contained in freq.

 fdiscond(sysr::FDFilterIF,SFDI,freq) -> (scond, β, γ)

Compute the detection and isolation sensitivity condition scond (and related quatities β and γ) for the binary structure vector SFDI associated to the stable transfer function matrix Rf(λ) of the transfer channel from the fault inputs to residuals of the fault detection filter internal form object sysr::FDFilterIF. If Rff(λ) is the transfer function matrix formed of those j-th columns of Rf(λ) for which SFDI[j] = 1, then: β - the H∞- index of Rff(λ), γ - the maximum of the columns norms of Rff(λ) and the fault detection sensitivity condition scond evaluated as scond := β/γ. If freq is a vector of real frequency values, then β and γ are evaluated over the frequencies contained in freq.

 fdiscond(sysr::FDFilterIF, SFDI, freq) -> (scond, β, γ)

Compute the detection and isolation sensitivity condition vector scond (and related quatities β and γ) for the q × mf binary structure matrix SFDI associated to the stable transfer function matrix Rf(λ) of the transfer channel from the fault inputs to residuals of the fault detection filter internal form object sysr::FDFilterIF. The filter sysr consists of N individual FDI filters sysr.sys[i], for i = 1, ..., N, where sysr.sys[i][:,sysr.faults] is the fault to residual channel of the i-th filter and Rfi(λ) is the corresponding transfer function matrix. The global transfer function matrix Rf(λ) is formed by row concatenation of the transfer function matrices of the N individual filters, i.e., Rf(λ) := [ Rf1(λ); Rf2(λ); ...; RfN(λ)]. The i-th element of the vectors scond, β and γ contain the quantities: β[i] - the H∞- index of the nonzero columns of Rffi(λ), γ - the maximum of the nonzero columns norms of Rfi(λ) and the correspomding fault detection sensitivity condition scond[i] evaluated as scond[i] := β[i]/γ[i], where Rffi(λ) is formed of those j-th columns of the i-th row of Rf(λ) for which S[i,j] = 1 and Rfi(λ) is the i-th row of Rf(λ). It is assumed that for each j such that SFDI[i,j] = 1, the j-th column of Rfi(λ) is nonzero and for each j such that SFDI[i,j] = 0, the j-th column of Rfi(λ) is zero. If freq is a vector of real frequency values, then β[i] and γ[i] are evaluated over the frequencies contained in freq.

 fdiscond(sysr::FDIFilterIF, SFDI, freq) -> (scond, β, γ)

Compute the detection and isolation sensitivity condition scond (and related quatities β and γ) for the N × mf binary structure matrix SFDI associated to the stable global transfer function matrix Rf(λ) of the transfer channel from the fault inputs to residuals of the fault detection and isolation filter internal form object sysr::FDIFilterIF. The filter sysr consists of N individual FDI filters sysr.sys[i], for i = 1, ..., N, where sysr.sys[i][:,sysr.faults] is the fault to residual channel of the i-th filter and Rfi(λ) is the corresponding transfer function matrix. The global transfer function matrix Rf(λ) is formed by row concatenation of the transfer function matrices of the N individual filters, i.e., Rf(λ) := [ Rf1(λ); Rf2(λ); ...; RfN(λ)]. It is assumed that for each j such that SFDI[i,j] = 1, the j-th column of Rfi(λ) is nonzero and for each j such that SFDI[i,j] = 0, the j-th column of Rfi(λ) is zero. The i-th element of the vectors scond, β and γ contain the quantities: β[i] - the H∞- index of the nonzero columns of Rfi(λ), γ - the maximum of the nonzero columns norms of Rf(λ) and the fault detection sensitivity condition scond evaluated as scond[i] := β[i]/γ[i]. If freq is a vector of real frequency values, then β[i] and γ[i] are evaluated over the frequencies contained in freq.

FaultDetectionTools.fdif2ngapFunction
 fdif2ngap(sysr::FDFilterIF, freq) -> (gap, β, γ)

Compute the fault-to-noise gap gap (and the related quantities β and γ) for the stable fault detection filter internal form object sysr::FDFilterIF. For the fault to residual channel of the filter sysr.sys[:,sysr.faults] with the corresponding transfer function matrix Rf(λ) and the noise to residual channel of the filter sysr.sys[:,sysr.noise] with the corresponding transfer function matrix Rw(λ), β is the H∞- index of Rf(λ), γ is the H∞-norm of Rw(λ) and gap is the fault-to-noise gap evaluated as gap := β/γ. If freq is a vector of real frequency values, then β and γ are evaluated over the frequencies contained in freq. gap = ∞ if there are no noise inputs and gap = 0 if there are no fault inputs.

 fdif2ngap(sysr::FDFilterIF, SFDI, freq) -> (gap, β, γ)

Compute the fault-to-noise gap gap (and the related quantities β and γ) for the stable fault detection filter internal form object sysr::FDFilterIF and the associated binary structure vector SFDI. sysr.sys[:,sysr.faults] is the fault to residual channel of the filter with the corresponding transfer function matrix Rf(λ) and sysr.sys[:,sysr.noise] is the noise to residual channel of the filter with the corresponding transfer function matrix Rw(λ). If Rff(λ) is the transfer function matrix formed of those j-th columns of Rf(λ) for which SFDI[j] = 1 and Rdf(λ) is the transfer function matrix formed of those j-th columns of Rf(λ) for which SFDI[j] = false, then: β is the H∞- index of Rff(λ), γ is the H∞-norm of [Rdf(λ) Rw(λ)] and gap is the fault-to-noise gap evaluated as gap := β/γ. If freq is a vector of real frequency values, then β and γ are evaluated over the frequencies contained in freq. gap = ∞ if [Rdf(λ) Rw(λ)] = 0 and gap = 0 if Rff(λ) = 0.

  fdif2ngap(sysr::FDFilterIF, SFDI, freq; atol = √ϵ) -> (gap, β, γ)

Compute the fault-to-noise gap gap (and the related quatities β and γ) for the q × mf binary structure matrix SFDI associated to the stable transfer function matrices Rf(λ) and Rw(λ) of the transfer channels from the fault inputs to residuals and noise inputs to residuals, respectively, of the fault detection filter internal form object sysr::FDFilterIF. The i-th element of the vectors gap, β and γ contain the quantities: β[i] - the H∞- index of Rffi(λ), γ[i] - the H∞-norm of [Rdfi(λ) Rwi(λ)] and the fault-to-noise gap gap evaluated as gap[i] := β[i]/γ[i], where Rffi(λ) is formed of those j-th columns of the i-th row of Rf(λ) for which S[i,j] = 1 and Rdfi(λ) is formed of those j-th columns of the i-th row of Rf(λ) for which S[i,j] = 0. gap[i] = ∞ if [Rdfi(λ) Rwi(λ)] = 0 and gap[i] = 0 if Rffi(λ) = 0. If freq is a vector of real frequency values, then β[i] and γ[i] are evaluated over the frequencies contained in freq. atol is an absolute tolerance for the norms Rwi(λ), such that norm values less than or equal to atol are considered zero (default: √ϵ, where ϵ is the working machine precision.)

 fdif2ngap(sysr::FDIFilterIF, SFDI, freq; atol = √ϵ) -> (gap, β, γ)

Compute the fault-to-noise gap gap (and the related quatities β and γ) for the N × mf binary structure matrix SFDI associated to the stable global transfer function matrices Rf(λ) and Rw(λ) of the transfer channels from the fault inputs to residuals and noise inputs to residuals, respectively, of the fault detection and isolation filter internal form object sysr::FDIFilterIF. The filter sysr consists of N individual FDI filters sysr.sys[i], for i = 1, ..., N, where sysr.sys[i][:,sysr.faults] is the fault to residual channel of the i-th filter with the corresponding transfer function matrix Rfi(λ) and sysr.sys[i][:,sysr.noise] is the noise to residual channel of the i-th filter with the corresponding transfer function matrix Rwi(λ). The global transfer function matrices Rf(λ) and Rw(λ) are formed by row concatenation of the transfer function matrices of the N individual filters, i.e., Rf(λ) := [ Rf1(λ); Rf2(λ); ...; RfN(λ)] and Rw(λ) := [ Rw1(λ); Rw2(λ); ...; RwN(λ)] Let Rffi(λ) be the transfer function matrix formed of those j-th columns of Rfi(λ) for which SFDI[i,j] = 1 and let Rdfi(λ) be the transfer function matrix formed of those j-th columns of Rfi(λ) for which SFDI[i,j] = 0. The i-th element of the vectors gap, β and γ contain the quantities: β[i] - the H∞- index of Rffi(λ), γ[i] - the H∞-norm of [Rdfi(λ) Rwi(λ)] and the fault-to-noise gap gap evaluated as gap[i] := β[i]/γ[i]. gap[i] = ∞ if [Rdfi(λ) Rwi(λ)] = 0 and gap[i] = 0 if Rffi(λ) = 0. If freq is a vector of real frequency values, then β[i] and γ[i] are evaluated over the frequencies contained in freq. atol is an absolute tolerance for the norms Rwi(λ), such that norm values less than or equal to atol are considered zero (default: √ϵ, where ϵ is the working machine precision.)

FaultDetectionTools.fdimmperfFunction
 γ = fdimmperf(sysr::FDFilterIF[, nrmflag])

Compute the model-matching performance γ of the fault detection filter internal form object sysr::FDFilterIF. If Rw(λ) is the transfer function matrix of the transfer channel from the noise inputs to residuals sysr.sys[:,sysr.noise], then γ is the H∞-norm of Rw(λ), if nrmflag = Inf (default) and the H2-norm of Rw(λ), if nrmflag = 2. The value of γ is infinite for an unstable filter or if nrmflag = 2 and the transfer function matrix Rw(λ) of a continuous-time system is not strictly proper.

 γ = fdimmperf(sysr::FDFilterIF, SFDI[, nrmflag])

Compute the model-matching performance γ of the fault detection filter internal form object sysr::FDFilterIF for a given binary structure vector SFDI. If Rf(λ) is the transfer function matrix of the transfer channel from the fault inputs to residuals sysr.sys[:,sysr.faults] and Rw(λ) is the transfer function matrix of the transfer channel from the noise inputs to residuals sysr.sys[:,sysr.noise], then γ is the H∞-norm of [Rdf(λ) Rw(λ)], if nrmflag = Inf (default) and the H2-norm of [Rdf(λ) Rw(λ)], if nrmflag = 2, where Rdf(λ) is the transfer function matrix formed by those j-th columns of Rf(λ) for which SFDI[j] = 0. The value of γ is infinite for an unstable filter or if nrmflag = 2 and the transfer function matrix [Rdf(λ) Rw(λ)] of a continuous-time system is not strictly proper.

 γ = fdimmperf(sysr::FDFilterIF, SFDI[, nrmflag])

Compute the model-matching performance γ of the fault detection filter internal form object sysr::FDFilterIF for a given binary structure matrix SFDI. If Rf(λ) is the transfer function matrix of the transfer channel from the fault inputs to residuals sysr.sys[:,sysr.faults] and Rw(λ) is the transfer function matrix of the transfer channel from the noise inputs to residuals sysr.sys[:,sysr.noise], then γ is the H∞-norm of [Rdf(λ) Rw(λ)], if nrmflag = Inf (default) and the H2-norm of [Rdf(λ) Rw(λ)], if nrmflag = 2, where Rdf(λ) = .!SFDI .* Rf(λ) (i.e., the element-wise product of .!SFDI and Rf(λ). The value of γ is infinite for an unstable filter or if nrmflag = 2 and the transfer function matrix [Rdf(λ) Rw(λ)] of a continuous-time system is not strictly proper.

 γ = fdimmperf(sysr::FDIFilterIF[, nrmflag])

Compute the model-matching performance γ of the stable fault detection and isolation filter internal form object sysr::FDIFilterIF. The filter sysr consists of N individual FDI filters sysr.sys[i], for i = 1, ..., N, where sysr.sys[i][:,sysr.noise] is the noise to residual channel of the i-th filter with the corresponding transfer function matrix Rwi(λ). Then, γ is an N-dimensional vector whose i-th component is the H∞-norm of Rwi(λ), if nrmflag = Inf (default) and the H2-norm of Rwi(λ), if nrmflag = 2. The i-th component of γ is infinite for an unstable filter or if nrmflag = 2 and the transfer function matrix Rwi(λ) of a continuous-time system is not strictly proper.

 γ = fdimmperf(sysr::FDIFilterIF, SFDI[, nrmflag])

Compute the model-matching performance γ of the stable fault detection and isolation filter internal form object sysr::FDIFilterIF and the associated binary structure matrix SFDI. The filter sysr consists of N individual FDI filters sysr.sys[i], for i = 1, ..., N, where sysr.sys[i][:,sysr.faults] is the fault to residual channel of the i-th filter with the corresponding transfer function matrix Rfi(λ) and sysr.sys[i][:,sysr.noise] is the noise to residual channel of the i-th filter with the corresponding transfer function matrix Rwi(λ). Then, γ is an N-dimensional vector whose i-th component is the H∞-norm of [Rfdi(λ) Rwi(λ)], if nrmflag = Inf (default) or the H2-norm of [Rfdi(λ) Rwi(λ)], if nrmflag = 2, where Rfdi(λ) is the transfer function matrix whose j-th column is zero if SFDI[i,j] = 1 and is equal to the j-th column of Rfi(λ) if SFDI[i,j] = 0. The i-th component of γ is infinite for an unstable filter or if nrmflag = 2 and the transfer function matrix [Rfdi(λ) Rwi(λ)] of a continuous-time system is not strictly proper.

 γ = fdimmperf(sysr::FDFilterIF, sysref::Union{FDFilterIF,FDIModel}[, nrmflag])

Compute the model-matching performance γ of the fault detection filter internal form object sysr::FDFilterIF with respect to the fault detection reference filter internal form sysref::FDFilterIF. If R(λ) is the transfer function matrix of the fault detection filter internal form sysr.sys and Mr(λ) is the transfer function matrix of the fault detection reference filter internal form sysref.sys, then γ is the H∞-norm of R(λ)-Mr(λ), if nrmflag = Inf (default) or the H2-norm of R(λ)-Mr(λ), if nrmflag = 2. The value of γ is infinite for an unstable difference R(λ)-Mr(λ) or if nrmflag = 2 and the transfer function matrix R(λ)-Mr(λ) of a continuous-time system is not strictly proper. In general, R(λ) and Mr(λ) are partitioned as R(λ) = [ Ru(λ) Rd(λ) Rf(λ) Rw(λ) Ra(λ) ] and Mr(λ) = [ Mru(λ) Mrd(λ) Mrf(λ) Mrw(λ) Mra(λ) ] in accordance with the partitioning of the inputs in control inputs, disturbance inputs, fault inputs, noise inputs and auxiliary inputs. Void components of Mr(λ) corresponding to non-void components in R(λ) are assumed to be zero.

 γ = fdimmperf(sysr::FDIFilterIF, sysref::FDIFilterIF[, nrmflag])

Compute the model-matching performance γ of the fault detection and isolation filter internal form object sysr::FDIFilterIF with respect to the fault detection and isolation reference filter internal form sysref::FDIFilterIF. If Ri(λ) is the transfer function matrix of the i-th fault detection and isolation filter internal form sysr.sys[i] and Mri(λ) is the transfer function matrix of the i-th fault detection and isolation reference filter internal form sysref.sys[i], then γ is a vector whose i-th component γ[i] is the H∞-norm of Ri(λ)-Mri(λ), if nrmflag = Inf (default) or the H2-norm of Ri(λ)-Mri(λ), if nrmflag = 2. The value of γ[i] is infinite for an unstable difference Ri(λ)-Mri(λ) or if nrmflag = 2 and the transfer function matrix Ri(λ)-Mri(λ) of a continuous-time system is not strictly proper. In general, Ri(λ) and Mri(λ) are partitioned as Ri(λ) = [ Rui(λ) Rdi(λ) Rfi(λ) Rwi(λ) Rai(λ) ] and Mri(λ) = [ Mrui(λ) Mrdi(λ) Mrfi(λ) Mrwi(λ) Mrai(λ) ] in accordance with the partitioning of the inputs in control inputs, disturbance inputs, fault inputs, noise inputs and auxiliary inputs. Void components of Mri(λ) corresponding to non-void components in Ri(λ) are assumed to be zero.