DOLFIN
DOLFIN C++ interface
Static Public Member Functions | List of all members
dolfin::CollisionPredicates Class Reference

#include <CollisionPredicates.h>

Static Public Member Functions

static bool collides (const MeshEntity &entity, const Point &point)
 
static bool collides (const MeshEntity &entity_0, const MeshEntity &entity_1)
 
static bool collides_segment_point (const Point &p0, const Point &p1, const Point &point, std::size_t gdim)
 Check whether segment p0-p1 collides with point.
 
static bool collides_segment_point_1d (double p0, double p1, double point)
 Check whether segment p0-p1 collides with point (1D version)
 
static bool collides_segment_point_2d (const Point &p0, const Point &p1, const Point &point)
 Check whether segment p0-p1 collides with point (2D version)
 
static bool collides_segment_point_3d (const Point &p0, const Point &p1, const Point &point)
 Check whether segment p0-p1 collides with point (3D version)
 
static bool collides_segment_segment (const Point &p0, const Point &p1, const Point &q0, const Point &q1, std::size_t gdim)
 Check whether segment p0-p1 collides with segment q0-q1.
 
static bool collides_segment_segment_1d (double p0, double p1, double q0, double q1)
 Check whether segment p0-p1 collides with segment q0-q1 (1D version)
 
static bool collides_segment_segment_2d (const Point &p0, const Point &p1, const Point &q0, const Point &q1)
 Check whether segment p0-p1 collides with segment q0-q1 (2D version)
 
static bool collides_segment_segment_3d (const Point &p0, const Point &p1, const Point &q0, const Point &q1)
 Check whether segment p0-p1 collides with segment q0-q1 (3D version)
 
static bool collides_triangle_point (const Point &p0, const Point &p1, const Point &p2, const Point &point, std::size_t gdim)
 Check whether triangle p0-p1-p2 collides with point.
 
static bool collides_triangle_point_2d (const Point &p0, const Point &p1, const Point &p2, const Point &point)
 Check whether triangle p0-p1-p2 collides with point (2D version)
 
static bool collides_triangle_point_3d (const Point &p0, const Point &p1, const Point &p2, const Point &point)
 Check whether triangle p0-p1-p2 collides with point (3D version)
 
static bool collides_triangle_segment (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, std::size_t gdim)
 Check whether triangle p0-p1-p2 collides with segment q0-q1.
 
static bool collides_triangle_segment_2d (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1)
 Check whether triangle p0-p1-p2 collides with segment q0-q1 (2D version)
 
static bool collides_triangle_segment_3d (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1)
 Check whether triangle p0-p1-p2 collides with segment q0-q1 (3D version)
 
static bool collides_triangle_triangle (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2, std::size_t gdim)
 Check whether triangle p0-p1-p2 collides with triangle q0-q1-q2.
 
static bool collides_triangle_triangle_2d (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2)
 Check whether triangle p0-p1-p2 collides with triangle q0-q1-q2 (2D version)
 
static bool collides_triangle_triangle_3d (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2)
 Check whether triangle p0-p1-p2 collides with triangle q0-q1-q2 (3D version)
 
static bool collides_tetrahedron_point_3d (const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &point)
 Check whether tetrahedron p0-p1-p2-p3 collides with point.
 
static bool collides_tetrahedron_segment_3d (const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1)
 Check whether tetrahedron p0-p1-p2-p3 collides with segment q0-q1.
 
static bool collides_tetrahedron_triangle_3d (const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1, const Point &q2)
 Check whether tetrahedron p0-p1-p2-p3 collides with triangle q0-q1-q2.
 
static bool collides_tetrahedron_tetrahedron_3d (const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1, const Point &q2, const Point &q3)
 Check whether tetrahedron p0-p1-p2-p3 collides with tetrahedron q0-q1-q2.
 

Detailed Description

This class implements algorithms for detecting pairwise collisions between mesh entities of varying dimensions.

Member Function Documentation

◆ collides() [1/2]

bool CollisionPredicates::collides ( const MeshEntity entity,
const Point point 
)
static

Check whether entity collides with point.

Arguments entity (MeshEntity) The entity. point (Point) The point.

Returns bool True iff entity collides with cell.

◆ collides() [2/2]

bool CollisionPredicates::collides ( const MeshEntity entity_0,
const MeshEntity entity_1 
)
static

Check whether two entities collide.

Arguments entity_0 (MeshEntity) The first entity. entity_1 (MeshEntity) The second entity.

Returns bool True iff entity collides with cell.


The documentation for this class was generated from the following files: