Geogram
Version 1.9.1-rc
A programming library of geometric algorithms
|
Geometric functions and utilities. More...
Functions | |
vec3 | barycenter (const vec3 &p1, const vec3 &p2) |
Computes the barycenter of two points in 3d. More... | |
vec3 | barycenter (const vec3 &p1, const vec3 &p2, const vec3 &p3) |
Computes the barycenter of three points in 3d. More... | |
double | cos_angle (const vec3 &a, const vec3 &b) |
Computes the cosine of the angle between two 3d vectors. More... | |
double | angle (const vec3 &a, const vec3 &b) |
Computes the angle between two 3d vectors. More... | |
double | det (const vec2 &a, const vec2 &b) |
Computes the determinant of two vectors. More... | |
vec3 | triangle_normal (const vec3 &p1, const vec3 &p2, const vec3 &p3) |
Computes the normal of a 3d triangle. More... | |
double | triangle_area_3d (const double *p1, const double *p2, const double *p3) |
Computes the area of a 3d triangle. More... | |
double | triangle_area (const vec3 &p1, const vec3 &p2, const vec3 &p3) |
Computes the area of a 3d triangle. More... | |
double | triangle_signed_area_2d (const double *p1, const double *p2, const double *p3) |
Computes the area of a 2d triangle. More... | |
double | triangle_signed_area (const vec2 &p1, const vec2 &p2, const vec2 &p3) |
Computes the area of a 2d triangle. More... | |
double | triangle_area_2d (const double *p1, const double *p2, const double *p3) |
Computes the area of a 2d triangle. More... | |
vec2 | triangle_circumcenter (const vec2 &p1, const vec2 &p2, const vec2 &p3) |
Computes the center of the circumscribed circle of a 2d triangle. More... | |
bool | has_nan (const vec3 &v) |
Tests whether a 3d vector has a NaN (not a number) coordinate. More... | |
vec3 | perpendicular (const vec3 &V) |
Computes a 3d vector orthogonal to another one. More... | |
double | tetra_signed_volume (const vec3 &p1, const vec3 &p2, const vec3 &p3, const vec3 &p4) |
Computes the signed volume of a 3d tetrahedron. More... | |
double | tetra_signed_volume (const double *p1, const double *p2, const double *p3, const double *p4) |
Computes the signed volume of a 3d tetrahedron. More... | |
double | tetra_volume (const vec3 &p1, const vec3 &p2, const vec3 &p3, const vec3 &p4) |
Computes the volume of a 3d tetrahedron. More... | |
vec3 | tetra_circum_center (const vec3 &p1, const vec3 &p2, const vec3 &p3, const vec3 &p4) |
Computes the center of the circumscribed sphere of 3d tetrahedron. More... | |
void | triangle_centroid (const vec3 &p, const vec3 &q, const vec3 &r, double a, double b, double c, vec3 &Vg, double &V) |
Computes the centroid of a 3d triangle with weighted points. More... | |
double | triangle_mass (const vec3 &p, const vec3 &q, const vec3 &r, double a, double b, double c) |
Computes the mass of a 3d triangle with weighted points. More... | |
vec3 | random_point_in_triangle (const vec3 &p1, const vec3 &p2, const vec3 &p3) |
Generates a random point in a 3d triangle. More... | |
template<class COORD_T > | |
double | distance2 (const COORD_T *p1, const COORD_T *p2, coord_index_t dim) |
Computes the squared distance between two nd points. More... | |
template<class COORD_T > | |
double | distance (const COORD_T *p1, const COORD_T *p2, coord_index_t dim) |
Computes the distance between two nd points. More... | |
template<class VEC > | |
double | distance2 (const VEC &p1, const VEC &p2) |
Computes the squared distance between two nd points. More... | |
template<class VEC > | |
double | distance (const VEC &p1, const VEC &p2) |
Computes the distance between two nd points. More... | |
template<class COORD_T > | |
double | triangle_area (const COORD_T *p1, const COORD_T *p2, const COORD_T *p3, coord_index_t dim) |
Computes the area of a nd triangle. More... | |
template<class COORD_T > | |
void | triangle_centroid (const COORD_T *p, const COORD_T *q, const COORD_T *r, COORD_T a, COORD_T b, COORD_T c, double *Vg, double &V, coord_index_t dim) |
Computes the centroid of a 3d triangle with weighted points. More... | |
template<class VEC > | |
double | triangle_area (const VEC &p1, const VEC &p2, const VEC &p3) |
Computes the area of a nd triangle. More... | |
template<class VEC > | |
double | triangle_mass (const VEC &p, const VEC &q, const VEC &r, double a, double b, double c) |
Computes the mass of a nd triangle with weighted points. More... | |
template<class POINT > | |
POINT | triangle_circumcenter (const POINT &Q1, const POINT &Q2, const POINT &Q3, double *denom=nullptr) |
Computes the center of the circumscribed circle of a nd triangle. More... | |
template<class VEC > | |
void | triangle_centroid (const VEC &p, const VEC &q, const VEC &r, double a, double b, double c, VEC &Vg, double &V) |
Computes the centroid of a nd triangle with weighted points. More... | |
template<class VEC > | |
VEC | random_point_in_triangle (const VEC &p1, const VEC &p2, const VEC &p3) |
Generates a random point in a nd triangle. More... | |
template<class VEC > | |
VEC | random_point_in_tetra (const VEC &p1, const VEC &p2, const VEC &p3, const VEC &p4) |
Generates a random point in a nd tetrahedron. More... | |
template<class VEC > | |
double | point_segment_squared_distance (const VEC &point, const VEC &V0, const VEC &V1, VEC &closest_point, double &lambda0, double &lambda1) |
Computes the point closest to a given point in a nd segment. More... | |
template<class VEC > | |
double | point_segment_squared_distance (const VEC &point, const VEC &V0, const VEC &V1) |
Computes the point closest to a given point in a nd segment. More... | |
template<class VEC > | |
double | point_triangle_squared_distance (const VEC &point, const VEC &V0, const VEC &V1, const VEC &V2, VEC &closest_point, double &lambda0, double &lambda1, double &lambda2) |
Computes the point closest to a given point in a nd triangle. More... | |
template<class VEC > | |
double | point_triangle_squared_distance (const VEC &p, const VEC &q1, const VEC &q2, const VEC &q3) |
Computes the squared distance between a point and a nd triangle. More... | |
double | tetra_volume_from_edge_lengths (double u, double U, double v, double V, double w, double W) |
Computes the volume of a tetrahedron from edge lengths. More... | |
template<class VEC > | |
double | tetra_volume (const VEC &p1, const VEC &p2, const VEC &p3, const VEC &p4) |
Computes the volume of a nd tetrahedron. More... | |
template<int DIM> | |
double | tetra_volume (const double *p1, const double *p2, const double *p3, const double *p4) |
Computes the volume of a nd tetrahedron. More... | |
template<> | |
double | tetra_volume< 3 > (const double *p1, const double *p2, const double *p3, const double *p4) |
Computes the volume of a 3d tetrahedron. More... | |
const vec3 & | mesh_vertex (const Mesh &M, index_t v) |
Gets a mesh vertex by its index. More... | |
const vec3 & | mesh_vertex_ref (const Mesh &M, index_t v) |
Gets a mesh vertex by its index. More... | |
vec3 & | mesh_vertex_ref (Mesh &M, index_t v) |
Gets a mesh vertex by its index. More... | |
const vec3 & | mesh_corner_vertex (const Mesh &M, index_t c) |
Gets a mesh vertex by an incident corner index. More... | |
vec3 & | mesh_corner_vertex_ref (Mesh &M, index_t c) |
Gets a mesh vertex by an incident corner index. More... | |
const vec3 & | mesh_vertex_normal (const Mesh &M, index_t v) |
Gets a mesh vertex normal by vertex index. More... | |
vec3 & | mesh_vertex_normal_ref (Mesh &M, index_t v) |
Gets a mesh vertex normal by vertex index. More... | |
const vec3 & | mesh_vertex_normal_ref (const Mesh &M, index_t v) |
Gets a mesh vertex normal by vertex index. More... | |
double | mesh_facet_area (const Mesh &M, index_t f, index_t dim=0) |
Computes the area of a facet. More... | |
vec3 | mesh_facet_normal (const Mesh &M, index_t f) |
Computes the normal to a mesh facet. More... | |
vec3 | mesh_facet_center (const Mesh &M, index_t f) |
Gets the centroid of the vertices of a facet in a mesh. More... | |
vec3 | mesh_cell_center (const Mesh &M, index_t c) |
Gets the centroid of the vertices of a cell in a mesh. More... | |
vec3 | mesh_tet_center (const Mesh &M, index_t t) |
Gets the centroid of a tetrahedron in a mesh. More... | |
vec3 | mesh_corner_vector (const Mesh &M, index_t c1) |
Gets a vector by a mesh corner. More... | |
double | mesh_normal_angle (const Mesh &M, index_t c) |
Computes the angle between the normal vectors of two mesh facets sharing an edge. More... | |
double | mesh_unsigned_normal_angle (const Mesh &M, index_t f1, index_t f2) |
Computes the angle between the normal vectors of two mesh facets sharing an edge. More... | |
double | mesh_area (const Mesh &M, index_t dim) |
Computes the total surface area of a mesh in arbitrary dimension. More... | |
double | mesh_area (const Mesh &M) |
Computes the total surface area of a mesh. More... | |
double | mesh_enclosed_volume (const Mesh &M) |
Computes the volume enclosed by a surfacic mesh. More... | |
const vec3 & | halfedge_vertex_from (const Mesh &M, const MeshHalfedges::Halfedge &H) |
Gets the origin point of a Halfedge. More... | |
const vec3 & | halfedge_vertex_to (const Mesh &M, const MeshHalfedges::Halfedge &H) |
Gets the arrow extremity point of a Halfedge. More... | |
vec3 | halfedge_vector (const Mesh &M, const MeshHalfedges::Halfedge &H) |
Gets a 3d vector that connects the origin with the arrow extremity of a Halfedge. More... | |
double | edge_length (const Mesh &M, const MeshHalfedges::Halfedge &H) |
Gets the length of a Halfedge. More... | |
index_t | colocate (const double *points, coord_index_t dim, index_t nb_points, vector< index_t > &old2new, double tolerance=0.0, index_t stride=0, const std::string &nn_algo="default") |
Finds sets of identical points in a point set. More... | |
index_t | colocate_by_lexico_sort (const double *points, coord_index_t dim, index_t nb_points, vector< index_t > &old2new, index_t stride) |
Finds sets of identical points in a point set. More... | |
Geometric functions and utilities.
Computes the angle between two 3d vectors.
Computes the angle between two 2d vectors.
[in] | a | first vector |
[in] | b | second vector |
a
and b
in radians, in the interval \( [ 0 \ldots \pi ] \).[in] | a | first vector |
[in] | b | second vector |
Definition at line 266 of file geometry.h.
Computes the barycenter of two points in 3d.
Computes the barycenter of two points in 2d.
[in] | p1 | first point |
[in] | p2 | second point |
p1
and p2
Definition at line 189 of file geometry.h.
Computes the barycenter of three points in 3d.
Computes the barycenter of three points in 2d.
[in] | p1 | first point |
[in] | p2 | second point |
[in] | p3 | third point |
p1
, p2
and p3
Definition at line 217 of file geometry.h.
index_t GEO::Geom::colocate | ( | const double * | points, |
coord_index_t | dim, | ||
index_t | nb_points, | ||
vector< index_t > & | old2new, | ||
double | tolerance = 0.0 , |
||
index_t | stride = 0 , |
||
const std::string & | nn_algo = "default" |
||
) |
Finds sets of identical points in a point set.
[in] | points | the point array |
[in] | dim | dimension of the points |
[in] | nb_points | number of points |
[out] | old2new | an array of size nb_points. |
[in] | tolerance | threshold for colocating points. |
[in] | stride | number of doubles between two consecutive points (set to dim if unspecified). |
[in] | nn_algo | factory name for nearest neighbor search. |
index_t GEO::Geom::colocate_by_lexico_sort | ( | const double * | points, |
coord_index_t | dim, | ||
index_t | nb_points, | ||
vector< index_t > & | old2new, | ||
index_t | stride | ||
) |
Finds sets of identical points in a point set.
This version uses a lexicographic sort. It does not have a 'tolerance' parameter (only points with exactly the same coordinates can be colocated).
[in] | points | the point array |
[in] | dim | dimension of the points |
[in] | nb_points | number of points |
[out] | old2new | an array of size nb_points. |
[in] | stride | number of doubles between two consecutive points (set to dim if unspecified). |
Computes the cosine of the angle between two 3d vectors.
Computes the cosine of the angle between two 2d vectors.
[in] | a | first vector |
[in] | b | second vector |
a
and b
Definition at line 249 of file geometry.h.
Computes the determinant of two vectors.
[in] | a | first vector |
[in] | b | second vector |
a
and b
Definition at line 292 of file geometry.h.
|
inline |
Computes the distance between two nd points.
[in] | p1 | a pointer to the coordinates of the first point |
[in] | p2 | a pointer to the coordinates of the second point |
[in] | dim | dimension (number of coordinates of the points) |
p1
and p2
COORD_T | the numeric type of the point coordinates |
Definition at line 84 of file geometry_nd.h.
|
inline |
Computes the distance between two nd points.
[in] | p1 | first point |
[in] | p2 | second point |
VEC | the class that represents the points. VEC needs to implement data(), that returns a pointer to the coordinates of the point. |
p1
and p2
Definition at line 120 of file geometry_nd.h.
|
inline |
Computes the squared distance between two nd points.
[in] | p1 | a pointer to the coordinates of the first point |
[in] | p2 | a pointer to the coordinates of the second point |
[in] | dim | dimension (number of coordinates of the points) |
p1
and p2
COORD_T | the numeric type of the point coordinates |
Definition at line 65 of file geometry_nd.h.
|
inline |
Computes the squared distance between two nd points.
[in] | p1 | first point |
[in] | p2 | second point |
VEC | the class that represents the points. VEC needs to implement data(), that returns a pointer to the coordinates of the point. |
p1
and p2
Definition at line 100 of file geometry_nd.h.
|
inline |
Gets the length of a Halfedge.
[in] | M | the Mesh |
[in] | H | the Halfedge |
H
Definition at line 358 of file mesh_halfedges.h.
|
inline |
Gets a 3d vector that connects the origin with the arrow extremity of a Halfedge.
[in] | M | the Mesh |
[in] | H | the Halfedge |
H
Definition at line 346 of file mesh_halfedges.h.
|
inline |
Gets the origin point of a Halfedge.
[in] | M | the mesh |
[in] | H | the Halfedge |
H
Definition at line 319 of file mesh_halfedges.h.
|
inline |
Gets the arrow extremity point of a Halfedge.
[in] | M | the mesh |
[in] | H | the Halfedge |
H
Definition at line 331 of file mesh_halfedges.h.
|
inline |
Tests whether a 3d vector has a NaN (not a number) coordinate.
Tests whether a 2d vector has a NaN (not a number) coordinate.
[in] | v | a 3d vector |
[in] | v | a 2d vector |
Definition at line 429 of file geometry.h.
|
inline |
Computes the total surface area of a mesh.
[in] | M | the mesh |
Definition at line 301 of file mesh_geometry.h.
Computes the total surface area of a mesh in arbitrary dimension.
[in] | M | the mesh |
[in] | dim | the dimension to be used for the computation |
M
computed in dim d
. Gets the centroid of the vertices of a cell in a mesh.
[in] | M | the mesh |
[in] | c | the index of the facet |
f
in M
Definition at line 219 of file mesh_geometry.h.
Gets a vector by a mesh corner.
[in] | M | a const reference to the mesh |
[in] | c1 | a corner index in M |
c1
and pointing at the next corner around the facet incident to c1
Definition at line 256 of file mesh_geometry.h.
Gets a mesh vertex by an incident corner index.
[in] | M | the mesh |
[in] | c | the index of a corner incident to the vertex |
vth
vertex of a mesh Definition at line 100 of file mesh_geometry.h.
Gets a mesh vertex by an incident corner index.
[in] | M | the mesh |
[in] | c | the index of a corner incident to the vertex |
vth
vertex of a mesh Definition at line 111 of file mesh_geometry.h.
double GEO::Geom::mesh_enclosed_volume | ( | const Mesh & | M | ) |
Computes the volume enclosed by a surfacic mesh.
[in] | M | a closed surfacic mesh. |
M
. Computes the area of a facet.
[in] | M | a const reference to the mesh |
[in] | f | index of the facet |
[in] | dim | dimension that will be used to compute the area |
dim
first coordinates of the vertices only Definition at line 159 of file mesh_geometry.h.
Gets the centroid of the vertices of a facet in a mesh.
[in] | M | the mesh |
[in] | f | the index of the facet |
f
in M
Definition at line 202 of file mesh_geometry.h.
Computes the normal to a mesh facet.
[in] | M | the mesh |
[in] | f | the facet index in M |
f
Computes the angle between the normal vectors of two mesh facets sharing an edge.
[in] | M | a const reference to the mesh |
[in] | c | a corner index in M |
Gets the centroid of a tetrahedron in a mesh.
[in] | M | the mesh |
[in] | t | the index of the tetrahedron |
t
in M
Definition at line 235 of file mesh_geometry.h.
Computes the angle between the normal vectors of two mesh facets sharing an edge.
[in] | M | a const reference to the mesh |
[in] | f1,f2 | two facets of the mesh |
f1
and f2
in radians Gets a mesh vertex by its index.
[in] | M | the mesh |
[in] | v | the index of the vertex |
vth
vertex of a mesh Definition at line 64 of file mesh_geometry.h.
Gets a mesh vertex normal by vertex index.
[in] | M | the mesh |
[in] | v | the index of the vertex |
v
Definition at line 122 of file mesh_geometry.h.
Gets a mesh vertex normal by vertex index.
[in] | M | the mesh |
[in] | v | the index of the vertex |
v
Definition at line 146 of file mesh_geometry.h.
Gets a mesh vertex normal by vertex index.
[in] | M | the mesh |
[in] | v | the index of the vertex |
v
Definition at line 134 of file mesh_geometry.h.
Gets a mesh vertex by its index.
[in] | M | the mesh |
[in] | v | the index of the vertex |
vth
vertex of a mesh Definition at line 76 of file mesh_geometry.h.
Gets a mesh vertex by its index.
[in] | M | the mesh |
[in] | v | the index of the vertex |
vth
vertex of a mesh Definition at line 88 of file mesh_geometry.h.
Computes a 3d vector orthogonal to another one.
[in] | V | a 3d vector |
V
|
inline |
Computes the point closest to a given point in a nd segment.
[in] | point | the query point |
[in] | V0 | first extremity of the segment |
[in] | V1 | second extremity of the segment |
VEC | the class that represents the points. |
V0
, V1
] Definition at line 441 of file geometry_nd.h.
|
inline |
Computes the point closest to a given point in a nd segment.
[in] | point | the query point |
[in] | V0 | first extremity of the segment |
[in] | V1 | second extremity of the segment |
[out] | closest_point | the point closest to point in the segment [V0 , V1 ] |
[out] | lambda0 | barycentric coordinate of the closest point relative to V0 |
[out] | lambda1 | barycentric coordinate of the closest point relative to V1 |
VEC | the class that represents the points. |
V0
, V1
] Definition at line 403 of file geometry_nd.h.
|
inline |
Computes the squared distance between a point and a nd triangle.
See http://www.geometrictools.com/LibMathematics/Distance/Distance.html
[in] | p | the query point |
[in] | q1 | first vertex of the triangle |
[in] | q2 | second vertex of the triangle |
[in] | q3 | third vertex of the triangle |
VEC | the class that represents the points. |
V0
, V1
, V2
) Definition at line 696 of file geometry_nd.h.
|
inline |
Computes the point closest to a given point in a nd triangle.
See http://www.geometrictools.com/LibMathematics/Distance/Distance.html
[in] | point | the query point |
[in] | V0 | first vertex of the triangle |
[in] | V1 | second vertex of the triangle |
[in] | V2 | third vertex of the triangle |
[out] | closest_point | the point closest to point in the triangle (V0 , V1 , V2 ) |
[out] | lambda0 | barycentric coordinate of the closest point relative to V0 |
[out] | lambda1 | barycentric coordinate of the closest point relative to V1 |
[out] | lambda2 | barycentric coordinate of the closest point relative to V2 |
VEC | the class that represents the points. |
V0
, V1
, V2
) Definition at line 475 of file geometry_nd.h.
|
inline |
Generates a random point in a nd tetrahedron.
Uses Greg Turk's second method (see article in Graphic Gems).
[in] | p1 | first vertex of the triangle |
[in] | p2 | second vertex of the triangle |
[in] | p3 | third vertex of the triangle |
[in] | p4 | fourth vertex of the triangle |
p1
, p2
, p3
, p4
) VEC | the class used to represent the vertices of the triangle |
Definition at line 364 of file geometry_nd.h.
|
inline |
Generates a random point in a nd triangle.
Uses Greg Turk's second method (see article in Graphic Gems).
[in] | p1 | first vertex of the triangle |
[in] | p2 | second vertex of the triangle |
[in] | p3 | third vertex of the triangle |
p1
, p2
, p3
) VEC | the class used to represent the vertices of the triangle |
Definition at line 338 of file geometry_nd.h.
|
inline |
Generates a random point in a 3d triangle.
Uses Greg Turk's second method. Reference: Greg Turk, Generating Random Points in Triangles, Graphics Gems, p. 24-28, code: p. 649-650.
[in] | p1 | first vertex of the triangle |
[in] | p2 | second vertex of the triangle |
[in] | p3 | third vertex of the triangle |
p1
, p2
, p3
) Definition at line 582 of file geometry.h.
vec3 GEO::Geom::tetra_circum_center | ( | const vec3 & | p1, |
const vec3 & | p2, | ||
const vec3 & | p3, | ||
const vec3 & | p4 | ||
) |
Computes the center of the circumscribed sphere of 3d tetrahedron.
[in] | p1 | first vertex of the tetrahedron |
[in] | p2 | second vertex of the tetrahedron |
[in] | p3 | third vertex of the tetrahedron |
[in] | p4 | fourth vertex of the tetrahedron |
p1
, p2
, p3
, p4
)
|
inline |
Computes the signed volume of a 3d tetrahedron.
[in] | p1 | first vertex of the tetrahedron |
[in] | p2 | second vertex of the tetrahedron |
[in] | p3 | third vertex of the tetrahedron |
[in] | p4 | fourth vertex of the tetrahedron |
p1
, p2
, p3
, p4
) Definition at line 479 of file geometry.h.
|
inline |
Computes the signed volume of a 3d tetrahedron.
[in] | p1 | first vertex of the tetrahedron |
[in] | p2 | second vertex of the tetrahedron |
[in] | p3 | third vertex of the tetrahedron |
[in] | p4 | fourth vertex of the tetrahedron |
p1
, p2
, p3
, p4
) Definition at line 463 of file geometry.h.
|
inline |
Computes the volume of a nd tetrahedron.
Uses a form of generalized Heron formula: W. Kahan, "What has the Volume of a Tetrahedron to do with Computer Programming Languages?"
[in] | p1 | first vertex of the tetrahedron |
[in] | p2 | second vertex of the tetrahedron |
[in] | p3 | third vertex of the tetrahedron |
[in] | p4 | fourth vertex of the tetrahedron |
DIM | dimension of the points |
p1
, p2
, p3
, p4
) Definition at line 783 of file geometry_nd.h.
|
inline |
Computes the volume of a nd tetrahedron.
Uses a form of generalized Heron formula: W. Kahan, "What has the Volume of a Tetrahedron to do with Computer Programming Languages?"
[in] | p1 | first vertex of the tetrahedron |
[in] | p2 | second vertex of the tetrahedron |
[in] | p3 | third vertex of the tetrahedron |
[in] | p4 | fourth vertex of the tetrahedron |
VEC | the class that represents the points |
p1
, p2
, p3
, p4
) Definition at line 757 of file geometry_nd.h.
|
inline |
Computes the volume of a 3d tetrahedron.
[in] | p1 | first vertex of the tetrahedron |
[in] | p2 | second vertex of the tetrahedron |
[in] | p3 | third vertex of the tetrahedron |
[in] | p4 | fourth vertex of the tetrahedron |
p1
, p2
, p3
, p4
) Definition at line 500 of file geometry.h.
|
inline |
Computes the volume of a 3d tetrahedron.
Partial specialization of tetra_volume() for DIM == 3. It uses the standard formula, much simpler than the N dimension version.
[in] | p1 | first vertex of the tetrahedron |
[in] | p2 | second vertex of the tetrahedron |
[in] | p3 | third vertex of the tetrahedron |
[in] | p4 | fourth vertex of the tetrahedron |
p1
, p2
, p3
, p4
) Definition at line 809 of file geometry_nd.h.
|
inline |
Computes the volume of a tetrahedron from edge lengths.
Uses a form of generalized Heron formula: W. Kahan, "What has the Volume of a Tetrahedron to do with Computer Programming Languages?"
[in] | u | distance between p1 and p4 |
[in] | U | distance between p2 and p3 |
[in] | v | distance between p2 and p4 |
[in] | V | distance between p3 and p1 |
[in] | w | distance between p3 and p4 |
[in] | W | distance between p1 and p2 |
Definition at line 720 of file geometry_nd.h.
|
inline |
Computes the area of a nd triangle.
Uses Heron formula (that computes the area from the lengths of the three edges).
[in] | p1 | a pointer to the coordinates of the first vertex of the triangle |
[in] | p2 | a pointer to the coordinates of the second vertex of the triangle |
[in] | p3 | a pointer to the coordinates of the third vertex of the triangle |
[in] | dim | dimension of the points |
COORD_T | the numeric type that represents the coordinates of the points |
p1
, p2
, p3
) Definition at line 143 of file geometry_nd.h.
|
inline |
Computes the area of a nd triangle.
Uses Heron formula (that compures the area from the lengths of the three edges).
[in] | p1 | first vertex of the triangle |
[in] | p2 | second vertex of the triangle |
[in] | p3 | third vertex of the triangle |
VEC | the class used to represent the vertices of the triangle |
p1
, p2
, p3
) Definition at line 212 of file geometry_nd.h.
Computes the area of a 3d triangle.
Computes the area of a 2d triangle.
[in] | p1,p2,p3 | the three vertices of the triangle |
p1
, p2
, p3
).[in] | p1 | first vertex of the triangle |
[in] | p2 | second vertex of the triangle |
[in] | p3 | third vertex of the triangle |
p1
, p2
, p3
). Definition at line 348 of file geometry.h.
|
inline |
Computes the area of a 2d triangle.
[in] | p1 | first vertex of the triangle |
[in] | p2 | second vertex of the triangle |
[in] | p3 | third vertex of the triangle |
p1
, p2
, p3
). Definition at line 406 of file geometry.h.
|
inline |
Computes the area of a 3d triangle.
[in] | p1,p2,p3 | the three vertices of the triangle |
p1
, p2
, p3
). Definition at line 326 of file geometry.h.
|
inline |
Computes the centroid of a 3d triangle with weighted points.
The integrated weight varies linearly in the triangle.
[in] | p | a pointer to the coordinates of the first vertex of the triangle |
[in] | q | a pointer to the coordinates of the second vertex of the triangle |
[in] | r | a pointer to the coordinates of the third vertex of the triangle |
[in] | a | the weight associated with vertex p |
[in] | b | the weight associated with vertex q |
[in] | c | the weight associated with vertex r |
[out] | Vg | the total weight times the centroid ( a pointer to a caller-allocated array of dim COORD_Ts) |
[out] | V | the total weight |
[in] | dim | the dimension of the vertices |
COORD_T | the numeric type that represents the coordinates of the points |
Definition at line 178 of file geometry_nd.h.
|
inline |
Computes the centroid of a nd triangle with weighted points.
The integrated weight varies linearly in the triangle.
[in] | p | first vertex of the triangle |
[in] | q | second vertex of the triangle |
[in] | r | third vertex of the triangle |
[in] | a | the weight associated with vertex p |
[in] | b | the weight associated with vertex q |
[in] | c | the weight associated with vertex r |
[out] | Vg | the total weight times the centroid |
[out] | V | the total weight |
VEC | the class used to represent the vertices of the triangle |
Definition at line 311 of file geometry_nd.h.
|
inline |
Computes the centroid of a 3d triangle with weighted points.
The integrated weight varies linearly in the triangle.
[in] | p | first vertex of the triangle |
[in] | q | second vertex of the triangle |
[in] | r | third vertex of the triangle |
[in] | a | the weight associated with vertex p |
[in] | b | the weight associated with vertex q |
[in] | c | the weight associated with vertex r |
[out] | Vg | the total weight times the centroid |
[out] | V | the total weight |
Definition at line 534 of file geometry.h.
POINT GEO::Geom::triangle_circumcenter | ( | const POINT & | Q1, |
const POINT & | Q2, | ||
const POINT & | Q3, | ||
double * | denom = nullptr |
||
) |
Computes the center of the circumscribed circle of a nd triangle.
[in] | Q1 | first vertex of the triangle |
[in] | Q2 | second vertex of the triangle |
[in] | Q3 | third vertex of the triangle |
[out] | denom | if the parameter is non null, it is set to the denominator of the barycentric coordinates of the circumcenter. |
POINT | the class used to represent the vertices of the triangle |
p1
, p2
, p3
). Definition at line 267 of file geometry_nd.h.
Computes the center of the circumscribed circle of a 2d triangle.
[in] | p1 | first vertex of the triangle |
[in] | p2 | second vertex of the triangle |
[in] | p3 | third vertex of the triangle |
p1
, p2
, p3
).
|
inline |
Computes the mass of a nd triangle with weighted points.
The integrated weight varies linearly in the triangle.
[in] | p | first vertex of the triangle |
[in] | q | second vertex of the triangle |
[in] | r | third vertex of the triangle |
[in] | a | the weight associated with vertex p |
[in] | b | the weight associated with vertex q |
[in] | c | the weight associated with vertex r |
VEC | the class used to represent the vertices of the triangle |
p
, a
), ( q
, b
), ( r
, c
) Definition at line 237 of file geometry_nd.h.
|
inline |
Computes the mass of a 3d triangle with weighted points.
The integrated weight varies linearly in the triangle.
[in] | p | first vertex of the triangle |
[in] | q | second vertex of the triangle |
[in] | r | third vertex of the triangle |
[in] | a | the weight associated with vertex p |
[in] | b | the weight associated with vertex q |
[in] | c | the weight associated with vertex r |
p
, a
), ( q
, b
), ( r
, c
) Definition at line 563 of file geometry.h.
Computes the normal of a 3d triangle.
[in] | p1,p2,p3 | the three vertices of the triangle |
p1
, p2
, p3
). Definition at line 315 of file geometry.h.
Computes the area of a 2d triangle.
[in] | p1 | first vertex of the triangle |
[in] | p2 | second vertex of the triangle |
[in] | p3 | third vertex of the triangle |
p1
, p2
, p3
), positive if the triangle is oriented clockwise, negative otherwise. Definition at line 380 of file geometry.h.
|
inline |
Computes the area of a 2d triangle.
[in] | p1 | first vertex of the triangle |
[in] | p2 | second vertex of the triangle |
[in] | p3 | third vertex of the triangle |
p1
, p2
, p3
), positive if the triangle is oriented clockwise, negative otherwise. Definition at line 362 of file geometry.h.