Skip to content

Commit

Permalink
Better implementations of DoMaybeGetFeasiblePoint for some subclasses…
Browse files Browse the repository at this point in the history
… of ConvexSet. (#19824)

fixup private doxygen
  • Loading branch information
cohnt authored Jul 23, 2023
1 parent 23e3325 commit 3fdd74c
Show file tree
Hide file tree
Showing 8 changed files with 61 additions and 0 deletions.
17 changes: 17 additions & 0 deletions geometry/optimization/cartesian_product.cc
Original file line number Diff line number Diff line change
Expand Up @@ -135,7 +135,24 @@ std::optional<VectorXd> CartesianProduct::DoMaybeGetPoint() const {
return std::nullopt;
}
}
return CartesianProduct::StackAndMaybeTransform(points);
}

std::optional<Eigen::VectorXd> CartesianProduct::DoMaybeGetFeasiblePoint()
const {
std::vector<VectorXd> points;
for (const auto& s : sets_) {
if (std::optional<VectorXd> point = s->MaybeGetFeasiblePoint()) {
points.push_back(std::move(*point));
} else {
return std::nullopt;
}
}
return CartesianProduct::StackAndMaybeTransform(points);
}

VectorXd CartesianProduct::StackAndMaybeTransform(
const std::vector<Eigen::VectorXd>& points) const {
// Stack the points.
const int y_ambient_dimension = A_ ? A_->rows() : ambient_dimension();
int start = 0;
Expand Down
9 changes: 9 additions & 0 deletions geometry/optimization/cartesian_product.h
Original file line number Diff line number Diff line change
Expand Up @@ -87,6 +87,15 @@ class CartesianProduct final : public ConvexSet {

std::optional<Eigen::VectorXd> DoMaybeGetPoint() const final;

std::optional<Eigen::VectorXd> DoMaybeGetFeasiblePoint() const final;

/* Given a list of vectors, one from each constituent set of this
CartesianProduct, this stacks them into a single vector. Then, if this
CartesianProduct has an associated transformation (in the form of an A_
matrix and b_ vector), it applies that transformation. */
Eigen::VectorXd StackAndMaybeTransform(
const std::vector<Eigen::VectorXd>& points) const;

bool DoPointInSet(const Eigen::Ref<const Eigen::VectorXd>& x,
double tol) const final;

Expand Down
4 changes: 4 additions & 0 deletions geometry/optimization/hyperellipsoid.cc
Original file line number Diff line number Diff line change
Expand Up @@ -193,6 +193,10 @@ bool Hyperellipsoid::DoIsEmpty() const {
return false;
}

std::optional<Eigen::VectorXd> Hyperellipsoid::DoMaybeGetFeasiblePoint() const {
return center_;
}

bool Hyperellipsoid::DoPointInSet(const Eigen::Ref<const VectorXd>& x,
double tol) const {
DRAKE_DEMAND(A_.cols() == x.size());
Expand Down
3 changes: 3 additions & 0 deletions geometry/optimization/hyperellipsoid.h
Original file line number Diff line number Diff line change
Expand Up @@ -117,6 +117,9 @@ class Hyperellipsoid final : public ConvexSet {

// N.B. No need to override DoMaybeGetPoint here.

/** Returns the center, which is always feasible. */
std::optional<Eigen::VectorXd> DoMaybeGetFeasiblePoint() const final;

bool DoPointInSet(const Eigen::Ref<const Eigen::VectorXd>& x,
double tol) const final;

Expand Down
16 changes: 16 additions & 0 deletions geometry/optimization/minkowski_sum.cc
Original file line number Diff line number Diff line change
Expand Up @@ -132,6 +132,22 @@ std::optional<VectorXd> MinkowskiSum::DoMaybeGetPoint() const {
return result;
}

std::optional<VectorXd> MinkowskiSum::DoMaybeGetFeasiblePoint() const {
std::optional<VectorXd> result;
for (const auto& s : sets_) {
if (std::optional<VectorXd> point = s->MaybeGetFeasiblePoint()) {
if (result.has_value()) {
*result += *point;
} else {
result = std::move(point);
}
} else {
return std::nullopt;
}
}
return result;
}

bool MinkowskiSum::DoPointInSet(const Eigen::Ref<const VectorXd>& x,
double) const {
// TODO(russt): Figure out if there is a general way to communicate tol
Expand Down
2 changes: 2 additions & 0 deletions geometry/optimization/minkowski_sum.h
Original file line number Diff line number Diff line change
Expand Up @@ -72,6 +72,8 @@ class MinkowskiSum final : public ConvexSet {

std::optional<Eigen::VectorXd> DoMaybeGetPoint() const final;

std::optional<Eigen::VectorXd> DoMaybeGetFeasiblePoint() const final;

bool DoPointInSet(const Eigen::Ref<const Eigen::VectorXd>& x,
double tol) const final;

Expand Down
8 changes: 8 additions & 0 deletions geometry/optimization/vpolytope.cc
Original file line number Diff line number Diff line change
Expand Up @@ -351,6 +351,14 @@ std::optional<VectorXd> VPolytope::DoMaybeGetPoint() const {
return std::nullopt;
}

std::optional<VectorXd> VPolytope::DoMaybeGetFeasiblePoint() const {
if (IsEmpty()) {
return std::nullopt;
} else {
return vertices_.col(0);
}
}

bool VPolytope::DoPointInSet(const Eigen::Ref<const VectorXd>& x,
double tol) const {
if (vertices_.cols() == 0) {
Expand Down
2 changes: 2 additions & 0 deletions geometry/optimization/vpolytope.h
Original file line number Diff line number Diff line change
Expand Up @@ -102,6 +102,8 @@ class VPolytope final : public ConvexSet {

std::optional<Eigen::VectorXd> DoMaybeGetPoint() const final;

std::optional<Eigen::VectorXd> DoMaybeGetFeasiblePoint() const final;

bool DoPointInSet(const Eigen::Ref<const Eigen::VectorXd>& x,
double tol) const final;

Expand Down

0 comments on commit 3fdd74c

Please sign in to comment.