xmsgrid  1.0
geoms.h
Go to the documentation of this file.
1 #pragma once
2 //------------------------------------------------------------------------------
8 //------------------------------------------------------------------------------
9 
10 //----- Included files ---------------------------------------------------------
11 
12 // 3. Standard Library Headers
13 #include <vector>
14 
15 // 4. External Library Headers
16 #include <xmscore/points/pt.h> // Pt3d, Pt2d, etc.
17 #include <xmscore/stl/vector.h>
18 
19 // 5. Shared Headers
20 
21 // 6. Non-shared Headers
22 
23 //----- Namespace declaration --------------------------------------------------
24 
25 namespace xms
26 {
28 enum Turn_enum {
33 };
34 
37  PT_ERROR = -1,
41 };
42 
43 //----- Forward declarations ---------------------------------------------------
44 
45 //----- Global functions -------------------------------------------------------
46 
47 bool gmPointInOrOnBox2d(const Pt3d& a_bMin, const Pt3d& a_bMax, const Pt3d& a_pt);
48 bool gmBoxesOverlap2d(const Pt3d& a_b1Min,
49  const Pt3d& a_b1Max,
50  const Pt3d& a_b2Min,
51  const Pt3d& a_b2Max);
53  const Pt3d& p2,
54  const Pt3d& p3,
55  double* a,
56  double* b,
57  double* c,
58  double* d);
60  const Pt3d* p2,
61  const Pt3d* p3,
62  double* a,
63  double* b,
64  double* c,
65  double* d);
66 double gmMiddleZ(const VecPt3d& a_points);
67 PtInOutOrOn_enum gmPtInCircumcircle(const Pt3d& pt, Pt3d circumcirclePts[3]);
68 double gmXyDistanceSquared(const Pt3d& pt1, const Pt3d& pt2);
69 bool gmCircumcircleWithTol(const Pt3d* pt1,
70  const Pt3d* pt2,
71  const Pt3d* pt3,
72  double* xc,
73  double* yc,
74  double* r2,
75  double tol);
76 int gmCartToBary(const Pt3d* cart, const Pt3d* orig, double coef[6], int dir, Pt3d* bary);
77 int gmBaryPrepare(const Pt3d* p1,
78  const Pt3d* p2,
79  const Pt3d* p3,
80  const Pt3d* norm,
81  Pt3d* orig,
82  double coef[6],
83  int* dir,
84  bool flag);
85 bool gmColinearWithTol(const Pt3d& p1, const Pt3d& p2, const Pt3d& p3, const double tol);
86 bool gmIntersectLineSegmentsWithTol(const Pt3d& one1,
87  const Pt3d& one2,
88  const Pt3d& two1,
89  const Pt3d& two2,
90  double* xi,
91  double* yi,
92  double* zi1,
93  double* zi2,
94  double tol);
95 bool gmCounterClockwiseTri(const Pt3d& vtx0, const Pt3d& vtx1, const Pt3d& vtx2);
96 double gmCross2D(const double& dx1, const double& dy1, const double& dx2, const double& dy2);
97 double gmCross2D(const Pt3d& a_origin, const Pt3d& a_A, const Pt3d& a_B);
98 double gmAngleBetween2DVectors(double dxp, double dyp, double dxn, double dyn);
99 double gmAngleBetween2DVectors(double dxp,
100  double dyp,
101  double dxn,
102  double dyn,
103  double a_magn,
104  double a_magp);
105 double gmAngleBetweenEdges(const Pt3d&, const Pt3d&, const Pt3d&);
106 double gmAngleBetweenEdges(const Pt2d& p1, const Pt2d& p2, const Pt2d& p3);
107 double gmComputeDeviationInDirection(const Pt3d&, const Pt3d&, const Pt3d&);
108 bool gmOnLineAndBetweenEndpointsWithTol(const Pt3d& a_pt1,
109  const Pt3d& a_pt2,
110  const double a_x,
111  const double a_y,
112  double a_tol);
113 bool gmOnLineWithTol(const Pt3d& vertex1,
114  const Pt3d& vertex2,
115  const double x,
116  const double y,
117  const double tol);
118 void gmAddToExtents(const Pt3d& a_pt, Pt3d& a_min, Pt3d& a_max);
119 void gmAddToExtents(const Pt3d& a_pt, Pt2d& a_min, Pt2d& a_max);
120 void gmAddToExtents(const Pt2d& a_pt, Pt3d& a_min, Pt3d& a_max);
121 double gmXyDistance(double x1, double y1, double x2, double y2);
122 double gmXyDistance(const Pt3d& pt1, const Pt3d& pt2);
123 Turn_enum gmTurn(const Pt3d& a_v1, const Pt3d& a_v2, const Pt3d& a_v3, double a_angtol = 0.0017453);
124 Pt3d gmComputeCentroid(const VecPt3d& a_points);
126 bool gmLinesIntersect(const Pt3d& one1, const Pt3d& one2, const Pt3d& two1, const Pt3d& two2);
127 bool gmLinesCross(const Pt3d& one1, const Pt3d& one2, const Pt3d& two1, const Pt3d& two2);
128 int gmPointInPolygon2D(const Pt3d*, size_t, double, double);
129 int gmPointInPolygon2D(const Pt3d*, size_t, double, double);
130 int gmPointInPolygon2D(const Pt3d*, size_t, Pt3d);
131 int gmPointInPolygon2D(const Pt2i*, size_t, Pt2d);
132 int gmPointInPolygon2D(const Pt2i*, size_t, Pt2i);
133 int gmPointInPolygon2D(const Pt2i*, size_t, Pt3d);
134 int gmPointInPolygon2D(const VecPt3d&, const Pt3d&);
135 template <typename T>
136 int gmPointInPolygon2D(const T* theverts,
137  const size_t numverts,
138  const double xpt,
139  const double ypt,
140  const double tol);
141 double gmComputeXyTol(const Pt3d& a_mn, const Pt3d& a_mx);
142 double gmXyTol(bool a_set = false, double a_value = 1e-9);
143 double gmZTol(bool a_set = false, double a_value = 1e-6);
144 bool gmEqualPointsXY(double x1, double y1, double x2, double y2);
145 bool gmEqualPointsXY(double x1, double y1, double x2, double y2, double tol);
146 bool gmEqualPointsXY(const Pt2i&, const Pt2i&);
147 bool gmEqualPointsXY(const Pt2d& a_pt1, const Pt2d& a_pt2, double tol);
148 bool gmEqualPointsXY(const Pt3d& a_pt1, const Pt3d& a_pt2, double tol);
149 bool gmEqualPointsXYZ(double x1, double y1, double z1, double x2, double y2, double z2);
150 bool gmEqualPointsXYZ(double x1, double y1, double z1, double x2, double y2, double z2, double tol);
151 bool gmEqualPointsXYZ(const Pt3d& pt1, const Pt3d& pt2, double tol);
152 bool gmPointInTriangleWithTol(const Pt3d* p1,
153  const Pt3d* p2,
154  const Pt3d* p3,
155  double x,
156  double y,
157  double tol);
158 bool gmInsideOrOnLineWithTol(const Pt3d* p1,
159  const Pt3d* p2,
160  const Pt3d* inpoint,
161  double x,
162  double y,
163  double tol,
164  double* dist = NULL);
165 double gmPolygonArea(const Pt3d* points, size_t npoints);
166 VecPt3d gmArrayToVecPt3d(double* a_array, int a_size);
167 void gmEnvelopeOfPts(const VecPt3d& a_pts, Pt3d& a_min, Pt3d& a_max);
168 void gmOrderPointsCounterclockwise(const VecPt3d& pts, VecInt& ccwOrder, int startindex = 0);
170 double gmFindClosestPtOnSegment(const Pt3d& pt1,
171  const Pt3d& pt2,
172  const Pt3d& pt,
173  Pt3d& newpt,
174  const double tol);
175 double gmPtDistanceAlongSegment(const Pt3d& pt1, const Pt3d& pt2, const Pt3d& pt, const double tol);
176 bool gmInsideOfLineWithTol(const Pt3d& a_vertex1,
177  const Pt3d& a_vertex2,
178  const Pt3d& a_oppositevertex,
179  const double a_x,
180  const double a_y,
181  const double a_tol);
182 void gmExtents2D(const VecPt3d& a_pts, Pt2d& a_min, Pt2d& a_max);
183 void gmExtents2D(const VecPt3d& a_points, Pt3d& a_min, Pt3d& a_max);
184 void gmExtents3D(const VecPt3d& a_points, Pt3d& a_min, Pt3d& a_max);
185 double gmPerpendicularAngle(const Pt3d& point1, const Pt3d& point2);
186 double gmBisectingAngle(const Pt3d& p1, const Pt3d& p2, const Pt3d& p3);
187 void gmComponentMagnitudes(double* x, double* y, double* mag, double* dir, bool tomagdir);
188 Pt3d gmCreateVector(const Pt3d& p1, const Pt3d& p2);
189 double gmConvertAngleToBetween0And360(double a_angle, bool a_InDegrees = true);
190 void gmCross3D(const Pt3d& a_vec1, const Pt3d& a_vec2, Pt3d* a_vec3);
191 inline double gmDot3D(const Pt3d& a_vec1, const Pt3d& a_vec2);
192 int gmIntersectTriangleAndLineSegment(const Pt3d& a_pt1,
193  const Pt3d& a_pt2,
194  const Pt3d& a_t0,
195  const Pt3d& a_t1,
196  const Pt3d& a_t2,
197  Pt3d& a_IntersectPt);
198 double gm2DDistanceToLineWithTol(const Pt3d* pt1, const Pt3d* pt2, double x, double y, double tol);
199 void gmGetConvexHull(const VecPt3d& a_pts, VecPt3d& a_hull, bool a_includeOn = false);
200 
201 } // namespace xms
double gmMiddleZ(const VecPt3d &a_points)
Returns the z value halfway between the max and min z. Different from the average z (where many point...
Definition: geoms.cpp:153
bool gmOnLineWithTol(const Pt3d &p1, const Pt3d &p2, const double x, const double y, const double tol)
Determines if a point (x,y) is on the line defined by p1 and p2. Assumes p1 and p2 aren&#39;t the same...
Definition: geoms.cpp:732
double gmZTol(bool a_set, double a_value)
Get or set (set first time) global z tolerance for float operations.
Definition: geoms.cpp:1292
bool gmInsideOfLineWithTol(const Pt3d &a_vertex1, const Pt3d &a_vertex2, const Pt3d &a_oppositevertex, const double a_x, const double a_y, const double a_tol)
Returns TRUE if the Point on the same side of the line (defined by vertex1 and vertex2) as oppositeve...
Definition: geoms.cpp:1782
double gmBisectingAngle(const Pt3d &a_p1, const Pt3d &a_p2, const Pt3d &a_p3)
Returns the angle (0-2pi) which bisects the edges p2-p1 and p2-p3 based on a ccw rotation from edge 1...
Definition: geoms.cpp:1901
std::vector< int > VecInt
void gmOrderPointsCounterclockwise(const VecPt3d &a_pts, VecInt &a_ccwOrder, int a_startindex)
Orders array of points counter clockwise. Given non-empty array of points. array of point indices ord...
Definition: geoms.cpp:1635
double gmXyDistance(double x1, double y1, double x2, double y2)
Compute the 2d distance between 2 points.
Definition: geoms.cpp:832
bool gmOnLineAndBetweenEndpointsWithTol(const Pt3d &a_pt1, const Pt3d &a_pt2, const double a_x, const double a_y, double a_tol)
determines if (x,y) is on the line passing through p1 & p2 and between p1 & p2.
Definition: geoms.cpp:767
bool gmPointInTriangleWithTol(const Pt3d *p1, const Pt3d *p2, const Pt3d *p3, double x, double y, double tol)
Returns true if (x,y) is in the tri formed by p1, p2, p3.
Definition: geoms.cpp:1423
PtInOutOrOn_enum gmPtInCircumcircle(const Pt3d &pt, Pt3d circumcirclePts[3])
Is a given point inside a circumcircle defined by three points?
Definition: geoms.cpp:172
double gmPerpendicularAngle(const Pt3d &a_pt1, const Pt3d &a_pt2)
Returns the angle in radians perpendicular to the two points.
Definition: geoms.cpp:1874
bool gmEqualPointsXY(double x1, double y1, double x2, double y2, double tolerance)
Returns true if the points are equal to within gmXyTol().
Definition: geoms.cpp:1308
Turn_enum gmTurn(const Pt3d &a_v1, const Pt3d &a_v2, const Pt3d &a_v3, double a_angtol)
Determine if angle a_v1, a_v2, a_v3 is a left turn, a right turn, colinear 180 degrees, or colinear 0 degrees.
Definition: geoms.cpp:863
double gmAngleBetween2DVectors(double dxp, double dyp, double dxn, double dyn)
Returns the angle (0-2PI) in radians between the edges p and n based on a ccw rotation from p to n...
Definition: geoms.cpp:603
bool gmIntersectLineSegmentsWithTol(const Pt3d &one1, const Pt3d &one2, const Pt3d &two1, const Pt3d &two2, double *xi, double *yi, double *zi1, double *zi2, double tol)
Intersects the plan projection of two line segments.
Definition: geoms.cpp:420
on
Definition: geoms.h:40
void gmGetConvexHull(const VecPt3d &a_pts, VecPt3d &a_hull, bool a_includeOn)
Calculate convex hull using Monotone chain aka Andrew&#39;s algorithm.
Definition: geoms.cpp:2211
turns back on same segment
Definition: geoms.h:32
void gmEnvelopeOfPts(const VecPt3d &a_pts, Pt3d &a_min, Pt3d &a_max)
Calculates the envelope of a vector of points.
Definition: geoms.cpp:1606
PtInOutOrOn_enum
point in, out, or on
Definition: geoms.h:36
Pt3d gmCreateVector(const Pt3d &a_p1, const Pt3d &a_p2)
creates a vector from a_p1 to a_p2
Definition: geoms.cpp:1982
bool gmInsideOrOnLineWithTol(const Pt3d *p1, const Pt3d *p2, const Pt3d *inpoint, const double x, const double y, const double tol, double *dist)
Returns true if the (x,y) is on the line segment (p1,p2) or on the same side of the line as "inpoint"...
Definition: geoms.cpp:1450
bool gmColinearWithTol(const Pt3d &p1, const Pt3d &p2, const Pt3d &p3, const double tol)
Returns true if (the three vertices are gmColinear. Result should be insensitive to the order of the ...
Definition: geoms.cpp:393
double gmPtDistanceAlongSegment(const Pt3d &a_pt1, const Pt3d &a_pt2, const Pt3d &a_pt, const double a_tol)
Finds the distance along a segment for the location closest to a_pt.
Definition: geoms.cpp:1749
turn right
Definition: geoms.h:30
double gmAngleBetweenEdges(const Pt3d &p1, const Pt3d &p2, const Pt3d &p3)
Returns the ccw angle (0-2pi) between p2-p1 and p2-p3.
Definition: geoms.cpp:655
VecPt3d gmArrayToVecPt3d(double *a_array, int a_size)
Useful in testing to create a VecPt3d from a C array of xy pairs.
Definition: geoms.cpp:1590
double gmPolygonArea(const Pt3d *pts, size_t npoints)
Compute 2d planview projection of area of polygon.
Definition: geoms.cpp:1537
void gmCalculateNormalizedPlaneCoefficients(const Pt3d &p1, const Pt3d &p2, const Pt3d &p3, double *a, double *b, double *c, double *d)
Calculates the plane coefficients for a triangle. Given point references calculate coefficents for pl...
Definition: geoms.cpp:105
bool gmLinesIntersect(const Pt3d &one1, const Pt3d &one2, const Pt3d &two1, const Pt3d &two2)
intersects the plan projection of two line segments (bad func name) segment 1 = one1,one2 = one1 + lambda(one2 - one1) segment 2 = two1,two2 = two1 + mu (two2 - two1)
Definition: geoms.cpp:966
void gmComponentMagnitudes(double *a_x, double *a_y, double *a_mag, double *a_dir, bool a_tomagdir)
converts the magnitude and angle to xy components or vice versa
Definition: geoms.cpp:1941
double gmComputeDeviationInDirection(const Pt3d &a_p0, const Pt3d &a_p1, const Pt3d &a_p2)
Computes the deviation in direction from one seg to next seg 1 is from a_p0 to a_p1 seg 2 is from a_p...
Definition: geoms.cpp:697
void gmCross3D(const Pt3d &a_vec1, const Pt3d &a_vec2, Pt3d *a_vec3)
Perform a cross product of Pt3d&#39;s.
Definition: geoms.cpp:2034
double gmComputeXyTol(const Pt3d &a_mn, const Pt3d &a_mx)
Given extents (min, max), compute a tolerance for the xy plane to be used with geometric functions...
Definition: geoms.cpp:1262
Pt3d gmComputeCentroid(const VecPt3d &a_points)
Computes the centroid of points (but not of a polygon). Shouldn&#39;t pass an empty vector (no checks are...
Definition: geoms.cpp:895
double gm2DDistanceToLineWithTol(const Pt3d *a_pt1, const Pt3d *a_pt2, double a_x, double a_y, double a_tol)
return the xy distance from (a_x,a_y) to line through (a_pt1, a_pt2)
Definition: geoms.cpp:2177
double gmXyTol(bool a_set, double a_value)
Get or set (set first time) global xy tolerance for float operations.
Definition: geoms.cpp:1279
double gmConvertAngleToBetween0And360(double a_angle, bool a_InDegrees)
given an angle, this function will return the corresponding angle that matches it in the range of 0 d...
Definition: geoms.cpp:1997
bool gmCircumcircleWithTol(const Pt3d *pt1, const Pt3d *pt2, const Pt3d *pt3, double *xc, double *yc, double *r2, double tol)
Computes center & radius squared for circumcircle of triangle made by the three points.
Definition: geoms.cpp:218
double gmCross2D(const double &dx1, const double &dy1, const double &dx2, const double &dy2)
Returns the cross product of two 2-d vectors.
Definition: geoms.cpp:579
error
Definition: geoms.h:37
double gmFindClosestPtOnSegment(const Pt3d &a_pt1, const Pt3d &a_pt2, const Pt3d &a_pt, Pt3d &a_newpt, const double a_tol)
Finds the closest point to another point on a segment.
Definition: geoms.cpp:1714
bool gmEqualPointsXYZ(double x1, double y1, double z1, double x2, double y2, double z2, double tolerance)
Returns true if the points are equal to within tolerance.
Definition: geoms.cpp:1376
XMS Namespace.
Definition: geoms.cpp:34
Pt2< int > Pt2i
out
Definition: geoms.h:39
turn left
Definition: geoms.h:29
bool gmLinesCross(const Pt3d &a_segment1Point1, const Pt3d &a_segment1Point2, const Pt3d &a_segment2Point1, const Pt3d &a_segment2Point2)
Determine whether 2 line segments cross. The segments may touch at the end points.
Definition: geoms.cpp:1013
int gmIntersectTriangleAndLineSegment(const Pt3d &a_pt1, const Pt3d &a_pt2, const Pt3d &a_t0, const Pt3d &a_t1, const Pt3d &a_t2, Pt3d &a_IntersectPt)
Determine if the line (described by a_pt1 and a_pt2) intersects the triangle (a_t0, a_t1, a_t2).
Definition: geoms.cpp:2075
Turn_enum
direction of turn between 3 points
Definition: geoms.h:28
bool gmPointInOrOnBox2d(const Pt3d &a_bMin, const Pt3d &a_bMax, const Pt3d &a_pt)
finds if a point is in or on a box in 2d
Definition: geoms.cpp:65
void gmAddToExtents(const Pt3d &a_pt, Pt3d &a_min, Pt3d &a_max)
Compares a_pt to a_min and a_max. If a_pt is less than a_min or greater than a_max, a_min and a_max are updated.
Definition: geoms.cpp:787
bool gmCounterClockwiseTri(const Pt3d &vtx0, const Pt3d &vtx1, const Pt3d &vtx2)
Returns true if the triangle is wrapped counter clockwise.
Definition: geoms.cpp:564
double gmDot3D(const Pt3d &a_vec1, const Pt3d &a_vec2)
Perform a dot product of Pt3d&#39;s.
Definition: geoms.cpp:2048
bool gmBoxesOverlap2d(const Pt3d &a_b1Min, const Pt3d &a_b1Max, const Pt3d &a_b2Min, const Pt3d &a_b2Max)
finds if 2 boxes overlap in 2d
Definition: geoms.cpp:79
continue onward
Definition: geoms.h:31
void gmExtents2D(const VecPt3d &a_points, Pt2d &a_min, Pt2d &a_max)
Get the 2D extents of a vector of points.
Definition: geoms.cpp:1821
int gmPointInPolygon2D(const T *a_verts, const size_t a_num, const double a_x, const double a_y, const double a_tol)
Determines whether (a_x, a_y) is inside=1, on=0, or outside=-1 the polygon defined by the given verti...
Definition: geoms.cpp:1057
int gmCartToBary(const Pt3d *cart, const Pt3d *orig, double coef[6], int dir, Pt3d *bary)
Compute Barycentric coords for point. Use gmBaryPrepare to get the coefficients and direction...
Definition: geoms.cpp:275
int gmBaryPrepare(const Pt3d *p1, const Pt3d *p2, const Pt3d *p3, const Pt3d *norm, Pt3d *orig, double coef[6], int *dir, bool flag)
For a tri - compute dir & coefs to compute Barycentric coords.
Definition: geoms.cpp:315
void gmExtents3D(const VecPt3d &a_points, Pt3d &a_min, Pt3d &a_max)
Get the 3D extents of a vector of points.
Definition: geoms.cpp:1856
in
Definition: geoms.h:38
Pt3d gmComputePolygonCentroid(const VecPt3d &pts)
Computes the plan view centroid of a non-self-intersecting polygon.
Definition: geoms.cpp:908
double gmXyDistanceSquared(const Pt3d &pt1, const Pt3d &pt2)
Calculate XY distance between two 3D points.
Definition: geoms.cpp:202
std::vector< Pt3d > VecPt3d