topic/diffentiators
Vincent Samy 2019-01-04 17:15:03 +09:00
rodzic 655db4a39a
commit 74f52009dc
9 zmienionych plików z 252 dodań i 55 usunięć

Wyświetl plik

@ -0,0 +1,4 @@
USE_MDFILE_AS_MAINPAGE = @CMAKE_SOURCE_DIR@/README.md
FILE_PATTERNS = *.h
EXTRA_PACKAGES = amsmath
EXTRACT_PRIVATE = YES

Wyświetl plik

@ -5,15 +5,40 @@
namespace fratio {
/*! \brief Transform an analog signal to a discrete signal and vice versa.
*
* \see https://en.wikipedia.org/wiki/Bilinear_transform
* \tparam T Floating (complex) types.
*/
template <typename T>
struct BilinearTransform {
using SubType = internal::complex_sub_type_t<T>;
using SubType = internal::complex_sub_type_t<T>; /*!< Sub-type of the complex if T is complex, T otherwise */
static_assert(std::is_floating_point<SubType>::value, "This struct can only accept floating point types (real and complex).");
/*! \brief Transformation from analog to discrete.
* \param fs Sampling frequency.
* \param sPlanePole Analog data.
* \param[out] zPlanePole Resulting discrete data.
*/
static void SToZ(SubType fs, const T& sPlanePole, T& zPlanePole);
static void SToZ(SubType fs, const vectX_t<T>& sPlanePoles, Eigen::Ref<vectX_t<T>>& zPlanePoles); // Can be optimized
/*! \brief Transformation from analog to discrete.
* \param fs Sampling frequency.
* \param sPlanePole Analog signal.
* \param[out] zPlanePole Resulting discrete signal.
*/
static void SToZ(SubType fs, const vectX_t<T>& sPlanePoles, Eigen::Ref<vectX_t<T>>& zPlanePoles); // Can be optimized maybe
/*! \brief Transformation from discrete to analog.
* \param fs Sampling frequency.
* \param zPlanePole Discrete data.
* \param[out] sPlanePole Resulting analog data.
*/
static void ZToS(SubType fs, const T& zPlanePole, T& sPlanePole);
static void ZToS(SubType fs, const vectX_t<T>& zPlanePoles, Eigen::Ref<vectX_t<T>>& sPlanePoles); // Can be optimized
/*! \brief Transformation from discrete to analog.
* \param fs Sampling frequency.
* \param zPlanePole Discrete signal.
* \param[out] sPlanePole Resulting analog signal.
*/
static void ZToS(SubType fs, const vectX_t<T>& zPlanePoles, Eigen::Ref<vectX_t<T>>& sPlanePoles); // Can be optimized maybe
};
template <typename T>

Wyświetl plik

@ -7,51 +7,127 @@
namespace fratio {
// https://www.dsprelated.com/showarticle/1119.php
// https://www.dsprelated.com/showarticle/1135.php
// https://www.dsprelated.com/showarticle/1128.php
// https://www.dsprelated.com/showarticle/1131.php
// https://www.mathworks.com/help/signal/ref/butter.html
/*! \brief Butterworth digital filter.
*
* This is a specialization of a digital filter in order to use a Butterworth filter.
* \see https://www.dsprelated.com/showarticle/1119.php
* \see https://www.dsprelated.com/showarticle/1135.php
* \see https://www.dsprelated.com/showarticle/1128.php
* \see https://www.dsprelated.com/showarticle/1131.php
* \see https://www.mathworks.com/help/signal/ref/butter.html
* \tparam Floating type.
*/
template <typename T>
class Butterworth : public DigitalFilter<T> {
public:
static T PI;
static T PI; /*!< pi depending on the type T */
public:
/*! \brief Type of butterworth filter0 */
enum class Type {
LowPass,
HighPass,
BandPass,
BandReject
LowPass, /*!< Define a low-pass filter */
HighPass, /*!< Define a high-pass filter */
BandPass, /*!< Define a band-pass filter */
BandReject /*!< Define a band-reject filter */
// AllPass
// CombPass
};
public:
// http://www.matheonics.com/Tutorials/Butterworth.html#Paragraph_3.2
// https://www.mathworks.com/help/signal/ref/buttord.html#d120e11079
/*! \brief Function to help you design a Butterworth filter.
*
* It finds optimal values of the order and cut-off frequency.
* \warning Works only for low-pass and high-pass filters.
* \see http://www.matheonics.com/Tutorials/Butterworth.html#Paragraph_3.2
* \see https://www.mathworks.com/help/signal/ref/buttord.html#d120e11079
* \param wPass Pass band edge.
* \param wStop Stop band edge.
* \param APass Maximum pass band attenuation.
* \param AStop Minimum stop band attenuation.
* \return A pair of { filter order, cut-off frequency }.
*/
static std::pair<int, T> findMinimumButter(T wPass, T wStop, T APass, T AStop); // Only for low and high pass
public:
/*! \brief Uninitialized constructor.
* \param type Filter type. Default is LowPass.
*/
Butterworth(Type type = Type::LowPass);
/*! \brief Constructor for both low-pass and high-pass filters.
* \param order Order of the filter.
* \param fc Cut-off frequency.
* \param fs Sampling frequency.
* \param type Filter type. Default is LowPass.
*/
Butterworth(int order, T fc, T fs, Type type = Type::LowPass);
/*! \brief Constructor for both band-pass and band-reject filters.
* \param order Order of the filter.
* \param fLower Lower bound frequency.
* \param fUpper Upper bound frequency.
* \param fs Sampling frequency.
* \param type Filter type. Default is BandPass.
*/
Butterworth(int order, T fLower, T fUpper, T fs, Type type = Type::BandPass);
/*! \brief Set filter set of parameters.
* \param order Order of the filter.
* \param fc Cut-off frequency.
* \param fs Sampling frequency.
*/
void setFilterParameters(int order, T fc, T fs);
/*! \brief Set filter set of parameters.
* \param order Order of the filter.
* \param fLower Lower bound frequency.
* \param fUpper Upper bound frequency.
* \param fs Sampling frequency.
*/
void setFilterParameters(int order, T fLower, T fUpper, T fs);
private:
void initialize(int order, T f1, T f2, T fs); // f1 = fc for LowPass/HighPass filter
/*! \brief Initialize the filter.
* \param order Order of the filter.
* \param f1 First frequency parameter.
* \param f2 Second frequency parameter.
* \param fs Sampling frequency.
*/
void initialize(int order, T f1, T f2, T fs);
/*! \brief Compute the digital filter representation for low-pass and high-pass.
* \param fc Cut-off frequency.
*/
void computeDigitalRep(T fc);
/*! \brief Compute the digital filter representation for band-pass and band-reject.
* \param fLower Lower bound frequency.
* \param fUpper Upper bound frequency.
*/
void computeBandDigitalRep(T fLower, T fUpper);
/*! \brief Generate an analog pole on the unit circle for low-pass and high-pass.
* \param k Step on the unit circle.
* \param fpw Continuous pre-warp cut-off frequency.
* \return Generated pole.
*/
std::complex<T> generateAnalogPole(int k, T fpw);
/*! \brief Generate an analog pole on the unit circle for band-pass and band-reject.
* \param k Step on the unit circle.
* \param fpw0 Continuous pre-warp frequency at geometric center.
* \param bw Bandwith.
* \return Pair of generated pole.
*/
std::pair<std::complex<T>, std::complex<T>> generateBandAnalogPole(int k, T fpw0, T bw);
vectXc_t<T> generateAnalogZeros(T f0 = T());
/*! \brief Generate all analog zeros.
* \param fpw0 Continuous pre-warp frequency at geometric center (Only use by the band-reject).
* \return Set of generated zeros.
*/
vectXc_t<T> generateAnalogZeros(T fpw0 = T());
/*! \brief Scale coefficients.
* \param aCoeff Unscaled poles.
* \param bCoeff Unscaled zeros.
* \param bpS Result of \f$H(fc)=1\f$ with \f$H(z)\f$ the discrete transfer function at the geometric center \f$fc\f$ (Only use by band-pass).
*/
void scaleAmplitude(const vectX_t<T>& aCoeff, Eigen::Ref<vectX_t<T>> bCoeff, const std::complex<T>& bpS = T());
private:
Type m_type;
int m_order;
T m_fs;
vectXc_t<T> m_poles;
Type m_type; /*!< Filter type */
int m_order; /*!< Filter order */
T m_fs; /*!< Filter sampling frequency */
};
} // namespace fratio

Wyświetl plik

@ -11,8 +11,8 @@ std::pair<int, T> Butterworth<T>::findMinimumButter(T wPass, T wStop, T APass, T
{
T num = std::log10((std::pow(T(10), T(0.1) * std::abs(AStop)) - 1) / (std::pow(T(10), T(0.1) * std::abs(APass)) - 1));
// pre-warp
T fwPass = std::tan(0.5 * PI * wPass);
T fwStop = std::tan(0.5 * PI * wStop);
T fwPass = std::tan(T(0.5) * PI * wPass);
T fwStop = std::tan(T(0.5) * PI * wStop);
T w;
if (wPass < wStop)
w = std::abs(fwStop / fwPass);
@ -65,6 +65,8 @@ void Butterworth<T>::setFilterParameters(int order, T fLower, T fUpper, T fs)
template <typename T>
void Butterworth<T>::initialize(int order, T f1, T f2, T fs)
{
// f1 = fc for LowPass/HighPass filter
// f1 = fLower, f2 = fUpper for BandPass/BandReject filter
if (order <= 0) {
m_status = FilterStatus::BAD_ORDER_SIZE;
return;

Wyświetl plik

@ -5,11 +5,20 @@
namespace fratio {
// https://www.mathworks.com/help/dsp/ug/how-is-moving-average-filter-different-from-an-fir-filter.html
/*! \brief Basic digital filter.
*
* This filter allows you to set any digital filter based on its coefficients.
* \tparam T Floating type.
*/
template <typename T>
class DigitalFilter : public GenericFilter<T> {
public:
/*! \brief Default uninitialized constructor. */
DigitalFilter() = default;
/*! \brief Constructor.
* \param aCoeff Denominator coefficients of the filter in decreasing order.
* \param bCoeff Numerator coefficients of the filter in decreasing order.
*/
DigitalFilter(const vectX_t<T>& aCoeff, const vectX_t<T>& bCoeff)
: GenericFilter<T>(aCoeff, bCoeff)
{

Wyświetl plik

@ -7,44 +7,104 @@
namespace fratio {
/*! \brief Low-level filter.
*
* It creates the basic and common functions of all linear filter that can written as a digital filter.
* This class can not be instantiated directly.
*
* \warning In Debug mode, all functions may throw if a filter is badly initialized.
* This not the case in Realese mode.
*
* \tparam T Floating type.
*/
template <typename T>
class GenericFilter {
static_assert(std::is_floating_point<T>::value && !std::is_const<T>::value, "Only accept non-complex floating point types.");
public:
/*! \brief Get the meaning of the filter status.
* \param status Filter status to get the meaning from.
* \return The meaning.
*/
static std::string filterStatus(FilterStatus status);
public:
// Careful: Only an assert check for the filter status
/*! \brief Filter a new data.
*
* This function is practical for online application that does not know the whole signal in advance.
* \param data New data to filter.
* \return Filtered data.
*/
T stepFilter(const T& data);
/*! \brief Filter a signal.
*
* Filter all data given by the signal.
* \param data Signal.
* \return Filtered signal.
*/
vectX_t<T> filter(const vectX_t<T>& data);
/*! \brief Filter a signal and store in a user-defined Eigen vector.
*
* Useful if the length of the signal is known in advance.
* \param[out] results Filtered signal.
* \param data Signal.
* \return False if vector's lengths do not match.
*/
bool getFilterResults(Eigen::Ref<vectX_t<T>> results, const vectX_t<T>& data);
/*! \brief Reset the data and filtered data. */
void resetFilter();
template <typename T2>
void setCoeffs(T2&& aCoeff, T2&& bCoeff);
/*! \brief Get digital filter coefficients.
*
* It will automatically resize the given vectors.
* \param[out] aCoeff Denominator coefficients of the filter in decreasing order.
* \param[out] bCoeff Numerator coefficients of the filter in decreasing order.
*/
void getCoeffs(vectX_t<T>& aCoeff, vectX_t<T>& bCoeff) const;
/*! \brief Return the current filter status. */
FilterStatus status() const noexcept { return m_status; }
/*! \brief Return the order the denominator polynome order of the filter. */
Eigen::Index aOrder() const noexcept { return m_aCoeff.size(); }
/*! \brief Return the order the numerator polynome order of the filter. */
Eigen::Index bOrder() const noexcept { return m_bCoeff.size(); }
protected:
/*! \brief Default uninitialized constructor. */
GenericFilter() = default;
/*! \brief Constructor.
* \param aCoeff Denominator coefficients of the filter in decreasing order.
* \param bCoeff Numerator coefficients of the filter in decreasing order.
*/
GenericFilter(const vectX_t<T>& aCoeff, const vectX_t<T>& bCoeff);
/*! \brief Default destructor. */
virtual ~GenericFilter() = default;
/*! \brief Set the new coefficients of the filters.
*
* It awaits a universal reference.
* \param aCoeff Denominator coefficients of the filter in decreasing order.
* \param bCoeff Numerator coefficients of the filter in decreasing order.
*/
template <typename T2>
void setCoeffs(T2&& aCoeff, T2&& bCoeff);
/*! \brief Normalized the filter coefficients such that aCoeff(0) = 1. */
void normalizeCoeffs();
/*! \brief Check for bad coefficients.
*
* Set the filter status to ready is everything is fine.
* \param aCoeff Denominator coefficients of the filter.
* \param bCoeff Numerator coefficients of the filter.
* \return True if the filter status is set on READY.
*/
bool checkCoeffs(const vectX_t<T>& aCoeff, const vectX_t<T>& bCoeff);
protected:
FilterStatus m_status;
FilterStatus m_status; /*!< Filter status */
private:
vectX_t<T> m_aCoeff;
vectX_t<T> m_bCoeff;
vectX_t<T> m_filteredData;
vectX_t<T> m_rawData;
vectX_t<T> m_aCoeff; /*!< Denominator coefficients of the filter */
vectX_t<T> m_bCoeff; /*!< Numerator coefficients of the filter */
vectX_t<T> m_filteredData; /*!< Last set of filtered data */
vectX_t<T> m_rawData; /*!< Last set of non-filtered data */
};
} // namespace fratio

Wyświetl plik

@ -5,15 +5,24 @@
namespace fratio {
/*! \brief Moving average digital filter.
*
* This is a specialization of a digital filter in order to use a moving average.
* \tparam T Floating type.
*/
template <typename T>
class MovingAverage : public DigitalFilter<T> {
public:
/*! \brief Default uninitialized constructor. */
MovingAverage() = default;
/*! \brief Constructor.
* \param windowSize Size of the moving average window.
*/
MovingAverage(int windowSize)
{
setWindowSize(windowSize);
}
/*! \brief Set the size of the moving average window. */
void setWindowSize(int windowSize)
{
if (windowSize <= 0) {
@ -23,6 +32,7 @@ public:
setCoeffs(vectX_t<T>::Constant(1, T(1)), vectX_t<T>::Constant(windowSize, T(1) / windowSize));
}
/*! \brief Get the size of the moving average window. */
int windowSize() const noexcept { return bOrder(); }
};

Wyświetl plik

@ -6,22 +6,31 @@
namespace fratio {
/*! \brief Compute polynome coefficients from roots.
*
* This is done through Vieta's algorithm: \see https://en.wikipedia.org/wiki/Vieta%27s_formulas
* \tparam T Floating type.
*/
template <typename T>
struct VietaAlgo {
static_assert(std::is_arithmetic<internal::complex_sub_type_t<T>>::value, "This struct can only accept arithmetic types or complex.");
// Vieta's computation: https://en.wikipedia.org/wiki/Vieta%27s_formulas
static vectX_t<T> polyCoeffFromRoot(const vectX_t<T>& poles);
/*! \brief Vieta's algorithm.
* \note The function return the coefficients in the decreasing order: \f$a_n X^n + a_{n-1}X^{n-1} + ... + a1X + a0\f$.
* \param roots Set of all roots of the polynome.
* \return Coefficients of the polynome.
*/
static vectX_t<T> polyCoeffFromRoot(const vectX_t<T>& roots);
};
template <typename T>
vectX_t<T> VietaAlgo<T>::polyCoeffFromRoot(const vectX_t<T>& poles)
vectX_t<T> VietaAlgo<T>::polyCoeffFromRoot(const vectX_t<T>& roots)
{
vectX_t<T> coeffs = vectX_t<T>::Zero(poles.size() + 1);
vectX_t<T> coeffs = vectX_t<T>::Zero(roots.size() + 1);
coeffs(0) = T(1);
for (Eigen::Index i = 0; i < poles.size(); ++i) {
for (Eigen::Index i = 0; i < roots.size(); ++i) {
for (Eigen::Index k = i + 1; k > 0; --k) {
coeffs(k) -= poles(i) * coeffs(k - 1);
coeffs(k) -= roots(i) * coeffs(k - 1);
}
}

Wyświetl plik

@ -5,24 +5,26 @@
namespace fratio {
template <typename T>
using vectX_t = Eigen::Matrix<T, Eigen::Dynamic, 1>;
template <typename T>
using vectXc_t = vectX_t<std::complex<T>>;
using vectX_t = Eigen::Matrix<T, Eigen::Dynamic, 1>; /*!< Eigen column-vector */
template <typename T>
using vectXc_t = vectX_t<std::complex<T>>; /*!< Eigen complex column-vector */
/*! \brief Filter status */
enum class FilterStatus {
// Generic filter
NONE,
READY,
BAD_ORDER_SIZE,
BAD_A_COEFF,
A_COEFF_MISSING,
B_COEFF_MISSING,
ALL_COEFF_MISSING = A_COEFF_MISSING | B_COEFF_MISSING,
NONE, /*!< Filter has not yet been initialized */
READY, /*!< Filter is ready to process data */
BAD_ORDER_SIZE, /*!< Order of the filter is bad */
BAD_A_COEFF, /*!< Denominator coefficients of the filter is bad */
A_COEFF_MISSING, /*!< Denominator coefficients have not been set */
B_COEFF_MISSING, /*!< Numerator coefficients have not been set */
ALL_COEFF_MISSING = A_COEFF_MISSING | B_COEFF_MISSING, /*!< Coefficients have not been set */
// Butterworth filter
BAD_FREQUENCY_VALUE,
BAD_CUTOFF_FREQUENCY,
BAD_BAND_FREQUENCY
BAD_FREQUENCY_VALUE, /*!< Given frequency is bad */
BAD_CUTOFF_FREQUENCY, /*!< Given ctu-off frequency is bad */
BAD_BAND_FREQUENCY /*!< Given band frequency is bad */
};
} // namespace fratio