pub struct ConvexSurfaceMesh2d { /* private fields */ }Expand description
The vertices and edges that make up a convex polygon.
ConvexPolytope::<2> and ConvexSurfaceMesh2d can both represent
2d convex polygons. The first is defined implicitly as the convex hull
of a set of points. It stores the given point set without any modification,
and can therefore be constructed quickly. The implicit convex hull is
formed by SupportMapping during intersection tests of
Convex(ConvexPolygon) with other Convex(_) types.
In contrast, ConvexSurfaceMesh2d explicitly computes the convex hull
on construction. After construction, the vertices of the shape include
only the points on the convex hull in a counter-clockwise order.
Using this representation, ConvexSurfaceMesh2d is able to provide
implementations of Volume, IsPointInside, and IntersectsAt.
The native ConvexSurfaceMesh2d–ConvexSurfaceMesh2d intersection test
is much faster than the generic Convex(ConvexPolygon) intersection test.
§Examples
Construction:
use hoomd_geometry::shape::ConvexSurfaceMesh2d;
let triangle = ConvexSurfaceMesh2d::from_point_set([
[1.0, -1.0].into(),
[-1.0, -1.0].into(),
[0.0, 1.0].into(),
])?;Intersection tests:
use hoomd_geometry::{Convex, IntersectsAt, shape::ConvexSurfaceMesh2d};
use hoomd_vector::{Angle, Cartesian};
use std::f64::consts::PI;
let rectangle = ConvexSurfaceMesh2d::from_point_set([
[-2.0, -1.0].into(),
[2.0, -1.0].into(),
[2.0, 1.0].into(),
[-2.0, 1.0].into(),
])?;
assert!(!rectangle.intersects_at(
&rectangle,
&[0.0, 2.1].into(),
&Angle::default()
));
assert!(rectangle.intersects_at(
&rectangle,
&[0.0, 2.1].into(),
&Angle::from(PI / 2.0)
));Implementations§
Source§impl ConvexSurfaceMesh2d
impl ConvexSurfaceMesh2d
Sourcepub fn from_point_set<I>(points: I) -> Result<Self, Error>where
I: IntoIterator<Item = Cartesian<2>>,
pub fn from_point_set<I>(points: I) -> Result<Self, Error>where
I: IntoIterator<Item = Cartesian<2>>,
Create an convex surface mesh from the convex hull of the given set of points.
The resulting shape contains a subset of the given points, including only the non-degenerate points on the convex hull. The result’s vertices are arranged in a counter-clockwise order.
§Errors
Error::DegeneratePolytopewhen there are fewer than 3 non-collinear points.
§Example
use hoomd_geometry::shape::ConvexSurfaceMesh2d;
let equilateral_triangle = ConvexSurfaceMesh2d::from_point_set([
[1.0, 0.0].into(),
[0.5, f64::sqrt(3.0) / 2.0].into(),
[-0.5, f64::sqrt(3.0) / 2.0].into(),
])?;Sourcepub fn construct_convex_hull(
points: Vec<Cartesian<2>>,
) -> Result<Vec<Cartesian<2>>, Error>
pub fn construct_convex_hull( points: Vec<Cartesian<2>>, ) -> Result<Vec<Cartesian<2>>, Error>
Compute the convex hull of points in two dimensions.
The resulting vector contains a subset of the given points, including only the non-degenerate points on the convex hull. The output vertices are arranged in a counter-clockwise order.
§Errors
Returns Error if the input vertices do not form a convex body with 3 or more points.
§Example
use hoomd_geometry::shape::ConvexSurfaceMesh2d;
let points = vec![
[1.0, 0.0].into(),
[0.5, f64::sqrt(3.0) / 2.0].into(),
[-0.5, f64::sqrt(3.0) / 2.0].into(),
];
let hull_vertices = ConvexSurfaceMesh2d::construct_convex_hull(points)?;Trait Implementations§
Source§impl BoundingSphereRadius for ConvexSurfaceMesh2d
impl BoundingSphereRadius for ConvexSurfaceMesh2d
Source§fn bounding_sphere_radius(&self) -> PositiveReal
fn bounding_sphere_radius(&self) -> PositiveReal
Radius of a circle that bounds the shape.
The circle has the same local origin as the shape self.
§Example
use approxim::assert_relative_eq;
use hoomd_geometry::{BoundingSphereRadius, shape::ConvexSurfaceMesh2d};
let triangle = ConvexSurfaceMesh2d::from_point_set([
[1.0, -1.0].into(),
[-1.0, -1.0].into(),
[0.0, 1.0].into(),
])?;
let bounding_radius = triangle.bounding_sphere_radius();
assert_relative_eq!(bounding_radius.get(), 2.0_f64.sqrt());Source§impl Clone for ConvexSurfaceMesh2d
impl Clone for ConvexSurfaceMesh2d
Source§fn clone(&self) -> ConvexSurfaceMesh2d
fn clone(&self) -> ConvexSurfaceMesh2d
1.0.0 · Source§fn clone_from(&mut self, source: &Self)
fn clone_from(&mut self, source: &Self)
source. Read moreSource§impl Debug for ConvexSurfaceMesh2d
impl Debug for ConvexSurfaceMesh2d
Source§impl<'de> Deserialize<'de> for ConvexSurfaceMesh2d
impl<'de> Deserialize<'de> for ConvexSurfaceMesh2d
Source§fn deserialize<__D>(__deserializer: __D) -> Result<Self, __D::Error>where
__D: Deserializer<'de>,
fn deserialize<__D>(__deserializer: __D) -> Result<Self, __D::Error>where
__D: Deserializer<'de>,
Source§impl<R> IntersectsAt<ConvexSurfaceMesh2d, Cartesian<2>, R> for ConvexSurfaceMesh2d
impl<R> IntersectsAt<ConvexSurfaceMesh2d, Cartesian<2>, R> for ConvexSurfaceMesh2d
Source§fn intersects_at(&self, other: &Self, v_ij: &Cartesian<2>, o_ij: &R) -> bool
fn intersects_at(&self, other: &Self, v_ij: &Cartesian<2>, o_ij: &R) -> bool
Test convex polygon intersections using the separating planes method.
When the number of vertices is small, separating planes is significantly
faster than the Xenocollide algorithm implemented for the Convex<ConvexPolygon>
type, even though separating planes is $O(n^2)$.
§Example
use hoomd_geometry::{IntersectsAt, shape::ConvexSurfaceMesh2d};
use hoomd_vector::{Angle, Cartesian};
use std::f64::consts::PI;
let rectangle = ConvexSurfaceMesh2d::from_point_set([
[-2.0, -1.0].into(),
[2.0, -1.0].into(),
[2.0, 1.0].into(),
[-2.0, 1.0].into(),
])?;
assert!(!rectangle.intersects_at(
&rectangle,
&[0.0, 2.1].into(),
&Angle::default()
));
assert!(rectangle.intersects_at(
&rectangle,
&[0.0, 2.1].into(),
&Angle::from(PI / 2.0)
));Source§fn approximate_separation_distance(
&self,
other: &S,
v_ij: &V,
o_ij: &R,
resolution: PositiveReal,
) -> f64where
V: InnerProduct,
fn approximate_separation_distance(
&self,
other: &S,
v_ij: &V,
o_ij: &R,
resolution: PositiveReal,
) -> f64where
V: InnerProduct,
Source§impl<R> IntersectsAtGlobal<ConvexSurfaceMesh2d, Cartesian<2>, R> for ConvexSurfaceMesh2d
impl<R> IntersectsAtGlobal<ConvexSurfaceMesh2d, Cartesian<2>, R> for ConvexSurfaceMesh2d
Source§impl IsPointInside<Cartesian<2>> for ConvexSurfaceMesh2d
impl IsPointInside<Cartesian<2>> for ConvexSurfaceMesh2d
Source§fn is_point_inside(&self, point: &Cartesian<2>) -> bool
fn is_point_inside(&self, point: &Cartesian<2>) -> bool
Check if a point is inside the convex polygon.
§Example
use hoomd_geometry::{IsPointInside, shape::ConvexSurfaceMesh2d};
let triangle = ConvexSurfaceMesh2d::from_point_set([
[1.0, -1.0].into(),
[-1.0, -1.0].into(),
[0.0, 1.0].into(),
])?;
assert!(triangle.is_point_inside(&[0.0, 0.0].into()));
assert!(triangle.is_point_inside(&[-0.9, -0.9].into()));
assert!(!triangle.is_point_inside(&[-1.5, 2.0].into()));Source§impl PartialEq for ConvexSurfaceMesh2d
impl PartialEq for ConvexSurfaceMesh2d
Source§impl Serialize for ConvexSurfaceMesh2d
impl Serialize for ConvexSurfaceMesh2d
Source§impl SupportMapping<Cartesian<2>> for ConvexSurfaceMesh2d
impl SupportMapping<Cartesian<2>> for ConvexSurfaceMesh2d
Source§fn support_mapping(&self, n: &Cartesian<2>) -> Cartesian<2>
fn support_mapping(&self, n: &Cartesian<2>) -> Cartesian<2>
Find the point on a shape that is the furthest in a given direction.
ConvexSurfaceMesh2d implements SupportMapping to enable
intersection tests between mixed convex types.
§Example
use hoomd_geometry::{
Convex, IntersectsAt,
shape::{Circle, ConvexSurfaceMesh2d, Sphero},
};
use hoomd_vector::{Angle, Cartesian};
use std::f64::consts::PI;
let circle = Convex(Circle {
radius: 0.5.try_into()?,
});
let rectangle = ConvexSurfaceMesh2d::from_point_set([
[-1.5, -1.0].into(),
[1.5, -1.0].into(),
[-1.5, 1.0].into(),
[1.5, 1.0].into(),
])?;
let rounded_rectangle = Convex(Sphero {
shape: rectangle,
rounding_radius: 0.5.try_into()?,
});
assert!(rounded_rectangle.intersects_at(
&circle,
&[2.4, 0.0].into(),
&Angle::default()
));
assert!(!rounded_rectangle.intersects_at(
&circle,
&[0.0, 2.4].into(),
&Angle::default()
));
assert!(circle.intersects_at(
&rounded_rectangle,
&[0.0, 2.4].into(),
&Angle::from(PI / 2.0)
));Source§impl<const MAX_VERTICES: usize> TryFrom<ConvexPolytope<2, MAX_VERTICES>> for ConvexSurfaceMesh2d
impl<const MAX_VERTICES: usize> TryFrom<ConvexPolytope<2, MAX_VERTICES>> for ConvexSurfaceMesh2d
Source§fn try_from(
v: ConvexPolytope<2, MAX_VERTICES>,
) -> Result<ConvexSurfaceMesh2d, Error>
fn try_from( v: ConvexPolytope<2, MAX_VERTICES>, ) -> Result<ConvexSurfaceMesh2d, Error>
Construct the convex hull of a ConvexPolytope<2>.
§Errors
Error::DegeneratePolytopewhen there are fewer than 3 non-collinear points.
§Example
Invalid conversion
use hoomd_geometry::shape::{ConvexPolygon, ConvexSurfaceMesh2d};
let equilateral_triangle = ConvexPolygon::regular(3);
let mesh = ConvexSurfaceMesh2d::try_from(equilateral_triangle)?;Source§impl Volume for ConvexSurfaceMesh2d
impl Volume for ConvexSurfaceMesh2d
Source§fn volume(&self) -> f64
fn volume(&self) -> f64
Compute the area of the convex polygon.
§Example
use approxim::assert_relative_eq;
use hoomd_geometry::{Volume, shape::ConvexSurfaceMesh2d};
let triangle = ConvexSurfaceMesh2d::from_point_set([
[-1.0, -1.0].into(),
[1.0, -1.0].into(),
[1.0, 1.0].into(),
])?;
let area = triangle.volume();
assert_relative_eq!(area, 2.0);impl StructuralPartialEq for ConvexSurfaceMesh2d
Auto Trait Implementations§
impl Freeze for ConvexSurfaceMesh2d
impl RefUnwindSafe for ConvexSurfaceMesh2d
impl Send for ConvexSurfaceMesh2d
impl Sync for ConvexSurfaceMesh2d
impl Unpin for ConvexSurfaceMesh2d
impl UnwindSafe for ConvexSurfaceMesh2d
Blanket Implementations§
Source§impl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere
T: ?Sized,
Source§fn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Source§impl<T> CloneToUninit for Twhere
T: Clone,
impl<T> CloneToUninit for Twhere
T: Clone,
Source§impl<T> IntoEither for T
impl<T> IntoEither for T
Source§fn into_either(self, into_left: bool) -> Either<Self, Self>
fn into_either(self, into_left: bool) -> Either<Self, Self>
self into a Left variant of Either<Self, Self>
if into_left is true.
Converts self into a Right variant of Either<Self, Self>
otherwise. Read moreSource§fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
self into a Left variant of Either<Self, Self>
if into_left(&self) returns true.
Converts self into a Right variant of Either<Self, Self>
otherwise. Read more