chebyshev_transform.hpp 6.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237
  1. // (C) Copyright Nick Thompson 2017.
  2. // Use, modification and distribution are subject to the
  3. // Boost Software License, Version 1.0. (See accompanying file
  4. // LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
  5. #ifndef BOOST_MATH_SPECIAL_CHEBYSHEV_TRANSFORM_HPP
  6. #define BOOST_MATH_SPECIAL_CHEBYSHEV_TRANSFORM_HPP
  7. #include <cmath>
  8. #include <type_traits>
  9. #include <boost/math/constants/constants.hpp>
  10. #include <boost/math/special_functions/chebyshev.hpp>
  11. #ifdef BOOST_HAS_FLOAT128
  12. #include <quadmath.h>
  13. #endif
  14. #ifdef __has_include
  15. # if __has_include(<fftw3.h>)
  16. # include <fftw3.h>
  17. # else
  18. # error "This feature is unavailable without fftw3 installed"
  19. #endif
  20. #endif
  21. namespace boost { namespace math {
  22. namespace detail{
  23. template <class T>
  24. struct fftw_cos_transform;
  25. template<>
  26. struct fftw_cos_transform<double>
  27. {
  28. fftw_cos_transform(int n, double* data1, double* data2)
  29. {
  30. plan = fftw_plan_r2r_1d(n, data1, data2, FFTW_REDFT10, FFTW_ESTIMATE);
  31. }
  32. ~fftw_cos_transform()
  33. {
  34. fftw_destroy_plan(plan);
  35. }
  36. void execute(double* data1, double* data2)
  37. {
  38. fftw_execute_r2r(plan, data1, data2);
  39. }
  40. static double cos(double x) { return std::cos(x); }
  41. static double fabs(double x) { return std::fabs(x); }
  42. private:
  43. fftw_plan plan;
  44. };
  45. template<>
  46. struct fftw_cos_transform<float>
  47. {
  48. fftw_cos_transform(int n, float* data1, float* data2)
  49. {
  50. plan = fftwf_plan_r2r_1d(n, data1, data2, FFTW_REDFT10, FFTW_ESTIMATE);
  51. }
  52. ~fftw_cos_transform()
  53. {
  54. fftwf_destroy_plan(plan);
  55. }
  56. void execute(float* data1, float* data2)
  57. {
  58. fftwf_execute_r2r(plan, data1, data2);
  59. }
  60. static float cos(float x) { return std::cos(x); }
  61. static float fabs(float x) { return std::fabs(x); }
  62. private:
  63. fftwf_plan plan;
  64. };
  65. template<>
  66. struct fftw_cos_transform<long double>
  67. {
  68. fftw_cos_transform(int n, long double* data1, long double* data2)
  69. {
  70. plan = fftwl_plan_r2r_1d(n, data1, data2, FFTW_REDFT10, FFTW_ESTIMATE);
  71. }
  72. ~fftw_cos_transform()
  73. {
  74. fftwl_destroy_plan(plan);
  75. }
  76. void execute(long double* data1, long double* data2)
  77. {
  78. fftwl_execute_r2r(plan, data1, data2);
  79. }
  80. static long double cos(long double x) { return std::cos(x); }
  81. static long double fabs(long double x) { return std::fabs(x); }
  82. private:
  83. fftwl_plan plan;
  84. };
  85. #ifdef BOOST_HAS_FLOAT128
  86. template<>
  87. struct fftw_cos_transform<__float128>
  88. {
  89. fftw_cos_transform(int n, __float128* data1, __float128* data2)
  90. {
  91. plan = fftwq_plan_r2r_1d(n, data1, data2, FFTW_REDFT10, FFTW_ESTIMATE);
  92. }
  93. ~fftw_cos_transform()
  94. {
  95. fftwq_destroy_plan(plan);
  96. }
  97. void execute(__float128* data1, __float128* data2)
  98. {
  99. fftwq_execute_r2r(plan, data1, data2);
  100. }
  101. static __float128 cos(__float128 x) { return cosq(x); }
  102. static __float128 fabs(__float128 x) { return fabsq(x); }
  103. private:
  104. fftwq_plan plan;
  105. };
  106. #endif
  107. }
  108. template<class Real>
  109. class chebyshev_transform
  110. {
  111. public:
  112. template<class F>
  113. chebyshev_transform(const F& f, Real a, Real b,
  114. Real tol = 500 * std::numeric_limits<Real>::epsilon(),
  115. size_t max_refinements = 16) : m_a(a), m_b(b)
  116. {
  117. if (a >= b)
  118. {
  119. throw std::domain_error("a < b is required.\n");
  120. }
  121. using boost::math::constants::half;
  122. using boost::math::constants::pi;
  123. using std::cos;
  124. using std::abs;
  125. Real bma = (b-a)*half<Real>();
  126. Real bpa = (b+a)*half<Real>();
  127. size_t n = 256;
  128. std::vector<Real> vf;
  129. size_t refinements = 0;
  130. while(refinements < max_refinements)
  131. {
  132. vf.resize(n);
  133. m_coeffs.resize(n);
  134. detail::fftw_cos_transform<Real> plan(static_cast<int>(n), vf.data(), m_coeffs.data());
  135. Real inv_n = 1/static_cast<Real>(n);
  136. for(size_t j = 0; j < n/2; ++j)
  137. {
  138. // Use symmetry cos((j+1/2)pi/n) = - cos((n-1-j+1/2)pi/n)
  139. Real y = detail::fftw_cos_transform<Real>::cos(pi<Real>()*(j+half<Real>())*inv_n);
  140. vf[j] = f(y*bma + bpa)*inv_n;
  141. vf[n-1-j]= f(bpa-y*bma)*inv_n;
  142. }
  143. plan.execute(vf.data(), m_coeffs.data());
  144. Real max_coeff = 0;
  145. for (auto const & coeff : m_coeffs)
  146. {
  147. if (detail::fftw_cos_transform<Real>::fabs(coeff) > max_coeff)
  148. {
  149. max_coeff = detail::fftw_cos_transform<Real>::fabs(coeff);
  150. }
  151. }
  152. size_t j = m_coeffs.size() - 1;
  153. while (abs(m_coeffs[j])/max_coeff < tol)
  154. {
  155. --j;
  156. }
  157. // If ten coefficients are eliminated, the we say we've done all
  158. // we need to do:
  159. if (n - j > 10)
  160. {
  161. m_coeffs.resize(j+1);
  162. return;
  163. }
  164. n *= 2;
  165. ++refinements;
  166. }
  167. }
  168. inline Real operator()(Real x) const
  169. {
  170. return chebyshev_clenshaw_recurrence(m_coeffs.data(), m_coeffs.size(), m_a, m_b, x);
  171. }
  172. // Integral over entire domain [a, b]
  173. Real integrate() const
  174. {
  175. Real Q = m_coeffs[0]/2;
  176. for(size_t j = 2; j < m_coeffs.size(); j += 2)
  177. {
  178. Q += -m_coeffs[j]/((j+1)*(j-1));
  179. }
  180. return (m_b - m_a)*Q;
  181. }
  182. const std::vector<Real>& coefficients() const
  183. {
  184. return m_coeffs;
  185. }
  186. Real prime(Real x) const
  187. {
  188. Real z = (2*x - m_a - m_b)/(m_b - m_a);
  189. Real dzdx = 2/(m_b - m_a);
  190. if (m_coeffs.size() < 2)
  191. {
  192. return 0;
  193. }
  194. Real b2 = 0;
  195. Real d2 = 0;
  196. Real b1 = m_coeffs[m_coeffs.size() -1];
  197. Real d1 = 0;
  198. for(size_t j = m_coeffs.size() - 2; j >= 1; --j)
  199. {
  200. Real tmp1 = 2*z*b1 - b2 + m_coeffs[j];
  201. Real tmp2 = 2*z*d1 - d2 + 2*b1;
  202. b2 = b1;
  203. b1 = tmp1;
  204. d2 = d1;
  205. d1 = tmp2;
  206. }
  207. return dzdx*(z*d1 - d2 + b1);
  208. }
  209. private:
  210. std::vector<Real> m_coeffs;
  211. Real m_a;
  212. Real m_b;
  213. };
  214. }}
  215. #endif