ConvexSurfaceMesh2d

Struct ConvexSurfaceMesh2d 

Source
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 ConvexSurfaceMesh2dConvexSurfaceMesh2d 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

Source

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
§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(),
])?;
Source

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)?;
Source

pub fn vertices(&self) -> &[Cartesian<2>]

The vertices of the convex polygon in counter-clockwise order.

Trait Implementations§

Source§

impl BoundingSphereRadius for ConvexSurfaceMesh2d

Source§

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

Source§

fn clone(&self) -> ConvexSurfaceMesh2d

Returns a duplicate of the value. Read more
1.0.0 · Source§

fn clone_from(&mut self, source: &Self)

Performs copy-assignment from source. Read more
Source§

impl Debug for ConvexSurfaceMesh2d

Source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result

Formats the value using the given formatter. Read more
Source§

impl<'de> Deserialize<'de> for ConvexSurfaceMesh2d

Source§

fn deserialize<__D>(__deserializer: __D) -> Result<Self, __D::Error>
where __D: Deserializer<'de>,

Deserialize this value from the given Serde deserializer. Read more
Source§

impl<R> IntersectsAt<ConvexSurfaceMesh2d, Cartesian<2>, R> for ConvexSurfaceMesh2d
where RotationMatrix<2>: From<R>, R: Copy,

Source§

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, ) -> f64
where V: InnerProduct,

Approximate the amount of overlap between two shapes. Read more
Source§

impl<R> IntersectsAtGlobal<ConvexSurfaceMesh2d, Cartesian<2>, R> for ConvexSurfaceMesh2d
where R: Rotation + Rotate<Cartesian<2>>, RotationMatrix<2>: From<R>,

Source§

fn intersects_at_global( &self, other: &Self, r_self: &Cartesian<2>, o_self: &R, r_other: &Cartesian<2>, o_other: &R, ) -> bool

Test whether the set of points in one shape intersects with the set of another (in the global frame). Read more
Source§

impl IsPointInside<Cartesian<2>> for ConvexSurfaceMesh2d

Source§

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

Source§

fn eq(&self, other: &ConvexSurfaceMesh2d) -> bool

Tests for self and other values to be equal, and is used by ==.
1.0.0 · Source§

fn ne(&self, other: &Rhs) -> bool

Tests for !=. The default implementation is almost always sufficient, and should not be overridden without very good reason.
Source§

impl Serialize for ConvexSurfaceMesh2d

Source§

fn serialize<__S>(&self, __serializer: __S) -> Result<__S::Ok, __S::Error>
where __S: Serializer,

Serialize this value into the given Serde serializer. Read more
Source§

impl SupportMapping<Cartesian<2>> for ConvexSurfaceMesh2d

Source§

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

Source§

fn try_from( v: ConvexPolytope<2, MAX_VERTICES>, ) -> Result<ConvexSurfaceMesh2d, Error>

Construct the convex hull of a ConvexPolytope<2>.

§Errors
§Example

Invalid conversion

use hoomd_geometry::shape::{ConvexPolygon, ConvexSurfaceMesh2d};

let equilateral_triangle = ConvexPolygon::regular(3);
let mesh = ConvexSurfaceMesh2d::try_from(equilateral_triangle)?;
Source§

type Error = Error

The type returned in the event of a conversion error.
Source§

impl Volume for ConvexSurfaceMesh2d

Source§

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);
Source§

impl StructuralPartialEq for ConvexSurfaceMesh2d

Auto Trait Implementations§

Blanket Implementations§

Source§

impl<T> Any for T
where T: 'static + ?Sized,

Source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
Source§

impl<T> Borrow<T> for T
where T: ?Sized,

Source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
Source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

Source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
Source§

impl<T> CloneToUninit for T
where T: Clone,

Source§

unsafe fn clone_to_uninit(&self, dest: *mut u8)

🔬This is a nightly-only experimental API. (clone_to_uninit)
Performs copy-assignment from self to dest. Read more
Source§

impl<T> From<T> for T

Source§

fn from(t: T) -> T

Returns the argument unchanged.

Source§

impl<T, U> Into<U> for T
where U: From<T>,

Source§

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

Source§

impl<T> IntoEither for T

Source§

fn into_either(self, into_left: bool) -> Either<Self, Self>

Converts 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 more
Source§

fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
where F: FnOnce(&Self) -> bool,

Converts 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
Source§

impl<T> ToOwned for T
where T: Clone,

Source§

type Owned = T

The resulting type after obtaining ownership.
Source§

fn to_owned(&self) -> T

Creates owned data from borrowed data, usually by cloning. Read more
Source§

fn clone_into(&self, target: &mut T)

Uses borrowed data to replace owned data, usually by cloning. Read more
Source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

Source§

type Error = Infallible

The type returned in the event of a conversion error.
Source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
Source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

Source§

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
Source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.
Source§

impl<T> DeserializeOwned for T
where T: for<'de> Deserialize<'de>,