Implements an ellipsoidal convex set represented by the quadratic form {x | (x-center)ᵀAᵀA(x-center) ≤ 1}
.
Note that A
need not be square; we require only that the matrix AᵀA is positive semi-definite.
Compare this with an alternative (very useful) parameterization of the ellipsoid: {Bu + center | |u|₂ ≤ 1}
, which is an affine scaling of the unit ball. This is related to the quadratic form by B = A⁻¹
, when A
is invertible, but the quadratic form can also represent unbounded sets. The affine scaling of the unit ball representation is available via the AffineBall class.
Note: the name Hyperellipsoid was taken here to avoid conflicting with geometry::Ellipsoid and to distinguish that this class supports N dimensions.
A hyperellipsoid can never be empty – it always contains its center. This includes the zero-dimensional case.
#include <drake/geometry/optimization/hyperellipsoid.h>
Public Member Functions | |
Hyperellipsoid () | |
Constructs a default (zero-dimensional, nonempty) set. More... | |
Hyperellipsoid (const Eigen::Ref< const Eigen::MatrixXd > &A, const Eigen::Ref< const Eigen::VectorXd > ¢er) | |
Constructs the ellipsoid. More... | |
Hyperellipsoid (const QueryObject< double > &query_object, GeometryId geometry_id, std::optional< FrameId > reference_frame=std::nullopt) | |
Constructs a Hyperellipsoid from a SceneGraph geometry and pose in the reference_frame frame, obtained via the QueryObject. More... | |
Hyperellipsoid (const AffineBall &ellipsoid) | |
Constructs a Hyperellipsoid from an AffineBall. More... | |
~Hyperellipsoid () final | |
const Eigen::MatrixXd & | A () const |
Returns the quadratic form matrix A. More... | |
const Eigen::VectorXd & | center () const |
Returns the center of the ellipsoid. More... | |
std::pair< double, Eigen::VectorXd > | MinimumUniformScalingToTouch (const ConvexSet &other) const |
Computes the smallest uniform scaling of this ellipsoid for which it still intersects other . More... | |
Hyperellipsoid | Scale (double scale) const |
Results a new Hyperellipsoid that is a scaled version of this about the center. More... | |
template<typename Archive > | |
void | Serialize (Archive *a) |
Passes this object to an Archive. More... | |
double | Volume () const |
Computes the volume for the hyperellipsoid set. More... | |
std::pair< std::unique_ptr< Shape >, math::RigidTransformd > | ToShapeWithPose () const |
Constructs a Ellipsoid shape description of this set. More... | |
bool | IsBounded (Parallelism parallelism=Parallelism::None()) const |
A Hyperellipsoid is bounded if the kernel of A is rank-zero. More... | |
Implements CopyConstructible, CopyAssignable, MoveConstructible, MoveAssignable | |
Hyperellipsoid (const Hyperellipsoid &)=default | |
Hyperellipsoid & | operator= (const Hyperellipsoid &)=default |
Hyperellipsoid (Hyperellipsoid &&)=default | |
Hyperellipsoid & | operator= (Hyperellipsoid &&)=default |
Public Member Functions inherited from ConvexSet | |
virtual | ~ConvexSet () |
std::unique_ptr< ConvexSet > | Clone () 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 Hyperellipsoid | MakeAxisAligned (const Eigen::Ref< const Eigen::VectorXd > &radius, const Eigen::Ref< const Eigen::VectorXd > ¢er) |
Constructs the an axis-aligned Hyperellipsoid with the implicit form (x₀-c₀)²/r₀² + (x₁-c₁)²/r₁² + ... More... | |
static Hyperellipsoid | MakeHypersphere (double radius, const Eigen::Ref< const Eigen::VectorXd > ¢er) |
Constructs a hypersphere with radius and center . More... | |
static Hyperellipsoid | MakeUnitBall (int dim) |
Constructs the L₂-norm unit ball in dim dimensions, {x | |x|₂ <= 1 }. More... | |
static Hyperellipsoid | MinimumVolumeCircumscribedEllipsoid (const Eigen::Ref< const Eigen::MatrixXd > &points, double rank_tol=1e-6) |
Constructs the minimum-volume ellipsoid which contains all of the points . 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 std::optional< Eigen::VectorXd > | DoMaybeGetPoint () const |
Non-virtual interface implementation for MaybeGetPoint(). More... | |
virtual bool | DoPointInSet (const Eigen::Ref< const Eigen::VectorXd > &x, double tol) const |
Non-virtual interface implementation for PointInSet(). More... | |
std::optional< symbolic::Variable > | HandleZeroAmbientDimensionConstraints (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 | |
ConvexSet & | operator= (const ConvexSet &)=default |
ConvexSet (ConvexSet &&)=default | |
ConvexSet & | operator= (ConvexSet &&)=default |
Static Protected Member Functions inherited from ConvexSet | |
static std::unique_ptr< ConvexSet > | AffineHullShortcut (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... | |
|
default |
|
default |
Hyperellipsoid | ( | ) |
Constructs a default (zero-dimensional, nonempty) set.
Hyperellipsoid | ( | const Eigen::Ref< const Eigen::MatrixXd > & | A, |
const Eigen::Ref< const Eigen::VectorXd > & | center | ||
) |
Constructs the ellipsoid.
Hyperellipsoid | ( | const QueryObject< double > & | query_object, |
GeometryId | geometry_id, | ||
std::optional< FrameId > | reference_frame = std::nullopt |
||
) |
Constructs a Hyperellipsoid from a SceneGraph geometry and pose in the reference_frame
frame, obtained via the QueryObject.
If reference_frame
frame is std::nullopt, then it will be expressed in the world frame.
std::exception | if geometry_id does not represent a shape that can be described as an Hyperellipsoid. |
|
explicit |
Constructs a Hyperellipsoid from an AffineBall.
|
final |
const Eigen::MatrixXd& A | ( | ) | const |
Returns the quadratic form matrix A.
const Eigen::VectorXd& center | ( | ) | const |
Returns the center of the ellipsoid.
bool IsBounded |
A Hyperellipsoid is bounded if the kernel of A is rank-zero.
parallelism | Ignored – no parallelization is used. |
|
static |
Constructs the an axis-aligned Hyperellipsoid with the implicit form (x₀-c₀)²/r₀² + (x₁-c₁)²/r₁² + ...
center
and r is shorthand for radius
.
|
static |
Constructs a hypersphere with radius
and center
.
|
static |
Constructs the L₂-norm unit ball in dim
dimensions, {x | |x|₂ <= 1 }.
Computes the smallest uniform scaling of this ellipsoid for which it still intersects other
.
√ minₓ (x-center)ᵀAᵀA(x-center) s.t. x ∈ other. Note that if center ∈ other, then we expect scaling = 0 and x = center (up to precision).
other
must have the same ambient_dimension as this. std::exception | if other is empty. |
std::exception | if ambient_dimension() == 0 |
|
static |
Constructs the minimum-volume ellipsoid which contains all of the points
.
This is commonly referred to as the outer Löwner-John ellipsoid.
points | is a d-by-n matrix, where d is the ambient dimension and each column represents one point. |
rank_tol | the singular values of the data matrix will be considered non-zero if they are strictly greater than rank_tol * max_singular_value . The default is 1e-6 to be compatible with common solver tolerances. This is used to detect if the data lies on a lower-dimensional affine space than the ambient dimension of the ellipsoid. If this is the case, then use AffineBall::MinimumVolumeCircumscribedEllipsoid instead. |
std::exception | if the MathematicalProgram fails to solve. If this were to happen (due to numerical issues), then increasing rank_tol should provide a mitigation. |
std::exception | if points includes NaNs or infinite values. |
std::exception | if the numerical data rank of points is less than d. |
|
default |
|
default |
Hyperellipsoid Scale | ( | double | scale | ) | const |
Results a new Hyperellipsoid that is a scaled version of this
about the center.
Any point on the boundary of the ellipsoid, x, is now translated to a new point, x*, such that ||x* - center|| = ||x - center|| * pow(scale, 1.0/ambient_dimension()). The volume of the resulting shape is scaled up by 'scale'.
scale
> 0. void Serialize | ( | Archive * | a | ) |
Passes this object to an Archive.
Refer to YAML Serialization for background.
std::pair<std::unique_ptr<Shape>, math::RigidTransformd> ToShapeWithPose |
Constructs a Ellipsoid shape description of this set.
Note that the choice of ellipsoid is not unique. This method chooses to order the Ellipsoid parameters a ≥ b ≥ c.
std::exception | if ambient_dimension() != 3 |
std::exception | if A is not invertible (has any eigenvalue less than sqrt(1e-12)). This tolerance is not carefully tuned (yet). We use Eigen's SelfAdjointEigenSolver to take the eigenvalues of AᵀA; this solver is listed as having accuracy "Good": https://eigen.tuxfamily.org/dox/group__TopicLinearAlgebraDecompositions.html but how does that translate into a numerical precision? |
double Volume | ( | ) | const |
Computes the volume for the hyperellipsoid set.