Drake
Drake C++ Documentation
GcsTrajectoryOptimization Class Referencefinal

Detailed Description

GcsTrajectoryOptimization implements a simplified motion planning optimization problem introduced in the paper "Motion Planning around Obstacles with Convex Optimization" by Tobia Marcucci, Mark Petersen, David von Wrangel, Russ Tedrake.

Instead of using the full time-scaling curve, this problem uses a single time-scaling variable for each region. This formulation yields continuous trajectories, which are not differentiable at the transition times between the regions since non-convex continuity constraints are not supported yet. However, it supports continuity on the path r(s) for arbitrary degree. The path r(s) can be reconstructed from the gcs solution q(t) with NormalizeSegmentTimes() and post-processed with e.g. Toppra to enforce acceleration bounds.

The ith piece of the composite trajectory is defined as q(t) = r((t - tᵢ) / hᵢ). r : [0, 1] → ℝⁿ is a the path, parametrized as a Bézier curve with order n. tᵢ and hᵢ are the initial time and duration of the ith sub-trajectory.

This class supports the notion of a Subgraph of regions. This has proven useful to facilitate multi-modal motion planning such as: Subgraph A: find a collision-free trajectory for the robot to a grasping posture, Subgraph B: find a collision-free trajectory for the robot with the object in its hand to a placing posture, etc.

Continuous Revolute Joints

This class also supports robots with continuous revolute joints (revolute joints that don't have any joint limits) and mobile bases. Adding or subtracting 2π to such a joint's angle leaves it unchanged; this logic is implemented behind the scenes. To use it, one should specify the joint indices that don't have limits, and ensure all sets satisfy the "convexity radius" property – their width along a dimension corresponding to a continuous revolute joint must be less than π. This can be enforced when constructing the convex sets, or after the fact with geometry::optimization::PartitionConvexSet. The GcsTrajectoryOptimization methods AddRegions and AddEdges will handle all of the intersection checks behind the scenes, including applying the appropriate logic to connect sets that "wrap around" 2π.

#include <drake/planning/trajectory_optimization/gcs_trajectory_optimization.h>

Classes

class  EdgesBetweenSubgraphs
 EdgesBetweenSubgraphs are defined as the connecting edges between two given subgraphs. More...
 
class  Subgraph
 A Subgraph is a subset of the larger graph. More...
 

Public Member Functions

 GcsTrajectoryOptimization (int num_positions, std::vector< int > continuous_revolute_joints=std::vector< int >())
 Constructs the motion planning problem. More...
 
 ~GcsTrajectoryOptimization ()
 
int num_positions () const
 Returns the number of position variables. More...
 
const std::vector< int > & continuous_revolute_joints ()
 Returns a list of indices corresponding to continuous revolute joints. More...
 
std::string GetGraphvizString (const std::optional< solvers::MathematicalProgramResult > &result=std::nullopt, bool show_slack=true, int precision=3, bool scientific=false) const
 Returns a Graphviz string describing the graph vertices and edges. More...
 
SubgraphAddRegions (const geometry::optimization::ConvexSets &regions, const std::vector< std::pair< int, int >> &edges_between_regions, int order, double h_min=0, double h_max=20, std::string name="", std::optional< const std::vector< Eigen::VectorXd >> edge_offsets=std::nullopt)
 Creates a Subgraph with the given regions and indices. More...
 
SubgraphAddRegions (const geometry::optimization::ConvexSets &regions, int order, double h_min=0, double h_max=20, std::string name="")
 Creates a Subgraph with the given regions. More...
 
void RemoveSubgraph (const Subgraph &subgraph)
 Remove a subgraph and all associated edges found in the subgraph and to and from other subgraphs. More...
 
EdgesBetweenSubgraphsAddEdges (const Subgraph &from_subgraph, const Subgraph &to_subgraph, const geometry::optimization::ConvexSet *subspace=nullptr)
 Connects two subgraphs with directed edges. More...
 
void AddTimeCost (double weight=1.0)
 Adds a minimum time cost to all regions in the whole graph. More...
 
void AddPathLengthCost (const Eigen::MatrixXd &weight_matrix)
 Adds multiple L2Norm Costs on the upper bound of the path length. More...
 
void AddPathLengthCost (double weight=1.0)
 Adds multiple L2Norm Costs on the upper bound of the path length. More...
 
void AddVelocityBounds (const Eigen::Ref< const Eigen::VectorXd > &lb, const Eigen::Ref< const Eigen::VectorXd > &ub)
 Adds a linear velocity constraint to the entire graph lb ≤ q̇(t) ≤ ub. More...
 
void AddPathContinuityConstraints (int continuity_order)
 Enforces derivative continuity constraints on the entire graph. More...
 
std::pair< trajectories::CompositeTrajectory< double >, solvers::MathematicalProgramResultSolvePath (const Subgraph &source, const Subgraph &target, const geometry::optimization::GraphOfConvexSetsOptions &options={})
 Formulates and solves the mixed-integer convex formulation of the shortest path problem on the whole graph. More...
 
std::pair< trajectories::CompositeTrajectory< double >, solvers::MathematicalProgramResultSolveConvexRestriction (const std::vector< const geometry::optimization::GraphOfConvexSets::Vertex * > &active_vertices, const geometry::optimization::GraphOfConvexSetsOptions &options={})
 Solves a trajectory optimization problem through specific vertices. More...
 
double EstimateComplexity () const
 Provide a heuristic estimate of the complexity of the underlying GCS mathematical program, for regression testing purposes. More...
 
std::vector< Subgraph * > GetSubgraphs () const
 Returns a vector of all subgraphs. More...
 
std::vector< EdgesBetweenSubgraphs * > GetEdgesBetweenSubgraphs () const
 Returns a vector of all edges between subgraphs. More...
 
const geometry::optimization::GraphOfConvexSetsgraph_of_convex_sets () const
 Getter for the underlying GraphOfConvexSets. More...
 
Does not allow copy, move, or assignment
 GcsTrajectoryOptimization (const GcsTrajectoryOptimization &)=delete
 
GcsTrajectoryOptimizationoperator= (const GcsTrajectoryOptimization &)=delete
 
 GcsTrajectoryOptimization (GcsTrajectoryOptimization &&)=delete
 
GcsTrajectoryOptimizationoperator= (GcsTrajectoryOptimization &&)=delete
 

Static Public Member Functions

static trajectories::CompositeTrajectory< doubleNormalizeSegmentTimes (const trajectories::CompositeTrajectory< double > &trajectory)
 Normalizes each trajectory segment to one second in duration. More...
 
static trajectories::CompositeTrajectory< doubleUnwrapToContinousTrajectory (const trajectories::CompositeTrajectory< double > &gcs_trajectory, std::vector< int > continuous_revolute_joints, std::optional< std::vector< int >> starting_rounds=std::nullopt)
 Unwraps a trajectory with continuous revolute joints into a continuous trajectory in the Euclidean space. More...
 

Constructor & Destructor Documentation

◆ GcsTrajectoryOptimization() [1/3]

◆ GcsTrajectoryOptimization() [2/3]

◆ GcsTrajectoryOptimization() [3/3]

GcsTrajectoryOptimization ( int  num_positions,
std::vector< int continuous_revolute_joints = std::vector< int >() 
)
explicit

Constructs the motion planning problem.

Parameters
num_positionsis the dimension of the configuration space.
continuous_revolute_jointsis a list of indices corresponding to continuous revolute joints, i.e., revolute joints which don't have any joint limits, and hence "wrap around" at 2π. Each entry in continuous_revolute_joints must be non-negative, less than num_positions, and unique. This feature is currently only supported within a single subgraph: continuous revolute joints won't be taken into account when constructing edges between subgraphs or checking if sets intersect through a subspace.

◆ ~GcsTrajectoryOptimization()

Member Function Documentation

◆ AddEdges()

EdgesBetweenSubgraphs& AddEdges ( const Subgraph from_subgraph,
const Subgraph to_subgraph,
const geometry::optimization::ConvexSet subspace = nullptr 
)

Connects two subgraphs with directed edges.

Parameters
from_subgraphis the subgraph to connect from. Must have been created from a call to AddRegions() on this object, not some other optimization program.
to_subgraphis the subgraph to connect to. Must have been created from a call to AddRegions() on this object, not some other optimization program.
subspaceis the subspace that the connecting control points must be in. Subspace is optional. Only edges that connect through the subspace will be added, and the subspace is added as a constraint on the connecting control points. Subspaces of type point or HPolyhedron are supported since other sets require constraints that are not yet supported by the GraphOfConvexSets::Edge constraint, e.g., set containment of a HyperEllipsoid is formulated via LorentzCone constraints. Workaround: Create a subgraph of zero order with the subspace as the region and connect it between the two subgraphs. This works because GraphOfConvexSet::Vertex , supports arbitrary instances of ConvexSets.

◆ AddPathContinuityConstraints()

void AddPathContinuityConstraints ( int  continuity_order)

Enforces derivative continuity constraints on the entire graph.

Parameters
continuity_orderis the order of the continuity constraint.

Note that the constraints are on the control points of the derivatives of r(s) and not q(t). This may result in discontinuities of the trajectory return by SolvePath() since the r(s) will get rescaled by the duration h to yield q(t). NormalizeSegmentTimes() will return r(s) with valid continuity.

Exceptions
std::exceptionif the continuity order is less than one since path continuity is enforced by default.

◆ AddPathLengthCost() [1/2]

void AddPathLengthCost ( const Eigen::MatrixXd &  weight_matrix)

Adds multiple L2Norm Costs on the upper bound of the path length.

Since we cannot directly compute the path length of a Bézier curve, we minimize the upper bound of the path integral by minimizing the sum of (weighted) distances between control points: ∑ |weight_matrix * (rᵢ₊₁ − rᵢ)|₂.

This cost will be added to the entire graph. Since the path length is only defined for Bézier curves that have two or more control points, this cost will only added to all subgraphs with order greater than zero. Note that this cost will be applied even to subgraphs added in the future.

Parameters
weight_matrixis the relative weight of each component for the cost. The diagonal of the matrix is the weight for each dimension. The off-diagonal elements are the weight for the cross terms, which can be used to penalize diagonal movement.
Precondition
weight_matrix must be of size num_positions() x num_positions().

◆ AddPathLengthCost() [2/2]

void AddPathLengthCost ( double  weight = 1.0)

Adds multiple L2Norm Costs on the upper bound of the path length.

Since we cannot directly compute the path length of a Bézier curve, we minimize the upper bound of the path integral by minimizing the sum of distances between control points. For Bézier curves, this is equivalent to the sum of the L2Norm of the derivative control points of the curve divided by the order.

This cost will be added to the entire graph. Since the path length is only defined for Bézier curves that have two or more control points, this cost will only added to all subgraphs with order greater than zero. Note that this cost will be applied even to subgraphs added in the future.

Parameters
weightis the relative weight of the cost.

◆ AddRegions() [1/2]

Subgraph& AddRegions ( const geometry::optimization::ConvexSets regions,
const std::vector< std::pair< int, int >> &  edges_between_regions,
int  order,
double  h_min = 0,
double  h_max = 20,
std::string  name = "",
std::optional< const std::vector< Eigen::VectorXd >>  edge_offsets = std::nullopt 
)

Creates a Subgraph with the given regions and indices.

Parameters
regionsrepresent the valid set a control point can be in. We retain a copy of the regions since other functions may access them. If any of the positions represent revolute joints without limits, each region has a maximum width of strictly less than π along dimensions corresponding to those joints.
edges_between_regionsis a list of pairs of indices into the regions vector. For each pair representing an edge between two regions, an edge is added within the subgraph. Note that the edges are directed so (i,j) will only add an edge from region i to region j.
orderis the order of the Bézier curve.
h_maxis the maximum duration to spend in a region (seconds). Some solvers struggle numerically with large values.
h_minis the minimum duration to spend in a region (seconds) if that region is visited on the optimal path. Some cost and constraints are only convex for h > 0. For example the perspective quadratic cost of the path energy ||ṙ(s)||² / h becomes non-convex for h = 0. Otherwise h_min can be set to 0.
nameis the name of the subgraph. If the passed name is an empty string, a default name will be provided.
edge_offsetsis an optional list of vectors. If defined, the list must contain the same number of entries as edges_between_regions. In other words, if defined, there must be one edge offset for each specified edge. For each pair of sets listed in edges_between_regions, the first set is translated (in configuration space) by the corresponding vector in edge_offsets before computing the constraints associated to that edge. This is used to add edges between sets that "wrap around" 2π along some dimension, due to, e.g., a continuous revolute joint. This edge offset corresponds to the translation component of the affine map τ_uv in equation (11) of "Non-Euclidean Motion Planning with Graphs of Geodesically-Convex Sets", and per the discussion in Subsection VI A, τ_uv has no rotation component.

◆ AddRegions() [2/2]

Subgraph& AddRegions ( const geometry::optimization::ConvexSets regions,
int  order,
double  h_min = 0,
double  h_max = 20,
std::string  name = "" 
)

Creates a Subgraph with the given regions.

This function will compute the edges between the regions based on the set intersections.

Parameters
regionsrepresent the valid set a control point can be in. We retain a copy of the regions since other functions may access them. If any of the positions represent continuous revolute joints, each region must have a maximum width of strictly less than π along dimensions corresponding to those joints.
orderis the order of the Bézier curve.
h_minis the minimum duration to spend in a region (seconds) if that region is visited on the optimal path. Some cost and constraints are only convex for h > 0. For example the perspective quadratic cost of the path energy ||ṙ(s)||² / h becomes non-convex for h = 0. Otherwise h_min can be set to 0.
h_maxis the maximum duration to spend in a region (seconds). Some solvers struggle numerically with large values.
nameis the name of the subgraph. A default name will be provided.
Exceptions
std::exceptionif any of the regions has a width of π or greater along dimensions corresponding to continuous revolute joints.

◆ AddTimeCost()

void AddTimeCost ( double  weight = 1.0)

Adds a minimum time cost to all regions in the whole graph.

The cost is the sum of the time scaling variables.

This cost will be added to the entire graph. Note that this cost will be applied even to subgraphs added in the future.

Parameters
weightis the relative weight of the cost.

◆ AddVelocityBounds()

void AddVelocityBounds ( const Eigen::Ref< const Eigen::VectorXd > &  lb,
const Eigen::Ref< const Eigen::VectorXd > &  ub 
)

Adds a linear velocity constraint to the entire graph lb ≤ q̇(t) ≤ ub.

Parameters
lbis the lower bound of the velocity.
ubis the upper bound of the velocity.

This constraint will be added to the entire graph. Since the velocity requires forming the derivative of the Bézier curve, this constraint will only added to all subgraphs with order greater than zero. Note that this constraint will be applied even to subgraphs added in the future.

Exceptions
std::exceptionif lb or ub are not of size num_positions().

◆ continuous_revolute_joints()

const std::vector<int>& continuous_revolute_joints ( )

Returns a list of indices corresponding to continuous revolute joints.

◆ EstimateComplexity()

double EstimateComplexity ( ) const

Provide a heuristic estimate of the complexity of the underlying GCS mathematical program, for regression testing purposes.

Here we sum the total number of variable appearances in our costs and constraints as a rough approximation of the complexity of the subproblems.

◆ GetEdgesBetweenSubgraphs()

std::vector<EdgesBetweenSubgraphs*> GetEdgesBetweenSubgraphs ( ) const

Returns a vector of all edges between subgraphs.

◆ GetGraphvizString()

std::string GetGraphvizString ( const std::optional< solvers::MathematicalProgramResult > &  result = std::nullopt,
bool  show_slack = true,
int  precision = 3,
bool  scientific = false 
) const

Returns a Graphviz string describing the graph vertices and edges.

If results is supplied, then the graph will be annotated with the solution values.

Parameters
show_slacksdetermines whether the values of the intermediate (slack) variables are also displayed in the graph.
precisionsets the floating point precision (how many digits are generated) of the annotations.
scientificsets the floating point formatting to scientific (if true) or fixed (if false).

◆ GetSubgraphs()

std::vector<Subgraph*> GetSubgraphs ( ) const

Returns a vector of all subgraphs.

◆ graph_of_convex_sets()

const geometry::optimization::GraphOfConvexSets& graph_of_convex_sets ( ) const

Getter for the underlying GraphOfConvexSets.

This is intended primarily for inspecting the resulting programs.

◆ NormalizeSegmentTimes()

static trajectories::CompositeTrajectory<double> NormalizeSegmentTimes ( const trajectories::CompositeTrajectory< double > &  trajectory)
static

Normalizes each trajectory segment to one second in duration.

Reconstructs the path r(s) from the solution trajectory q(t) of SolvePath() s.t. each segment of the resulting trajectory will be one second long. The start time will match the original start time.

Parameters
trajectoryThe solution trajectory returned by SolvePath().
Exceptions
std::exceptionif not all trajectory segments of the CompositeTrajectory are of type BezierCurve<double>

◆ num_positions()

int num_positions ( ) const

Returns the number of position variables.

◆ operator=() [1/2]

◆ operator=() [2/2]

GcsTrajectoryOptimization& operator= ( const GcsTrajectoryOptimization )
delete

◆ RemoveSubgraph()

void RemoveSubgraph ( const Subgraph subgraph)

Remove a subgraph and all associated edges found in the subgraph and to and from other subgraphs.

Precondition
The subgraph must have been created from a call to AddRegions() on this object.

◆ SolveConvexRestriction()

std::pair<trajectories::CompositeTrajectory<double>, solvers::MathematicalProgramResult> SolveConvexRestriction ( const std::vector< const geometry::optimization::GraphOfConvexSets::Vertex * > &  active_vertices,
const geometry::optimization::GraphOfConvexSetsOptions options = {} 
)

Solves a trajectory optimization problem through specific vertices.

This method allows for targeted optimization by considering only selected active vertices, reducing the problem's complexity. See geometry::optimization::GraphOfConvexSets::SolveConvexRestriction(). This API prefers a sequence of vertices over edges, as a user may know which regions the solution should pass through. GcsTrajectoryOptimization::AddRegions() automatically manages edge creation and intersection checks, which makes passing a sequence of edges less convenient.

Parameters
active_verticesA sequence of ordered vertices of subgraphs to be included in the problem.
optionsinclude all settings for solving the shortest path problem.
Precondition
There must be at least two vertices in active_vertices.
Exceptions
std::exceptionif the vertices are not connected.
std::exceptionif two vertices are connected by multiple edges. This may happen if one connects two graphs through multiple subspaces, which is currently not supported with this method.
std::exceptionif the program cannot be written as a convex optimization consumable by one of the standard solvers.

◆ SolvePath()

std::pair<trajectories::CompositeTrajectory<double>, solvers::MathematicalProgramResult> SolvePath ( const Subgraph source,
const Subgraph target,
const geometry::optimization::GraphOfConvexSetsOptions options = {} 
)

Formulates and solves the mixed-integer convex formulation of the shortest path problem on the whole graph.

See also
geometry::optimization::GraphOfConvexSets::SolveShortestPath() for further details.
Parameters
sourcespecifies the source subgraph. Must have been created from a call to AddRegions() on this object, not some other optimization program. If the source is a subgraph with more than one region, an empty set will be added and optimizer will choose the best region to start in. To start in a particular point, consider adding a subgraph of order zero with a single region of type Point.
targetspecifies the target subgraph. Must have been created from a call to AddRegions() on this object, not some other optimization program. If the target is a subgraph with more than one region, an empty set will be added and optimizer will choose the best region to end in. To end in a particular point, consider adding a subgraph of order zero with a single region of type Point.
optionsinclude all settings for solving the shortest path problem. The following default options will be used if they are not provided in options:
  • options.convex_relaxation = true,
  • options.max_rounded_paths = 5,
  • options.preprocessing = true.
See also
geometry::optimization::GraphOfConvexSetsOptions for further details.

◆ UnwrapToContinousTrajectory()

static trajectories::CompositeTrajectory<double> UnwrapToContinousTrajectory ( const trajectories::CompositeTrajectory< double > &  gcs_trajectory,
std::vector< int continuous_revolute_joints,
std::optional< std::vector< int >>  starting_rounds = std::nullopt 
)
static

Unwraps a trajectory with continuous revolute joints into a continuous trajectory in the Euclidean space.

Trajectories produced by GcsTrajectoryOptimization for robotic systems with continuous revolute joints may include apparent discontinuities, where a multiple of 2π is instantaneously added to a joint value at the boundary between two adjacent segments of the trajectory. This function removes such discontinuities by adding or subtracting the appropriate multiple of 2π, "unwrapping" the trajectory into a continuous representation suitable for downstream tasks that do not take the joint wraparound into account.

Parameters
gcs_trajectoryThe trajectory to unwrap.
continuous_revolute_jointsThe indices of the continuous revolute joints.
starting_roundsA vector of integers that sets the starting rounds for each continuous revolute joint. Given integer k for the starting_round of a joint, its initial position will be wrapped into [2πk , 2π(k+1)). If the starting rounds are not provided, the initial position of trajectory will be unchanged.
Returns
an unwrapped (continous in Euclidean space) CompositeTrajectory.
Exceptions
std::exceptionif starting_rounds.size()!=continuous_revolute_joints.size().
std::exceptionif continuous_revolute_joints contain repeated indices and/or indices outside the range [0, gcs_trajectory.rows()).
std::exceptionif the gcs_trajectory is not continuous on the manifold defined by the continuous_revolute_joints, i.e., the shift between two consecutive segments is not an integer multiple of 2π (within a tolerance of 1e-10 radians).
std::exceptionif all the segments are not of type BezierCurve. Other types are not supported yet. Note that currently the output of GcsTrajectoryOptimization::SolvePath() is a CompositeTrajectory of BezierCurves.

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