Skip to content

Commit 04bd67a

Browse files
committed
Polynomial roots via eigenvalues of the companion matrix
1 parent a53b013 commit 04bd67a

File tree

5 files changed

+190
-3
lines changed

5 files changed

+190
-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

+67
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,11 @@
2930
#include <type_traits>
3031
#include <iterator>
3132

33+
#if __has_include(<Eigen/Eigenvalues>)
34+
#include <complex> // roots are complex numbers.
35+
#include <Eigen/Eigenvalues>
36+
#endif
37+
3238
namespace boost{ namespace math{ namespace tools{
3339

3440
template <class T>
@@ -575,6 +581,67 @@ class polynomial
575581
m_data.erase(std::find_if(m_data.rbegin(), m_data.rend(), [](const T& x)->bool { return x != T(0); }).base(), m_data.end());
576582
}
577583

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

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