Local Discriminant Basis
WaveletsExt.LDB.AsymmetricRelativeEntropyWaveletsExt.LDB.BasisDiscriminantMeasureWaveletsExt.LDB.DiscriminantMeasureWaveletsExt.LDB.DiscriminantPowerWaveletsExt.LDB.EarthMoverDistanceWaveletsExt.LDB.EnergyMapWaveletsExt.LDB.FishersClassSeparabilityWaveletsExt.LDB.HellingerDistanceWaveletsExt.LDB.LocalDiscriminantBasisWaveletsExt.LDB.LpDistanceWaveletsExt.LDB.ProbabilityDensityWaveletsExt.LDB.ProbabilityDensityDMWaveletsExt.LDB.RobustFishersClassSeparabilityWaveletsExt.LDB.SignaturesWaveletsExt.LDB.SignaturesDMWaveletsExt.LDB.SymmetricRelativeEntropyWaveletsExt.LDB.TimeFrequencyWaveletsExt.LDB.change_nfeaturesWaveletsExt.LDB.discriminant_measureWaveletsExt.LDB.discriminant_powerWaveletsExt.LDB.energy_mapWaveletsExt.LDB.fit!WaveletsExt.LDB.fit_transformWaveletsExt.LDB.fitdec!WaveletsExt.LDB.inverse_transformWaveletsExt.LDB.pairwise_discriminant_measureWaveletsExt.LDB.transform
Energy Maps
Public API
WaveletsExt.LDB.EnergyMap — TypeEnergyMapEnergy map for Local Discriminant Basis. Current available types are:
WaveletsExt.LDB.TimeFrequency — TypeTimeFrequency <: EnergyMapAn energy map based on time frequencies, a measure based on the differences of derived quantities from projection $Z_i$, such as mean class energies or cumulants.
See also: EnergyMap, ProbabilityDensity, Signatures
WaveletsExt.LDB.ProbabilityDensity — TypeProbabilityDensity <: EnergyMapAn energy map based on probability density, a measure based on the differences among the pdfs of $Z_i$. Since we do not know the true density functions of the coefficients, the PDFs are estimated using the Average Shifted Histogram (ASH).
See also: EnergyMap, TimeFrequency, Signatures
WaveletsExt.LDB.Signatures — TypeSignatures <: EnergyMapAn energy map based on signatures, a measure that uses the Earth Mover's Distance (EMD) to compute the discriminating measure of a coordinate. Signatures provide us with a fully data-driven representation, which can be efficiently used with EMD. This representation is more efficient than a histogram and is able to represent complex data structure with fewer samples.
Here, a signature for the coefficients in the $j$-th level, $k$-th node, $l$-th index of class $c$ is defined as
$s_{j,k,l}^{(c)} = \{(\alpha_{i;j,k,l}^{(c)}, w_{i;j,k,l}^{(c)})\}_{i=1}^{N_c}$
where $\alpha_{i;j,k,l}^{(c)}$ and $w_{i;j,k,l}^{(c)}$ are the expansion coefficients and weights at location $(j,k,l)$ for signal $i$ of class $c$ respectively. Currently, the two valid types of weights are :equal and :pdf.
Parameters
weight::Symbol: (Default::equal) Type of weight to be used to compute $w_{i;j,k,l}^{(c)}$. Available methods are:equaland:pdf.
See also: EnergyMap, TimeFrequency, ProbabilityDensity
WaveletsExt.LDB.energy_map — Functionenergy_map(Xw, y, method)Computes the energy map based on the decomposed input signals Xw and corresponding labels y.
Arguments
Xw::AbstractArray{S} where S<:AbstractFloat: Decomposed 1D or 2D signals.y::AbstractVector{T} where T: Corresponding labels to the signals inXw.method::EnergyMap: Type of energy map to compute. Supported types are:TimeFrequency()ProbabilityDensity()Signatures()
Returns
Γ: Computed energy map. Depending on the inputmethod, the data structure ofΓis:TimeFrequency()$\Longrightarrow$Array{S,3}for 1D signals orArray{S,4}for 2D signals.ProbabilityDensity()$\Longrightarrow$Array{S,4}for 1D signals orArray{S,5}for 2D signals.Signatures(weight = :equal): $\Longrightarrow$ Vector{NamedTuple{(:coef, :weight), Tuple{Array{S}, S}}}Signatures(weight = :pdf: $\Longrightarrow$ Vector{NamedTuple{(:coef, :weight), Tuple{Array{S}, Array{S}}}}
When setting method = Signatures() of any type of weight, the output is a vector of named tuples, where each named tuple contains the coefficients and weights of all the signals of the same label.
Examples
using Wavelets, WaveletsExt
X, y = generateclassdata(ClassData(:tri, 5, 5, 5))
Xw = wpdall(X, wavelet(WT.haar))
energy_map(Xw, y, TimeFrequency())
energy_map(Xw, y, ProbabilityDensity())
energy_map(Xw, y, Signatures())See also: EnergyMap. TimeFrequency, ProbabilityDensity
Discriminant Measures
Public API
WaveletsExt.LDB.DiscriminantMeasure — TypeDiscriminant measure for Local Discriminant Basis. Current available subtypes are:
WaveletsExt.LDB.ProbabilityDensityDM — TypeDiscriminant measure for Probability Density and Time Frequency based energy maps. Current available measures are:
WaveletsExt.LDB.SignaturesDM — TypeDiscriminant measure for Signatures based energy maps. Current available measures are:
WaveletsExt.LDB.AsymmetricRelativeEntropy — TypeAsymmetricRelativeEntropy <: ProbabilityDensityDMAsymmetric Relative Entropy discriminant measure for the Probability Density and Time Frequency based energy maps. This measure is also known as cross entropy and Kullback-Leibler divergence.
Equation:
$D(p,q) = \sum p(x) \log \frac{p(x)}{q(x)}$
where $\int_{-\infty}^{\infty} p(t) dt = \int_{-\infty}^{\infty} q(t) dt = 1$ ($p(t)$ and $q(t)$ are probability density functions.)
WaveletsExt.LDB.SymmetricRelativeEntropy — TypeSymmetricRelativeEntropy <: ProbabilityDensityDMSymmetric Relative Entropy discriminant measure for the Probability Density and Time Frequency energy maps. Similar idea to the Asymmetric Relative Entropy, but this aims to make the measure more symmetric.
Equation: Denote the Asymmetric Relative Entropy as $D_A(p,q)$, then
$D(p,q) = D_A(p,q) + D_A(q,p) = \sum p(x) \log \frac{p(x)}{q(x)} + q(x) \log \frac{q(x)}{p(x)}$
where $\int_{-\infty}^{\infty} p(t) dt = \int_{-\infty}^{\infty} q(t) dt = 1$ ($p(t)$ and $q(t)$ are probability density functions.)
See also: AsymmetricRelativeEntropy
WaveletsExt.LDB.HellingerDistance — TypeHellingerDistance <: ProbabilityDensityDMHellinger Distance discriminant measure for the Probability Density energy map.
Equation:
$H(p,q) = \sum_{i=1}^n (\sqrt{p_i} - \sqrt{q_i})^2$
WaveletsExt.LDB.LpDistance — TypeLpDistance <: ProbabilityDensityDM$\ell^p$ Distance discriminant measure for the Probability Density and Time Frequency based energy maps. The default $p$ value is set to 2.
Equation:
$W(q,r) = ||q-r||_p^p = \sum_{i=1}^n (q_i - r_i)^p$
WaveletsExt.LDB.EarthMoverDistance — TypeEarthMoverDistance <: SignaturesDMEarth Mover Distance discriminant measure for the Signatures energy map.
Equation:
$E(P,Q) = \frac{\sum_{k=1}^{m+n+1} |\hat p_k - \hat q_k| (r_{k+1} - r_k)}{w_\Sigma}$
where
- $r_1, r_2, \ldots, r_{m+n}$ is the sorted list of $p_1, \ldots, p_m, q_1, \ldots, q_n$
- $\hat p_k = \sum_{p_i \leq r_k} w_{p_i}$
- $\hat q_k = \sum_{q_i \leq r_k} w_{q_i}$
WaveletsExt.LDB.discriminant_measure — Functiondiscriminant_measure(Γ, dm)Computes the discriminant measure of each subspace calculated from the energy maps.
Arguments
Γ: Energy map computed fromenergy_mapfunction. The data structures ofΓ, depending on their corresponding energy map, should be:TimeFrequency():AbstractArray{T,3}for 1D signals orAbstractArray{T,4}for 2D signals.ProbabilityDensity():AbstractArray{T,4}for 1D signals orAbstractArray{T,5}for 2D signals.Signatures():AbstractVector{NamedTuple{(:coef, :weight), Tuple{S₁,S₂}}}.
dm::DiscriminantMeasure: Type of Discriminant Measure. The type ofdmmust match the type ofΓ.
Returns
D::Array{T}: Discriminant measure at each coefficient of the decomposed signals.
Examples
using Wavelets, WaveletsExt
X, y = generateclassdata(ClassData(:tri, 5, 5, 5))
Xw = wpdall(X, wavelet(WT.haar))
Γ = energy_map(Xw, y, TimeFrequency()); discriminant_measure(Γ, AsymmetricRelativeEntropy())
Γ = energy_map(Xw, y, ProbabilityDensity()); discriminant_measure(Γ, LpDistance())
Γ = energy_map(Xw, y, Signatures()); discriminant_measure(Γ, EarthMoverDistance())See also: pairwise_discriminant_measure
Private API
WaveletsExt.LDB.pairwise_discriminant_measure — Functionpairwise_discriminant_measure(Γ₁, Γ₂, dm)Computes the discriminant measure between 2 classes based on their energy maps.
Arguments
Γ₁::AbstractArray{T} where T<:AbstractFloat: Energy map for class 1.Γ₂::AbstractArray{T} where T<:AbstractFloat: Energy map for class 2.dm::DiscriminantMeasure: Type of discriminant measure.
Returns
::Array{T}: Discriminant measure betweenΓ₁andΓ₂.
pairwise_discriminant_measure(p, q, dm)Computes the discriminant measure between 2 classes at an index $i$, where Γ₁[i] = p and Γ₂[i] = q.
Arguments
p::T where T<:AbstractFloat: Coefficient at index $i$ fromΓ₁.q::T where T<:AbstractFloat: Coefficient at index $i$ fromΓ₂.dm::DiscriminantMeasure: Type of discriminant measure.
Returns
::T: Discriminant measure betweenpandq.
Computation of Discriminant Powers
Public API
WaveletsExt.LDB.DiscriminantPower — TypeDiscriminant Power measure for the Local Discriminant Basis. Current available measures are
WaveletsExt.LDB.BasisDiscriminantMeasure — TypeBasisDiscriminantMeasure <: DiscriminantPowerThis is the discriminant measure of a single basis function computed in a previous step to construct the energy maps.
WaveletsExt.LDB.FishersClassSeparability — TypeFishersClassSeparability <: DiscriminantPowerThe Fisher's class separability of the expansion coefficients in the basis function.
Equation: $\frac{\sum_{c=1}^C \pi_c\{{\rm mean}_i(\alpha_{\lambda,i}^{(c)}) - {\rm mean}_c \cdot {\rm mean}_i(\alpha_{\lambda,i}^{(c)})\}^2}{\sum_{c=1}^C \pi_c {\rm var}_i(\alpha_{\lambda,i}^{(c)})}$
WaveletsExt.LDB.RobustFishersClassSeparability — TypeRobustFishersClassSeparability <: DiscriminantPowerThe robust version of Fisher's class separability of the expansion coefficients in the basis function.
Equation: $\frac{\sum_{c=1}^C \pi_c\{{\rm med}_i(\alpha_{\lambda,i}^{(c)}) - {\rm med}_c \cdot {\rm med}_i(\alpha_{\lambda,i}^{(c)})\}^2}{\sum_{c=1}^C \pi_c {\rm mad}_i(\alpha_{\lambda,i}^{(c)})}$
WaveletsExt.LDB.discriminant_power — Functiondiscriminant_power(D, tree, dp)
discriminant_power(coefs, y, dp)Returns the discriminant power of each leaf from the local discriminant basis (LDB) tree.
Arguments
D::AbstractArray{T} where T<:AbstractFloat: Discriminant measures.tree::BitVector: Best basis tree for selecting coefficients with largest discriminant measures.coefs::AbstractArray{T} where T<:AbstractFloat: Best basis coefficients for the input signals.y::AbstractVector{S} where S: Labels corresponding to each signal incoefs.dp::DiscriminantPower: The measure for discriminant power.
discriminant_power(D, tree, dp) only works for dp = BasisDiscriminantMeasure(), whereas discriminant_power(coefs, y, dp) works for dp = FishersClassSeparability() and dp = RobustFishersClassSeparability().
Returns
power::Array{T}: The discriminant power at each index ofDorcoefs.order::Vector{T}: The order of discriminant power in descending order.
Feature Extraction and Transformation
Public API
WaveletsExt.LDB.LocalDiscriminantBasis — TypeLocalDiscriminantBasisClass type for the Local Discriminant Basis (LDB), a feature selection algorithm developed by N. Saito and R. Coifman in "Local Discriminant Bases and Their Applications" in the Journal of Mathematical Imaging and Vision, Vol 5, 337-358 (1995). This struct contains the following field values:
Parameters and Attributes:
wt::DiscreteWavelet: (Default:wavelet(WT.haar)) A discrete wavelet for transform purposes.max_dec_level::Union{Integer, Nothing}: (Default:nothing) Max level of wavelet packet decomposition to be computed.dm::DiscriminantMeasure: (Default:AsymmetricRelativeEntropy()) the discriminant measure for the LDB algorithm. Supported measures are theAsymmetricRelativeEntropy(),LpDistance(),SymmetricRelativeEntropy(), andHellingerDistance().en::EnergyMap: (Default:TimeFrequency()) the type of energy map used. Supported maps areTimeFrequency(),ProbabilityDensity(), andSignatures().dp::DiscriminantPower(): (Default:BasisDiscriminantMeasure()) the measure of discriminant power among expansion coefficients. Supported measures areBasisDiscriminantMeasure(),FishersClassSeparability(), andRobustFishersClassSeparability().top_k::Union{Integer, Nothing}: (Default:nothing) the top-k coefficients used in each node to determine the discriminant measure.n_features::Union{Integer, Nothing}: (Default:nothing) the dimension of output after undergoing feature selection and transformation.sz::Union{Vector{T} where T<:Integer, Nothing}: (Default:nothing) Size of signalΓ::Union{AbstractArray{T} where T<:AbstractFloat, AbstractArray{NamedTuple{(:coef, :weight), Tuple{S₁, S₂}}} where {S₁<:Array{T} where T<:AbstractFloat, S2<:Union{T, Array{T}} where T<:AbstractFloat}, Nothing}: (Default:nothing) computed energy mapDM::Union{AbstractArray{<:AbstractFloat}, Nothing}: (Default:nothing) computed discriminant measurecost::Union{AbstractVector{<:AbstractFloat}, Nothing}: (Default:nothing) computed wavelet packet decomposition (WPD) tree cost based on the discriminant measureDM.tree::Union{BitVector, Nothing}: (Default:nothing) computed best WPD tree based on the discriminant measureDM.DP::Union{AbstractVector{<:AbstractFloat}, Nothing}: (Default:nothing) computed discriminant powerorder::Union{AbstractVector{Integer}, Nothing}: (Default:nothing) ordering ofDPby descending order.
WaveletsExt.LDB.fit! — Functionfit!(f, X, y)Fits the Local Discriminant Basis feature selection algorithm f onto the signals X with labels y.
Arguments
f::LocalDiscriminantBasis: Local discriminant basis object.X::AbstractArray{S} where S<:AbstractFloat: Group of signals of the same size, arranged in the 2nd (1D signals) or 3rd (2D signals) dimension.y::AbstractVector{T} where T: Labels corresponding to each signal in X.
Returns
f::LocalDiscriminantBasis: The same LDB objectfwith updated fields.
Examples
using Wavelets, WaveletsExt
X, y = generateclassdata(ClassData(:tri, 5, 5, 5))
f = LocalDiscriminantBasis()
fit!(f, X, y)See also: LocalDiscriminantBasis, fit_transform, transform, inverse_transform, change_nfeatures
WaveletsExt.LDB.fitdec! — Functionfitdec!(f, Xw, y)Fits the Local Discriminant Basis feature selection algorithm f onto the decomposed signals Xw with labels y.
Arguments
f::LocalDiscriminantBasis: Local discriminant basis object.Xw::AbstractArray{S} where S<:AbstractFloat: Group of decomposed signals of the same size, arranged in the 3rd (1D signals) or 4th (2D signals) dimension.y::AbstractVector{T} where T: Labels corresponding to each signal in X.
Returns
f::LocalDiscriminantBasis: The same LDB objectfwith updated fields.
Examples
using Wavelets, WaveletsExt
X, y = generateclassdata(ClassData(:tri, 5, 5, 5))
f = LocalDiscriminantBasis()
Xw = wpdall(X, f.wt)
fitdec!(f, Xw, y)See also: LocalDiscriminantBasis, fit!, fit_transform, transform, inverse_transform, change_nfeatures
WaveletsExt.LDB.fit_transform — Functionfit_transform(f, X, y)Fit and transform the signals X with labels y based on the LDB class f. This function essentially combines both fit! and transform into one, and is the preferred method when trying to extract the LDB coefficients of the training dataset. For the test dataset, one needs to first run fit! on the training dataset before running transform on the test dataset.
Arguments
f::LocalDiscriminantBasis: Local discriminant basis object.Xw::AbstractArray{S} where S<:AbstractFloat: Group of decomposed signals of the same size, arranged in the 3rd (1D signals) or 4th (2D signals) dimension.y::AbstractVector{T} where T: Labels corresponding to each signal in X.
Returns
Xc::Matrix{S}: A matrix of the extracted LDB features fromX, where each column corresponds to the coefficients extracted from the corresponding signal.
Examples
using Wavelets, WaveletsExt
X, y = generateclassdata(ClassData(:tri, 5, 5, 5))
f = LocalDiscriminantBasis()
Xc = fit_transform(f, X, y)See also: LocalDiscriminantBasis, fit!, transform, inverse_transform, change_nfeatures
WaveletsExt.LDB.transform — Functiontransform(f, X)Extract the LDB features from signals X based on the LDB tree constructed from previously fitted signals.
Arguments
f::LocalDiscriminantBasis: Local discriminant basis object.X::AbstractArray{S} where S<:AbstractFloat: Group of signals of the same size, arranged in the 2nd (1D signals) or 3rd (2D signals) dimension.
Returns
Xc::Matrix{S}: A matrix of the extracted LDB features fromX, where each column corresponds to the coefficients extracted from the corresponding signal.
Examples
using Wavelets, WaveletsExt
X, y = generateclassdata(ClassData(:tri, 5, 5, 5))
f = LocalDiscriminantBasis()
fit!(f, X, y)
Xc = transform(f, X)See also: LocalDiscriminantBasis, fit!, fit_transform, inverse_transform, change_nfeatures
WaveletsExt.LDB.inverse_transform — Functioninverse_transform(f, X)Compute the inverse transform on the feature matrix X to form the original signal based on the LDB class f. This function can be used to reconstruct the signals and compare the reconstructed features between different classes.
Arguments
f::LocalDiscriminantBasis: Local discriminant basis object.X::Matrix{S} where S<:AbstractFloat: A matrix of the extracted LDB features from original signals, where each column corresponds to the coefficients extracted from the corresponding signal.
Returns
Xₜ::Array{S}: Reconstructed signals.
Examples
X, y = generateclassdata(ClassData(:tri, 5, 5, 5))
f = LocalDiscriminantBasis()
Xc = fit_transform(f, X, y)
Xₜ = inverse_transform(f, Xc)See also: LocalDiscriminantBasis, fit!, fit_transform, transform, change_nfeatures
WaveletsExt.LDB.change_nfeatures — Functionchange_nfeatures(f, x, n_features)Change the number of features extracted from the signals from f.n_features to n_features.
!!! warning If the input n_features is larger than f.n_features, it results in the regeneration of signals based on the current f.n_features before reselecting the features. This will cause additional features to be less accurate and effective. Often times, the "additional" features tend to be zeros.
Arguments
f::LocalDiscriminantBasis: Local discriminant basis object.x::AbstractMatrix{S} where S<:AbstractFloat: A matrix of the extracted LDB features from original signals, where each column corresponds to the coefficients extracted from the corresponding signal.n_features::Integer: Desired number of features in selected by LDB.
Returns
xₜ::Matrix{S}: Matrix of extracted LDB features, where the number of rows is nown_featuresinstead of the previousf.n_features.
See also: LocalDiscriminantBasis, fit!, fit_transform, transform, inverse_transform