Points

Name

Points -- point object and related functions.

Synopsis


#include <gts.h>


#define     GTS_POINT_CLASS                 (klass)
#define     GTS_POINT                       (obj)
#define     GTS_IS_POINT                    (obj)
struct      GtsPointClass;
struct      GtsPoint;

GtsPointClass* gts_point_class              (void);
GtsPoint*   gts_point_new                   (GtsPointClass *klass,
                                             gdouble x,
                                             gdouble y,
                                             gdouble z);
void        gts_point_set                   (GtsPoint *p,
                                             gdouble x,
                                             gdouble y,
                                             gdouble z);
#define     gts_point_is_in_rectangle       (p, p1, p2)
GtsPoint*   gts_segment_triangle_intersection
                                            (GtsSegment *s,
                                             GtsTriangle *t,
                                             gboolean boundary,
                                             GtsPointClass *klass);
void        gts_point_transform             (GtsPoint *p,
                                             GtsMatrix *m);
gdouble     gts_point_distance              (GtsPoint *p1,
                                             GtsPoint *p2);
gdouble     gts_point_distance2             (GtsPoint *p1,
                                             GtsPoint *p2);
gdouble     gts_point_orientation_3d        (GtsPoint *p1,
                                             GtsPoint *p2,
                                             GtsPoint *p3,
                                             GtsPoint *p4);
gint        gts_point_orientation_3d_sos    (GtsPoint *p1,
                                             GtsPoint *p2,
                                             GtsPoint *p3,
                                             GtsPoint *p4);
enum        GtsIntersect;
gdouble     gts_point_in_circle             (GtsPoint *p,
                                             GtsPoint *p1,
                                             GtsPoint *p2,
                                             GtsPoint *p3);
gdouble     gts_point_in_triangle_circle    (GtsPoint *p,
                                             GtsTriangle *t);
GtsIntersect gts_point_is_in_triangle       (GtsPoint *p,
                                             GtsTriangle *t);
gdouble     gts_point_orientation           (GtsPoint *p1,
                                             GtsPoint *p2,
                                             GtsPoint *p3);
gint        gts_point_orientation_sos       (GtsPoint *p1,
                                             GtsPoint *p2,
                                             GtsPoint *p3);
gdouble     gts_point_segment_distance2     (GtsPoint *p,
                                             GtsSegment *s);
gdouble     gts_point_segment_distance      (GtsPoint *p,
                                             GtsSegment *s);
void        gts_point_segment_closest       (GtsPoint *p,
                                             GtsSegment *s,
                                             GtsPoint *closest);
gdouble     gts_point_triangle_distance     (GtsPoint *p,
                                             GtsTriangle *t);
void        gts_point_triangle_closest      (GtsPoint *p,
                                             GtsTriangle *t,
                                             GtsPoint *closest);
gdouble     gts_point_triangle_distance2    (GtsPoint *p,
                                             GtsTriangle *t);
gboolean    gts_point_is_inside_surface     (GtsPoint *p,
                                             GNode *tree,
                                             gboolean is_open);

Description

Three-dimensional points and a few geometrical functions.

Details

GTS_POINT_CLASS()

#define     GTS_POINT_CLASS(klass)

Casts klass to GtsPointClass.

klass :

a descendant of the GtsPointClass.


GTS_POINT()

#define     GTS_POINT(obj)

Casts obj to GtsPoint.

obj :

a descendant of the GtsPoint.


GTS_IS_POINT()

#define     GTS_IS_POINT(obj)

Evaluates to TRUE if obj is a descendant of a GtsPoint, FALSE otherwise.

obj :

a pointer to test.


struct GtsPointClass

struct GtsPointClass {

  GtsObjectClass parent_class;
  gboolean binary;
};

The class for GtsPoint. No virtual functions are associated.


struct GtsPoint

struct GtsPoint {

  GtsObject object;

  gdouble x, y, z; /* must be contiguous (cast to robust functions) */
};

The point object.

GtsObject object

The parent object.

gdouble x

x coordinate.

gdouble y

y coordinate.

gdouble z

z coordinate.


gts_point_class ()

GtsPointClass* gts_point_class              (void);

Returns :

the GtsPointClass.


gts_point_new ()

GtsPoint*   gts_point_new                   (GtsPointClass *klass,
                                             gdouble x,
                                             gdouble y,
                                             gdouble z);

klass :

a GtsPointClass.

x :

the x-coordinate.

y :

the y-coordinate.

z :

the z-coordinate.

Returns :

a new GtsPoint.


gts_point_set ()

void        gts_point_set                   (GtsPoint *p,
                                             gdouble x,
                                             gdouble y,
                                             gdouble z);

Sets the coordinates of p.

p :

a GtsPoint.

x :

the x-coordinate.

y :

the y-coordinate.

z :

the z-coordinate.


gts_point_is_in_rectangle()

#define     gts_point_is_in_rectangle(p, p1, p2)

Evaluates to TRUE if p is contained in the box (boundary included) defined by its two corners p1 and p2, FALSE otherwise.

p :

a GtsPoint.

p1 :

the lower-left front corner of the box.

p2 :

the upper-right back corner of the box.


gts_segment_triangle_intersection ()

GtsPoint*   gts_segment_triangle_intersection
                                            (GtsSegment *s,
                                             GtsTriangle *t,
                                             gboolean boundary,
                                             GtsPointClass *klass);

Checks if s intersects t. If this is the case, creates a new point pi intersection of s with t.

This function is geometrically robust in the sense that it will not return a point if s and t do not intersect and will return a point if s and t do intersect. However, the point coordinates are subject to round-off errors.

Note that this function will not return any point if s is contained in the plane defined by t.

s :

a GtsSegment.

t :

a GtsTriangle.

boundary :

if TRUE, the boundary of t is taken into account.

klass :

a GtsPointClass to be used for the new point.

Returns :

a summit of t (if boundary is set to TRUE), one of the endpoints of s or a new GtsPoint, intersection of s with t or NULL if s and t don't intersect.


gts_point_transform ()

void        gts_point_transform             (GtsPoint *p,
                                             GtsMatrix *m);

Transform the coordinates of p according to m. (p[] becomes m[][].p[]).

p :

a GtsPoint.

m :

the GtsMatrix representing the transformation to apply to the coordinates of p.


gts_point_distance ()

gdouble     gts_point_distance              (GtsPoint *p1,
                                             GtsPoint *p2);

p1 :

a GtsPoint.

p2 :

another GtsPoint.

Returns :

the Euclidean distance between p1 and p2.


gts_point_distance2 ()

gdouble     gts_point_distance2             (GtsPoint *p1,
                                             GtsPoint *p2);

p1 :

a GtsPoint.

p2 :

another GtsPoint.

Returns :

the square of the Euclidean distance between p1 and p2.


gts_point_orientation_3d ()

gdouble     gts_point_orientation_3d        (GtsPoint *p1,
                                             GtsPoint *p2,
                                             GtsPoint *p3,
                                             GtsPoint *p4);

Checks if p4 lies above, below or on the plane passing through the points p1, p2 and p3. Below is defined so that p1, p2 and p3 appear in counterclockwise order when viewed from above the plane. The returned value is an approximation of six times the signed volume of the tetrahedron defined by the four points. This function uses adaptive floating point arithmetic and is consequently geometrically robust.

p1 :

a GtsPoint.

p2 :

a GtsPoint.

p3 :

a GtsPoint.

p4 :

a GtsPoint.

Returns :

a positive value if p4 lies below, a negative value if p4 lies above the plane, zero if the four points are coplanar.


gts_point_orientation_3d_sos ()

gint        gts_point_orientation_3d_sos    (GtsPoint *p1,
                                             GtsPoint *p2,
                                             GtsPoint *p3,
                                             GtsPoint *p4);

Checks if p4 lies above or below the plane passing through the points p1, p2 and p3. Below is defined so that p1, p2 and p3 appear in counterclockwise order when viewed from above the plane. This function uses adaptive floating point arithmetic and is consequently geometrically robust.

Simulation of Simplicity (SoS) is used to break ties when the orientation is degenerate (i.e. p4 lies on the plane defined by p1, p2 and p3).

p1 :

a GtsPoint.

p2 :

a GtsPoint.

p3 :

a GtsPoint.

p4 :

a GtsPoint.

Returns :

+1 if p4 lies below, -1 if p4 lies above the plane.


enum GtsIntersect

typedef enum 
{ 
  GTS_OUT = -1,
  GTS_ON = 0,
  GTS_IN = 1
} GtsIntersect;


gts_point_in_circle ()

gdouble     gts_point_in_circle             (GtsPoint *p,
                                             GtsPoint *p1,
                                             GtsPoint *p2,
                                             GtsPoint *p3);

Tests if the planar projection (x, y) of p is inside or outside the circle defined by the planar projection of p1, p2 and p3.

p :

a GtsPoint.

p1 :

a GtsPoint.

p2 :

a GtsPoint.

p3 :

a GtsPoint.

Returns :

a positive number if p lies inside, a negative number if p lies outside and zero if p lies on the circle.


gts_point_in_triangle_circle ()

gdouble     gts_point_in_triangle_circle    (GtsPoint *p,
                                             GtsTriangle *t);

Tests if the planar projection (x, y) of p is inside or outside the circumcircle of the planar projection of t. This function is geometrically robust.

p :

a GtsPoint.

t :

a GtsTriangle.

Returns :

a positive number if p lies inside, a negative number if p lies outside and zero if p lies on the circumcircle of t.


gts_point_is_in_triangle ()

GtsIntersect gts_point_is_in_triangle       (GtsPoint *p,
                                             GtsTriangle *t);

Tests if the planar projection (x, y) of p is inside, outside or on the boundary of the planar projection of t. This function is geometrically robust.

p :

a GtsPoint.

t :

a GtsTriangle.

Returns :

GTS_IN if p is inside t, GTS_ON if p is on the boundary of t, GTS_OUT otherwise.


gts_point_orientation ()

gdouble     gts_point_orientation           (GtsPoint *p1,
                                             GtsPoint *p2,
                                             GtsPoint *p3);

Checks for orientation of the projection of three points on the (x,y) plane. The result is also an approximation of twice the signed area of the triangle defined by the three points. This function uses adaptive floating point arithmetic and is consequently geometrically robust.

p1 :

a GtsPoint.

p2 :

a GtsPoint.

p3 :

a GtsPoint.

Returns :

a positive value if p1, p2 and p3 appear in counterclockwise order, a negative value if they appear in clockwise order and zero if they are colinear.


gts_point_orientation_sos ()

gint        gts_point_orientation_sos       (GtsPoint *p1,
                                             GtsPoint *p2,
                                             GtsPoint *p3);

Checks for orientation of the projection of three points on the (x,y) plane.

Simulation of Simplicity (SoS) is used to break ties when the orientation is degenerate (i.e. p3 lies on the line defined by p1 and p2).

p1 :

a GtsPoint.

p2 :

a GtsPoint.

p3 :

a GtsPoint.

Returns :

a positive value if p1, p2 and p3 appear in counterclockwise order or a negative value if they appear in clockwise order.


gts_point_segment_distance2 ()

gdouble     gts_point_segment_distance2     (GtsPoint *p,
                                             GtsSegment *s);

p :

a GtsPoint.

s :

a GtsSegment.

Returns :

the square of the minimun Euclidean distance between p and s.


gts_point_segment_distance ()

gdouble     gts_point_segment_distance      (GtsPoint *p,
                                             GtsSegment *s);

p :

a GtsPoint.

s :

a GtsSegment.

Returns :

the minimun Euclidean distance between p and s.


gts_point_segment_closest ()

void        gts_point_segment_closest       (GtsPoint *p,
                                             GtsSegment *s,
                                             GtsPoint *closest);

Set the coordinates of closest to the coordinates of the point belonging to s closest to p.

p :

a GtsPoint.

s :

a GtsSegment.

closest :

a GtsPoint.


gts_point_triangle_distance ()

gdouble     gts_point_triangle_distance     (GtsPoint *p,
                                             GtsTriangle *t);

p :

a GtsPoint.

t :

a GtsTriangle.

Returns :

the minimun Euclidean distance between p and t.


gts_point_triangle_closest ()

void        gts_point_triangle_closest      (GtsPoint *p,
                                             GtsTriangle *t,
                                             GtsPoint *closest);

Set the coordinates of closest to those of the point belonging to t and closest to p.

p :

a GtsPoint.

t :

a GtsTriangle.

closest :

a GtsPoint.


gts_point_triangle_distance2 ()

gdouble     gts_point_triangle_distance2    (GtsPoint *p,
                                             GtsTriangle *t);

p :

a GtsPoint.

t :

a GtsTriangle.

Returns :

the square of the minimun Euclidean distance between p and t.


gts_point_is_inside_surface ()

gboolean    gts_point_is_inside_surface     (GtsPoint *p,
                                             GNode *tree,
                                             gboolean is_open);

p :

a GtsPoint.

tree :

a bounding box tree of the faces of a closed, orientable surface (see gts_bb_tree_surface()).

is_open :

TRUE if the surface defined by tree is "open" i.e. its volume is negative, FALSE otherwise.

Returns :

TRUE if p is inside the surface defined by tree, FALSE otherwise.