Drake
Drake C++ Documentation
Hyperrectangle Class Referencefinal

Detailed Description

Axis-aligned hyperrectangle in Rᵈ defined by its lower bounds and upper bounds as {x| lb ≤ x ≤ ub}.

#include <drake/geometry/optimization/hyperrectangle.h>

Public Member Functions

 Hyperrectangle ()
 Constructs a default (zero-dimensional, nonempty) hyperrectangle. More...
 
 Hyperrectangle (const Eigen::Ref< const Eigen::VectorXd > &lb, const Eigen::Ref< const Eigen::VectorXd > &ub)
 Constructs a hyperrectangle from its lower and upper bounds. More...
 
const Eigen::VectorXd & lb () const
 Get the lower bounds of the hyperrectangle. More...
 
const Eigen::VectorXd & ub () const
 Get the upper bounds of the hyperrectangle. More...
 
Eigen::VectorXd UniformSample (RandomGenerator *generator) const
 Draws a uniform sample from the set. More...
 
Eigen::VectorXd Center () const
 Get the center of the hyperrectangle. More...
 
HPolyhedron MakeHPolyhedron () const
 Helper to convert this hyperrectangle to an HPolyhedron. More...
 
std::optional< HyperrectangleMaybeGetIntersection (const Hyperrectangle &other) const
 Constructs the intersection of two Hyperrectangle by taking the pointwise maximum of the lower bounds and the pointwise minimums of the upper bounds. More...
 
template<typename Archive >
void Serialize (Archive *a)
 Passes this object to an Archive. More...
 
Implements CopyConstructible, CopyAssignable, MoveConstructible, MoveAssignable
 Hyperrectangle (const Hyperrectangle &)=default
 
Hyperrectangleoperator= (const Hyperrectangle &)=default
 
 Hyperrectangle (Hyperrectangle &&)=default
 
Hyperrectangleoperator= (Hyperrectangle &&)=default
 
- Public Member Functions inherited from ConvexSet
virtual ~ConvexSet ()
 
std::unique_ptr< ConvexSetClone () const
 Creates a unique deep copy of this set. More...
 
int ambient_dimension () const
 Returns the dimension of the vector space in which the elements of this set are evaluated. More...
 
bool IntersectsWith (const ConvexSet &other) const
 Returns true iff the intersection between this and other is non-empty. More...
 
bool IsBounded () const
 Returns true iff the set is bounded, e.g., there exists an element-wise finite lower and upper bound for the set. More...
 
bool IsEmpty () const
 Returns true iff the set is empty. More...
 
std::optional< Eigen::VectorXd > MaybeGetPoint () const
 If this set trivially contains exactly one point, returns the value of that point. More...
 
std::optional< Eigen::VectorXd > MaybeGetFeasiblePoint () const
 Returns a feasible point within this convex set if it is nonempty, and nullopt otherwise. More...
 
bool PointInSet (const Eigen::Ref< const Eigen::VectorXd > &x, double tol=0) const
 Returns true iff the point x is contained in the set. More...
 
std::pair< VectorX< symbolic::Variable >, std::vector< solvers::Binding< solvers::Constraint > > > AddPointInSetConstraints (solvers::MathematicalProgram *prog, const Eigen::Ref< const solvers::VectorXDecisionVariable > &vars) const
 Adds a constraint to an existing MathematicalProgram enforcing that the point defined by vars is inside the set. More...
 
std::vector< solvers::Binding< solvers::Constraint > > AddPointInNonnegativeScalingConstraints (solvers::MathematicalProgram *prog, const Eigen::Ref< const solvers::VectorXDecisionVariable > &x, const symbolic::Variable &t) const
 Let S be this convex set. More...
 
std::vector< solvers::Binding< solvers::Constraint > > AddPointInNonnegativeScalingConstraints (solvers::MathematicalProgram *prog, const Eigen::Ref< const Eigen::MatrixXd > &A, const Eigen::Ref< const Eigen::VectorXd > &b, const Eigen::Ref< const Eigen::VectorXd > &c, double d, const Eigen::Ref< const solvers::VectorXDecisionVariable > &x, const Eigen::Ref< const solvers::VectorXDecisionVariable > &t) const
 Let S be this convex set. More...
 
std::pair< std::unique_ptr< Shape >, math::RigidTransformd > ToShapeWithPose () const
 Constructs a Shape and a pose of the set in the world frame for use in the SceneGraph geometry ecosystem. More...
 
double CalcVolume () const
 Computes the exact volume for the convex set. More...
 
SampledVolume CalcVolumeViaSampling (RandomGenerator *generator, const double desired_rel_accuracy=1e-2, const int max_num_samples=1e4) const
 Calculates an estimate of the volume of the convex set using sampling and performing Monte Carlo integration. More...
 
std::optional< std::pair< std::vector< double >, Eigen::MatrixXd > > Projection (const Eigen::Ref< const Eigen::MatrixXd > &points) const
 Computes in the L₂ norm the distance and the nearest point in this convex set to every column of points. More...
 
bool has_exact_volume () const
 Returns true if the exact volume can be computed for this convex set instance. More...
 

Static Public Member Functions

static std::optional< HyperrectangleMaybeCalcAxisAlignedBoundingBox (const ConvexSet &set)
 Returns the minimum axis-aligned bounding box of a convex set, for sets with finite volume. More...
 

Additional Inherited Members

- Protected Member Functions inherited from ConvexSet
 ConvexSet (int ambient_dimension, bool has_exact_volume)
 For use by derived classes to construct a ConvexSet. More...
 
template<typename Archive >
void Serialize (Archive *a)
 Implements non-virtual base class serialization. More...
 
virtual std::vector< std::optional< double > > DoProjectionShortcut (const Eigen::Ref< const Eigen::MatrixXd > &points, EigenPtr< Eigen::MatrixXd > projected_points) const
 Non-virtual interface implementation for DoProjectionShortcut(). More...
 
virtual bool DoPointInSet (const Eigen::Ref< const Eigen::VectorXd > &x, double tol) const
 Non-virtual interface implementation for PointInSet(). More...
 
std::optional< symbolic::VariableHandleZeroAmbientDimensionConstraints (solvers::MathematicalProgram *prog, const ConvexSet &set, std::vector< solvers::Binding< solvers::Constraint >> *constraints) const
 Instances of subclasses such as CartesianProduct and MinkowskiSum can have constituent sets with zero ambient dimension, which much be handled in a special manner when calling methods such as DoAddPointInSetConstraints. More...
 
 ConvexSet (const ConvexSet &)=default
 
ConvexSetoperator= (const ConvexSet &)=default
 
 ConvexSet (ConvexSet &&)=default
 
ConvexSetoperator= (ConvexSet &&)=default
 

Constructor & Destructor Documentation

◆ Hyperrectangle() [1/4]

Hyperrectangle ( const Hyperrectangle )
default

◆ Hyperrectangle() [2/4]

Hyperrectangle ( Hyperrectangle &&  )
default

◆ Hyperrectangle() [3/4]

Constructs a default (zero-dimensional, nonempty) hyperrectangle.

◆ Hyperrectangle() [4/4]

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

Constructs a hyperrectangle from its lower and upper bounds.

Precondition
lb.size() == ub.size()
lb and ub are finite.
lb(i) <= ub(i) for all i

Member Function Documentation

◆ Center()

Eigen::VectorXd Center ( ) const

Get the center of the hyperrectangle.

◆ lb()

const Eigen::VectorXd& lb ( ) const

Get the lower bounds of the hyperrectangle.

◆ MakeHPolyhedron()

HPolyhedron MakeHPolyhedron ( ) const

Helper to convert this hyperrectangle to an HPolyhedron.

◆ MaybeCalcAxisAlignedBoundingBox()

static std::optional<Hyperrectangle> MaybeCalcAxisAlignedBoundingBox ( const ConvexSet set)
static

Returns the minimum axis-aligned bounding box of a convex set, for sets with finite volume.

(std::nullopt otherwise).

◆ MaybeGetIntersection()

std::optional<Hyperrectangle> MaybeGetIntersection ( const Hyperrectangle other) const

Constructs the intersection of two Hyperrectangle by taking the pointwise maximum of the lower bounds and the pointwise minimums of the upper bounds.

Returns std::nullopt if the intersection is empty.

Precondition
this and other need to have the same ambient dimension.

◆ operator=() [1/2]

Hyperrectangle& operator= ( const Hyperrectangle )
default

◆ operator=() [2/2]

Hyperrectangle& operator= ( Hyperrectangle &&  )
default

◆ Serialize()

void Serialize ( Archive *  a)

Passes this object to an Archive.

Refer to YAML Serialization for background.

◆ ub()

const Eigen::VectorXd& ub ( ) const

Get the upper bounds of the hyperrectangle.

◆ UniformSample()

Eigen::VectorXd UniformSample ( RandomGenerator generator) const

Draws a uniform sample from the set.


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