daubechies_coefficients.cpp 8.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221
  1. /*
  2. * Copyright Nick Thompson, 2018
  3. * Use, modification and distribution are subject to the
  4. * Boost Software License, Version 1.0. (See accompanying file
  5. * LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
  6. */
  7. #include <iostream>
  8. #include <vector>
  9. #include <string>
  10. #include <complex>
  11. #include <bitset>
  12. #include <boost/assert.hpp>
  13. #include <boost/multiprecision/cpp_bin_float.hpp>
  14. #include <boost/math/constants/constants.hpp>
  15. #include <boost/math/tools/polynomial.hpp>
  16. #include <boost/math/tools/roots.hpp>
  17. #include <boost/math/special_functions/binomial.hpp>
  18. #include <boost/multiprecision/cpp_complex.hpp>
  19. #include <boost/multiprecision/complex128.hpp>
  20. #include <boost/math/quadrature/gauss_kronrod.hpp>
  21. using std::string;
  22. using boost::math::tools::polynomial;
  23. using boost::math::binomial_coefficient;
  24. using boost::math::tools::schroder_iterate;
  25. using boost::math::tools::halley_iterate;
  26. using boost::math::tools::newton_raphson_iterate;
  27. using boost::math::tools::complex_newton;
  28. using boost::math::constants::half;
  29. using boost::math::constants::root_two;
  30. using boost::math::constants::pi;
  31. using boost::math::quadrature::gauss_kronrod;
  32. using boost::multiprecision::cpp_bin_float_100;
  33. using boost::multiprecision::cpp_complex_100;
  34. template<class Complex>
  35. std::vector<std::pair<Complex, Complex>> find_roots(size_t p)
  36. {
  37. // Initialize the polynomial; see Mallat, A Wavelet Tour of Signal Processing, equation 7.96
  38. BOOST_ASSERT(p>0);
  39. typedef typename Complex::value_type Real;
  40. std::vector<Complex> coeffs(p);
  41. for (size_t k = 0; k < coeffs.size(); ++k)
  42. {
  43. coeffs[k] = Complex(binomial_coefficient<Real>(p-1+k, k), 0);
  44. }
  45. polynomial<Complex> P(std::move(coeffs));
  46. polynomial<Complex> Pcopy = P;
  47. polynomial<Complex> Pcopy_prime = P.prime();
  48. auto orig = [&](Complex z) { return std::make_pair<Complex, Complex>(Pcopy(z), Pcopy_prime(z)); };
  49. polynomial<Complex> P_prime = P.prime();
  50. // Polynomial is of degree p-1.
  51. std::vector<Complex> roots(p-1, {std::numeric_limits<Real>::quiet_NaN(),std::numeric_limits<Real>::quiet_NaN()});
  52. size_t i = 0;
  53. while(P.size() > 1)
  54. {
  55. Complex guess = {0.0, 1.0};
  56. std::cout << std::setprecision(std::numeric_limits<Real>::digits10+3);
  57. auto f = [&](Complex x)->std::pair<Complex, Complex>
  58. {
  59. return std::make_pair<Complex, Complex>(P(x), P_prime(x));
  60. };
  61. Complex r = complex_newton(f, guess);
  62. using std::isnan;
  63. if(isnan(r.real()))
  64. {
  65. int i = 50;
  66. do {
  67. // Try a different guess
  68. guess *= Complex(1.0,-1.0);
  69. r = complex_newton(f, guess);
  70. std::cout << "New guess: " << guess << ", result? " << r << std::endl;
  71. } while (isnan(r.real()) && i-- > 0);
  72. if (isnan(r.real()))
  73. {
  74. std::cout << "Polynomial that killed the process: " << P << std::endl;
  75. throw std::logic_error("Newton iteration did not converge");
  76. }
  77. }
  78. // Refine r with the original function.
  79. // We only use the polynomial division to ensure we don't get the same root over and over.
  80. // However, the division induces error which can grow quickly-or slowly! See Numerical Recipes, section 9.5.1.
  81. r = complex_newton(orig, r);
  82. if (isnan(r.real()))
  83. {
  84. throw std::logic_error("Found a root for the deflated polynomial which is not a root for the original. Indicative of catastrophic numerical error.");
  85. }
  86. // Test the root:
  87. using std::sqrt;
  88. Real tol = sqrt(sqrt(std::numeric_limits<Real>::epsilon()));
  89. if (norm(Pcopy(r)) > tol)
  90. {
  91. std::cout << "This is a bad root: P" << r << " = " << Pcopy(r) << std::endl;
  92. std::cout << "Reduced polynomial leading to bad root: " << P << std::endl;
  93. throw std::logic_error("Donezo.");
  94. }
  95. BOOST_ASSERT(i < roots.size());
  96. roots[i] = r;
  97. ++i;
  98. polynomial<Complex> q{-r, {1,0}};
  99. // This optimization breaks at p = 11. I have no clue why.
  100. // Unfortunate, because I expect it to be considerably more stable than
  101. // repeatedly dividing by the complex root.
  102. /*polynomial<Complex> q;
  103. if (r.imag() > sqrt(std::numeric_limits<Real>::epsilon()))
  104. {
  105. // Then the complex conjugate is also a root:
  106. using std::conj;
  107. using std::norm;
  108. BOOST_ASSERT(i < roots.size());
  109. roots[i] = conj(r);
  110. ++i;
  111. q = polynomial<Complex>({{norm(r), 0}, {-2*r.real(),0}, {1,0}});
  112. }
  113. else
  114. {
  115. // The imaginary part is numerical noise:
  116. r.imag() = 0;
  117. q = polynomial<Complex>({-r, {1,0}});
  118. }*/
  119. auto PR = quotient_remainder(P, q);
  120. // I should validate that the remainder is small, but . . .
  121. //std::cout << "Remainder = " << PR.second<< std::endl;
  122. P = PR.first;
  123. P_prime = P.prime();
  124. }
  125. std::vector<std::pair<Complex, Complex>> Qroots(p-1);
  126. for (size_t i = 0; i < Qroots.size(); ++i)
  127. {
  128. Complex y = roots[i];
  129. Complex z1 = static_cast<Complex>(1) - static_cast<Complex>(2)*y + static_cast<Complex>(2)*sqrt(y*(y-static_cast<Complex>(1)));
  130. Complex z2 = static_cast<Complex>(1) - static_cast<Complex>(2)*y - static_cast<Complex>(2)*sqrt(y*(y-static_cast<Complex>(1)));
  131. Qroots[i] = {z1, z2};
  132. }
  133. return Qroots;
  134. }
  135. template<class Complex>
  136. std::vector<typename Complex::value_type> daubechies_coefficients(std::vector<std::pair<Complex, Complex>> const & Qroots)
  137. {
  138. typedef typename Complex::value_type Real;
  139. size_t p = Qroots.size() + 1;
  140. // Choose the minimum abs root; see Mallat, discussion just after equation 7.98
  141. std::vector<Complex> chosen_roots(p-1);
  142. for (size_t i = 0; i < p - 1; ++i)
  143. {
  144. if(norm(Qroots[i].first) <= 1)
  145. {
  146. chosen_roots[i] = Qroots[i].first;
  147. }
  148. else
  149. {
  150. BOOST_ASSERT(norm(Qroots[i].second) <= 1);
  151. chosen_roots[i] = Qroots[i].second;
  152. }
  153. }
  154. polynomial<Complex> R{1};
  155. for (size_t i = 0; i < p-1; ++i)
  156. {
  157. Complex ak = chosen_roots[i];
  158. R *= polynomial<Complex>({-ak/(static_cast<Complex>(1)-ak), static_cast<Complex>(1)/(static_cast<Complex>(1)-ak)});
  159. }
  160. polynomial<Complex> a{{half<Real>(), 0}, {half<Real>(),0}};
  161. polynomial<Complex> poly = root_two<Real>()*pow(a, p)*R;
  162. std::vector<Complex> result = poly.data();
  163. // If we reverse, we get the Numerical Recipes and Daubechies convention.
  164. // If we don't reverse, we get the Pywavelets and Mallat convention.
  165. // I believe this is because of the sign convention on the DFT, which differs between Daubechies and Mallat.
  166. // You implement a dot product in Daubechies/NR convention, and a convolution in PyWavelets/Mallat convention.
  167. // I won't reverse so I can spot check against Pywavelets: http://wavelets.pybytes.com/wavelet/
  168. //std::reverse(result.begin(), result.end());
  169. std::vector<Real> h(result.size());
  170. for (size_t i = 0; i < result.size(); ++i)
  171. {
  172. Complex r = result[i];
  173. BOOST_ASSERT(r.imag() < sqrt(std::numeric_limits<Real>::epsilon()));
  174. h[i] = r.real();
  175. }
  176. // Quick sanity check: We could check all vanishing moments, but that sum is horribly ill-conditioned too!
  177. Real sum = 0;
  178. Real scale = 0;
  179. for (size_t i = 0; i < h.size(); ++i)
  180. {
  181. sum += h[i];
  182. scale += h[i]*h[i];
  183. }
  184. BOOST_ASSERT(abs(scale -1) < sqrt(std::numeric_limits<Real>::epsilon()));
  185. BOOST_ASSERT(abs(sum - root_two<Real>()) < sqrt(std::numeric_limits<Real>::epsilon()));
  186. return h;
  187. }
  188. int main()
  189. {
  190. typedef boost::multiprecision::cpp_complex<100> Complex;
  191. for(size_t p = 1; p < 200; ++p)
  192. {
  193. auto roots = find_roots<Complex>(p);
  194. auto h = daubechies_coefficients(roots);
  195. std::cout << "h_" << p << "[] = {";
  196. for (auto& x : h) {
  197. std::cout << x << ", ";
  198. }
  199. std::cout << "} // = h_" << p << "\n\n\n\n";
  200. }
  201. }