Drake
Drake C++ Documentation
Loading...
Searching...
No Matches
GenericPolynomial< BasisElement > Class Template Reference

Detailed Description

template<typename BasisElement>
class drake::symbolic::GenericPolynomial< BasisElement >

Represents symbolic generic polynomials using a given basis (for example, monomial basis, Chebyshev basis, etc).

A generic symbolic polynomial keeps a mapping from a basis element of indeterminates to its coefficient in a symbolic expression. A generic polynomial p has to satisfy an invariant such that p.decision_variables() ∩ p.indeterminates() = ∅. We have CheckInvariant() method to check the invariant. For polynomials using different basis, you could refer to section 3.1.5 of Semidefinite Optimization and Convex Algebraic Geometry on the pros/cons of each basis.

We provide two instantiations of this template

#include <drake/common/symbolic/generic_polynomial.h>

Public Types

using MapType = std::map<BasisElement, Expression>
 Type of mapping from basis element to coefficient.

Public Member Functions

 GenericPolynomial ()=default
 Constructs a zero polynomial.
 GenericPolynomial (std::nullptr_t)
 Constructs a default value.
 GenericPolynomial (MapType init)
 Constructs a generic polynomial from a map, basis_element → coefficient.
 GenericPolynomial (const BasisElement &m)
 Constructs a generic polynomial from a single basis element m.
 GenericPolynomial (const Expression &e)
 Constructs a polynomial from an expression e.
 GenericPolynomial (const Expression &e, Variables indeterminates)
 Constructs a polynomial from an expression e by decomposing it with respect to indeterminates.
const Variablesindeterminates () const
 Returns the indeterminates of this generic polynomial.
const Variablesdecision_variables () const
 Returns the decision variables of this generic polynomial.
void SetIndeterminates (const Variables &new_indeterminates)
 Sets the indeterminates to new_indeterminates.
const MapTypebasis_element_to_coefficient_map () const
 Returns the map from each basis element to its coefficient.
int Degree (const Variable &v) const
 Returns the highest degree of this generic polynomial in an indeterminate v.
int TotalDegree () const
 Returns the total degree of this generic polynomial.
Expression ToExpression () const
 Returns an equivalent symbolic expression of this generic polynomial.
GenericPolynomial< BasisElement > Differentiate (const Variable &x) const
 Differentiates this generic polynomial with respect to the variable x.
template<typename Derived>
Eigen::Matrix< GenericPolynomial< BasisElement >, 1, Derived::RowsAtCompileTime > Jacobian (const Eigen::MatrixBase< Derived > &vars) const
 Computes the Jacobian matrix J of the generic polynomial with respect to vars.
double Evaluate (const Environment &env) const
 Evaluates this generic polynomial under a given environment env.
GenericPolynomial< BasisElement > EvaluatePartial (const Environment &env) const
 Partially evaluates this generic polynomial using an environment env.
GenericPolynomial< BasisElement > EvaluatePartial (const Variable &var, double c) const
 Partially evaluates this generic polynomial by substituting var with c.
GenericPolynomial< BasisElement > & AddProduct (const Expression &coeff, const BasisElement &m)
 Adds coeff * m to this generic polynomial.
GenericPolynomial< BasisElement > RemoveTermsWithSmallCoefficients (double coefficient_tol) const
 Removes the terms whose absolute value of the coefficients are smaller than or equal to coefficient_tol.
GenericPolynomial< BasisElement > & operator+= (const GenericPolynomial< BasisElement > &p)
GenericPolynomial< BasisElement > & operator+= (const BasisElement &m)
GenericPolynomial< BasisElement > & operator+= (double c)
GenericPolynomial< BasisElement > & operator+= (const Variable &v)
GenericPolynomial< BasisElement > & operator-= (const GenericPolynomial< BasisElement > &p)
GenericPolynomial< BasisElement > & operator-= (const BasisElement &m)
GenericPolynomial< BasisElement > & operator-= (double c)
GenericPolynomial< BasisElement > & operator-= (const Variable &v)
GenericPolynomial< BasisElement > & operator*= (const GenericPolynomial< BasisElement > &p)
GenericPolynomial< BasisElement > & operator*= (const BasisElement &m)
GenericPolynomial< BasisElement > & operator*= (double c)
GenericPolynomial< BasisElement > & operator*= (const Variable &v)
GenericPolynomial< BasisElement > & operator/= (double c)
bool EqualTo (const GenericPolynomial< BasisElement > &p) const
 Returns true if this and p are structurally equal.
bool EqualToAfterExpansion (const GenericPolynomial< BasisElement > &p) const
 Returns true if this generic polynomial and p are equal after expanding the coefficients.
bool CoefficientsAlmostEqual (const GenericPolynomial< BasisElement > &p, double tol) const
 Returns true if this polynomial and p are almost equal (the difference in the corresponding coefficients are all less than tol), after expanding the coefficients.
Formula operator== (const GenericPolynomial< BasisElement > &p) const
 Returns a symbolic formula representing the condition where this polynomial and p are the same.
Formula operator!= (const GenericPolynomial< BasisElement > &p) const
 Returns a symbolic formula representing the condition where this polynomial and p are not the same.
Implements CopyConstructible, CopyAssignable, MoveConstructible, MoveAssignable
 GenericPolynomial (const GenericPolynomial &)=default
GenericPolynomialoperator= (const GenericPolynomial &)=default
 GenericPolynomial (GenericPolynomial &&)=default
GenericPolynomialoperator= (GenericPolynomial &&)=default

Friends

template<class HashAlgorithm>
void hash_append (HashAlgorithm &hasher, const GenericPolynomial< BasisElement > &item) noexcept
 Implements the hash_append generic hashing concept.

Member Typedef Documentation

◆ MapType

template<typename BasisElement>
using MapType = std::map<BasisElement, Expression>

Type of mapping from basis element to coefficient.

Constructor & Destructor Documentation

◆ GenericPolynomial() [1/8]

template<typename BasisElement>
GenericPolynomial ( )
default

Constructs a zero polynomial.

◆ GenericPolynomial() [2/8]

template<typename BasisElement>
GenericPolynomial ( const GenericPolynomial< BasisElement > & )
default

◆ GenericPolynomial() [3/8]

template<typename BasisElement>
GenericPolynomial ( GenericPolynomial< BasisElement > && )
default

◆ GenericPolynomial() [4/8]

template<typename BasisElement>
GenericPolynomial ( std::nullptr_t )
explicit

Constructs a default value.

This overload is used by Eigen when EIGEN_INITIALIZE_MATRICES_BY_ZERO is enabled.

◆ GenericPolynomial() [5/8]

template<typename BasisElement>
GenericPolynomial ( MapType init)
explicit

Constructs a generic polynomial from a map, basis_element → coefficient.

For example

{{MonomialBasisElement(x, 2), a}, {MonomialBasisElement(x, 3), a+b}})
GenericPolynomial()=default
Constructs a zero polynomial.
MonomialBasisElement represents a monomial, a product of powers of variables with non-negative intege...
Definition monomial_basis_element.h:26
x
Definition light_parameter.h:147

constructs a polynomial ax²+(a+b)x³.

◆ GenericPolynomial() [6/8]

template<typename BasisElement>
GenericPolynomial ( const BasisElement & m)

Constructs a generic polynomial from a single basis element m.

Note
that all variables in m are considered as indeterminates. Namely the constructed generic polynomial contains the map with a single key m, with the coefficient being 1.

◆ GenericPolynomial() [7/8]

template<typename BasisElement>
GenericPolynomial ( const Expression & e)
explicit

Constructs a polynomial from an expression e.

Note that all variables in e are considered as indeterminates.

Exceptions
std::exceptionif e is not a polynomial.

◆ GenericPolynomial() [8/8]

template<typename BasisElement>
GenericPolynomial ( const Expression & e,
Variables indeterminates )

Constructs a polynomial from an expression e by decomposing it with respect to indeterminates.

Note
The indeterminates for the polynomial are indeterminates. Even if a variable in indeterminates does not show up in e, that variable is still registered as an indeterminate in this polynomial, as this->indeterminates() be the same as indeterminates.
Exceptions
std::exceptionif e is not a polynomial in indeterminates.

Member Function Documentation

◆ AddProduct()

template<typename BasisElement>
GenericPolynomial< BasisElement > & AddProduct ( const Expression & coeff,
const BasisElement & m )

Adds coeff * m to this generic polynomial.

◆ basis_element_to_coefficient_map()

template<typename BasisElement>
const MapType & basis_element_to_coefficient_map ( ) const
nodiscard

Returns the map from each basis element to its coefficient.

◆ CoefficientsAlmostEqual()

template<typename BasisElement>
bool CoefficientsAlmostEqual ( const GenericPolynomial< BasisElement > & p,
double tol ) const

Returns true if this polynomial and p are almost equal (the difference in the corresponding coefficients are all less than tol), after expanding the coefficients.

◆ decision_variables()

template<typename BasisElement>
const Variables & decision_variables ( ) const
nodiscard

Returns the decision variables of this generic polynomial.

◆ Degree()

template<typename BasisElement>
int Degree ( const Variable & v) const
nodiscard

Returns the highest degree of this generic polynomial in an indeterminate v.

◆ Differentiate()

template<typename BasisElement>
GenericPolynomial< BasisElement > Differentiate ( const Variable & x) const
nodiscard

Differentiates this generic polynomial with respect to the variable x.

Note that a variable x can be either a decision variable or an indeterminate.

◆ EqualTo()

template<typename BasisElement>
bool EqualTo ( const GenericPolynomial< BasisElement > & p) const
nodiscard

Returns true if this and p are structurally equal.

◆ EqualToAfterExpansion()

template<typename BasisElement>
bool EqualToAfterExpansion ( const GenericPolynomial< BasisElement > & p) const
nodiscard

Returns true if this generic polynomial and p are equal after expanding the coefficients.

◆ Evaluate()

template<typename BasisElement>
double Evaluate ( const Environment & env) const
nodiscard

Evaluates this generic polynomial under a given environment env.

Exceptions
std::exceptionif there is a variable in this generic polynomial whose assignment is not provided by env.

◆ EvaluatePartial() [1/2]

template<typename BasisElement>
GenericPolynomial< BasisElement > EvaluatePartial ( const Environment & env) const
nodiscard

Partially evaluates this generic polynomial using an environment env.

Exceptions
std::exceptionif NaN is detected during evaluation.

◆ EvaluatePartial() [2/2]

template<typename BasisElement>
GenericPolynomial< BasisElement > EvaluatePartial ( const Variable & var,
double c ) const
nodiscard

Partially evaluates this generic polynomial by substituting var with c.

Exceptions
std::exceptionif NaN is detected at any point during evaluation.

◆ indeterminates()

template<typename BasisElement>
const Variables & indeterminates ( ) const
nodiscard

Returns the indeterminates of this generic polynomial.

◆ Jacobian()

template<typename BasisElement>
template<typename Derived>
Eigen::Matrix< GenericPolynomial< BasisElement >, 1, Derived::RowsAtCompileTime > Jacobian ( const Eigen::MatrixBase< Derived > & vars) const

Computes the Jacobian matrix J of the generic polynomial with respect to vars.

J(0,i) contains ∂f/∂vars(i). vars should be an Eigen column vector of symbolic variables.

◆ operator!=()

template<typename BasisElement>
Formula operator!= ( const GenericPolynomial< BasisElement > & p) const

Returns a symbolic formula representing the condition where this polynomial and p are not the same.

◆ operator*=() [1/4]

template<typename BasisElement>
GenericPolynomial< BasisElement > & operator*= ( const BasisElement & m)

◆ operator*=() [2/4]

template<typename BasisElement>
GenericPolynomial< BasisElement > & operator*= ( const GenericPolynomial< BasisElement > & p)

◆ operator*=() [3/4]

template<typename BasisElement>
GenericPolynomial< BasisElement > & operator*= ( const Variable & v)

◆ operator*=() [4/4]

template<typename BasisElement>
GenericPolynomial< BasisElement > & operator*= ( double c)

◆ operator+=() [1/4]

template<typename BasisElement>
GenericPolynomial< BasisElement > & operator+= ( const BasisElement & m)

◆ operator+=() [2/4]

template<typename BasisElement>
GenericPolynomial< BasisElement > & operator+= ( const GenericPolynomial< BasisElement > & p)

◆ operator+=() [3/4]

template<typename BasisElement>
GenericPolynomial< BasisElement > & operator+= ( const Variable & v)

◆ operator+=() [4/4]

template<typename BasisElement>
GenericPolynomial< BasisElement > & operator+= ( double c)

◆ operator-=() [1/4]

template<typename BasisElement>
GenericPolynomial< BasisElement > & operator-= ( const BasisElement & m)

◆ operator-=() [2/4]

template<typename BasisElement>
GenericPolynomial< BasisElement > & operator-= ( const GenericPolynomial< BasisElement > & p)

◆ operator-=() [3/4]

template<typename BasisElement>
GenericPolynomial< BasisElement > & operator-= ( const Variable & v)

◆ operator-=() [4/4]

template<typename BasisElement>
GenericPolynomial< BasisElement > & operator-= ( double c)

◆ operator/=()

template<typename BasisElement>
GenericPolynomial< BasisElement > & operator/= ( double c)

◆ operator=() [1/2]

template<typename BasisElement>
GenericPolynomial & operator= ( const GenericPolynomial< BasisElement > & )
default

◆ operator=() [2/2]

template<typename BasisElement>
GenericPolynomial & operator= ( GenericPolynomial< BasisElement > && )
default

◆ operator==()

template<typename BasisElement>
Formula operator== ( const GenericPolynomial< BasisElement > & p) const

Returns a symbolic formula representing the condition where this polynomial and p are the same.

◆ RemoveTermsWithSmallCoefficients()

template<typename BasisElement>
GenericPolynomial< BasisElement > RemoveTermsWithSmallCoefficients ( double coefficient_tol) const
nodiscard

Removes the terms whose absolute value of the coefficients are smaller than or equal to coefficient_tol.

For example, if the generic polynomial is 2x² + 3xy + 10⁻⁴x - 10⁻⁵, then after calling RemoveTermsWithSmallCoefficients(1e-3), the returned polynomial becomes 2x² + 3xy.

Parameters
coefficient_tolA positive scalar.
Return values
polynomial_cleanedA generic polynomial whose terms with small coefficients are removed.

◆ SetIndeterminates()

template<typename BasisElement>
void SetIndeterminates ( const Variables & new_indeterminates)

Sets the indeterminates to new_indeterminates.

Changing the indeterminates will change basis_element_to_coefficient_map(), and also potentially the degree of the polynomial. Here is an example.

// p is a quadratic polynomial with x being the only indeterminate.
{x});
// p.basis_element_to_coefficient_map() contains {1: c, x: b, x*x:a}.
std::cout << p.TotalDegree(); // prints 2.
// Now set (a, b, c) to the indeterminates. p becomes a linear
// polynomial of a, b, c.
p.SetIndeterminates({a, b, c});
// p.basis_element_to_coefficient_map() now is {a: x * x, b: x, c: 1}.
std::cout << p.TotalDegree(); // prints 1.
Represents symbolic generic polynomials using a given basis (for example, monomial basis,...
Definition generic_polynomial.h:35
int TotalDegree() const
Returns the total degree of this generic polynomial.

This function can be expensive, as it potentially reconstructs the polynomial (using the new indeterminates) from the expression.

◆ ToExpression()

template<typename BasisElement>
Expression ToExpression ( ) const
nodiscard

Returns an equivalent symbolic expression of this generic polynomial.

◆ TotalDegree()

template<typename BasisElement>
int TotalDegree ( ) const
nodiscard

Returns the total degree of this generic polynomial.

◆ hash_append

template<typename BasisElement>
template<class HashAlgorithm>
void hash_append ( HashAlgorithm & hasher,
const GenericPolynomial< BasisElement > & item )
friend

Implements the hash_append generic hashing concept.


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