Skip to content

Commit 5df7bca

Browse files
committed
Polynomial roots via eigenvalues of the companion matrix
1 parent a53b013 commit 5df7bca

File tree

5 files changed

+202
-3
lines changed

5 files changed

+202
-3
lines changed

doc/internals/polynomial.qbk

+12-2
Original file line numberDiff line numberDiff line change
@@ -48,7 +48,8 @@
4848

4949
T operator()(T z) const;
5050

51-
51+
// Only available if Eigen library is available:
52+
std::vector<std::complex<T>> roots() const;
5253

5354
// modify:
5455
void set_zero();
@@ -183,11 +184,20 @@ for output
183184

184185
The source code is at [@../../example/polynomial_arithmetic.cpp polynomial_arithmetic.cpp]
185186

186-
[endsect] [/section:polynomials Polynomials]
187+
[h4 Roots]
188+
189+
Polynomial roots are recovered by the eigenvalues of the companion matrix.
190+
This requires that the Eigen C++ library to be available; otherwise calling `.roots()` is a compile-time error.
191+
In addition, the polynomial coefficients must be of floating point type, or a static assertion will fail.
192+
193+
[endsect]
194+
195+
[/section:polynomials Polynomials]
187196

188197
[/
189198
Copyright 2006 John Maddock and Paul A. Bristow.
190199
Copyright 2015 Jeremy William Murphy.
200+
Copyright 2024 Nick Thompson.
191201

192202
Distributed under the Boost Software License, Version 1.0.
193203
(See accompanying file LICENSE_1_0.txt or copy at

include/boost/math/tools/polynomial.hpp

+70
Original file line numberDiff line numberDiff line change
@@ -1,5 +1,6 @@
11
// (C) Copyright John Maddock 2006.
22
// (C) Copyright Jeremy William Murphy 2015.
3+
// (C) Copyright Nick Thompson 2024.
34

45

56
// Use, modification and distribution are subject to the
@@ -29,6 +30,15 @@
2930
#include <type_traits>
3031
#include <iterator>
3132

33+
#if __has_include(<Eigen/Eigenvalues>)
34+
#include <Eigen/Eigenvalues>
35+
#include <complex> // roots are complex numbers.
36+
#if __has_include(<boost/multiprecision/eigen.hpp>)
37+
#include <boost/multiprecision/eigen.hpp>
38+
#include <boost/math/cstdfloat/cstdfloat_complex_std.hpp>
39+
#endif
40+
#endif
41+
3242
namespace boost{ namespace math{ namespace tools{
3343

3444
template <class T>
@@ -575,6 +585,66 @@ class polynomial
575585
m_data.erase(std::find_if(m_data.rbegin(), m_data.rend(), [](const T& x)->bool { return x != T(0); }).base(), m_data.end());
576586
}
577587

588+
#if __has_include(<Eigen/Eigenvalues>)
589+
/*
590+
* Polynomial root recovery by the eigenvalues of the companion matrix.
591+
* N.B.: This algorithm is not the state of the art; a faster algorithm is
592+
* "Fast and backward stable computation of roots of polynomials" by Aurentz et al.
593+
*/
594+
[[nodiscard]] auto roots() const {
595+
// At least as of Eigen 3.4.0, we cannot provide the eigensolver with complex numbers.
596+
using std::polar;
597+
using std::complex;
598+
using Complex = decltype(polar(T(), T()));
599+
if (m_data.size() == 1) {
600+
return std::vector<Complex>();
601+
}
602+
// There is a temptation to split off the linear and quadratic case, since
603+
// they are so easy. Resist the temptation! Your best unit tests will become
604+
// tautological.
605+
std::size_t n = m_data.size() - 1;
606+
Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic> C(n, n);
607+
C << Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic>::Zero(n,n);
608+
for (std::size_t i = 0; i < n; ++i) {
609+
// Remember the class invariant m_data.back() != 0 from the normalize() call?
610+
// Reaping blessings right here y'all:
611+
C(i, n - 1) = -m_data[i] / m_data.back();
612+
}
613+
for (std::size_t i = 0; i < n - 1; ++i) {
614+
C(i + 1, i) = 1;
615+
}
616+
Eigen::EigenSolver<decltype(C)> es;
617+
es.compute(C, /*computeEigenvectors=*/ false);
618+
auto info = es.info();
619+
if (info != Eigen::ComputationInfo::Success) {
620+
std::ostringstream oss;
621+
oss << __FILE__ << ":" << __LINE__ << ":" << __func__ << ": Eigen's eigensolver did not succeed.";
622+
switch (info) {
623+
case Eigen::ComputationInfo::NumericalIssue:
624+
oss << " Problem: numerical issue.";
625+
break;
626+
case Eigen::ComputationInfo::NoConvergence:
627+
oss << " Problem: no convergence.";
628+
break;
629+
case Eigen::ComputationInfo::InvalidInput:
630+
oss << " Problem: Invalid input.";
631+
break;
632+
default:
633+
oss << " Problem: Unknown.";
634+
}
635+
BOOST_MATH_THROW_EXCEPTION(std::runtime_error(oss.str()));
636+
}
637+
// Don't want to expose Eigen types to the rest of the world;
638+
// Eigen is a detail of this algorithm, so big sad copy:
639+
auto eigen_zeros = es.eigenvalues();
640+
std::vector<Complex> zeros(eigen_zeros.size());
641+
for (std::size_t i = 0; i < zeros.size(); ++i) {
642+
zeros[i] = eigen_zeros[i];
643+
}
644+
return zeros;
645+
}
646+
#endif
647+
578648
private:
579649
template <class U, class R>
580650
polynomial& addition(const U& value, R op)

test/Jamfile.v2

+1
Original file line numberDiff line numberDiff line change
@@ -1068,6 +1068,7 @@ test-suite interpolators :
10681068
[ run jso_test.cpp : : : [ requires cxx17_if_constexpr cxx17_std_apply ] ]
10691069
[ run random_search_test.cpp : : : [ requires cxx17_if_constexpr cxx17_std_apply ] ]
10701070
[ run cma_es_test.cpp : : : [ requires cxx17_if_constexpr cxx17_std_apply ] [ check-target-builds ../../multiprecision/config//has_eigen : : <build>no ] ]
1071+
[ run polynomial_roots_test.cpp : : : [ requires cxx17_if_constexpr cxx17_std_apply ] [ check-target-builds ../../multiprecision/config//has_eigen : : <build>no ] ]
10711072
[ compile compile_test/random_search_incl_test.cpp : [ requires cxx17_if_constexpr cxx17_std_apply ] ]
10721073
[ compile compile_test/differential_evolution_incl_test.cpp : [ requires cxx17_if_constexpr cxx17_std_apply ] ]
10731074
[ compile compile_test/jso_incl_test.cpp : [ requires cxx17_if_constexpr cxx17_std_apply ] ]

test/polynomial_roots_test.cpp

+118
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,118 @@
1+
/*
2+
* Copyright Nick Thompson, 2024
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 <vector>
8+
#include <iostream>
9+
#include <list>
10+
#include <random>
11+
#include <cmath>
12+
#include <complex>
13+
#include <utility>
14+
#include <limits>
15+
#include <algorithm>
16+
#include <boost/math/tools/polynomial.hpp>
17+
using boost::math::tools::polynomial;
18+
#ifdef BOOST_HAS_FLOAT128
19+
#include <boost/multiprecision/float128.hpp>
20+
using boost::multiprecision::float128;
21+
#endif
22+
#include <boost/multiprecision/cpp_bin_float.hpp>
23+
#include "math_unit_test.hpp"
24+
25+
#if __has_include(<Eigen/Eigenvalues>)
26+
27+
void test_random_coefficients() {
28+
std::random_device rd;
29+
uint32_t seed = rd();
30+
std::mt19937_64 mt(seed);
31+
std::uniform_real_distribution<double> unif(-1, 1);
32+
std::size_t n = seed % 3 + 3;
33+
std::vector<double> coeffs(n, std::numeric_limits<double>::quiet_NaN());
34+
for (std::size_t i = 0; i < coeffs.size(); ++i) {
35+
coeffs[i] = unif(mt);
36+
}
37+
coeffs[coeffs.size() -1] = 1.0;
38+
auto p = polynomial(std::move(coeffs));
39+
auto roots = p.roots();
40+
CHECK_EQUAL(roots.size(), p.size() - 1);
41+
std::complex<double> root_product = -1;
42+
std::complex<double> root_sum = 0.0;
43+
for (auto const & root : roots) {
44+
root_product *= static_cast<std::complex<double>>(root);
45+
root_sum += static_cast<std::complex<double>>(root);
46+
}
47+
if (p.size() & 1) {
48+
root_product *= -1;
49+
}
50+
CHECK_ULP_CLOSE(root_product.real(), p[0], 1000);
51+
CHECK_LE(root_product.imag(), 1e-6);
52+
53+
CHECK_LE(root_sum.imag(), 1e-7);
54+
CHECK_ULP_CLOSE(root_sum.real(), -p[p.size() - 2], 1000);
55+
}
56+
57+
58+
59+
void test_wilkinson_polynomial() {
60+
// CoefficientList[Expand[(x - 1)*(x - 2)*(x - 3)*(x - 4)*(x - 5)*(x - 6)*(x - 7)*(x - 8)*(x - 9)*(x - 10)], x]
61+
std::vector<float> coeffs{3628800.0, -10628640.0, 12753576.0, -8409500.0, 3416930.0, -902055.0, 157773.0, -18150.0, 1320.0, -55.0 ,1.0};
62+
auto p = polynomial(std::move(coeffs));
63+
auto roots = p.roots();
64+
CHECK_EQUAL(roots.size(), p.size() - 1);
65+
std::complex<double> root_product = -1;
66+
std::complex<double> root_sum = 0.0;
67+
for (auto const & root : roots) {
68+
root_product *= static_cast<std::complex<double>>(root);
69+
root_sum += static_cast<std::complex<double>>(root);
70+
}
71+
if (p.size() & 1) {
72+
root_product *= -1;
73+
}
74+
CHECK_ABSOLUTE_ERROR(root_product.real(), double(p[0]), double(1e-3*p[0]));
75+
CHECK_LE(root_product.imag(), 1e-8);
76+
77+
CHECK_LE(root_sum.imag(), 1e-8);
78+
CHECK_ABSOLUTE_ERROR(root_sum.real(), -double(p[p.size() - 2]), 1e-5);
79+
80+
std::complex<double> c = 0.0;
81+
for (std::size_t i = 0; i < roots.size(); ++i) {
82+
auto ri = static_cast<std::complex<double>>(roots[i]);
83+
for (std::size_t j = i + 1; j < roots.size(); ++j) {
84+
c += ri*static_cast<std::complex<double>>(roots[j]);
85+
}
86+
}
87+
CHECK_ULP_CLOSE(p[p.size()-3], static_cast<float>(c.real()), 10);
88+
CHECK_ABSOLUTE_ERROR(0.0, c.imag(), 1e-8);
89+
90+
}
91+
92+
template<typename T>
93+
void test_singular_companion()
94+
{
95+
std::vector<T> coeffs{0.0, 0.0, 1.0};
96+
auto p = polynomial(std::move(coeffs));
97+
auto roots = p.roots();
98+
CHECK_EQUAL(roots.size(), p.size() - 1);
99+
for (std::size_t i = 0; i < roots.size() - 1; ++i) {
100+
CHECK_ABSOLUTE_ERROR(T(0), roots[i].real(), std::numeric_limits<T>::epsilon());
101+
CHECK_ABSOLUTE_ERROR(T(0), roots[i].imag(), std::numeric_limits<T>::epsilon());
102+
}
103+
}
104+
105+
106+
int main()
107+
{
108+
test_random_coefficients();
109+
test_singular_companion<float>();
110+
test_singular_companion<double>();
111+
#if BOOST_HAS_FLOAT128
112+
test_singular_companion<float128>();
113+
#endif
114+
test_singular_companion<boost::multiprecision::cpp_bin_float_100>();
115+
test_wilkinson_polynomial();
116+
return boost::math::test::report_errors();
117+
}
118+
#endif

test/test_autodiff_5.cpp

+1-1
Original file line numberDiff line numberDiff line change
@@ -58,7 +58,7 @@ BOOST_AUTO_TEST_CASE_TEMPLATE(chebyshev_hpp, T, all_float_types) {
5858
auto x = x_sampler.next();
5959
BOOST_CHECK_CLOSE(
6060
boost::math::chebyshev_t(n, make_fvar<T, m>(x)).derivative(0u),
61-
boost::math::chebyshev_t(n, x), 40 * test_constants::pct_epsilon());
61+
boost::math::chebyshev_t(n, x), 80 * test_constants::pct_epsilon());
6262
// Lower accuracy with clang/apple:
6363
BOOST_CHECK_CLOSE(
6464
boost::math::chebyshev_u(n, make_fvar<T, m>(x)).derivative(0u),

0 commit comments

Comments
 (0)