Drake
Drake C++ Documentation

Detailed Description

Implements an ellipsoidal convex set represented as an affine scaling of the unit ball {Bu + center | |u|₂ ≤ 1}.

B must be a square matrix.

Compare this with an alternative parametrization of the ellipsoid: {x | (x-center)ᵀAᵀA(x-center) ≤ 1}, which utilizes a quadratic form. The two representations are related by B = A⁻¹ if A and B are invertible.

The quadratic form parametrization is implemented in Hyperellipsoid. It can represent unbounded sets, but not sets along a lower-dimensional affine subspace. The AffineBall parametrization can represent sets along a lower-dimensional affine subspace, but not unbounded sets.

An AffineBall can never be empty – it always contains its center. This includes the zero-dimensional case.

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

Public Member Functions

 AffineBall ()
 Constructs a default (zero-dimensional, nonempty) set. More...
 
 AffineBall (const Eigen::Ref< const Eigen::MatrixXd > &B, const Eigen::Ref< const Eigen::VectorXd > &center)
 Constructs the ellipsoid from a transformation matrix B and translation center. More...
 
 AffineBall (const Hyperellipsoid &ellipsoid)
 Constructs an AffineBall from a Hyperellipsoid. More...
 
 ~AffineBall () final
 
const Eigen::MatrixXd & B () const
 Returns the affine transformation matrix B. More...
 
const Eigen::VectorXd & center () const
 Returns the center of the ellipsoid. More...
 
template<typename Archive >
void Serialize (Archive *a)
 Passes this object to an Archive. More...
 
bool IsBounded (Parallelism parallelism=Parallelism::None()) const
 Every AffineBall is bounded by construction. More...
 
Implements CopyConstructible, CopyAssignable, MoveConstructible, MoveAssignable
 AffineBall (const AffineBall &)=default
 
AffineBalloperator= (const AffineBall &)=default
 
 AffineBall (AffineBall &&)=default
 
AffineBalloperator= (AffineBall &&)=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 (Parallelism parallelism=Parallelism::None()) 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 AffineBall MakeAxisAligned (const Eigen::Ref< const Eigen::VectorXd > &radius, const Eigen::Ref< const Eigen::VectorXd > &center)
 Constructs an axis-aligned AffineBall with the implicit form (x₀-c₀)²/r₀² + (x₁-c₁)²/r₁² + ... More...
 
static AffineBall MakeHypersphere (double radius, const Eigen::Ref< const Eigen::VectorXd > &center)
 Constructs a hypersphere with radius and center. More...
 
static AffineBall MakeUnitBall (int dim)
 Constructs the L₂-norm unit ball in dim dimensions, {x | |x|₂ <= 1 }. More...
 
static AffineBall MinimumVolumeCircumscribedEllipsoid (const Eigen::Ref< const Eigen::MatrixXd > &points, double rank_tol=1e-6)
 Constructs the minimum-volume ellipsoid which contains all of the given points. More...
 
static AffineBall MakeAffineBallFromLineSegment (const Eigen::Ref< const Eigen::VectorXd > &x_1, const Eigen::Ref< const Eigen::VectorXd > &x_2, const double epsilon=1e-3)
 Constructs an affine ball such that its main diameter is the line segment from x_1 to x_2, and the length of all other diameters are 2 * epsilon. 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::optional< bool > DoIsBoundedShortcutParallel (Parallelism) const
 Non-virtual interface implementation for DoIsBoundedShortcutParallel(). 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
 
- Static Protected Member Functions inherited from ConvexSet
static std::unique_ptr< ConvexSetAffineHullShortcut (const ConvexSet &self, std::optional< double > tol)
 When there is a more efficient strategy to compute the affine hull of this set, returns affine hull as an AffineSubspace. More...
 

Constructor & Destructor Documentation

◆ AffineBall() [1/5]

AffineBall ( const AffineBall )
default

◆ AffineBall() [2/5]

AffineBall ( AffineBall &&  )
default

◆ AffineBall() [3/5]

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

◆ AffineBall() [4/5]

AffineBall ( const Eigen::Ref< const Eigen::MatrixXd > &  B,
const Eigen::Ref< const Eigen::VectorXd > &  center 
)

Constructs the ellipsoid from a transformation matrix B and translation center.

B describes the linear transformation that is applied to the unit ball in order to produce the ellipsoid, and center describes the translation of the center of the ellipsoid from the origin.

Precondition
B.rows() == B.cols().
B.cols() == center.size().

◆ AffineBall() [5/5]

AffineBall ( const Hyperellipsoid ellipsoid)
explicit

Constructs an AffineBall from a Hyperellipsoid.

Precondition
ellipsoid.IsBounded().

◆ ~AffineBall()

~AffineBall ( )
final

Member Function Documentation

◆ B()

const Eigen::MatrixXd& B ( ) const

Returns the affine transformation matrix B.

◆ center()

const Eigen::VectorXd& center ( ) const

Returns the center of the ellipsoid.

◆ IsBounded()

bool IsBounded

Every AffineBall is bounded by construction.

Parameters
parallelismIgnored – no parallelization is used.
Note
See parent class's documentation for more details.

◆ MakeAffineBallFromLineSegment()

static AffineBall MakeAffineBallFromLineSegment ( const Eigen::Ref< const Eigen::VectorXd > &  x_1,
const Eigen::Ref< const Eigen::VectorXd > &  x_2,
const double  epsilon = 1e-3 
)
static

Constructs an affine ball such that its main diameter is the line segment from x_1 to x_2, and the length of all other diameters are 2 * epsilon.

Precondition
x_1.size() == x_2.size().
epsilon >= 0.
Exceptions
std::runtime_errorif ‖x_1 - x_2‖₂ is less than 1e-9.

◆ MakeAxisAligned()

static AffineBall MakeAxisAligned ( const Eigen::Ref< const Eigen::VectorXd > &  radius,
const Eigen::Ref< const Eigen::VectorXd > &  center 
)
static

Constructs an axis-aligned AffineBall with the implicit form (x₀-c₀)²/r₀² + (x₁-c₁)²/r₁² + ...

  • (x_N - c_N)²/r_N² ≤ 1, where c is shorthand for center and r is shorthand for radius.
    Precondition
    radius.size() == center.size().
    radius[i] >= 0, for all i.

◆ MakeHypersphere()

static AffineBall MakeHypersphere ( double  radius,
const Eigen::Ref< const Eigen::VectorXd > &  center 
)
static

Constructs a hypersphere with radius and center.

Precondition
radius >= 0.

◆ MakeUnitBall()

static AffineBall MakeUnitBall ( int  dim)
static

Constructs the L₂-norm unit ball in dim dimensions, {x | |x|₂ <= 1 }.

Precondition
dim >= 0.

◆ MinimumVolumeCircumscribedEllipsoid()

static AffineBall MinimumVolumeCircumscribedEllipsoid ( const Eigen::Ref< const Eigen::MatrixXd > &  points,
double  rank_tol = 1e-6 
)
static

Constructs the minimum-volume ellipsoid which contains all of the given points.

This is commonly referred to as the outer Löwner-John ellipsoid.

If all of the points lie along a proper affine subspace, this method instead computes the minimum-k-volume ellipsoid, where k is the affine dimension of the convex hull of points.

Parameters
pointsis a d-by-n matrix, where d is the ambient dimension and each column represents one point.
rank_tolthe tolerance used to detect if the data lies in an affine subspace. The affine ball is computed in the subspace spanned by the left singular vectors of the data matrix whose associated singular values are larger than rank_tol * max_singular_value. The default is 1e-6 to be compatible with common solver tolerances.
Exceptions
std::exceptionif the MathematicalProgram fails to solve. This can happen due to numerical issues caused by rank_tol being set too low.
std::exceptionif points includes NaNs or infinite values.
Precondition
points.rows() >= 1.
points.cols() >= 1.

◆ operator=() [1/2]

AffineBall& operator= ( AffineBall &&  )
default

◆ operator=() [2/2]

AffineBall& operator= ( const AffineBall )
default

◆ Serialize()

void Serialize ( Archive *  a)

Passes this object to an Archive.

Refer to YAML Serialization for background.


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