pydrake.symbolic

Symbolic variable, variables, monomial, expression, polynomial, and formula

pydrake.symbolic.abs(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.acos(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.arccos(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.arcsin(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.arctan(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.arctan2(y: pydrake.symbolic.Expression, x: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.asin(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.atan(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.atan2(y: pydrake.symbolic.Expression, x: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.CalcMonomialBasisOrderUpToOne(x: pydrake.symbolic.Variables, sort_monomial: bool = False) numpy.ndarray[object[m, 1]]

Generates all the monomials of x, such that the degree for x(i) is no larger than 1 for every x(i) in x.

Parameter x:

The variables whose monomials are generated.

Parameter sort_monomial:

If true, the returned monomials are sorted in the graded reverse lexicographic order. For example if x = (x₀, x₁) with x₀< x₁, then this function returns [x₀x₁, x₁, x₀, 1]. If sort_monomial=false, then we return the monomials in an arbitrary order.

pydrake.symbolic.CalcPolynomialWLowerTriangularPart(*args, **kwargs)

Overloaded function.

  1. CalcPolynomialWLowerTriangularPart(monomial_basis: numpy.ndarray[object[m, 1]], gram_lower: numpy.ndarray[numpy.float64[m, 1]]) -> pydrake.symbolic.Polynomial

Returns the polynomial m(x)ᵀ * Q * m(x), where m(x) is the monomial basis, and Q is the Gram matrix.

Parameter monomial_basis:

m(x) in the documentation. A vector of monomials.

Parameter gram_lower:

The lower triangular entries in Q, stacked columnwise into a vector.

  1. CalcPolynomialWLowerTriangularPart(monomial_basis: numpy.ndarray[object[m, 1]], gram_lower: numpy.ndarray[object[m, 1]]) -> pydrake.symbolic.Polynomial

Returns the polynomial m(x)ᵀ * Q * m(x), where m(x) is the monomial basis, and Q is the Gram matrix.

Parameter monomial_basis:

m(x) in the documentation. A vector of monomials.

Parameter gram_lower:

The lower triangular entries in Q, stacked columnwise into a vector.

  1. CalcPolynomialWLowerTriangularPart(monomial_basis: numpy.ndarray[object[m, 1]], gram_lower: numpy.ndarray[object[m, 1]]) -> pydrake.symbolic.Polynomial

Returns the polynomial m(x)ᵀ * Q * m(x), where m(x) is the monomial basis, and Q is the Gram matrix.

Parameter monomial_basis:

m(x) in the documentation. A vector of monomials.

Parameter gram_lower:

The lower triangular entries in Q, stacked columnwise into a vector.

pydrake.symbolic.ceil(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.cos(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.cosh(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.DecomposeAffineExpression(e: pydrake.symbolic.Expression, map_var_to_index: Dict[int, int]) Tuple[numpy.ndarray[numpy.float64[1, n]], float]

Decomposes an affine combination e = c0 + c1 * v1 + … cn * vn into the following:

constant term : c0 coefficient vector : [c1, …, cn] variable vector : [v1, …, vn]

Then, it extracts the coefficient and the constant term. A map from variable ID to int, map_var_to_index, is used to decide a variable’s index in a linear combination.

Precondition: 1. coeffs is a row vector of double, whose length matches with the size of map_var_to_index. 2. e.is_polynomial() is true. 3. e is an affine expression. 4. all values in map_var_to_index should be in the range [0, map_var_to_index.size())

$Parameter e:

The symbolic affine expression

Parameter map_var_to_index:

A mapping from variable ID to variable index, such that map_var_to_index[vi.get_ID()] = i.

Parameter coeffs:

A row vector. coeffs(i) = ci.

Parameter constant_term:

c0 in the equation above.

Returns

num_variable. Number of variables in the expression. 2 * x(0) + 3 has 1 variable, 2 * x(0) + 3 * x(1) - 2 * x(0) has 1 variable.

pydrake.symbolic.DecomposeAffineExpressions(*args, **kwargs)

Overloaded function.

  1. DecomposeAffineExpressions(expressions: numpy.ndarray[object[m, 1]], vars: numpy.ndarray[object[m, 1]]) -> Tuple[numpy.ndarray[numpy.float64[m, n]], numpy.ndarray[numpy.float64[m, 1]]]

Decomposes expressions into M * vars + v.

Raises

RuntimeError if expressions is not affine in vars.

Precondition:

M.rows() == expressions.rows() && M.cols() == vars.rows().

Precondition:

v.rows() == expressions.rows().

  1. DecomposeAffineExpressions(v: numpy.ndarray[object[m, 1]]) -> Tuple[numpy.ndarray[numpy.float64[m, n]], numpy.ndarray[numpy.float64[m, 1]], numpy.ndarray[object[m, 1]]]

Given a vector of affine expressions v, decompose it to \(v = A vars + b\)

Parameter v:

A vector of affine expressions

Parameter A:

The matrix containing the linear coefficients.

Parameter b:

The vector containing all the constant terms.

Parameter vars:

All variables.

pydrake.symbolic.DecomposeLinearExpressions(expressions: numpy.ndarray[object[m, 1]], vars: numpy.ndarray[object[m, 1]]) numpy.ndarray[numpy.float64[m, n]]

Decomposes expressions into M * vars.

Raises

RuntimeError if expressions is not linear in vars.

Precondition:

M.rows() == expressions.rows() && M.cols() == vars.rows().

pydrake.symbolic.DecomposeLumpedParameters(f: numpy.ndarray[object[m, 1]], parameters: numpy.ndarray[object[m, 1]]) Tuple[numpy.ndarray[object[m, n]], numpy.ndarray[object[m, 1]], numpy.ndarray[object[m, 1]]]

Given a vector of Expressions f and a list of parameters we define all additional variables in f to be a vector of “non-parameter variables”, n. This method returns a factorization of f into an equivalent “data matrix”, W, which depends only on the non-parameter variables, and a “lumped parameter vector”, α, which depends only on parameters: f = W(n)*α(parameters) + w0(n).

Note

The current implementation makes some simple attempts to minimize the number of lumped parameters, but more simplification could be implemented relatively easily. Optimal simplification, however, involves the complexity of comparing two arbitrary Expressions (see Expression::EqualTo for more details).

Raises
  • RuntimeError if f is not decomposable in this way (cells

  • containing parameters may only be added or multiplied with

  • cells containing non-parameter variables)

Returns

W(n), α(parameters), and w0(n).

pydrake.symbolic.DecomposeQuadraticPolynomial(poly: pydrake.symbolic.Polynomial, map_var_to_index: Dict[int, int]) Tuple[numpy.ndarray[numpy.float64[m, n]], numpy.ndarray[numpy.float64[m, 1]], float]

Given a quadratic polynomial poly, decomposes it into the form 0.5 * x’ Q * x + b’ * x + c

Parameter poly:

Quadratic polynomial to decompose.

Parameter map_var_to_index:

maps variables in poly.GetVariables() to the index in the vector x.

Parameter Q:

[out] The Hessian of the quadratic expression.

Precondition:

The size of Q should be num_variables x num_variables. Q is a symmetric matrix.

Parameter b:

[out] The linear term of the quadratic expression.

Precondition:

The size of b should be num_variables.

Parameter c:

[out] The constant term of the quadratic expression.

pydrake.symbolic.Evaluate(*args, **kwargs)

Overloaded function.

  1. Evaluate(m: numpy.ndarray[object[m, n]], env: Dict[pydrake.symbolic.Variable, float] = {}, generator: pydrake.common.RandomGenerator = None) -> numpy.ndarray[numpy.float64[m, n]]

Evaluates a symbolic matrix m using env and random_generator.

If there is a random variable in m which is unassigned in env, this function uses random_generator to sample a value and use the value to substitute all occurrences of the random variable in m.

Returns

a matrix of double whose size is the size of m.

Raises
  • RuntimeError if NaN is detected during evaluation.

  • RuntimeError if m includes unassigned random variables but

  • random_generator` is nullptr

  1. Evaluate(m: numpy.ndarray[object[m, n]], env: Dict[pydrake.symbolic.Variable, float]) -> numpy.ndarray[numpy.float64[m, n]]

Evaluates a matrix m of symbolic polynomials using env.

Returns

a matrix of double whose size is the size of m.

Raises

RuntimeError if NaN is detected during evaluation.

pydrake.symbolic.EvenDegreeMonomialBasis(vars: pydrake.symbolic.Variables, degree: int) numpy.ndarray[object[m, 1]]

Returns all even degree monomials up to a given degree under the graded reverse lexicographic order. A monomial has an even degree if its total degree is even. So xy is an even degree monomial (degree 2) while x²y is not (degree 3). Note that graded reverse lexicographic order uses the total order among Variable which is based on a variable’s unique ID. For example, for a given variable ordering x > y > z, EvenDegreeMonomialBasis({x, y, z}, 2) returns a column vector [x², xy, y², xz, yz, z², 1].

Precondition:

vars is a non-empty set.

Precondition:

degree is a non-negative integer.

pydrake.symbolic.exp(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
class pydrake.symbolic.Expression

Represents a symbolic form of an expression.

Its syntax tree is as follows:

Click to expand C++ code...
E := Var | Constant | E + ... + E | E * ... * E | E / E | log(E)
| abs(E) | exp(E) | sqrt(E) | pow(E, E) | sin(E) | cos(E) | tan(E)
| asin(E) | acos(E) | atan(E) | atan2(E, E) | sinh(E) | cosh(E) | tanh(E)
| min(E, E) | max(E, E) | ceil(E) | floor(E) | if_then_else(F, E, E)
| NaN | uninterpreted_function(name, {v_1, ..., v_n})

In the implementation, Expression directly stores Constant values inline, but in all other cases stores a shared pointer to a const ExpressionCell class that is a super-class of different kinds of symbolic expressions (i.e., ExpressionAdd, ExpressionMul, ExpressionLog, ExpressionSin), which makes it efficient to copy, move, and assign to an Expression.

Note

-E is represented as -1 * E internally.

Note

A subtraction E1 - E2 is represented as E1 + (-1 * E2) internally.

The following simple simplifications are implemented:

Click to expand C++ code...
E + 0             ->  E
0 + E             ->  E
E - 0             ->  E
E - E             ->  0
E * 1             ->  E
1 * E             ->  E
E * 0             ->  0
0 * E             ->  0
E / 1             ->  E
E / E             ->  1
pow(E, 0)         ->  1
pow(E, 1)         ->  E
E * E             ->  E^2 (= pow(E, 2))
sqrt(E * E)       ->  |E| (= abs(E))
sqrt(E) * sqrt(E) -> E

Constant folding is implemented:

Click to expand C++ code...
E(c1) + E(c2)  ->  E(c1 + c2)    // c1, c2 are constants
E(c1) - E(c2)  ->  E(c1 - c2)
E(c1) * E(c2)  ->  E(c1 * c2)
E(c1) / E(c2)  ->  E(c1 / c2)
f(E(c))        ->  E(f(c))       // c is a constant, f is a math function

For the math functions which are only defined over a restricted domain (namely, log, sqrt, pow, asin, acos), we check the domain of argument(s), and throw ValueError exception if a function is not well-defined for a given argument(s).

Relational operators over expressions (==, !=, <, >, <=, >=) return symbolic::Formula instead of bool. Those operations are declared in formula.h file. To check structural equality between two expressions a separate function, Expression::EqualTo, is provided.

Regarding the arithmetic of an Expression when operating on NaNs, we have the following rules: 1. NaN values are extremely rare during typical computations. Because they are difficult to handle symbolically, we will round that up to “must never occur”. We allow the user to form ExpressionNaN cells in a symbolic tree. For example, the user can initialize an Expression to NaN and then overwrite it later. However, evaluating a tree that has NaN in its evaluated sub-trees is an error (see rule (3) below). 2. It’s still valid for code to check isnan in order to fail-fast. So we provide isnan(const Expression&) for the common case of non-NaN value returning False. This way, code can fail-fast with double yet still compile with Expression. 3. If there are expressions that embed separate cases (if_then_else), some of the sub-expressions may be not used in evaluation when they are in the not-taken case (for NaN reasons or any other reason). Bad values within those not-taken branches does not cause exceptions. 4. The isnan check is different than if_then_else. In the latter, the ExpressionNaN is within a dead sub-expression branch. In the former, it appears in an evaluated trunk. That goes against rule (1) where a NaN anywhere in a computation (other than dead code) is an error.

note for Drake developers: under the hood of Expression, we have an internal::BoxedCell helper class that uses NaN for pointer tagging; that’s a distinct concept from the Expression::NaN() rules enumerated just above.

symbolic::Expression can be used as a scalar type of Eigen types.

__init__(*args, **kwargs)

Overloaded function.

  1. __init__(self: pydrake.symbolic.Expression) -> None

Default constructor. It constructs Zero().

  1. __init__(self: pydrake.symbolic.Expression, constant: float) -> None

Constructs a constant.

  1. __init__(self: pydrake.symbolic.Expression, var: pydrake.symbolic.Variable) -> None

Constructs an expression from var.

Precondition:

var is not a BOOLEAN variable.

abs(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
acos(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
arccos(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
arcsin(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
arctan(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
arctan2(self: pydrake.symbolic.Expression, x: pydrake.symbolic.Expression) pydrake.symbolic.Expression

Uses self for y in arctan2(y, x).

asin(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
atan(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
atan2(self: pydrake.symbolic.Expression, x: pydrake.symbolic.Expression) pydrake.symbolic.Expression

Uses self for y in atan2(y, x).

ceil(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
cos(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
cosh(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
Differentiate(self: pydrake.symbolic.Expression, x: pydrake.symbolic.Variable) pydrake.symbolic.Expression

Differentiates this symbolic expression with respect to the variable var.

Raises

RuntimeError if it is not differentiable.

EqualTo(self: pydrake.symbolic.Expression, arg0: pydrake.symbolic.Expression) bool

Checks structural equality.

Two expressions e1 and e2 are structurally equal when they have the same internal AST(abstract-syntax tree) representation. Please note that we can have two computationally (or extensionally) equivalent expressions which are not structurally equal. For example, consider:

e1 = 2 * (x + y) e2 = 2x + 2y

Obviously, we know that e1 and e2 are evaluated to the same value for all assignments to x and y. However, e1 and e2 are not structurally equal by the definition. Note that e1 is a multiplication expression (is_multiplication(e1) is true) while e2 is an addition expression (is_addition(e2) is true).

One main reason we use structural equality in EqualTo is due to Richardson’s Theorem. It states that checking ∀x. E(x) = F(x) is undecidable when we allow sin, asin, log, exp in E and F. Read https://en.wikipedia.org/wiki/Richardson%27s_theorem for details.

Note that for polynomial cases, you can use Expand method and check if two polynomial expressions p1 and p2 are computationally equal. To do so, you check the following:

p1.Expand().EqualTo(p2.Expand())

Evaluate(*args, **kwargs)

Overloaded function.

  1. Evaluate(self: pydrake.symbolic.Expression, env: Dict[pydrake.symbolic.Variable, float] = {}, generator: pydrake.common.RandomGenerator = None) -> float

Evaluates using a given environment (by default, an empty environment) and a random number generator. If there is a random variable in this expression which is unassigned in env, this method uses random_generator to sample a value and use the value to substitute all occurrences of the variable in this expression.

Raises
  • RuntimeError if there exists a non-random variable in this

  • expression whose assignment is not provided by env.

  • RuntimeError if an unassigned random variable is detected while

  • random_generator` is nullptr

  • RuntimeError if NaN is detected during evaluation.

  1. Evaluate(self: pydrake.symbolic.Expression, generator: pydrake.common.RandomGenerator) -> float

Evaluates using an empty environment and a random number generator. It uses random_generator to sample values for the random variables in this expression.

See the above overload for the exceptions that it might throw.

EvaluatePartial(self: pydrake.symbolic.Expression, env: Dict[pydrake.symbolic.Variable, float]) pydrake.symbolic.Expression

Partially evaluates this expression using an environment env. Internally, this method promotes env into a substitution (Variable → Expression) and call Evaluate::Substitute with it.

Raises

RuntimeError if NaN is detected during evaluation.

exp(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
Expand(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression

Expands out products and positive integer powers in expression. For example, (x + 1) * (x - 1) is expanded to x^2 - 1 and (x + y)^2 is expanded to x^2 + 2xy + y^2. Note that Expand applies recursively to sub-expressions. For instance, sin(2 * (x + y)) is expanded to sin(2x + 2y). It also simplifies “division by constant” cases. See “drake/common/test/symbolic_expansion_test.cc” to find the examples.

Raises

RuntimeError if NaN is detected during expansion.

floor(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
get_kind(self: pydrake.symbolic.Expression) pydrake.symbolic.ExpressionKind

Returns expression kind.

GetVariables(self: pydrake.symbolic.Expression) pydrake.symbolic.Variables

Collects variables in expression.

is_polynomial(self: pydrake.symbolic.Expression) bool

Checks if this symbolic expression is convertible to Polynomial.

Jacobian(self: pydrake.symbolic.Expression, vars: numpy.ndarray[object[m, 1]]) numpy.ndarray[object[1, n]]

Let f be this Expression, computes a row vector of derivatives, [∂f/∂vars(0), ... , ∂f/∂vars(n-1)] with respect to the variables vars.

log(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
max(self: pydrake.symbolic.Expression, arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
min(self: pydrake.symbolic.Expression, arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pow(*args, **kwargs)

Overloaded function.

  1. pow(self: pydrake.symbolic.Expression, arg0: float) -> pydrake.symbolic.Expression

  2. pow(self: pydrake.symbolic.Expression, arg0: pydrake.symbolic.Expression) -> pydrake.symbolic.Expression

sin(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
sinh(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
sqrt(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
Substitute(*args, **kwargs)

Overloaded function.

  1. Substitute(self: pydrake.symbolic.Expression, var: pydrake.symbolic.Variable, e: pydrake.symbolic.Expression) -> pydrake.symbolic.Expression

Returns a copy of this expression replacing all occurrences of var with e.

Raises

RuntimeError if NaN is detected during substitution.

  1. Substitute(self: pydrake.symbolic.Expression, s: Dict[pydrake.symbolic.Variable, pydrake.symbolic.Expression]) -> pydrake.symbolic.Expression

Returns a copy of this expression replacing all occurrences of the variables in s with corresponding expressions in s. Note that the substitutions occur simultaneously. For example, (x / y).Substitute({{x, y}, {y, x}}) gets (y / x).

Raises

RuntimeError if NaN is detected during substitution.

tan(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
tanh(self: pydrake.symbolic.Expression) pydrake.symbolic.Expression
to_string(self: pydrake.symbolic.Expression) str

Returns string representation of Expression.

Unapply(self: pydrake.symbolic.Expression) object

Given an expression, returns a tuple (ctor, args) that would re-create an equivalent expression when called as ctor(*args). This is a useful way to unpack the contents of a compound expression, e.g., to obtain the terms of an addition. In many cases (all arithmetic operators, trig functions, max, min, sqrt, etc.), the returned args will all be either Expressions or floats; in other cases (e.g., if_then_else) some args will other types (e.g., a Formula). To check the form (i.e., kind) of an expression, use e.get_kind(); do not try to infer it from the returned ctor’s identity.

class pydrake.symbolic.ExpressionKind

Kinds of symbolic expressions. The constants here are carefully chosen to support nanboxing. For all elements except Constant, the bit pattern must have have 0x7FF0 bits set but must not be exactly 0x7FF0 nor 0xFFF0 (reserved for ±infinity). Refer to the details in boxed_cell.h for more information.

Members:

Constant : constant (double)

Var : variable

Add : addition (+)

Mul : multiplication (*)

Div : division (/)

Log : logarithms

Abs : absolute value function

Exp : exponentiation

Sqrt : square root

Pow : power function

Sin : sine

Cos : cosine

Tan : tangent

Asin : arcsine

Acos : arccosine

Atan : arctangent

Atan2 : arctangent2 (atan2(y,x) = atan(y/x))

Sinh : hyperbolic sine

Cosh : hyperbolic cosine

Tanh : hyperbolic tangent

Min : min

Max : max

Ceil : ceil

Floor : floor

IfThenElse : if then else

NaN : NaN

UninterpretedFunction : Uninterpreted function

__init__(self: pydrake.symbolic.ExpressionKind, value: int) None
Abs = <ExpressionKind.Abs: 32758>
Acos = <ExpressionKind.Acos: 32766>
Add = <ExpressionKind.Add: 32754>
Asin = <ExpressionKind.Asin: 32765>
Atan = <ExpressionKind.Atan: 32767>
Atan2 = <ExpressionKind.Atan2: 65521>
Ceil = <ExpressionKind.Ceil: 65527>
Constant = <ExpressionKind.Constant: 0>
Cos = <ExpressionKind.Cos: 32763>
Cosh = <ExpressionKind.Cosh: 65523>
Div = <ExpressionKind.Div: 32756>
Exp = <ExpressionKind.Exp: 32759>
Floor = <ExpressionKind.Floor: 65528>
IfThenElse = <ExpressionKind.IfThenElse: 65529>
Log = <ExpressionKind.Log: 32757>
Max = <ExpressionKind.Max: 65526>
Min = <ExpressionKind.Min: 65525>
Mul = <ExpressionKind.Mul: 32755>
property name
NaN = <ExpressionKind.NaN: 65530>
Pow = <ExpressionKind.Pow: 32761>
Sin = <ExpressionKind.Sin: 32762>
Sinh = <ExpressionKind.Sinh: 65522>
Sqrt = <ExpressionKind.Sqrt: 32760>
Tan = <ExpressionKind.Tan: 32764>
Tanh = <ExpressionKind.Tanh: 65524>
UninterpretedFunction = <ExpressionKind.UninterpretedFunction: 65531>
property value
Var = <ExpressionKind.Var: 32753>
pydrake.symbolic.ExtractVariablesFromExpression(*args, **kwargs)

Overloaded function.

  1. ExtractVariablesFromExpression(e: pydrake.symbolic.Expression) -> Tuple[numpy.ndarray[object[m, 1]], Dict[int, int]]

Given an expression e, extracts all variables inside e.

Parameter e:

A symbolic expression.

Returns pair:

pair.first is the variables in e. pair.second is the mapping from the variable ID to the index in pair.first, such that pair.second[pair.first(i).get_id()] = i

  1. ExtractVariablesFromExpression(expressions: numpy.ndarray[object[m, 1]]) -> Tuple[numpy.ndarray[object[m, 1]], Dict[int, int]]

Overloads ExtractVariablesFromExpression but with a vector of expressions.

pydrake.symbolic.floor(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.forall(vars: pydrake.symbolic.Variables, f: pydrake.symbolic.Formula) pydrake.symbolic.Formula

Returns a formula f, universally quantified by variables vars.

class pydrake.symbolic.Formula

Represents a symbolic form of a first-order logic formula.

It has the following grammar:

Click to expand C++ code...
F := ⊥ | ⊤ | Var | E = E | E ≠ E | E > E | E ≥ E | E < E | E ≤ E
| E ∧ ... ∧ E | E ∨ ... ∨ E | ¬F | ∀ x₁, ..., xn. F

In the implementation, Formula is a simple wrapper including a shared pointer to FormulaCell class which is a super-class of different kinds of symbolic formulas (i.e. FormulaAnd, FormulaOr, FormulaEq). Note that it includes a shared pointer, not a unique pointer, to allow sharing sub-expressions.

Note

The sharing of sub-expressions is not yet implemented.

The following simple simplifications are implemented:

Click to expand C++ code...
E1 = E2        ->  True    (if E1 and E2 are structurally equal)
E1 ≠ E2        ->  False   (if E1 and E2 are structurally equal)
E1 > E2        ->  False   (if E1 and E2 are structurally equal)
E1 ≥ E2        ->  True    (if E1 and E2 are structurally equal)
E1 < E2        ->  False   (if E1 and E2 are structurally equal)
E1 ≤ E2        ->  True    (if E1 and E2 are structurally equal)
F1 ∧ F2        ->  False   (if either F1 or F2 is False)
F1 ∨ F2        ->  True    (if either F1 or F2 is True)
¬(¬(F))        ->  F

We flatten nested conjunctions (or disjunctions) at the construction. A conjunction (resp. disjunction) takes a set of conjuncts (resp. disjuncts). Note that any duplicated conjunct/disjunct is removed. For example, both of f1 && (f2 && f1) and (f1 && f2) && f1 are flattened to f1 && f2 && f1 and simplified into f1 && f2. As a result, the two are identified as the same formula.

Note

Formula class has an explicit conversion operator to bool. It evaluates a symbolic formula under an empty environment. If a symbolic formula includes variables, the conversion operator throws an exception. This operator is only intended for third-party code doing things like (imag(SymbolicExpression(0)) == SymbolicExpression(0)) { ... }; that we found in Eigen3 codebase. In general, a user of this class should explicitly call Evaluate from within Drake for readability.

__init__(*args, **kwargs)

Overloaded function.

  1. __init__(self: pydrake.symbolic.Formula) -> None

Default constructor. Sets the value to Formula::False, to be consistent with value-initialized `bool`s.

  1. __init__(self: pydrake.symbolic.Formula, value: bool) -> None

Constructs from a bool. This overload is also used by Eigen when EIGEN_INITIALIZE_MATRICES_BY_ZERO is enabled.

  1. __init__(self: pydrake.symbolic.Formula, var: pydrake.symbolic.Variable) -> None

Constructs a formula from var.

Precondition:

var is of BOOLEAN type.

EqualTo(self: pydrake.symbolic.Formula, arg0: pydrake.symbolic.Formula) bool

Checks structural equality.

Evaluate(self: pydrake.symbolic.Formula, env: Dict[pydrake.symbolic.Variable, float] = {}) bool

Evaluates using a given environment (by default, an empty environment) and a random number generator. If there is a random variable in this formula which is unassigned in env, it uses random_generator to sample a value and use it to substitute all occurrences of the random variable in this formula.

Raises
  • RuntimeError if a variable v is needed for an evaluation but

  • not provided by env.

  • RuntimeError if an unassigned random variable is detected while

  • random_generator` is nullptr

static False_() pydrake.symbolic.Formula

Symbolic formula representing false.

get_kind(self: pydrake.symbolic.Formula) pydrake.symbolic.FormulaKind
GetFreeVariables(self: pydrake.symbolic.Formula) pydrake.symbolic.Variables

Gets free variables (unquantified variables).

Substitute(*args, **kwargs)

Overloaded function.

  1. Substitute(self: pydrake.symbolic.Formula, var: pydrake.symbolic.Variable, e: pydrake.symbolic.Expression) -> pydrake.symbolic.Formula

Returns a copy of this formula replacing all occurrences of var with e.

Raises

RuntimeError if NaN is detected during substitution.

  1. Substitute(self: pydrake.symbolic.Formula, var: pydrake.symbolic.Variable, e: pydrake.symbolic.Variable) -> pydrake.symbolic.Formula

Returns a copy of this formula replacing all occurrences of var with e.

Raises

RuntimeError if NaN is detected during substitution.

  1. Substitute(self: pydrake.symbolic.Formula, var: pydrake.symbolic.Variable, e: float) -> pydrake.symbolic.Formula

Returns a copy of this formula replacing all occurrences of var with e.

Raises

RuntimeError if NaN is detected during substitution.

  1. Substitute(self: pydrake.symbolic.Formula, s: Dict[pydrake.symbolic.Variable, pydrake.symbolic.Expression]) -> pydrake.symbolic.Formula

Returns a copy of this formula replacing all occurrences of the variables in s with corresponding expressions in s. Note that the substitutions occur simultaneously. For example, (x / y > 0).Substitute({{x, y}, {y, x}}) gets (y / x > 0).

Raises

RuntimeError if NaN is detected during substitution.

to_string(self: pydrake.symbolic.Formula) str

Returns string representation of Formula.

static True_() pydrake.symbolic.Formula

Symbolic formula representing true.

Unapply(self: pydrake.symbolic.Formula) object

Given a formula, returns a tuple (ctor, args) that would re-create an equivalent formula when called as ctor(*args). This is a useful way to unpack the contents of a compound formula, e.g., to obtain the terms of a comparison. For relational formulae (==, <, >, etc.) the returned args will both be of type Expression. For compound formulae (and, or, not) the returned args will be of type Formula. To check the form (i.e., kind) of a formula, use f.get_kind(); do not try to infer it from he returned ctor’s identity.

class pydrake.symbolic.FormulaKind

Kinds of symbolic formulas.

Members:

False_ : ⊥

True_ : ⊤

Var : Boolean Variable

Eq : =

Neq : !=

Gt : >

Geq : >=

Lt : <

Leq : <=

And : Conjunction (∧)

Or : Disjunction (∨)

Not : Negation (¬)

Forall : Universal quantification (∀)

Isnan : NaN check predicate

PositiveSemidefinite : Positive semidefinite matrix

__init__(self: pydrake.symbolic.FormulaKind, value: int) None
And = <FormulaKind.And: 9>
Eq = <FormulaKind.Eq: 3>
False_ = <FormulaKind.False_: 0>
Forall = <FormulaKind.Forall: 12>
Geq = <FormulaKind.Geq: 6>
Gt = <FormulaKind.Gt: 5>
Isnan = <FormulaKind.Isnan: 13>
Leq = <FormulaKind.Leq: 8>
Lt = <FormulaKind.Lt: 7>
property name
Neq = <FormulaKind.Neq: 4>
Not = <FormulaKind.Not: 11>
Or = <FormulaKind.Or: 10>
PositiveSemidefinite = <FormulaKind.PositiveSemidefinite: 14>
True_ = <FormulaKind.True_: 1>
property value
Var = <FormulaKind.Var: 2>
pydrake.symbolic.from_sympy(x: Union[float, int, bool, sympy.Expr], *, memo: Dict = None) Union[float, int, bool, pydrake.symbolic.Variable, pydrake.symbolic.Expression, pydrake.symbolic.Formula]

Converts a SymPy Expr to the corresponding pydrake object.

Certain expressions are not supported and will raise NotImplementedError.

This function aims to support the latest contemporaneous version of SymPy as of Drake’s release date.

Parameters
  • x – The SymPy object to be converted.

  • memo – (Optional) Mapping between SymPy variables and Drake variables. Converting a SymPy variable to a pydrake.symbolic.Variable will look up the SymPy variable as a key in this dictionary. If a value is found, then it will be used as the Drake variable. If no value is found, then (for now) raises an exception. In the future, we might support automatically creating variables, but that is not yet implemented.

See also pydrake.symbolic.to_sympy().

pydrake.symbolic.GetVariableVector(expressions: numpy.ndarray[object[m, 1]]) numpy.ndarray[object[m, 1]]

Constructs a vector of variables from the vector of variable expressions.

Raises
  • RuntimeError if there is an expression in vec which is not a

  • variable.

pydrake.symbolic.if_then_else(f_cond: drake::symbolic::Formula, e_then: pydrake.symbolic.Expression, e_else: pydrake.symbolic.Expression) pydrake.symbolic.Expression

Constructs if-then-else expression.

Click to expand C++ code...
if_then_else(cond, expr_then, expr_else)

The value returned by the above if-then-else expression is expr_then if cond is evaluated to true. Otherwise, it returns expr_else.

The semantics is similar to the C++’s conditional expression constructed by its ternary operator, @c ?:. However, there is a key difference between the C++’s conditional expression and our if_then_else expression in a way the arguments are evaluated during the construction.

  • In case of the C++’s conditional expression, `` cond ? expr_then :

expr_else``, the then expression expr_then (respectively, the else expression expr_else) is only evaluated when the conditional expression cond is evaluated to true (respectively, when cond is evaluated to false).

  • In case of the symbolic expression, ``if_then_else(cond, expr_then,

expr_else)``, however, both arguments expr_then and expr_else are evaluated first and then passed to the if_then_else function.

Note

This function returns an expression and it is different from the C++’s if-then-else statement.

Note

While it is still possible to define `` min, max, abs`` math functions using if_then_else expression, it is highly recommended to use the provided native definitions for them because it allows solvers to detect specific math functions and to have a room for special optimizations.

Note

More information about the C++’s conditional expression and ternary operator is available at http://en.cppreference.com/w/cpp/language/operator_other#Conditional_operator.

pydrake.symbolic.intersect(vars1: pydrake.symbolic.Variables, vars2: pydrake.symbolic.Variables) pydrake.symbolic.Variables

Returns the intersection of vars1 and vars2.

This function has a time complexity of O(N₁ + N₂) where N₁ and N₂ are the size of vars1 and vars2 respectively.

pydrake.symbolic.inv(arg0: numpy.ndarray[object[m, n]]) numpy.ndarray[object[m, n]]
pydrake.symbolic.IsAffine(*args, **kwargs)

Overloaded function.

  1. IsAffine(m: numpy.ndarray[object[m, n], flags.f_contiguous], vars: pydrake.symbolic.Variables) -> bool

Checks if every element in m is affine in vars.

Note

If m is an empty matrix, it returns true.

  1. IsAffine(m: numpy.ndarray[object[m, n], flags.f_contiguous]) -> bool

Checks if every element in m is affine.

Note

If m is an empty matrix, it returns true.

pydrake.symbolic.isnan(*args, **kwargs)

Overloaded function.

  1. isnan(x: pydrake.symbolic.Expression) -> drake::symbolic::Formula

  2. isnan(e: pydrake.symbolic.Expression) -> pydrake.symbolic.Formula

Returns a Formula for the predicate isnan(e) to the given expression. This serves as the argument-dependent lookup related to std::isnan(double).

When this formula is evaluated, there are two possible outcomes: - Returns false if the e.Evaluate() is not NaN. - Throws RuntimeError if NaN is detected during evaluation. Note that the evaluation of isnan(e) never returns true.

pydrake.symbolic.Jacobian(*args, **kwargs)

Overloaded function.

  1. Jacobian(f: numpy.ndarray[object[m, 1]], vars: numpy.ndarray[object[m, 1]]) -> numpy.ndarray[object[m, n]]

Computes the Jacobian matrix J of the vector function f with respect to vars. J(i,j) contains ∂f(i)/∂vars(j).

For example, Jacobian([x * cos(y), x * sin(y), x^2], {x, y}) returns the following 3x2 matrix:

Click to expand C++ code...
= |cos(y)   -x * sin(y)|
   |sin(y)    x * cos(y)|
   | 2 * x             0|
Precondition:

{vars is non-empty}.

  1. Jacobian(f: numpy.ndarray[object[m, 1]], vars: numpy.ndarray[object[m, 1]]) -> numpy.ndarray[object[m, n]]

Computes the Jacobian matrix J of the vector function f with respect to vars. J(i,j) contains ∂f(i)/∂vars(j).

Precondition:

vars is non-empty.

pydrake.symbolic.log(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.logical_and(*formulas)
pydrake.symbolic.logical_not(arg0: pydrake.symbolic.Formula) pydrake.symbolic.Formula
pydrake.symbolic.logical_or(*formulas)
pydrake.symbolic.MakeMatrixBinaryVariable(rows: int, cols: int, name: str) numpy.ndarray[object[m, n]]

Creates a dynamically-sized Eigen matrix of symbolic binary variables.

Parameter rows:

The number of rows in the new matrix.

Parameter cols:

The number of cols in the new matrix.

Parameter name:

The common prefix for variables. The (i, j)-th element will be named as name(i, j).

pydrake.symbolic.MakeMatrixBooleanVariable(rows: int, cols: int, name: str) numpy.ndarray[object[m, n]]

Creates a dynamically-sized Eigen matrix of symbolic Boolean variables.

Parameter rows:

The number of rows in the new matrix.

Parameter cols:

The number of cols in the new matrix.

Parameter name:

The common prefix for variables. The (i, j)-th element will be named as name(i, j).

pydrake.symbolic.MakeMatrixContinuousVariable(rows: int, cols: int, name: str) numpy.ndarray[object[m, n]]

Creates a dynamically-sized Eigen matrix of symbolic continuous variables.

Parameter rows:

The number of rows in the new matrix.

Parameter cols:

The number of cols in the new matrix.

Parameter name:

The common prefix for variables. The (i, j)-th element will be named as name(i, j).

pydrake.symbolic.MakeMatrixVariable(rows: int, cols: int, name: str, type: pydrake.symbolic.Variable.Type = <Type.CONTINUOUS: 0>) numpy.ndarray[object[m, n]]

Creates a dynamically-sized Eigen matrix of symbolic variables.

Parameter rows:

The number of rows in the new matrix.

Parameter cols:

The number of cols in the new matrix.

Parameter name:

The common prefix for variables. The (i, j)-th element will be named as name(i, j).

Parameter type:

The type of variables in the matrix.

pydrake.symbolic.MakeVectorBinaryVariable(rows: int, name: str) numpy.ndarray[object[m, 1]]

Creates a dynamically-sized Eigen vector of symbolic binary variables.

Parameter rows:

The size of vector.

Parameter name:

The common prefix for variables. The i-th element will be named as name(i).

pydrake.symbolic.MakeVectorBooleanVariable(rows: int, name: str) numpy.ndarray[object[m, 1]]

Creates a dynamically-sized Eigen vector of symbolic Boolean variables.

Parameter rows:

The size of vector.

Parameter name:

The common prefix for variables. The i-th element will be named as name(i).

pydrake.symbolic.MakeVectorContinuousVariable(rows: int, name: str) numpy.ndarray[object[m, 1]]

Creates a dynamically-sized Eigen vector of symbolic continuous variables.

Parameter rows:

The size of vector.

Parameter name:

The common prefix for variables. The i-th element will be named as name(i).

pydrake.symbolic.MakeVectorVariable(rows: int, name: str, type: pydrake.symbolic.Variable.Type = <Type.CONTINUOUS: 0>) numpy.ndarray[object[m, 1]]

Creates a dynamically-sized Eigen vector of symbolic variables.

Parameter rows:

The size of vector.

Parameter name:

The common prefix for variables. The i-th element will be named as name(i).

Parameter type:

The type of variables in the vector.

pydrake.symbolic.max(arg0: pydrake.symbolic.Expression, arg1: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.min(arg0: pydrake.symbolic.Expression, arg1: pydrake.symbolic.Expression) pydrake.symbolic.Expression
class pydrake.symbolic.Monomial

Represents a monomial, a product of powers of variables with non-negative integer exponents. Note that it does not include the coefficient part of a monomial.

__init__(*args, **kwargs)

Overloaded function.

  1. __init__(self: pydrake.symbolic.Monomial) -> None

Constructs a monomial equal to 1. Namely the total degree is zero.

  1. __init__(self: pydrake.symbolic.Monomial, var: pydrake.symbolic.Variable) -> None

Constructs a monomial from var.

  1. __init__(self: pydrake.symbolic.Monomial, var: pydrake.symbolic.Variable, exponent: int) -> None

Constructs a monomial from var and exponent.

  1. __init__(self: pydrake.symbolic.Monomial, powers: Dict[pydrake.symbolic.Variable, int]) -> None

Constructs a monomial from powers.

Raises

RuntimeError if powers includes a negative exponent.

  1. __init__(self: pydrake.symbolic.Monomial, vars: numpy.ndarray[object[m, 1]], exponents: numpy.ndarray[numpy.int32[m, 1]]) -> None

Constructs a monomial from a vector of variables vars and their corresponding integer exponents exponents. For example, %Monomial([x, y, z], [2, 0, 1]) constructs a monomial x²z.

Precondition:

The size of vars should be the same as the size of exponents.

Raises

RuntimeError if exponents includes a negative integer.

degree(self: pydrake.symbolic.Monomial, v: pydrake.symbolic.Variable) int

Returns the degree of this monomial in a variable v.

EqualTo(self: pydrake.symbolic.Monomial, arg0: pydrake.symbolic.Monomial) bool
Evaluate(*args, **kwargs)

Overloaded function.

  1. Evaluate(self: pydrake.symbolic.Monomial, env: Dict[pydrake.symbolic.Variable, float]) -> float

Evaluates under a given environment env.

Raises
  • RuntimeError if there is a variable in this monomial whose

  • assignment is not provided by env.

  1. Evaluate(self: pydrake.symbolic.Monomial, vars: numpy.ndarray[object[m, 1]], vars_values: numpy.ndarray[numpy.float64[m, n], flags.f_contiguous]) -> numpy.ndarray[numpy.float64[m, 1]]

Evaluates the monomial for a batch of data. We return monomial_vals such that monomial_vals(j) is obtained by substituting vars(i) with vars_values(i, j), note that vars_values.rows() == vars.rows() and vars_values.cols() == monomial_vals.rows().

Parameter vars:

The variables whose value will be substituted. vars must contain all variables in this->GetVariables(). Also vars cannot contain any duplicate variables, namely vars(i) != vars(j) if i != j.

Parameter vars_values:

The i’th column of vars_values is the i’th data for vars.

Raises
  • RuntimeError if vars doesn't contain all the variables in

  • this->GetVariables()`

EvaluatePartial(self: pydrake.symbolic.Monomial, env: Dict[pydrake.symbolic.Variable, float]) Tuple[float, pydrake.symbolic.Monomial]

Partially evaluates using a given environment env. The evaluation result is of type pair<double, Monomial>. The first component (a double) represents the coefficient part while the second component represents the remaining parts of the monomial which was not evaluated.

Example 1. Evaluate with a fully-specified environment (x³*y²).EvaluatePartial({{x, 2}, {y, 3}}) = (2³ * 3² = 8 * 9 = 72, Monomial{} = 1).

Example 2. Evaluate with a partial environment (x³*y²).EvaluatePartial({{x, 2}}) = (2³ = 8, y²).

get_powers(self: pydrake.symbolic.Monomial) Dict[pydrake.symbolic.Variable, int]

Returns the internal representation of Monomial, the map from a base (Variable) to its exponent (int).

GetVariables(self: pydrake.symbolic.Monomial) pydrake.symbolic.Variables

Returns the set of variables in this monomial.

pow_in_place(self: pydrake.symbolic.Monomial, p: int) pydrake.symbolic.Monomial

Returns this monomial raised to p.

Raises

RuntimeError if p is negative.

ToExpression(self: pydrake.symbolic.Monomial) pydrake.symbolic.Expression

Returns a symbolic expression representing this monomial.

total_degree(self: pydrake.symbolic.Monomial) int

Returns the total degree of this monomial.

pydrake.symbolic.MonomialBasis(*args, **kwargs)

Overloaded function.

  1. MonomialBasis(vars: numpy.ndarray[object[m, 1]], degree: int) -> numpy.ndarray[object[m, 1]]

Returns all monomials up to a given degree under the graded reverse lexicographic order. Note that graded reverse lexicographic order uses the total order among Variable which is based on a variable’s unique ID. For example, for a given variable ordering x > y > z, MonomialBasis({x, y, z}, 2) returns a column vector [x^2, xy, y^2, xz, yz, z^2, x, y, z, 1].

Precondition:

vars is a non-empty set.

Precondition:

degree is a non-negative integer.

  1. MonomialBasis(vars: pydrake.symbolic.Variables, degree: int) -> numpy.ndarray[object[m, 1]]

Returns all monomials up to a given degree under the graded reverse lexicographic order. Note that graded reverse lexicographic order uses the total order among Variable which is based on a variable’s unique ID. For example, for a given variable ordering x > y > z, MonomialBasis({x, y, z}, 2) returns a column vector [x^2, xy, y^2, xz, yz, z^2, x, y, z, 1].

Precondition:

vars is a non-empty set.

Precondition:

degree is a non-negative integer.

  1. MonomialBasis(vars_degree: Dict[pydrake.symbolic.Variables, int]) -> numpy.ndarray[object[m, 1]]

Returns all the monomials (in graded reverse lexicographic order) such that the total degree for each set of variables is no larger than a specific degree. For example if x_set = {x₀, x₁} and y_set = {y₀, y₁}, then MonomialBasis({{x_set, 2}, {y_set, 1}}) will include all the monomials, whose total degree of x_set is no larger than 2, and the total degree of y_set is no larger than 1. Hence it can include monomials such as x₀x₁y₀, but not x₀y₀y₁ because the total degree for y_set is 2. So it would return the following set of monomials (ignoring the ordering) {x₀²y₀, x₀²y₁, x₀x₁y₀, x₀x₁y₁, x₁²y₀, x₁²y₀, x₀y₀, x₀y₁, x₁y₀, x₁y₁, x₀², x₀x₁, x₁², x₀, x₁, y₀, y₁, 1}.

Parameter variables_degree:

(vars, degree) maps each set of variables vars to the maximal degree of these variables in the monomial. Namely the summation of the degree of each variable in vars is no larger than degree.

Precondition:

The variables in variables_degree don’t overlap.

Precondition:

The degree in variables_degree are non-negative.

pydrake.symbolic.OddDegreeMonomialBasis(vars: pydrake.symbolic.Variables, degree: int) numpy.ndarray[object[m, 1]]

Returns all odd degree monomials up to a given degree under the graded reverse lexicographic order. A monomial has an odd degree if its total degree is odd. So x²y is an odd degree monomial (degree 3) while xy is not (degree 2). Note that graded reverse lexicographic order uses the total order among Variable which is based on a variable’s unique ID. For example, for a given variable ordering x > y > z, OddDegreeMonomialBasis({x, y, z}, 3) returns a column vector [x³, x²y, xy², y³, x²z, xyz, y²z, xz², yz², z³, x, y, z]

Precondition:

vars is a non-empty set.

Precondition:

degree is a non-negative integer.

class pydrake.symbolic.Polynomial

Represents symbolic polynomials. A symbolic polynomial keeps a mapping from a monomial of indeterminates to its coefficient in a symbolic expression.

A polynomial p has to satisfy an invariant such that p.decision_variables() p.indeterminates() = . We have CheckInvariant() method to check the invariant.

** Operation between a Polynomial and a Variable **

Note that for arithmetic operations ⊕ (where ⊕ can be +,-,*) between a Polynomial p and a Variable v, if the variable v is an indeterminate of the polynomial p, then we regard v as a monomial with indeterminate v with coefficient 1; on the other hand, if the variable v is not an indeterminate of p, then we regard v as a monomial with 0 degree and coefficient v, and v will be appended to decision_variables() of the resulted polynomial. For example, if p = ax²+y where (x,y) are indeterminates and a is a decision variable, then in the operation p + y, the result stores the monomial-to-coefficient mapping as {(x² -> a), (y -> 2)}; on the other hand, in the operation p + b, b is regarded as a 0-degree monomial with coefficient b, and p + b stores the monomial-to-coefficient mapping as {(x² -> a), (y -> 1), (1 -> b)}, with (p + b).decision_variables() being {a, b}. If you want to append the variable v to the indeterminates of the p⊕v, then explicitly convert it to a monomial as p ⊕ symbolic::Monomial(v).

__init__(*args, **kwargs)

Overloaded function.

  1. __init__(self: pydrake.symbolic.Polynomial) -> None

Constructs a zero polynomial.

  1. __init__(self: pydrake.symbolic.Polynomial, map: Dict[pydrake.symbolic.Monomial, pydrake.symbolic.Expression]) -> None

Constructs a polynomial from a map, Monomial → Expression.

  1. __init__(self: pydrake.symbolic.Polynomial, m: pydrake.symbolic.Monomial) -> None

Constructs a polynomial from a monomial m. Note that all variables in m are considered as indeterminates. Note that this implicit conversion is desirable to have a dot product of two Eigen::Vector<Monomial>s return a Polynomial.

  1. __init__(self: pydrake.symbolic.Polynomial, e: pydrake.symbolic.Expression) -> None

Constructs a polynomial from an expression e. Note that all variables in e are considered as indeterminates.

Raises

RuntimeError if e is not a polynomial.

  1. __init__(self: pydrake.symbolic.Polynomial, e: pydrake.symbolic.Expression, indeterminates: pydrake.symbolic.Variables) -> None

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

Note

It collects the intersection of the variables appeared in e and the provided indeterminates.

Raises

RuntimeError if e is not a polynomial in indeterminates.

  1. __init__(self: pydrake.symbolic.Polynomial, e: pydrake.symbolic.Expression, indeterminates: numpy.ndarray[object[m, 1]]) -> None

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

Note

It collects the intersection of the variables appeared in e and the provided indeterminates.

Raises

RuntimeError if e is not a polynomial in indeterminates.

AddProduct(self: pydrake.symbolic.Polynomial, coeff: pydrake.symbolic.Expression, m: pydrake.symbolic.Monomial) pydrake.symbolic.Polynomial

Adds coeff * m to this polynomial.

CoefficientsAlmostEqual(self: pydrake.symbolic.Polynomial, p: pydrake.symbolic.Polynomial, tolerance: float) bool

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

decision_variables(self: pydrake.symbolic.Polynomial) pydrake.symbolic.Variables

Returns the decision variables of this polynomial.

Degree(self: pydrake.symbolic.Polynomial, v: pydrake.symbolic.Variable) int

Returns the highest degree of this polynomial in a variable v.

Differentiate(self: pydrake.symbolic.Polynomial, x: pydrake.symbolic.Variable) pydrake.symbolic.Polynomial

Differentiates this polynomial with respect to the variable x. Note that a variable x can be either a decision variable or an indeterminate.

EqualTo(self: pydrake.symbolic.Polynomial, arg0: pydrake.symbolic.Polynomial) bool

Returns true if this polynomial and p are structurally equal.

Evaluate(self: pydrake.symbolic.Polynomial, env: Dict[pydrake.symbolic.Variable, float]) float

Evaluates this polynomial under a given environment env.

Raises
  • RuntimeError if there is a variable in this polynomial whose

  • assignment is not provided by env.

EvaluateIndeterminates(self: pydrake.symbolic.Polynomial, indeterminates: numpy.ndarray[object[m, 1]], indeterminates_values: numpy.ndarray[numpy.float64[m, n], flags.f_contiguous]) numpy.ndarray[numpy.float64[m, 1]]

Evaluates the polynomial at a batch of indeterminates values.

Parameter indeterminates:

Must include all this->indeterminates()

Parameter indeterminates_values:

Each column of indeterminates_values stores one specific value of indeterminates; indeterminates_values.rows() == indeterminates.rows().

Returns

polynomial_values polynomial_values(j) is obtained by substituting indeterminates(i) in this polynomial with indeterminates_values(i, j) for all i.

Raises
  • RuntimeError if any coefficient in this polynomial is not a

  • constant.

EvaluatePartial(*args, **kwargs)

Overloaded function.

  1. EvaluatePartial(self: pydrake.symbolic.Polynomial, env: Dict[pydrake.symbolic.Variable, float]) -> pydrake.symbolic.Polynomial

Partially evaluates this polynomial using an environment env.

Raises

RuntimeError if NaN is detected during evaluation.

  1. EvaluatePartial(self: pydrake.symbolic.Polynomial, var: pydrake.symbolic.Variable, c: float) -> pydrake.symbolic.Polynomial

Partially evaluates this polynomial by substituting var with c.

Raises

RuntimeError if NaN is detected at any point during evaluation.

EvaluateWithAffineCoefficients(self: pydrake.symbolic.Polynomial, indeterminates: numpy.ndarray[object[m, 1]], indeterminates_values: numpy.ndarray[numpy.float64[m, n], flags.f_contiguous]) Tuple[numpy.ndarray[numpy.float64[m, n]], numpy.ndarray[object[m, 1]], numpy.ndarray[numpy.float64[m, 1]]]

Evaluates the polynomial at a batch of indeterminate values. For a polynomial whose coefficients are affine expressions of decision variables, we evaluate this polynomial on a batch of indeterminate values, and return the matrix representation of the evaluated affine expressions. For example if p(x) = (a+1)x² + b*x where a, b are decision variables, if we evaluate this polynomial on x = 1 and x = 2, then p(x) = a+b+1 and 4a+2b+4 respectively. We return the evaluation result as A * decision_variables + b, where A.row(i) * decision_variables + b(i) is the evaluation of the polynomial on indeterminates_values.col(i).

Parameter indeterminates:

Must include all this->indeterminates()

Parameter indeterminates_values:

A matrix representing a batch of values. Each column of indeterminates_values stores one specific value of indeterminates, where indeterminates_values.rows() == indeterminates.rows().

Parameter A:

The coefficient of the evaluation results.

Parameter decision_variables:

The decision variables in the evaluation results.

Parameter b:

The constant terms in the evaluation results.

Raises
  • RuntimeError if the coefficients of this polynomial is not an

  • affine expression of its decision variables. For example, the

  • polynomial (2+sin(a)) * x² + 1 (where a is a decision variable

  • and x is a indeterminate) doesn't have affine expression as

  • its coefficient 2+sin(a)

Expand(self: pydrake.symbolic.Polynomial) pydrake.symbolic.Polynomial

Expands each coefficient expression and returns the expanded polynomial. If any coefficient is equal to 0 after expansion, then remove that term from the returned polynomial.

indeterminates(self: pydrake.symbolic.Polynomial) pydrake.symbolic.Variables

Returns the indeterminates of this polynomial.

Integrate(*args, **kwargs)

Overloaded function.

  1. Integrate(self: pydrake.symbolic.Polynomial, x: pydrake.symbolic.Variable) -> pydrake.symbolic.Polynomial

Integrates this polynomial with respect to an indeterminate x. Integration with respect to decision variables is not supported yet. If x is not an indeterminate nor decision variable, then it will be added to the list of indeterminates.

Raises

RuntimeError if x is a decision variable.

  1. Integrate(self: pydrake.symbolic.Polynomial, x: pydrake.symbolic.Variable, a: float, b: float) -> pydrake.symbolic.Polynomial

Computes the definite integrate of this polynomial with respect to the indeterminate x over the domain [a, b]. Integration with respect to decision variables is not supported yet.

Raises

RuntimeError if x is a decision variable.

IsEven(self: pydrake.symbolic.Polynomial) bool

Returns true if the polynomial is even, namely p(x) = p(-x). Meaning that the coefficient for all odd-degree monomials are 0. Returns false otherwise. Note that this is different from the p.TotalDegree() being an even number.

IsOdd(self: pydrake.symbolic.Polynomial) bool

Returns true if the polynomial is odd, namely p(x) = -p(-x). Meaning that the coefficient for all even-degree monomials are 0. Returns false otherwise. Note that this is different from the p.TotalDegree() being an odd number.

Jacobian(self: pydrake.symbolic.Polynomial, vars: numpy.ndarray[object[m, 1]]) numpy.ndarray[object[1, n]]

Computes the Jacobian matrix J of the polynomial with respect to vars. J(0,i) contains ∂f/∂vars(i).

monomial_to_coefficient_map(self: pydrake.symbolic.Polynomial) Dict[pydrake.symbolic.Monomial, pydrake.symbolic.Expression]

Returns the mapping from a Monomial to its corresponding coefficient of this polynomial. We maintain the invariance that for any [monomial, coeff] pair in monomial_to_coefficient_map(), symbolic:is_zero(coeff) is false.

RemoveTermsWithSmallCoefficients(self: pydrake.symbolic.Polynomial, coefficient_tol: float) pydrake.symbolic.Polynomial

Removes the terms whose absolute value of the coefficients are smaller than or equal to coefficient_tol. For example, if the polynomial is 2x² + 3xy + 10⁻⁴x - 10⁻⁵, then after calling RemoveTermsWithSmallCoefficients(1e-3), the returned polynomial becomes 2x² + 3xy.

Parameter coefficient_tol:

A positive scalar.

Returns polynomial_cleaned:

A polynomial whose terms with small coefficients are removed.

Roots(self: pydrake.symbolic.Polynomial) numpy.ndarray[numpy.complex128[m, 1]]

Returns the roots of a univariate polynomial with constant coefficients as a column vector. There is no specific guarantee on the order of the returned roots.

Raises
  • RuntimeError if this is not univariate with constant

  • coefficients.

SetIndeterminates(self: pydrake.symbolic.Polynomial, new_indeterminates: pydrake.symbolic.Variables) None

Sets the indeterminates to new_indeterminates. Changing the indeterminates would change monomial_to_coefficient_map(), and also potentially the degree of the polynomial. Here is an example.

Click to expand C++ code...
// p is a quadratic polynomial with x being the indeterminates.
symbolic::Polynomial p(a * x * x + b * x + c, {x});
// p.monomial_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.monomial_to_coefficient_map() now is {a: x * x, b: x, c: 1}.
std::cout << p.TotalDegree(); // prints 1.
SubstituteAndExpand(self: pydrake.symbolic.Polynomial, indeterminate_substitution: Dict[pydrake.symbolic.Variable, pydrake.symbolic.Polynomial], substitutions_cached_data: drake::symbolic::Polynomial::SubstituteAndExpandCacheData = None) pydrake.symbolic.Polynomial

Substitutes the monomials of this polynomial with new polynomial expressions and expand the polynomial to the monomial basis. For example, consider the substitution x = a(1-y) into the polynomial x¹⁴ + (x−1)². Repeatedly expanding the powers of x can take a long time using factory methods, so we store intermediate computations in the substitution map to avoid recomputing very high powers.

Parameter indeterminate_substitution:

The substitutions of every indeterminate with the new desired expression. This map must contain each element of indeterminates(). For performance reasons, it is recommended that this map contains Expanded polynomials as its values, but this is not necessary.

Parameter substitutions_cached_data:

A container caching the higher order expansions of the indeterminate_substitutions. Typically, the first time an indeterminate_substitution is performed, this will be empty. If the same indeterminate_substitutions is used for multiple polynomials, passing this value will enable the user to re-use the expansions across multiple calls. For example, suppose we wish to perform the substitution x = a(1-y) into the polynomials p1 = x¹⁴ + (x−1)² and p2 = x⁷. A user may call p1.SubstituteAndExpand({x : a(1-y), substitutions_cached_data}) where substitutions_cached_data is a pointer to an empty container. As part of computing the expansion of p1, the expansion of x⁷ may get computed and stored in substitutions_cached_data, and so a subsequent call of p2.SubstituteAndExpand({x : a(1-y)}, substitutions_cached_data) would be very fast.

Never reuse substitutions_cached_data if indeterminate_substitutions changes as this function will then compute an incorrect result.

Note that this function is NOT responsible for ensuring that substitutions_cached_data is consistent i.e. this method will not throw an error if substitutions_cached_data contains the inconsistent substitutions {x: y, x²: 2y}. To ensure correct results, ensure that the passed substitutions_cached_data object is consistent with indeterminate_substitutions. The easiest way to do this is to pass a pointer to an empty substitutions_cached_data or nullopt to this function.

ToExpression(self: pydrake.symbolic.Polynomial) pydrake.symbolic.Expression

Returns an equivalent symbolic expression of this polynomial.

TotalDegree(self: pydrake.symbolic.Polynomial) int

Returns the total degree of this polynomial.

pydrake.symbolic.positive_semidefinite(m: numpy.ndarray[object[m, n], flags.f_contiguous]) pydrake.symbolic.Formula

Returns a symbolic formula constraining m to be a positive-semidefinite matrix. By definition, a symmetric matrix m is positive-semidefinte if xᵀ m x ≥ 0 for all vector x ∈ ℝⁿ.

Raises

RuntimeError if m is not symmetric.

Note

This method checks if m is symmetric, which can be costly. If you want to avoid it, please consider using positive_semidefinite(m.triangularView<Eigen::Lower>()) or positive_semidefinite(m.triangularView<Eigen::Upper>()) instead of positive_semidefinite(m).

pydrake.symbolic.pow(*args, **kwargs)

Overloaded function.

  1. pow(arg0: pydrake.symbolic.Expression, arg1: float) -> pydrake.symbolic.Expression

  2. pow(arg0: pydrake.symbolic.Expression, arg1: pydrake.symbolic.Expression) -> pydrake.symbolic.Expression

class pydrake.symbolic.RationalFunction

Represents symbolic rational function. A function f(x) is a rational function, if f(x) = p(x) / q(x), where both p(x) and q(x) are polynomials of x. Note that rational functions are closed under (+, -, x, /). One application of rational function is in polynomial optimization, where we represent (or approximate) functions using rational functions, and then convert the constraint f(x) = h(x) (where h(x) is a polynomial) to a polynomial constraint p(x) - q(x) * h(x) = 0, or convert the inequality constraint f(x) >= h(x) as p(x) - q(x) * h(x) >= 0 if we know q(x) > 0.

This class represents a special subset of the symbolic::Expression. While a symbolic::Expression can represent a rational function, extracting the numerator and denominator, generally, is quite difficult; for instance, from p1(x) / q1(x) + p2(x) / q2(x) + … + pn(x) / qn(x). This class’s explicit structure facilitates this decomposition.

__init__(*args, **kwargs)

Overloaded function.

  1. __init__(self: pydrake.symbolic.RationalFunction) -> None

Constructs a zero rational function 0 / 1.

  1. __init__(self: pydrake.symbolic.RationalFunction, numerator: pydrake.symbolic.Polynomial, denominator: pydrake.symbolic.Polynomial) -> None

Constructs the rational function: numerator / denominator.

Parameter numerator:

The numerator of the fraction.

Parameter denominator:

The denominator of the fraction.

Precondition:

denominator is not equal to zero.

Precondition:

None of the indeterminates in the numerator can be decision variables in the denominator; similarly none of the indeterminates in the denominator can be decision variables in the numerator.

  1. __init__(self: pydrake.symbolic.RationalFunction, p: pydrake.symbolic.Polynomial) -> None

Constructs the rational function: p / 1. Note that we use 1 as the denominator.

Parameter p:

The numerator of the rational function.

  1. __init__(self: pydrake.symbolic.RationalFunction, m: pydrake.symbolic.Monomial) -> None

Constructs the rational function: m / 1 for any type which can be cast to a monomial

Parameter m:

The numerator of the rational function.

  1. __init__(self: pydrake.symbolic.RationalFunction, c: float) -> None

Constructs the rational function: c / 1. Note that we use 1 as the denominator.

Parameter c:

The numerator of the rational function.

  1. __init__(self: pydrake.symbolic.RationalFunction) -> None

Constructs a zero rational function 0 / 1.

denominator(self: pydrake.symbolic.RationalFunction) pydrake.symbolic.Polynomial

Getter for the denominator.

EqualTo(self: pydrake.symbolic.RationalFunction, f: pydrake.symbolic.RationalFunction) bool

Returns true if this rational function and f are structurally equal.

Evaluate(self: pydrake.symbolic.RationalFunction, env: Dict[pydrake.symbolic.Variable, float]) float

Evaluates this rational function under a given environment env.

Raises
  • RuntimeError if there is a variable in this rational function

  • whose assignment is not provided by env.

numerator(self: pydrake.symbolic.RationalFunction) pydrake.symbolic.Polynomial

Getter for the numerator.

SetIndeterminates(self: pydrake.symbolic.RationalFunction, new_indeterminates: pydrake.symbolic.Variables) None

Sets the indeterminates of the numerator and denominator polynomials

ToExpression(self: pydrake.symbolic.RationalFunction) pydrake.symbolic.Expression

Returns an equivalent symbolic expression of this rational function.

pydrake.symbolic.ReplaceBilinearTerms(e: pydrake.symbolic.Expression, x: numpy.ndarray[object[m, 1]], y: numpy.ndarray[object[m, 1]], W: numpy.ndarray[object[m, n], flags.f_contiguous]) pydrake.symbolic.Expression

Replaces all the bilinear product terms in the expression e, with the corresponding terms in W, where W represents the matrix x * yᵀ, such that after replacement, e does not have bilinear terms involving x and y. For example, if e = x(0)*y(0) + 2 * x(0)*y(1) + x(1) * y(1) + 3 * x(1), e has bilinear terms x(0)*y(0), x(0) * y(1) and x(2) * y(1), if we call ReplaceBilinearTerms(e, x, y, W) where W(i, j) represent the term x(i) * y(j), then this function returns W(0, 0) + 2 * W(0, 1) + W(1, 1) + 3 * x(1).

Parameter e:

An expression potentially contains bilinear products between x and y.

Parameter x:

The bilinear product between x and y will be replaced by the corresponding term in W.

Raises

RuntimeError if x contains duplicate entries.

Parameter y:

The bilinear product between x and y will be replaced by the corresponding term in ``W.

Raises

RuntimeError if y` contains duplicate entries.

Parameter W:

Bilinear product term x(i) * y(j) will be replaced by W(i, j). If W(i,j) is not a single variable, but an expression, then this expression cannot contain a variable in either x or y.

Raises
  • RuntimeError, if W(i, j) is not a single variable, and also

  • contains a variable in x or y.

Precondition:

W.rows() == x.rows() and W.cols() == y.rows().

Returns

The symbolic expression after replacing x(i) * y(j) with W(i, j).

pydrake.symbolic.sin(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
class pydrake.symbolic.SinCos

Represents a pair of Variables corresponding to sin(q) and cos(q).

__init__(self: pydrake.symbolic.SinCos, s: pydrake.symbolic.Variable, c: pydrake.symbolic.Variable, type: pydrake.symbolic.SinCosSubstitutionType = <SinCosSubstitutionType.kAngle: 0>) None
property c

cos variable.

property s

sin variable.

property type

Allows a user to specify non-default substitutions, such as using half-angle formulas.

class pydrake.symbolic.SinCosSubstitutionType

Members:

kAngle : Substitutes s <=> sin(q), c <=> cos(q).

kHalfAnglePreferSin : Substitutes s <=> sin(q/2), c <=> cos(q/2), and prefers sin when the

choice is ambiguous; e.g. cos(q) => 1 - 2s².

kHalfAnglePreferCos : Substitutes s <=> sin(q/2), c <=> cos(q/2), and prefers cos when the

choice is ambiguous; e.g. cos(q) => 2c² - 1.

__init__(self: pydrake.symbolic.SinCosSubstitutionType, value: int) None
kAngle = <SinCosSubstitutionType.kAngle: 0>
kHalfAnglePreferCos = <SinCosSubstitutionType.kHalfAnglePreferCos: 2>
kHalfAnglePreferSin = <SinCosSubstitutionType.kHalfAnglePreferSin: 1>
property name
property value
pydrake.symbolic.sinh(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.sqrt(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.Substitute(*args, **kwargs)

Overloaded function.

  1. Substitute(m: numpy.ndarray[object[m, n]], subst: Dict[pydrake.symbolic.Variable, pydrake.symbolic.Expression]) -> numpy.ndarray[object[m, n]]

Substitutes a symbolic matrix m using a given substitution subst.

Returns

a matrix of symbolic expressions whose size is the size of m.

Raises

RuntimeError if NaN is detected during substitution.

  1. Substitute(m: numpy.ndarray[object[m, n]], var: pydrake.symbolic.Variable, e: pydrake.symbolic.Expression) -> numpy.ndarray[object[m, n]]

Substitutes var with e in a symbolic matrix m.

Returns

a matrix of symbolic expressions whose size is the size of m.

Raises

RuntimeError if NaN is detected during substitution.

  1. Substitute(e: pydrake.symbolic.Expression, subs: Dict[pydrake.symbolic.Variable, pydrake.symbolic.SinCos]) -> pydrake.symbolic.Expression

Given a substitution map q => {s, c}, substitutes instances of sin(q) and cos(q) in e with s and c, with partial support for trigonometric expansions. For instance,

Click to expand C++ code...
Variable x{"x"}, y{"y"};
Variable sx{"sx"}, cx{"cx"}, sy{"sy"}, cy{"cy"};
SinCosSubstitution subs;
subs.emplace(x, SinCos(sx, cx));
subs.emplace(y, SinCos(sy, cy));
Expression e = Substitute(x * sin(x + y), subs);

will result in the expression x * (sx*cy + cx*sy).

Parameter half_angle:

If true, then the same workflow replaces instances of sin(q/2) and cos(q/2) in e will be replaced with s, and c.

Default: false.

The half-angle representation is more natural in many analysis computations for robots, for instance: https://underactuated.csail.mit.edu/lyapunov.html#trig_quadratic

Raises
  • RuntimeError if a trigonometric function is not a trigonometric

  • polynomial in q or if the e requires a trigonometric

  • expansion that not supported yet.

  1. Substitute(m: numpy.ndarray[object[m, n]], subs: Dict[pydrake.symbolic.Variable, pydrake.symbolic.SinCos]) -> numpy.ndarray[object[m, n]]

Matrix version of sin/cos substitution.

class pydrake.symbolic.SubstituteAndExpandCacheData

An encapsulated data type for use with the method SubstituteAndExpand.

__init__(self: pydrake.symbolic.SubstituteAndExpandCacheData) None
get_data(self: pydrake.symbolic.SubstituteAndExpandCacheData) Dict[pydrake.symbolic.Monomial, pydrake.symbolic.Polynomial]
pydrake.symbolic.SubstituteStereographicProjection(e: drake::symbolic::Polynomial, sin_cos: List[pydrake.symbolic.SinCos], t: numpy.ndarray[object[m, 1]]) drake::symbolic::RationalFunction

Substitutes the variables representing sine and cosine functions with their stereographic projection. We replace cosθᵢ with (1-tᵢ²)/(1+tᵢ²), and sinθᵢ with 2tᵢ/(1+tᵢ²), and get a rational polynomial. The indeterminates of this rational polynomial are t together with the indeterminates in e that are not cosθ or sinθ. If the input expression doesn’t contain the sine and cosine functions, then the returned rational has denominator being 1. Notice that the indeterminates of e can include variables other than cosθ and sinθ, and we impose no requirements on these variables that are not cosθ or sinθ.

Parameter e:

The symbolic polynomial to be substituted.

Parameter sin_cos:

sin_cos(i) is the pair of variables (sᵢ, cᵢ), (where sᵢ=sinθᵢ, cᵢ=cosθᵢ) as documented above.

Parameter t:

New variables to express cos and sin as rationals of t. tᵢ = tan(θᵢ/2).

Precondition:

t.rows() == sin_cos.size()

Returns

e_rational The rational polynomial of e after replacement. The indeterminates of the polynomials are t together with the indeterminates in e that are not cosθ or sinθ. Example

Click to expand C++ code...
std::vector<SinCos> sin_cos;
sin_cos.emplace_back(symbolic::Variable("s0"), symbolic::Variable("c0"));
sin_cos.emplace_back(symbolic::Variable("s1"), symbolic::Variable("c1"));
Vector2<symbolic::Variable> t(symbolic::Variable("t0"),
                              symbolic::Variable("t1"));
const auto e_rational =
SubstituteStereographicProjection(t(0) * sin_cos[0].s*sin_cos[1].c + 1,
                                  sin_cos, t);
// e_rational should be
// (2*t0*t0*(1-t1*t1) + (1+t0*t0)*(1+t1*t1))
// --------------------------------------------
//        ((1+t0*t0)*(1+t1*t1))
pydrake.symbolic.tan(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.tanh(arg0: pydrake.symbolic.Expression) pydrake.symbolic.Expression
pydrake.symbolic.TaylorExpand(f: drake::symbolic::Expression, a: Dict[pydrake.symbolic.Variable, float], order: int) drake::symbolic::Expression

Returns the Taylor series expansion of f around a of order order.

Parameter f:

Symbolic expression to approximate using Taylor series expansion.

Parameter a:

Symbolic environment which specifies the point of approximation. If a partial environment is provided, the unspecified variables are treated as symbolic variables (e.g. decision variable).

Parameter order:

Positive integer which specifies the maximum order of the resulting polynomial approximating f around a.

pydrake.symbolic.to_sympy(x: Union[float, int, bool, pydrake.symbolic.Variable, pydrake.symbolic.Expression, pydrake.symbolic.Formula], *, memo: Dict = None) Union[float, int, bool, sympy.Expr]

Converts a pydrake object to the corresponding SymPy Expr.

Certain expressions are not supported and will raise NotImplementedError. (Most acutely, note that int is not yet supported.)

This function aims to support the latest contemporaneous version of SymPy as of Drake’s release date.

Parameters
  • x – The pydrake object to be converted.

  • memo – (Optional) Mapping between Drake variables and SymPy variables. Converting a pydrake.symbolic.Variable to SymPy will look up the Drake variable’s v.get_id() as a key in this dictionary. If a value is found, it will be used as the SymPy atom for that Drake variable. Otherwise, a new SymPy variable will be created and used, and both mappings {drake_var.get_id(): sympy_var} and {sympy_var: drake_var} will be inserted into memo.

See also pydrake.symbolic.from_sympy().

pydrake.symbolic.ToLatex(*args, **kwargs)

Overloaded function.

  1. ToLatex(e: pydrake.symbolic.Expression, precision: int = 3) -> str

Generates a LaTeX string representation of e with floating point coefficients displayed using precision.

  1. ToLatex(f: pydrake.symbolic.Formula, precision: int = 3) -> str

Generates a LaTeX string representation of f with floating point coefficients displayed using precision.

  1. ToLatex(M: numpy.ndarray[object[m, n]], precision: int = 3) -> str

Generates a LaTeX string representation of M with floating point coefficients displayed using precision.

  1. ToLatex(M: numpy.ndarray[numpy.float64[m, n]], precision: int = 3) -> str

Generates a LaTeX string representation of M with floating point coefficients displayed using precision.

pydrake.symbolic.uninterpreted_function(name: str, arguments: List[pydrake.symbolic.Expression]) pydrake.symbolic.Expression

Constructs an uninterpreted-function expression with name and arguments. An uninterpreted function is an opaque function that has no other property than its name and a list of its arguments. This is useful to applications where it is good enough to provide abstract information of a function without exposing full details. Declaring sparsity of a system is a typical example.

class pydrake.symbolic.Variable

Represents a symbolic variable.

Note

Expression::Evaluate and Formula::Evaluate methods take a symbolic environment (Variable → double) and a random number generator. When an expression or a formula includes random variables, Evaluate methods use the random number generator to draw a number for a random variable from the given distribution. Then this numeric value is used to substitute all the occurrences of the corresponding random variable in an expression or a formula.

__init__(self: pydrake.symbolic.Variable, name: str, type: pydrake.symbolic.Variable.Type = <Type.CONTINUOUS: 0>) None

Constructs a variable with a string. If not specified, it has CONTINUOUS type by default.

abs(self: pydrake.symbolic.Variable) drake::symbolic::Expression
acos(self: pydrake.symbolic.Variable) drake::symbolic::Expression
arccos(self: pydrake.symbolic.Variable) drake::symbolic::Expression
arcsin(self: pydrake.symbolic.Variable) drake::symbolic::Expression
arctan(self: pydrake.symbolic.Variable) drake::symbolic::Expression
arctan2(self: pydrake.symbolic.Variable, x: pydrake.symbolic.Variable) drake::symbolic::Expression

Uses self for y in arctan2(y, x).

asin(self: pydrake.symbolic.Variable) drake::symbolic::Expression
atan(self: pydrake.symbolic.Variable) drake::symbolic::Expression
atan2(self: pydrake.symbolic.Variable, x: pydrake.symbolic.Variable) drake::symbolic::Expression

Uses self for y in atan2(y, x).

ceil(self: pydrake.symbolic.Variable) drake::symbolic::Expression
cos(self: pydrake.symbolic.Variable) drake::symbolic::Expression
cosh(self: pydrake.symbolic.Variable) drake::symbolic::Expression
EqualTo(self: pydrake.symbolic.Variable, arg0: pydrake.symbolic.Variable) bool

Checks the equality of two variables based on their ID values.

exp(self: pydrake.symbolic.Variable) drake::symbolic::Expression
floor(self: pydrake.symbolic.Variable) drake::symbolic::Expression
get_id(self: pydrake.symbolic.Variable) int
get_name(self: pydrake.symbolic.Variable) str
get_type(self: pydrake.symbolic.Variable) pydrake.symbolic.Variable.Type
is_dummy(self: pydrake.symbolic.Variable) bool

Checks if this is the variable created by the default constructor.

log(self: pydrake.symbolic.Variable) drake::symbolic::Expression
max(self: pydrake.symbolic.Variable, arg0: pydrake.symbolic.Variable) drake::symbolic::Expression
min(self: pydrake.symbolic.Variable, arg0: pydrake.symbolic.Variable) drake::symbolic::Expression
pow(*args, **kwargs)

Overloaded function.

  1. pow(self: pydrake.symbolic.Variable, arg0: float) -> drake::symbolic::Expression

  2. pow(self: pydrake.symbolic.Variable, arg0: drake::symbolic::Expression) -> drake::symbolic::Expression

sin(self: pydrake.symbolic.Variable) drake::symbolic::Expression
sinh(self: pydrake.symbolic.Variable) drake::symbolic::Expression
sqrt(self: pydrake.symbolic.Variable) drake::symbolic::Expression
tan(self: pydrake.symbolic.Variable) drake::symbolic::Expression
tanh(self: pydrake.symbolic.Variable) drake::symbolic::Expression
class Type

Members:

CONTINUOUS : A CONTINUOUS variable takes a double value.

INTEGER : An INTEGER variable takes an int value.

BINARY : A BINARY variable takes an integer value from {0, 1}.

BOOLEAN : A BOOLEAN variable takes a bool value.

RANDOM_UNIFORM : A random variable whose value will be drawn from uniform real

distributed ∈ [0,1).

RANDOM_GAUSSIAN : A random variable whose value will be drawn from mean-zero,

unit-variance normal.

RANDOM_EXPONENTIAL : A random variable whose value will be drawn from exponential

distribution with λ=1.

__init__(self: pydrake.symbolic.Variable.Type, value: int) None
BINARY = <Type.BINARY: 2>
BOOLEAN = <Type.BOOLEAN: 3>
CONTINUOUS = <Type.CONTINUOUS: 0>
INTEGER = <Type.INTEGER: 1>
property name
RANDOM_EXPONENTIAL = <Type.RANDOM_EXPONENTIAL: 6>
RANDOM_GAUSSIAN = <Type.RANDOM_GAUSSIAN: 5>
RANDOM_UNIFORM = <Type.RANDOM_UNIFORM: 4>
property value
class pydrake.symbolic.Variables

Represents a set of variables.

This class is based on std::set<Variable>. The intent is to add things that we need including set-union (Variables::insert, operator+, operator+=), set-minus (Variables::erase, operator-, operator-=), and subset/superset checking functions (Variables::IsSubsetOf, Variables::IsSupersetOf, Variables::IsStrictSubsetOf, Variables::IsStrictSupersetOf).

__init__(*args, **kwargs)

Overloaded function.

  1. __init__(self: pydrake.symbolic.Variables) -> None

Default constructor.

  1. __init__(self: pydrake.symbolic.Variables, arg0: numpy.ndarray[object[m, 1]]) -> None

Constructs from an Eigen vector of variables.

empty(self: pydrake.symbolic.Variables) bool

Checks if this set is empty or not.

EqualTo(self: pydrake.symbolic.Variables, arg0: pydrake.symbolic.Variables) bool
erase(*args, **kwargs)

Overloaded function.

  1. erase(self: pydrake.symbolic.Variables, key: pydrake.symbolic.Variable) -> int

Erases key from a set. Return number of erased elements (0 or 1).

  1. erase(self: pydrake.symbolic.Variables, vars: pydrake.symbolic.Variables) -> int

Erases variables in vars from a set. Return number of erased elements ([0, vars.size()]).

include(self: pydrake.symbolic.Variables, key: pydrake.symbolic.Variable) bool

Return true if key is included in the Variables.

insert(*args, **kwargs)

Overloaded function.

  1. insert(self: pydrake.symbolic.Variables, var: pydrake.symbolic.Variable) -> None

Inserts a variable var into a set.

  1. insert(self: pydrake.symbolic.Variables, vars: pydrake.symbolic.Variables) -> None

Inserts variables in vars into a set.

IsStrictSubsetOf(self: pydrake.symbolic.Variables, vars: pydrake.symbolic.Variables) bool

Return true if vars is a strict subset of the Variables.

IsStrictSupersetOf(self: pydrake.symbolic.Variables, vars: pydrake.symbolic.Variables) bool

Return true if vars is a strict superset of the Variables.

IsSubsetOf(self: pydrake.symbolic.Variables, vars: pydrake.symbolic.Variables) bool

Return true if vars is a subset of the Variables.

IsSupersetOf(self: pydrake.symbolic.Variables, vars: pydrake.symbolic.Variables) bool

Return true if vars is a superset of the Variables.

size(self: pydrake.symbolic.Variables) int

Returns the number of elements.

to_string(self: pydrake.symbolic.Variables) str

Returns string representation of Variables.