pydrake.trajectories

class pydrake.trajectories.PiecewisePolynomial

Bases: pydrake.trajectories.PiecewiseTrajectory

A scalar multi-variate piecewise polynomial.

PiecewisePolynomial represents a list of contiguous segments in a scalar independent variable (typically corresponding to time) with Polynomials defined at each segment. We call the output from evaluating the PiecewisePolynomial at the scalar independent variable “the output”, and that output can be either a Eigen MatrixX<T> (if evaluated using value()) or a scalar (if evaluated using scalar_value()).

An example of a piecewise polynomial is a function of m segments in time, where a different polynomial is defined for each segment. For a specific example, consider the absolute value function over the interval [-1, 1]. We can define a PiecewisePolynomial over this interval using breaks at t = { -1.0, 0.0, 1.0 }, and “knots” of abs(t).

// Construct the PiecewisePolynomial.
const std::vector<double> breaks = { -1.0, 0.0, 1.0 };
std::vector<Eigen::MatrixXd> knots(3);
for (int i = 0; i < static_cast<int>(breaks.size()); ++i) {
  knots[i].resize(1, 1);
  knots[i](0, 0) = std::abs(breaks[i]);
}
const auto pp = PiecewisePolynomial<double>::FirstOrderHold(breaks, knots);
const int row = 0, col = 0;

// Evaluate the PiecewisePolynomial at some values.
std::cout << pp.value(-.5)(row, col) << std::endl;    // Outputs 0.5.
std::cout << pp.value(0.0)(row, col) << std::endl;    // Outputs 0.0;

// Show how we can evaluate the first derivative (outputs -1.0).
std::cout << pp.derivative(1).value(-.5)(row, col) << std::endl;

PiecewisePolynomial objects can be added, subtracted, and multiplied. They cannot be divided because Polynomials are not closed under division.

Warning

PiecewisePolynomial silently clips input evaluations outside of the range defined by the breaks. So pp.value(-2.0, row, col) in the example above would evaluate to -1.0. See value().

Template parameter T:
is a scalar type.

Explicit instantiations are provided for:

  • double
ConcatenateInTime(self: pydrake.trajectories.PiecewisePolynomial, other: pydrake.trajectories.PiecewisePolynomial) → None

Concatenates other to the end of this.

Warning

The resulting PiecewisePolynomial will only be continuous to the degree that the first Polynomial of other is continuous with the last Polynomial of this. See warning about evaluating discontinuous derivatives at breaks in derivative().

Parameter other:
PiecewisePolynomial instance to concatenate.
Raises:
  • RuntimeError if trajectories’ dimensions do not match each other
  • (either rows() or cols() does not match between this and
  • other).
Raises:
  • RuntimeError if this->end_time() and other->start_time()
  • are not within PiecewiseTrajectory<T>::kEpsilonTime from each
  • other.
Cubic(*args, **kwargs)

Overloaded function.

  1. Cubic(breaks: numpy.ndarray[float64[m, 1]], knots: numpy.ndarray[float64[m, n], flags.f_contiguous], knots_dot_start: numpy.ndarray[float64[m, 1]], knots_dot_end: numpy.ndarray[float64[m, 1]]) -> pydrake.trajectories.PiecewisePolynomial

Version of Cubic(breaks, knots, knots_dot_start, knots_dot_end) that uses vector knots and Eigen VectorXd / MatrixX<T> inputs. Each column of knots represents a knot point.

Precondition:
knots.cols() == breaks.size().
Raises:
  • RuntimeError under the conditions specified under
  • coefficient_construction_methods.
  1. Cubic(breaks: numpy.ndarray[float64[m, 1]], knots: numpy.ndarray[float64[m, n], flags.f_contiguous], knots_dot: numpy.ndarray[float64[m, n], flags.f_contiguous]) -> pydrake.trajectories.PiecewisePolynomial

Version of Cubic(breaks, knots, knots_dot) that uses vector knots and Eigen VectorXd / MatrixX<T> inputs. Corresponding columns of knots and knots_dot are used as the knot point and independent variable derivative, respectively.

Precondition:
knots.cols() == knots_dot.cols() == breaks.size().
  1. Cubic(breaks: numpy.ndarray[float64[m, 1]], knots: numpy.ndarray[float64[m, n], flags.f_contiguous], periodic_end: bool) -> pydrake.trajectories.PiecewisePolynomial

Version of Cubic(breaks, knots) that uses vector knots and Eigen VectorXd / MatrixX<T> inputs. Each column of knots represents a knot point.

Precondition:
knots.cols() == breaks.size().
FirstOrderHold(arg0: numpy.ndarray[float64[m, 1]], arg1: numpy.ndarray[float64[m, n], flags.f_contiguous]) → pydrake.trajectories.PiecewisePolynomial

Version of FirstOrderHold(breaks, knots) that uses vector knots and Eigen VectorXd / MatrixX<T> inputs. Each column of knots represents a knot point.

Precondition:
knots.cols() == breaks.size()
Raises:
  • RuntimeError under the conditions specified under
  • coefficient_construction_methods.
Pchip(breaks: numpy.ndarray[float64[m, 1]], knots: numpy.ndarray[float64[m, n], flags.f_contiguous], zero_end_point_derivatives: bool = False) → pydrake.trajectories.PiecewisePolynomial

Version of Pchip(breaks, knots, zero_end_point_derivatives) that uses vector knots and Eigen VectorXd / MatrixX<T> inputs. Each column of knots represents a knot point.

Precondition:
knots.cols() == breaks.size().
Raises:
  • RuntimeError under the conditions specified under
  • coefficient_construction_methods.
ZeroOrderHold(arg0: numpy.ndarray[float64[m, 1]], arg1: numpy.ndarray[float64[m, n], flags.f_contiguous]) → pydrake.trajectories.PiecewisePolynomial

Version of ZeroOrderHold(breaks, knots) that uses vector knots and Eigen VectorXd/MatrixX<T> inputs. Each column of knots represents a knot point.

Precondition:
knots.cols() == breaks.size()
Raises:
  • RuntimeError under the conditions specified under
  • coefficient_construction_methods.
cols(self: pydrake.trajectories.PiecewisePolynomial) → int

Returns the column count of the output matrices.

Raises:RuntimeError if empty().
derivative(self: pydrake.trajectories.PiecewisePolynomial, derivative_order: int = 1) → pydrake.trajectories.PiecewisePolynomial

Returns a PiecewisePolynomial where each segment is the specified derivative of the corresponding segment in this. Any rules or limitations of Polynomial::derivative() also apply to this function.

Derivatives evaluated at non-differentiable points return the value at the left hand side of the interval.

Parameter derivative_order:
The order of the derivative, namely, if derivative_order = n, the n’th derivative of the polynomial will be returned.

Warning

In the event of discontinuous derivatives evaluated at breaks, it is not defined which polynomial (i.e., to the left or right of the break) will be the one that is evaluated at the break.

getPolynomial(self: pydrake.trajectories.PiecewisePolynomial, segment_index: int, row: int = 0, col: int = 0) → pydrake.polynomial.Polynomial

Gets the Polynomial with the given matrix row and column index that corresponds to the given segment index. Equivalent to getPolynomialMatrix(segment_index)(row, col).

Note

Calls PiecewiseTrajectory<T>::segment_number_range_check() to validate segment_index.

getPolynomialMatrix(self: pydrake.trajectories.PiecewisePolynomial, segment_index: int) → numpy.ndarray[object[m, n]]

Gets the matrix of Polynomials corresponding to the given segment index.

Warning

segment_index is not checked for validity.

getSegmentPolynomialDegree(self: pydrake.trajectories.PiecewisePolynomial, segment_index: int, row: int = 0, col: int = 0) → int

Gets the degree of the Polynomial with the given matrix row and column index that corresponds to the given segment index. Equivalent to getPolynomial(segment_index, row, col).GetDegree().

isApprox(self: pydrake.trajectories.PiecewisePolynomial, other: pydrake.trajectories.PiecewisePolynomial, tol: float) → bool

Checks whether a PiecewisePolynomial is approximately equal to this one.

Checks that every coefficient of other is within tol of the corresponding coefficient of this PiecewisePolynomial.

Raises:
  • RuntimeError if any Polynomial in either PiecewisePolynomial is
  • not univariate.
rows(self: pydrake.trajectories.PiecewisePolynomial) → int

Returns the row count of the output matrices.

Raises:RuntimeError if empty().
setPolynomialMatrixBlock(self: pydrake.trajectories.PiecewisePolynomial, replacement: numpy.ndarray[object[m, n]], segment_index: int, row_start: int = 0, col_start: int = 0) → None

Replaces the specified block of the PolynomialMatrix at the given segment index.

Note

Calls PiecewiseTrajectory<T>::segment_number_range_check() to validate segment_index.

Warning

This code relies upon Eigen to verify that the replacement block is not too large.

shiftRight(self: pydrake.trajectories.PiecewisePolynomial, offset: float) → None

Adds offset to all of the breaks. offset need not be a non-negative number.

Note

has no effect if empty().

slice(self: pydrake.trajectories.PiecewisePolynomial, start_segment_index: int, num_segments: int) → pydrake.trajectories.PiecewisePolynomial

Returns the PiecewisePolynomial comprising the num_segments segments starting at the specified start_segment_index.

Note

Calls PiecewiseTrajectory<T>::segment_number_range_check() to validate segment_index.

value(self: pydrake.trajectories.PiecewisePolynomial, t: float) → numpy.ndarray[float64[m, n]]

Evaluates the PiecewisePolynomial at the given time t.

Parameter t:
The time at which to evaluate the PiecewisePolynomial.
Returns:The matrix of evaluated values.

Warning

If t does not lie in the range that the polynomial is defined over, the polynomial will silently be evaluated at the closest point to t. For example, value(-1) will return value(0) for a polynomial defined over [0, 1].

Warning

See warning in polynomial_construction_warning.

class pydrake.trajectories.PiecewiseTrajectory

Bases: pydrake.trajectories.Trajectory

Abstract class that implements the basic logic of maintaining consequent segments of time (delimited by breaks) to implement a trajectory that is represented by simpler logic in each segment or “piece”.

Template parameter T:
is the scalar type.

Instantiated templates for the following kinds of T’s are provided:

  • double
duration(self: pydrake.trajectories.PiecewiseTrajectory, segment_index: int) → float
end_time(*args, **kwargs)

Overloaded function.

  1. end_time(self: pydrake.trajectories.PiecewiseTrajectory, segment_index: int) -> float
  2. end_time(self: pydrake.trajectories.PiecewiseTrajectory) -> float
get_number_of_segments(self: pydrake.trajectories.PiecewiseTrajectory) → int
get_segment_index(self: pydrake.trajectories.PiecewiseTrajectory, t: float) → int
get_segment_times(self: pydrake.trajectories.PiecewiseTrajectory) → List[float]
is_time_in_range(self: pydrake.trajectories.PiecewiseTrajectory, t: float) → bool

Returns true iff t >= getStartTime() && t <= getEndTime().

start_time(*args, **kwargs)

Overloaded function.

  1. start_time(self: pydrake.trajectories.PiecewiseTrajectory, segment_index: int) -> float
  2. start_time(self: pydrake.trajectories.PiecewiseTrajectory) -> float
class pydrake.trajectories.Trajectory

A Trajectory represents a time-varying matrix, indexed by a single scalar double time.

Template parameter T:
is a Scalar type for the data that is returned.