Eigen::PolynomialSolver< Scalar_, Deg_ > Class Template Reference

A polynomial solver. More...

+ Inheritance diagram for Eigen::PolynomialSolver< Scalar_, Deg_ >:

Public Types

typedef Matrix< Scalar, Deg_, Deg_ > CompanionMatrixType
 
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, Scalar, std::complex< Scalar > > ComplexScalar
 
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, ComplexEigenSolver< CompanionMatrixType >, EigenSolver< CompanionMatrixType > > EigenSolverType
 
typedef PolynomialSolverBase< Scalar_, Deg_ > PS_Base
 
- Public Types inherited from Eigen::PolynomialSolverBase< Scalar_, Deg_ >
typedef DenseIndex Index
 
typedef NumTraits< Scalar >::Real RealScalar
 
typedef Matrix< RootType, Deg_, 1 > RootsType
 
typedef std::complex< RealScalarRootType
 
typedef Scalar_ Scalar
 

Public Member Functions

template<typename OtherPolynomial >
void compute (const OtherPolynomial &poly)
 
 PolynomialSolver ()
 
template<typename OtherPolynomial >
 PolynomialSolver (const OtherPolynomial &poly)
 
- Public Member Functions inherited from Eigen::PolynomialSolverBase< Scalar_, Deg_ >
const RealScalarabsGreatestRealRoot (bool &hasArealRoot, const RealScalar &absImaginaryThreshold=NumTraits< Scalar >::dummy_precision()) const
 
const RealScalarabsSmallestRealRoot (bool &hasArealRoot, const RealScalar &absImaginaryThreshold=NumTraits< Scalar >::dummy_precision()) const
 
const RealScalargreatestRealRoot (bool &hasArealRoot, const RealScalar &absImaginaryThreshold=NumTraits< Scalar >::dummy_precision()) const
 
const RootTypegreatestRoot () const
 
 PolynomialSolverBase ()
 
template<typename OtherPolynomial >
 PolynomialSolverBase (const OtherPolynomial &poly)
 
template<typename Stl_back_insertion_sequence >
void realRoots (Stl_back_insertion_sequence &bi_seq, const RealScalar &absImaginaryThreshold=NumTraits< Scalar >::dummy_precision()) const
 
const RootsTyperoots () const
 
const RealScalarsmallestRealRoot (bool &hasArealRoot, const RealScalar &absImaginaryThreshold=NumTraits< Scalar >::dummy_precision()) const
 
const RootTypesmallestRoot () const
 

Protected Attributes

EigenSolverType m_eigenSolver
 
RootsType m_roots
 
- Protected Attributes inherited from Eigen::PolynomialSolverBase< Scalar_, Deg_ >
RootsType m_roots
 

Additional Inherited Members

- Protected Member Functions inherited from Eigen::PolynomialSolverBase< Scalar_, Deg_ >
template<typename squaredNormBinaryPredicate >
const RootTypeselectComplexRoot_withRespectToNorm (squaredNormBinaryPredicate &pred) const
 
template<typename squaredRealPartBinaryPredicate >
const RealScalarselectRealRoot_withRespectToAbsRealPart (squaredRealPartBinaryPredicate &pred, bool &hasArealRoot, const RealScalar &absImaginaryThreshold=NumTraits< Scalar >::dummy_precision()) const
 
template<typename RealPartBinaryPredicate >
const RealScalarselectRealRoot_withRespectToRealPart (RealPartBinaryPredicate &pred, bool &hasArealRoot, const RealScalar &absImaginaryThreshold=NumTraits< Scalar >::dummy_precision()) const
 
template<typename OtherPolynomial >
void setPolynomial (const OtherPolynomial &poly)
 

Detailed Description

template<typename Scalar_, int Deg_>
class Eigen::PolynomialSolver< Scalar_, Deg_ >

A polynomial solver.

Computes the complex roots of a real polynomial.

Parameters
Scalar_the scalar type, i.e., the type of the polynomial coefficients
Deg_the degree of the polynomial, can be a compile time value or Dynamic. Notice that the number of polynomial coefficients is Deg_+1.

This class implements a polynomial solver and provides convenient methods such as

  • real roots,
  • greatest, smallest complex roots,
  • real roots with greatest, smallest absolute real value.
  • greatest, smallest real roots.

WARNING: this polynomial solver is experimental, part of the unsupported Eigen modules.

Currently a QR algorithm is used to compute the eigenvalues of the companion matrix of the polynomial to compute its roots. This supposes that the complex moduli of the roots are all distinct: e.g. there should be no multiple roots or conjugate roots for instance. With 32bit (float) floating types this problem shows up frequently. However, almost always, correct accuracy is reached even in these cases for 64bit (double) floating types and small polynomial degree (<20).

Definition at line 333 of file PolynomialSolver.h.

Member Typedef Documentation

◆ CompanionMatrixType

template<typename Scalar_ , int Deg_>
typedef Matrix<Scalar,Deg_,Deg_> Eigen::PolynomialSolver< Scalar_, Deg_ >::CompanionMatrixType

Definition at line 341 of file PolynomialSolver.h.

◆ ComplexScalar

template<typename Scalar_ , int Deg_>
typedef std::conditional_t<NumTraits<Scalar>::IsComplex, Scalar, std::complex<Scalar> > Eigen::PolynomialSolver< Scalar_, Deg_ >::ComplexScalar

Definition at line 345 of file PolynomialSolver.h.

◆ EigenSolverType

template<typename Scalar_ , int Deg_>
typedef std::conditional_t<NumTraits<Scalar>::IsComplex, ComplexEigenSolver<CompanionMatrixType>, EigenSolver<CompanionMatrixType> > Eigen::PolynomialSolver< Scalar_, Deg_ >::EigenSolverType

Definition at line 344 of file PolynomialSolver.h.

◆ PS_Base

template<typename Scalar_ , int Deg_>
typedef PolynomialSolverBase<Scalar_,Deg_> Eigen::PolynomialSolver< Scalar_, Deg_ >::PS_Base

Definition at line 338 of file PolynomialSolver.h.

Constructor & Destructor Documentation

◆ PolynomialSolver() [1/2]

template<typename Scalar_ , int Deg_>
template<typename OtherPolynomial >
Eigen::PolynomialSolver< Scalar_, Deg_ >::PolynomialSolver ( const OtherPolynomial &  poly)
inline

Definition at line 389 of file PolynomialSolver.h.

389  {
390  compute( poly ); }
void compute(const OtherPolynomial &poly)

◆ PolynomialSolver() [2/2]

template<typename Scalar_ , int Deg_>
Eigen::PolynomialSolver< Scalar_, Deg_ >::PolynomialSolver ( )
inline

Definition at line 392 of file PolynomialSolver.h.

392 {}

Member Function Documentation

◆ compute()

template<typename Scalar_ , int Deg_>
template<typename OtherPolynomial >
void Eigen::PolynomialSolver< Scalar_, Deg_ >::compute ( const OtherPolynomial &  poly)
inline

Computes the complex roots of a new polynomial.

Definition at line 350 of file PolynomialSolver.h.

351  {
352  eigen_assert( Scalar(0) != poly[poly.size()-1] );
353  eigen_assert( poly.size() > 1 );
354  if(poly.size() > 2 )
355  {
356  internal::companion<Scalar,Deg_> companion( poly );
357  companion.balance();
358  m_eigenSolver.compute( companion.denseMatrix() );
359  m_roots = m_eigenSolver.eigenvalues();
360  // cleanup noise in imaginary part of real roots:
361  // if the imaginary part is rather small compared to the real part
362  // and that cancelling the imaginary part yield a smaller evaluation,
363  // then it's safe to keep the real part only.
364  RealScalar coarse_prec = RealScalar(std::pow(4,poly.size()+1))*NumTraits<RealScalar>::epsilon();
365  for(Index i = 0; i<m_roots.size(); ++i)
366  {
369  coarse_prec) )
370  {
371  ComplexScalar as_real_root = ComplexScalar(numext::real(m_roots[i]));
372  if( numext::abs(poly_eval(poly, as_real_root))
373  <= numext::abs(poly_eval(poly, m_roots[i])))
374  {
375  m_roots[i] = as_real_root;
376  }
377  }
378  }
379  }
380  else if(poly.size () == 2)
381  {
382  m_roots.resize(1);
383  m_roots[0] = -poly[0]/poly[1];
384  }
385  }
int i
#define eigen_assert(x)
constexpr void resize(Index rows, Index cols)
NumTraits< Scalar >::Real RealScalar
EigenSolverType m_eigenSolver
std::conditional_t< NumTraits< Scalar >::IsComplex, Scalar, std::complex< Scalar > > ComplexScalar
T poly_eval(const Polynomials &poly, const T &x)
bool isMuchSmallerThan(const Scalar &x, const OtherScalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_precision())
EIGEN_ALWAYS_INLINE std::enable_if_t< NumTraits< T >::IsSigned||NumTraits< T >::IsComplex, typename NumTraits< T >::Real > abs(const T &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
Eigen::AutoDiffScalar< EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE(internal::remove_all_t< DerType >, typename internal::traits< internal::remove_all_t< DerType >>::Scalar, product) > pow(const Eigen::AutoDiffScalar< DerType > &x, const typename internal::traits< internal::remove_all_t< DerType >>::Scalar &y)
const adouble & real(const adouble &x)
Definition: AdolcForward:72
adouble imag(const adouble &)
Definition: AdolcForward:73

Member Data Documentation

◆ m_eigenSolver

template<typename Scalar_ , int Deg_>
EigenSolverType Eigen::PolynomialSolver< Scalar_, Deg_ >::m_eigenSolver
protected

Definition at line 396 of file PolynomialSolver.h.

◆ m_roots

template<typename Scalar_ , int Deg_>
RootsType Eigen::PolynomialSolverBase< Scalar_, Deg_ >::m_roots
protected

Definition at line 292 of file PolynomialSolver.h.


The documentation for this class was generated from the following file: