Page MenuHomec4science

low_pass_filter.hh
No OneTemporary

File Metadata

Created
Thu, Jul 25, 09:19

low_pass_filter.hh

/**
* @file low_pass_filter.hh
*
* @author Guillaume Anciaux <guillaume.anciaux@epfl.ch>
*
* @date Mon Jan 07 16:36:30 2013
*
* @brief This implements a low pass filter
*
* @section LICENSE
*
* Copyright (©) 2010-2011 EPFL (Ecole Polytechnique Fédérale de Lausanne)
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* LibMultiScale is free software: you can redistribute it and/or modify it
* under the
* terms of the GNU Lesser General Public License as published by the Free
* Software Foundation, either version 3 of the License, or (at your option) any
* later version.
*
* LibMultiScale is distributed in the hope that it will be useful, but
* WITHOUT ANY
* WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
* A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more
* details.
*
* You should have received a copy of the GNU Lesser General Public License
* along with LibMultiScale. If not, see <http://www.gnu.org/licenses/>.
*
*/
#ifndef __LIBMULTISCALE_LOW_PASS_FILTER_HH__
#define __LIBMULTISCALE_LOW_PASS_FILTER_HH__
/* -------------------------------------------------------------------------- */
#include "function_interface.hh"
#include <gsl/gsl_sf_bessel.h>
/* -------------------------------------------------------------------------- */
__BEGIN_LIBMULTISCALE__
/** Class LowPassFilter
*
* This class implements access to the function
* @f$ F(x) = f_c \cdot \frac{J_1(2\pi f_c x)}{x} @f$
*
* where @f$J_1@f$ is the Bessel function of the first kind and
* @f$f_c = \frac{1}{\lambda}@f$.
*
* Also to define it numerically the limit when @f$x \to 0@f$ is:
*
* @f$F(x) \to f_c^2 \pi@f$
*/
template <UInt Dim> class LowPassFilter : public FunctionInterface {
/* ------------------------------------------------------------------------ */
/* Constructors/Destructors */
/* ------------------------------------------------------------------------ */
public:
LowPassFilter(Real cutoff_wavelength) {
fc = 1. / cutoff_wavelength;
phase = 2 * M_PI * fc;
};
virtual ~LowPassFilter(){};
/* ------------------------------------------------------------------------ */
/* Methods */
/* ------------------------------------------------------------------------ */
public:
template <UInt derivation_order> inline Real compute(Real x);
/* ------------------------------------------------------------------------ */
/* Accessors */
/* ------------------------------------------------------------------------ */
public:
/* ------------------------------------------------------------------------ */
/* Class Members */
/* ------------------------------------------------------------------------ */
private:
Real fc;
Real phase;
};
/* -------------------------------------------------------------------------- */
template <UInt Dim>
template <UInt derivation_order>
inline Real LowPassFilter<Dim>::compute(Real x) {
LM_FATAL("order " << derivation_order << "not defined yet");
return 0.0;
}
/* -------------------------------------------------------------------------- */
/** @f[
* F(x) = f_c \frac{J_1(2\pi f_c x)}{x}
* @f]
* which tends at the limit where @f$x \to 0@f$ to @f$ F(0) = f_c^2 \pi @f$
*/
template <> template <> inline Real LowPassFilter<1>::compute<0>(Real x) {
if (x) {
return sin(phase * x) / (M_PI * x);
}
return 2. * fc;
}
/* -------------------------------------------------------------------------- */
/** @f[
* F(x) = f_c \frac{J_1(2\pi f_c x)}{x}
* @f]
* which tends at the limit where @f$x \to 0@f$ to @f$ F(0) = f_c^2 \pi @f$
*/
template <> template <> inline Real LowPassFilter<2>::compute<0>(Real x) {
if (x) {
Real J1 = gsl_sf_bessel_Jn(1, phase * x);
return fc * J1 / x;
}
return fc * fc * M_PI;
}
/* -------------------------------------------------------------------------- */
/** @f[
* F'(x) = 2 f_c \left( \frac{\pi f_c}{x} J_0(2\pi f_c x)
* - \frac{J_1(2\pi f_c x)}{x^2} \right)
* @f]
*/
template <> template <> inline Real LowPassFilter<2>::compute<1>(Real x) {
if (x) {
Real J1 = gsl_sf_bessel_Jn(1, phase * x);
Real J0 = gsl_sf_bessel_Jn(0, phase * x);
return 2 * fc * (M_PI * fc / x * J0 - 1. / (x * x) * J1);
}
return 0.0;
}
/* -------------------------------------------------------------------------- */
/** @f[
* F''(x) = - \frac{2f_c}{x^3}
* \left(
* (2 \pi^2 f_c^2 x^2 - 3 )J_1(2\pi f_c x)
* + 3 J_0(2\pi f_c x) \pi f_c x
* \right)
*
* @f]
*/
template <> template <> inline Real LowPassFilter<2>::compute<2>(Real x) {
if (x) {
Real J1 = gsl_sf_bessel_Jn(1, phase * x);
Real J0 = gsl_sf_bessel_Jn(0, phase * x);
Real _1_x = 1. / x;
Real _1_x2 = _1_x * _1_x;
Real _1_x3 = _1_x2 * _1_x;
return -2. * fc * (J1 * (2. * M_PI * M_PI * fc * fc * _1_x - 3. * _1_x3) +
J0 * 3. * M_PI * fc * _1_x2);
}
return -1. * fc * fc * fc * fc * M_PI * M_PI * M_PI;
}
/* -------------------------------------------------------------------------- */
/** @f[
* F^{(3)}(x) =
* -4 f_c (2 Pi^3 f_c^3 x^3 J[0]-5 J[1] Pi^2 f_c^2 x^2+6 J[1]-6 J[0] Pi f_c
* x)/x^4
*
* @f]
*/
template <> template <> inline Real LowPassFilter<2>::compute<3>(Real x) {
if (x) {
Real J1 = gsl_sf_bessel_Jn(1, phase * x);
Real J0 = gsl_sf_bessel_Jn(0, phase * x);
return -4 * fc * (+J0 * (2 * M_PI * M_PI * M_PI * fc * fc * fc * x * x * x -
6 * M_PI * fc * x) +
J1 * (6 - 5 * M_PI * M_PI * fc * fc * x * x)) /
(x * x * x * x);
}
return 0;
}
/* -------------------------------------------------------------------------- */
__END_LIBMULTISCALE__
#endif /* __LIBMULTISCALE_LOW_PASS_FILTER_HH__ */

Event Timeline