An affine subspace (also known as a "flat", a "linear variety", or a "linear manifold") is a vector subspace of some Euclidean space, potentially translated so as to not pass through the origin.
Examples include points, lines, and planes (not necessarily through the origin).
An affine subspace is described by a basis of its corresponding vector subspace, plus a translation. This description is not unique as any point in the affine subspace can be used as a translation, and any basis of the corresponding vector subspace is valid.
An affine subspace can never be empty, because a vector subspace can never be empty. Thus, the translation will always be contained in the flat. An affine subspace is bounded if it is a point, which is when the basis has zero columns.
#include <drake/geometry/optimization/affine_subspace.h>
Public Member Functions | |
AffineSubspace () | |
Constructs a default (zero-dimensional, nonempty) affine subspace. More... | |
AffineSubspace (const Eigen::Ref< const Eigen::MatrixXd > &basis, const Eigen::Ref< const Eigen::VectorXd > &translation) | |
Constructs the affine subspace from an n-by-m matrix describing the basis, where n is the ambient dimension, and m is the dimension of the subspace, and from an n-dimensional vector describing the translation. More... | |
AffineSubspace (const ConvexSet &set, std::optional< double > tol=std::nullopt) | |
Constructs an affine subspace as the affine hull of another convex set. More... | |
~AffineSubspace () final | |
const Eigen::MatrixXd & | basis () const |
Returns the basis in an n-by-m matrix, where n is the ambient dimension, and m is the number of vectors in the basis. More... | |
const Eigen::VectorXd & | translation () const |
Returns the translation as a length n vector. More... | |
template<typename Archive > | |
void | Serialize (Archive *a) |
Passes this object to an Archive. More... | |
int | AffineDimension () const |
Returns the affine dimension of this set. More... | |
Eigen::MatrixXd | ToLocalCoordinates (const Eigen::Ref< const Eigen::MatrixXd > &x) const |
Given a point x in the standard basis of the ambient space, returns the coordinates of x in the basis of the AffineSubspace, with the zero point at translation_. More... | |
Eigen::MatrixXd | ToGlobalCoordinates (const Eigen::Ref< const Eigen::MatrixXd > &y) const |
Given a point y in the basis of the AffineSubspace, with the zero point at translation_, returns the coordinates of y in the standard basis of the ambient space. More... | |
bool | ContainedIn (const AffineSubspace &other, double tol=1e-15) const |
Returns true if this AffineSubspace is contained in other . More... | |
bool | IsNearlyEqualTo (const AffineSubspace &other, double tol=1e-15) const |
Returns true if the two AffineSubspaces describe the same set, by checking that each set is contained in the other. More... | |
Eigen::MatrixXd | OrthogonalComplementBasis () const |
Returns an orthonormal basis of the vector subspace which is orthogonal to this AffineSubspace. More... | |
bool | IsBounded (Parallelism parallelism=Parallelism::None()) const |
An AffineSubspace is bounded if and only if it is zero-dimensional (i.e., a point). More... | |
Implements CopyConstructible, CopyAssignable, MoveConstructible, MoveAssignable | |
AffineSubspace (const AffineSubspace &)=default | |
AffineSubspace & | operator= (const AffineSubspace &)=default |
AffineSubspace (AffineSubspace &&)=default | |
AffineSubspace & | operator= (AffineSubspace &&)=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... | |
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 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 |
AffineSubspace | ( | ) |
Constructs a default (zero-dimensional, nonempty) affine subspace.
|
explicit |
Constructs the affine subspace from an n-by-m matrix describing the basis, where n is the ambient dimension, and m is the dimension of the subspace, and from an n-dimensional vector describing the translation.
The set is {x | x = translation + basis*y, y ∈ Rᵐ} The columns must be linearly independent.
|
explicit |
Constructs an affine subspace as the affine hull of another convex set.
The generic approach is to find a feasible point in the set, and then iteratively compute feasible vectors until we have a basis that spans the set. If you pass in a convex set whose points are matrix-valued (e.g. a Spectrahedron), then the affine subspace will work over a flattened representation of those coordinates. (So a Spectrahedron with n-by-n matrices will output an AffineSubspace with ambient dimension (n * (n+1)) / 2.)
tol
sets the numerical precision of the computation. For each dimension, a pair of feasible points are constructed, so as to maximize the displacement in that dimension. If their displacement along that dimension is larger than tol, then the vector connecting the points is added as a basis vector.
std::exception | if set is empty. |
std::exception | if tol < 0 . |
For several subclasses of ConvexSet, there is a closed-form computation (or more efficient numerical computation) that is preferred.
tol
is used as the numerical tolerance for the rank of the matrix. Pass std::nullopt
for tol
to use Eigen's automatic tolerance computation.tol
is ignored.tol
is propagated to the constituent calls. (This is not done if the Cartesian product has an associated affine transformation.)tol
is ignored.tol
has the same meaning as in the generic affine hull computation.tol
is ignored. This also encompasses sets which are obviously a singleton point, as determined via MaybeGetPoint.tol
is used as the numerical tolerance for the rank of the matrix. Pass std::nullopt
for tol
to use Eigen's automatic tolerance computation.
|
final |
int AffineDimension | ( | ) | const |
Returns the affine dimension of this set.
For an affine subspace, this is simply the number of columns in the basis_ matrix. A point will have affine dimension zero.
const Eigen::MatrixXd& basis | ( | ) | const |
Returns the basis in an n-by-m matrix, where n is the ambient dimension, and m is the number of vectors in the basis.
bool ContainedIn | ( | const AffineSubspace & | other, |
double | tol = 1e-15 |
||
) | const |
Returns true
if this
AffineSubspace is contained in other
.
This is computed by checking if translation()
is in other
and then checking if each basis vector is in the span of the basis of other
. The latter step requires finding a least-squares solution, so a nonzero tolerance (tol
) is almost always necessary. (You may have to adjust the default tolerance depending on the dimension of your space and the scale of your basis vectors.)
bool IsBounded |
An AffineSubspace is bounded if and only if it is zero-dimensional (i.e., a point).
parallelism | Ignored – no parallelization is used. |
bool IsNearlyEqualTo | ( | const AffineSubspace & | other, |
double | tol = 1e-15 |
||
) | const |
Returns true if the two AffineSubspaces describe the same set, by checking that each set is contained in the other.
|
default |
|
default |
Eigen::MatrixXd OrthogonalComplementBasis | ( | ) | const |
Returns an orthonormal basis of the vector subspace which is orthogonal to this AffineSubspace.
void Serialize | ( | Archive * | a | ) |
Passes this object to an Archive.
Refer to YAML Serialization for background.
Eigen::MatrixXd ToGlobalCoordinates | ( | const Eigen::Ref< const Eigen::MatrixXd > & | y | ) | const |
Given a point y in the basis of the AffineSubspace, with the zero point at translation_, returns the coordinates of y in the standard basis of the ambient space.
If the AffineSubspace is a point, it has an empty basis, so the only possible local coordinates are also empty (and should be passed in as a length-zero vector). Each column of the input should be a vector in the affine subspace, represented in its local coordinates, and the corresponding column of the output will be its representation in the coordinate system of the ambient space.
Eigen::MatrixXd ToLocalCoordinates | ( | const Eigen::Ref< const Eigen::MatrixXd > & | x | ) | const |
Given a point x in the standard basis of the ambient space, returns the coordinates of x in the basis of the AffineSubspace, with the zero point at translation_.
The component of x that is orthogonal to the AffineSubspace (if it exists) is discarded, so ToGlobalCoordinates(ToLocalCoordinates(x)) is equivalent to Project(x). Note that if the AffineSubspace is a point, the basis is empty, so the local coordinates will also be empty (and returned as a length-zero vector). Each column of the input should be a vector in the ambient space, and the corresponding column of the output will be its representation in the local coordinates of the affine subspace.
const Eigen::VectorXd& translation | ( | ) | const |
Returns the translation as a length n vector.