xmsgrid
1.0
|
XMS Namespace. More...
Namespaces |
Classes | |
struct | BarycentricVals |
Structure for Barycentric coordinate calculations. More... | |
class | DaStreamReader |
class | DaStreamWriter |
struct | edgerecord |
Defines an edge that intersects a breakline. More... | |
class | ETestMessagingState |
class | fSatisfies |
a class used with the boost::geometry::index::satisfies function More... | |
class | GmExtents2d |
2D geometric extents (min/max). More... | |
class | GmExtents3d |
3D geometric extents (min/max). More... | |
class | GmMultiPolyIntersectionSorter |
Base class for sorting intersections from GmMultiPolyIntersector. More... | |
class | GmMultiPolyIntersectionSorterTerse |
Class for sorting intersections from GmMultiPolyIntersector in a terse way (with no duplicate info). More... | |
class | GmMultiPolyIntersector |
See GmMultiPolyIntersectorImpl comments. More... | |
struct | GmMultiPolyIntersectorData |
Struct used by GmMultiPolyIntersector. More... | |
class | GmMultiPolyIntersectorImpl |
Used to intersect a line segment with any number of polygons in 2D. Returns the intersected polygons in order along with t values. May be an alternative to SurfaceIter and feiTraverseLineSegment2. More... | |
class | GmPolygon |
Interface to a boost::geometry::polygon class. More... | |
class | GmPolygonImpl |
Wraps a boost polygon for point in poly, min distance from point to poly etc. More... | |
class | GmPolyLinePtRedistributer |
Redistributes the point locations on a polyline based on a size. More... | |
class | GmPolyLinePtRedistributerImpl |
Implementation of GmPolyLinePtRedistributer. More... | |
class | GmPtSearch |
Spatial index for searching points. More... | |
class | GmPtSearchImpl |
Implementation of GmPtSearch. Generic class for searching location data. Uses a boost R-tree to query a set of locations. Works for 2D and 3D. More... | |
class | GmTriSearch |
Spatial index for searching triangles. More... | |
class | GmTriSearchImpl |
Implementation of GmTriSearch. More... | |
class | idx_pt |
class for indexing the points More... | |
class | idx_tri |
provides indexing for the spatial index More... | |
class | ix |
An intersection point of a line with a polygon. More... | |
class | ltPt2 |
class | ltPt3 |
class | ltPt3_2D |
class | Observer |
class | Progress |
class | ProgressListener |
class | Pt2 |
class | Pt3 |
class | Pt4 |
class | SharedSingleton |
class | Singleton |
class | StCommaNumpunct |
class | StTemp2DigitExponents |
class | TrAutoFixFourTrianglePts |
class | TrAutoFixFourTrianglePtsImpl |
Used to delete points that are connected to 4 triangles and then retriangulate the void. More... | |
class | TrBreaklineAdder |
Adds breaklines to a tin. More... | |
class | TrBreaklineAdderImpl |
Adds breaklines to a triangulation. More... | |
class | TrOuterTriangleDeleter |
class | TrOuterTriangleDeleterImpl |
Used to delete tin triangles that are outside given polygons. The polygons may include holes - polygons inside polygons. More... | |
class | TrTin |
class | TrTinImpl |
Class to encapsulate a tin made simply of arrays of points, triangles and adjacency information. Also methods to manipulate it. More... | |
class | TrTriangulator |
Base class used to derive a class to triangulate points. More... | |
class | TrTriangulatorPoints |
Class to triangulate simple points. More... | |
class | XmEdge |
Two integer values representing an edge of an XmUGrid. By default has a direction. Can be sorted to have minimum index first. More... | |
class | XmLog |
class | XmUGrid |
Geometry for an unstructured grid. An XmUGrid is defined as a vector of 3d points and a stream of cells. Throughout this interface, we use these terms: More... | |
Typedefs | |
typedef boost::geometry::model::d2::point_xy< double > | GmBstPt2d |
Boost geometry 2d point. | |
typedef boost::geometry::model::polygon< GmBstPt2d > | GmBstPoly2d |
Boost polygon 2d. | |
typedef boost::geometry::model::polygon< Pt3d > | GmBstPoly3d |
Boost polygon 3d. | |
typedef boost::geometry::model::box< Pt3d > | GmBstBox3d |
Boost box 3d. | |
typedef boost::geometry::model::linestring< Pt3d > | GmBstLine3d |
Boost line 3d. | |
typedef boost::geometry::ring_type< GmBstPoly3d >::type | GmBstRing3d |
Boost ring 3d. | |
typedef std::pair< GmBstBox3d, int > | ValueBox |
Pair used in rtree. | |
typedef bgi::rtree< ValueBox, bgi::quadratic< 8 > > | RtreeBox |
Rtree typedef. | |
typedef bg::model::box< bPt > | box |
box | |
typedef size_t | value |
value | |
typedef bgi::quadratic< 8 > | qRtree |
qRtree | |
typedef std::vector< edgerecord > | VecEdge |
Vector of edgerecord. | |
Enumerations | |
enum | Turn_enum { TURN_LEFT, TURN_RIGHT, TURN_COLINEAR_180, TURN_COLINEAR_0 } |
direction of turn between 3 points More... | |
enum | PtInOutOrOn_enum { PT_ERROR = -1, PT_IN, PT_OUT, PT_ON } |
point in, out, or on More... | |
enum | GmMultiPolyIntersectorQueryEnum { GMMPIQ_COVEREDBY, GMMPIQ_INTERSECTS } |
Type of query. | |
enum | XmUGridCellType { XMU_INVALID_CELL_TYPE = -1, XMU_EMPTY_CELL = 0, XMU_VERTEX = 1, XMU_POLY_VERTEX = 2, XMU_LINE = 3, XMU_POLY_LINE = 4, XMU_TRIANGLE = 5, XMU_TRIANGLE_STRIP = 6, XMU_POLYGON = 7, XMU_PIXEL = 8, XMU_QUAD = 9, XMU_TETRA = 10, XMU_VOXEL = 11, XMU_HEXAHEDRON = 12, XMU_WEDGE = 13, XMU_PYRAMID = 14, XMU_PENTAGONAL_PRISM = 15, XMU_HEXAGONAL_PRISM = 16, XMU_QUADRATIC_EDGE = 21, XMU_QUADRATIC_TRIANGLE = 22, XMU_QUADRATIC_QUAD = 23, XMU_QUADRATIC_POLYGON = 36, XMU_QUADRATIC_TETRA = 24, XMU_QUADRATIC_HEXAHEDRON = 25, XMU_QUADRATIC_WEDGE = 26, XMU_QUADRATIC_PYRAMID = 27, XMU_BIQUADRATIC_QUAD = 28, XMU_TRIQUADRATIC_HEXAHEDRON = 29, XMU_QUADRATIC_LINEAR_QUAD = 30, XMU_QUADRATIC_LINEAR_WEDGE = 31, XMU_BIQUADRATIC_QUADRATIC_WEDGE = 32, XMU_BIQUADRATIC_QUADRATIC_HEXAHEDRON = 33, XMU_BIQUADRATIC_TRIANGLE = 34, XMU_CUBIC_LINE = 35, XMU_CONVEX_POINT_SET = 41, XMU_POLYHEDRON = 42, XMU_PARAMETRIC_CURVE = 51, XMU_PARAMETRIC_SURFACE = 52, XMU_PARAMETRIC_TRI_SURFACE = 53, XMU_PARAMETRIC_QUAD_SURFACE = 54, XMU_PARAMETRIC_TETRA_REGION = 55, XMU_PARAMETRIC_HEX_REGION = 56, XMU_HIGHER_ORDER_EDGE = 60, XMU_HIGHER_ORDER_TRIANGLE = 61, XMU_HIGHER_ORDER_QUAD = 62, XMU_HIGHER_ORDER_POLYGON = 63, XMU_HIGHER_ORDER_TETRAHEDRON = 64, XMU_HIGHER_ORDER_WEDGE = 65, XMU_HIGHER_ORDER_PYRAMID = 66, XMU_HIGHER_ORDER_HEXAHEDRON = 67, XMU_NUMBER_OF_CELL_TYPES } |
Matches cell types from VTK (see vtkCellType.h) | |
enum | XmUGridFaceOrientation { XMU_ORIENTATION_UNKNOWN = -1, XMU_ORIENTATION_SIDE = 0, XMU_ORIENTATION_TOP = 1, XMU_ORIENTATION_BOTTOM = 2, XMU_ORIENTATION_NUMBER } |
The orientation of a 3D face must be one of these. | |
Functions | |
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 More... | |
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 More... | |
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 plane (ax+by+cz+d=0). More... | |
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 pointers calculate coefficents for plane (ax+by+cz+d=0). More... | |
double | gmMiddleZ (const VecPt3d &a_points) |
Returns the z value halfway between the max and min z. Different from the average z (where many points can skew the average). More... | |
PtInOutOrOn_enum | gmPtInCircumcircle (const Pt3d &pt, Pt3d circumcirclePts[3]) |
Is a given point inside a circumcircle defined by three points? More... | |
double | gmXyDistanceSquared (const Pt3d &pt1, const Pt3d &pt2) |
Calculate XY distance between two 3D points. More... | |
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. More... | |
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. More... | |
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. More... | |
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 vertices. More... | |
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. More... | |
bool | gmCounterClockwiseTri (const Pt3d &vtx0, const Pt3d &vtx1, const Pt3d &vtx2) |
Returns true if the triangle is wrapped counter clockwise. More... | |
double | gmCross2D (const double &dx1, const double &dy1, const double &dx2, const double &dy2) |
Returns the cross product of two 2-d vectors. More... | |
double | gmCross2D (const Pt3d &a_origin, const Pt3d &a_A, const Pt3d &a_B) |
2D cross product of two points More... | |
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. More... | |
double | gmAngleBetween2DVectors (double dxp, double dyp, double dxn, double dyn, double a_magn, double a_magp) |
Returns the angle (0-2PI) in radians between the edges p and n based on a ccw rotation from p to n. More... | |
double | gmAngleBetweenEdges (const Pt3d &p1, const Pt3d &p2, const Pt3d &p3) |
Returns the ccw angle (0-2pi) between p2-p1 and p2-p3. More... | |
double | gmAngleBetweenEdges (const Pt2d &p1, const Pt2d &p2, const Pt2d &p3) |
Returns the ccw angle (0-2pi) between p2-p1 and p2-p3. More... | |
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_p1 to a_p2. More... | |
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't the same. More... | |
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. More... | |
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. More... | |
void | gmAddToExtents (const Pt3d &a_pt, Pt2d &a_min, Pt2d &a_max) |
void | gmAddToExtents (const Pt2d &a_pt, Pt3d &a_min, Pt3d &a_max) |
double | gmXyDistance (double x1, double y1, double x2, double y2) |
Compute the 2d distance between 2 points. More... | |
double | gmXyDistance (const Pt3d &pt1, const Pt3d &pt2) |
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. More... | |
Pt3d | gmComputeCentroid (const VecPt3d &a_points) |
Computes the centroid of points (but not of a polygon). Shouldn't pass an empty vector (no checks are performed). More... | |
Pt3d | gmComputePolygonCentroid (const VecPt3d &pts) |
Computes the plan view centroid of a non-self-intersecting polygon. More... | |
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) More... | |
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. More... | |
template<typename T > | |
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 vertices. More... | |
int | gmPointInPolygon2D (const Pt3d *a_verts, size_t a_num, double a_x, double a_y) |
int | gmPointInPolygon2D (const Pt3d *a_verts, size_t a_num, Pt3d a_pt) |
int | gmPointInPolygon2D (const Pt2i *a_verts, size_t a_num, Pt2d a_pt) |
int | gmPointInPolygon2D (const Pt2i *a_verts, size_t a_num, Pt3d a_pt) |
int | gmPointInPolygon2D (const Pt2i *a_verts, size_t a_num, Pt2i a_pt) |
int | gmPointInPolygon2D (const VecPt3d &a_verts, const Pt3d &a_pt) |
template int | gmPointInPolygon2D< Pt3d > (const Pt3d *a_verts, const size_t a_num, const double a_x, const double a_y, const double a_tol) |
template int | gmPointInPolygon2D< Pt2d > (const Pt2d *a_verts, const size_t a_num, const double a_x, const double a_y, const double a_tol) |
template int | gmPointInPolygon2D< Pt2i > (const Pt2i *a_verts, const size_t a_num, const double a_x, const double a_y, const double a_tol) |
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. More... | |
double | gmXyTol (bool a_set, double a_value) |
Get or set (set first time) global xy tolerance for float operations. More... | |
double | gmZTol (bool a_set, double a_value) |
Get or set (set first time) global z tolerance for float operations. More... | |
bool | gmEqualPointsXY (double x1, double y1, double x2, double y2, double tolerance) |
Returns true if the points are equal to within gmXyTol(). More... | |
bool | gmEqualPointsXY (double x1, double y1, double x2, double y2) |
bool | gmEqualPointsXY (const Pt2d &a_pt1, const Pt2d &a_pt2, double tol) |
bool | gmEqualPointsXY (const Pt3d &a_pt1, const Pt3d &a_pt2, double tol) |
bool | gmEqualPointsXY (const Pt2i &point1, const Pt2i &point2) |
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. More... | |
bool | gmEqualPointsXYZ (double x1, double y1, double z1, double x2, double y2, double z2) |
bool | gmEqualPointsXYZ (const Pt3d &pt1, const Pt3d &pt2, double tol) |
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. More... | |
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". ASSERTs in debug if "inpoint" is on the line (within tol). More... | |
double | gmPolygonArea (const Pt3d *pts, size_t npoints) |
Compute 2d planview projection of area of polygon. More... | |
VecPt3d | gmArrayToVecPt3d (double *a_array, int a_size) |
Useful in testing to create a VecPt3d from a C array of xy pairs. More... | |
void | gmEnvelopeOfPts (const VecPt3d &a_pts, Pt3d &a_min, Pt3d &a_max) |
Calculates the envelope of a vector of points. More... | |
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 ordered counter clockwise based on the angle from the centroid where angle starts at point at startindex. More... | |
void | gmOrderPointsCounterclockwise (VecPt3d &a_pts) |
Overload to gmOrderPointsCounterclockwise. More... | |
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. More... | |
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. More... | |
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 oppositevertex. More... | |
void | gmExtents2D (const VecPt3d &a_points, Pt2d &a_min, Pt2d &a_max) |
Get the 2D extents of a vector of points. More... | |
void | gmExtents2D (const VecPt3d &a_points, Pt3d &a_min, Pt3d &a_max) |
void | gmExtents3D (const VecPt3d &a_points, Pt3d &a_min, Pt3d &a_max) |
Get the 3D extents of a vector of points. More... | |
double | gmPerpendicularAngle (const Pt3d &a_pt1, const Pt3d &a_pt2) |
Returns the angle in radians perpendicular to the two points. More... | |
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 to edge 2. More... | |
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 More... | |
Pt3d | gmCreateVector (const Pt3d &a_p1, const Pt3d &a_p2) |
creates a vector from a_p1 to a_p2 More... | |
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 deg to 360 deg More... | |
void | gmCross3D (const Pt3d &a_vec1, const Pt3d &a_vec2, Pt3d *a_vec3) |
Perform a cross product of Pt3d's. More... | |
double | gmDot3D (const Pt3d &a_vec1, const Pt3d &a_vec2) |
Perform a dot product of Pt3d's. More... | |
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). More... | |
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) More... | |
void | gmGetConvexHull (const VecPt3d &a_pts, VecPt3d &a_hull, bool a_includeOn) |
Calculate convex hull using Monotone chain aka Andrew's algorithm. More... | |
static double | iDistanceToRing (const GmBstRing3d &a_ring, const Pt3d &a_pt) |
Return the minimum distance from the point to a line on the ring. More... | |
void | iCartToBary (const Pt3d &a_pt1, const Pt3d &a_pt2, const Pt3d &a_pt3, BarycentricVals &a_b) |
Converts from Cartesian coordinates to Barycentric coordinates (local triangle coords). More... | |
void | iGetBarycentricCoords (const Pt3d &a_p, BarycentricVals &a_b, Pt3d &weights) |
Calculates weights at the 3 points of a triangle given the location "a_p" and the Barycentric coords of the triangle. More... | |
int | mxLUDecomp (double **mat, int n, int *indx, double *d) |
Decompose an n x n matrix into Upper & Lower trianglar (in place), from "Numerical Recipes in C" p 43. More... | |
int | mxLUBcksub (double **mat, int n, const int *indx, double *b) |
solve [mat]{x} = {b} by back substitution (mat = LU of mat), from "Numerical Recipes in C" p 44. More... | |
bool | mxiLudcmp3 (double mat[3][3], int indx[3], double *d) |
Decomposes a 3 x 3 matrix into Upper & Lower trianglar(in place) More... | |
void | mxiLubksb3 (const double mat[3][3], const int indx[3], double b[3]) |
solve [mat]{x} = {b} by back substitution (mat = LU of mat) More... | |
bool | mxSolveNxN (double **A, double *x, double *b, int n) |
Solves the matrix equation [A]*{x}={b} for {x} where [A] is NxN and {x} and {b} are Nx1 vectors. The matrix A is replaced with the LU decomposition of A. More... | |
bool | mxSolveBandedEquations (double **a, double *x, double *b, int numeqs, int halfbandwidth) |
Solves the matrix equation [a][x]=[b] using gauss elimination. More... | |
bool | mxSolve3x3 (double A[3][3], double x[3], double b[3]) |
Solves the matrix equation [A]*{x}={b} for {x} where [A] is 3x3 and {x} and {b} are 3x1 vectors. More... | |
int | mxInvert4x4 (const double matrix[4][4], double inv[4][4]) |
Compute the inverse of a transformation matrix. Checks for special case of Orthogonal 4x4 and does easy inverse otherwise uses Gauss-Jordan Elimination with full pivoting and scaling with partial unraveling of loops for speed. See Numerical Recipes in C pages 36-37. More... | |
void | mxPointMult (Pt3d *pt, const double ctm[4][4]) |
Multiplies a point by a transformation matrix. More... | |
void | mxCopy4x4 (double copy[4][4], const double matrix[4][4]) |
Copy a transformation matrix. Simply uses memcpy. More... | |
void | mxIdent4x4 (double matrix[4][4]) |
Sets a transformation matrix to identity. More... | |
void | mxMult4x4 (double product[4][4], const double matrix1[4][4], const double matrix2[4][4]) |
Multiplies two transformation matrices. More... | |
void | mxRotate4x4 (double xrot, double yrot, double zrot, double rMatt[4][4]) |
return a rotation matrix with specified x, y, z rotations z - y - x order More... | |
void | mxTranslate4x4 (const Pt3d &a_translation, double a_mx[4][4]) |
return a translation matrix with specified translation z - y - x order More... | |
bool | trTriangulateIt (TrTriangulator &a_Client) |
Do delaunay divide-and-conquer triangulation. More... | |
double | trArea (const Pt3d &a_pt1, const Pt3d &a_pt2, const Pt3d &a_pt3) |
Return the signed planar area of the triangle (CCW Positive). More... | |
void | trBuildGridTrianglePolys (int rows, int cols, VecPt3d &a_points, VecInt2d &a_polys) |
Create something like this: More... | |
void | trBuildGridPolys (int rows, int cols, VecPt3d &pts, VecInt2d &polys) |
Build something like this: More... | |
int | trIncrementIndex (int i) |
Faster than a % operation and we do this a lot. More... | |
int | trDecrementIndex (int i) |
Faster than a % operation and we do this a lot. More... | |
void | trRenumberOnDelete (const SetInt &a_deleting, VecInt &a_vec) |
Boost serialize function. More... | |
VecPt3d | ConvexHull (const std::vector< Pt3< double >> &a_points) |
Returns the convex hull of a set of points (Deprecated). More... | |
VecInt | ConvexHullWithIndices (const std::vector< int > &a_points, std::shared_ptr< XmUGrid > a_ugrid) |
Returns the convex hull of a set of points (Deprecated). More... | |
bool | DoLineSegmentsCross (const std::pair< Pt3< double >, Pt3< double >> &a_segment1, const std::pair< Pt3< double >, Pt3< double >> &a_segment2) |
Determine whether 2 line segments intersect (Deprecated). More... | |
bool | DoLineSegmentsCross (const Pt3< double > &a_segment1Point1, const Pt3< double > &a_segment1Point2, const Pt3< double > &a_segment2Point1, const Pt3< double > &a_segment2Point2) |
Determine whether 2 line segments cross (Deprecated). More... | |
bool | XmEdgesEquivalent (const XmEdge &a_edge1, const XmEdge &a_edge2) |
Test if two edges are the same ignoring direction. More... | |
std::shared_ptr< XmUGrid > | TEST_XmUGrid1Left90Tri () |
Builds a 1 cell (left 90 degree triangle) 2D XmUGrid for testing. More... | |
std::shared_ptr< XmUGrid > | TEST_XmUGridSimpleQuad () |
Builds a 2 cell (Quad) 2D XmUGrid for testing. More... | |
std::shared_ptr< XmUGrid > | TEST_XmUGrid2dLinear () |
Builds an XmUGrid with supported 1D and 2D linear cells for testing. More... | |
std::shared_ptr< XmUGrid > | TEST_XmUGrid3dLinear () |
Builds an XmUGrid with supported 3D linear cells for testing. More... | |
std::shared_ptr< XmUGrid > | TEST_XmUGridHexagonalPolyhedron () |
Builds a 1 cell hexagon with cell type polyhedron. More... | |
std::shared_ptr< xms::XmUGrid > | TEST_XmUBuildQuadUGrid (int a_rows, int a_cols) |
Builds a UGrid of Quads at 1 spacing for rows & cols specified. More... | |
std::shared_ptr< xms::XmUGrid > | TEST_XmUBuildQuadUGrid (int a_rows, int a_cols, const xms::Pt3d &a_origin) |
Builds a UGrid of Quads at 1 spacing for rows & cols specified. More... | |
std::shared_ptr< xms::XmUGrid > | TEST_XmUBuildHexahedronUgrid (int a_rows, int a_cols, int a_lays) |
Builds a UGrid of Quads at 1 spacing for rows & cols specified. More... | |
std::shared_ptr< xms::XmUGrid > | TEST_XmUBuildHexahedronUgrid (int a_rows, int a_cols, int a_lays, const xms::Pt3d &a_origin) |
Builds a UGrid of Quads at 1 spacing for rows & cols specified. More... | |
std::shared_ptr< xms::XmUGrid > | TEST_XmUBuildPolyhedronUgrid (int a_rows, int a_cols, int a_lays) |
Builds a UGrid of Quads at 1 spacing for rows & cols specified. More... | |
std::shared_ptr< xms::XmUGrid > | TEST_XmUBuildPolyhedronUgrid (int a_rows, int a_cols, int a_lays, const xms::Pt3d &a_origin) |
Builds a UGrid of Quads at 1 spacing for rows & cols specified. More... | |
std::shared_ptr< xms::XmUGrid > | TEST_XmUBuild3DChevronUgrid () |
Builds a UGrid with one 3D Chevron polyhedron. More... | |
std::shared_ptr< XmUGrid > | XmReadUGridFromAsciiFile (const std::string &a_filePath) |
Read XmUGrid from an ASCII file. More... | |
std::shared_ptr< XmUGrid > | XmReadUGridFromStream (std::istream &a_inStream) |
Read an XmUGrid from ASCII text from an input stream. More... | |
void | XmWriteUGridToAsciiFile (std::shared_ptr< XmUGrid > a_ugrid, const std::string &a_filePath) |
Write an XmUGrid to an ASCII file. More... | |
void | XmWriteUGridToStream (std::shared_ptr< XmUGrid > a_ugrid, std::ostream &a_outStream) |
Save an XmUGrid ASCII text to output stream. More... | |
void | XmWriteUGridToStream (const XmUGrid &a_ugrid, std::ostream &a_outStream, bool a_binary) |
Save an XmUGrid ASCII text to output stream. More... | |
XMS Namespace.
enum xms::Turn_enum |
Returns the convex hull of a set of points (Deprecated).
[in] | a_points | The list of points |
Definition at line 55 of file XmGeometry.cpp.
References gmGetConvexHull().
Referenced by ConvexHullWithIndices().
VecInt xms::ConvexHullWithIndices | ( | const std::vector< int > & | a_points, |
std::shared_ptr< XmUGrid > | a_ugrid | ||
) |
Returns the convex hull of a set of points (Deprecated).
[in] | a_points | The list of indices to points |
[in] | a_ugrid | The ugrid that contains the points. |
Definition at line 70 of file XmGeometry.cpp.
References ConvexHull().
bool xms::DoLineSegmentsCross | ( | const std::pair< Pt3< double >, Pt3< double >> & | a_segment1, |
const std::pair< Pt3< double >, Pt3< double >> & | a_segment2 | ||
) |
Determine whether 2 line segments intersect (Deprecated).
[in] | a_segment1 | The first line segment |
[in] | a_segment2 | The second line segment |
Definition at line 99 of file XmGeometry.cpp.
References gmLinesCross().
bool xms::DoLineSegmentsCross | ( | const Pt3< double > & | a_segment1Point1, |
const Pt3< double > & | a_segment1Point2, | ||
const Pt3< double > & | a_segment2Point1, | ||
const Pt3< double > & | a_segment2Point2 | ||
) |
Determine whether 2 line segments cross (Deprecated).
[in] | a_segment1Point1 | First point 3d of line segment 1 |
[in] | a_segment1Point2 | Second point 3d of line segment 1 |
[in] | a_segment2Point1 | First point 3d of line segment 2 |
[in] | a_segment2Point2 | Second point 3d of line segment 2 |
Definition at line 115 of file XmGeometry.cpp.
References gmLinesCross().
double xms::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)
[in] | a_pt1 | first point of line segment |
[in] | a_pt2 | second point of line segment |
[in] | a_x | x location of point used to compute distance to line |
[in] | a_y | y location of point used to compute distance to line |
[in] | a_tol | tolerance used in geometric computations |
Definition at line 2177 of file geoms.cpp.
References sqr().
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.
[in] | a_pt | A point. |
[in,out] | a_min | Minimum. |
[in,out] | a_max | Maximum. |
Definition at line 787 of file geoms.cpp.
Referenced by xms::GmMultiPolyIntersectorImpl::BuildRTree(), xms::GmMultiPolyIntersectorImpl::CalculateBuffer(), xms::TrTinImpl::GetExtents(), xms::XmUGrid::Impl::GetExtentsFromPoints(), gmExtents2D(), gmExtents3D(), and xms::TrTriangulatorPoints::UpdateAreaTolerance().
double xms::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.
[in] | dxp | x of vector p. |
[in] | dyp | y of vector p. |
[in] | dxn | x of vector n. |
[in] | dyn | y of vector n. |
Definition at line 603 of file geoms.cpp.
References sqr().
Referenced by gmAngleBetweenEdges().
double xms::gmAngleBetween2DVectors | ( | double | dxp, |
double | dyp, | ||
double | dxn, | ||
double | dyn, | ||
double | a_magn, | ||
double | a_magp | ||
) |
Returns the angle (0-2PI) in radians between the edges p and n based on a ccw rotation from p to n.
[in] | dxp | x of vector p. |
[in] | dyp | y of vector p. |
[in] | dxn | x of vector n. |
[in] | dyn | y of vector n. |
[in] | a_magn | Magnitude of n. |
[in] | a_magp | Magnitude of p. |
Definition at line 622 of file geoms.cpp.
References gmCross2D(), and XM_PI.
Returns the ccw angle (0-2pi) between p2-p1 and p2-p3.
[in] | p1 | Point 1. |
[in] | p2 | Point 2. |
[in] | p3 | Point 3. |
Definition at line 655 of file geoms.cpp.
References gmAngleBetween2DVectors().
Referenced by xms::TrTinImpl::SwapEdge().
Returns the ccw angle (0-2pi) between p2-p1 and p2-p3.
[in] | p1 | Point 1. |
[in] | p2 | Point 2. |
[in] | p3 | Point 3. |
Definition at line 672 of file geoms.cpp.
References gmAngleBetween2DVectors().
VecPt3d xms::gmArrayToVecPt3d | ( | double * | a_array, |
int | a_size | ||
) |
Useful in testing to create a VecPt3d from a C array of xy pairs.
[in] | a_array | Array of xy pairs ([x][y][x][y]...). |
[in] | a_size | Array size. |
Definition at line 1590 of file geoms.cpp.
Referenced by GmPolygonUnitTests::SetUpPolyWithHole().
int xms::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.
p1 | = 1st point of triangle. |
p2 | = 2nd point of triangle. |
p3 | = 3rd point of triangle. |
norm | = triangle normal. |
orig | = triangle origin (min xyz of p1, p2 and p3). |
coef | = ? |
dir | = 0 if norm points to x, 1 if to y, 2 if to z. |
flag | = if true, origin will be computed from the triangle. If false, just calculate the coefficients and the origin is zero. |
Definition at line 315 of file geoms.cpp.
References Mmin(), and XM_SUCCESS.
Referenced by iCartToBary().
Returns the angle (0-2pi) which bisects the edges p2-p1 and p2-p3 based on a ccw rotation from edge 1 to edge 2.
[in] | a_p1 | The first point defining the edge (a_p1:a_p2) |
[in] | a_p2 | The second point defining two edges (a_p1:a_p2, a_p2:a_p3) |
[in] | a_p3 | The third point defining the edge (a_p2:a_p3) |
Definition at line 1901 of file geoms.cpp.
References gmCross2D(), sqr(), and XM_PI.
void xms::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 plane (ax+by+cz+d=0).
[in] | p1 | First point. |
[in] | p2 | Second point. |
[in] | p3 | Third point. |
[out] | a | Coefficient a. |
[out] | b | Coefficient b. |
[out] | c | Coefficient c. |
[out] | d | Coefficient d. |
void xms::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 pointers calculate coefficents for plane (ax+by+cz+d=0).
[in] | p1 | First point. |
[in] | p2 | Second point. |
[in] | p3 | Third point. |
[out] | a | Coefficient a. |
[out] | b | Coefficient b. |
[out] | c | Coefficient c. |
[out] | d | Coefficient d. |
int xms::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.
[in] | cart | cart? |
[in] | orig | orig? |
[in] | coef | 6 coefficients. |
[in] | dir | Direction. |
[out] | bary | point? |
Definition at line 275 of file geoms.cpp.
References XM_SUCCESS.
Referenced by iGetBarycentricCoords().
bool xms::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.
pt1 | first of 3 locations defining a triangle |
pt2 | second of 3 locations defining a triangle |
pt3 | third of 3 locations defining a triangle |
xc | circumcenter x coord |
yc | circumcenter y coord |
r2 | radius squared of circumcircle |
tol | tolerance for geometric comparison |
Definition at line 218 of file geoms.cpp.
References sqr().
Referenced by gmPtInCircumcircle().
Returns true if (the three vertices are gmColinear. Result should be insensitive to the order of the vertices.
[in] | p1 | Point 1. |
[in] | p2 | Point 2. |
[in] | p3 | Point 3. |
[in] | tol | Tolerance |
Definition at line 393 of file geoms.cpp.
References gmOnLineWithTol().
Referenced by gmIntersectLineSegmentsWithTol(), and xms::TrTriangulatorPoints::ReceiveTriangle().
void xms::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
[in,out] | a_x | vector x component either specified or calculated from a_mag, a_dir |
[in,out] | a_y | vector y component either specified or calculated from a_mag, a_dir |
[in,out] | a_mag | vector magnitude either specified or calculated from a_x, a_y |
[in,out] | a_dir | vector direction (degrees) either specified or calculated from a_x, a_y |
[in] | a_tomagdir | flag, when true a_x, a_y are used to calculate a_mag, a_dir; when false a_mag, a_dir are used to calculate a_x, a_y |
Definition at line 1941 of file geoms.cpp.
References sqr(), XM_PI, and XM_ZERO_TOL.
Computes the centroid of points (but not of a polygon). Shouldn't pass an empty vector (no checks are performed).
[in] | a_points | Points. |
Definition at line 895 of file geoms.cpp.
Referenced by xms::GmMultiPolyIntersectionSorterTerse::RemoveDuplicateEdges(), GeomsXmsngUnitTests::test_gmComputeCentroid(), and trBuildGridTrianglePolys().
double xms::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_p1 to a_p2.
[in] | a_p0 | 1st point on 1st segment. |
[in] | a_p1 | 2nd point on 1st segment and 1st point on second segment. |
[in] | a_p2 | 2nd point on 2nd segment. |
Computes the plan view centroid of a non-self-intersecting polygon.
[in] | pts | Points. |
Definition at line 908 of file geoms.cpp.
Referenced by xms::XmUGrid::Impl::GetCellCentroid(), and GeomsXmsngUnitTests::test_gmComputePolygonCentroid().
Given extents (min, max), compute a tolerance for the xy plane to be used with geometric functions.
a_mn | Minimum. |
a_mx | Maximum. |
Definition at line 1262 of file geoms.cpp.
References gmXyDistance().
Referenced by xms::GmMultiPolyIntersectorImpl::CalculateBuffer(), xms::TrBreaklineAdderImpl::ComputeTolerance(), and xms::TrTriangulatorPoints::UpdateAreaTolerance().
double xms::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 deg to 360 deg
[in] | a_angle | The angle to convert to 0 to 360 |
[in] | a_InDegrees | Flag to tell if the angle is in degrees vs radians |
Definition at line 1997 of file geoms.cpp.
References GTEQ_TOL(), LT_TOL(), and XM_PI.
Returns true if the triangle is wrapped counter clockwise.
[in] | vtx0 | Vertex 0. |
[in] | vtx1 | Vertex 1. |
[in] | vtx2 | Vertex 2. |
Definition at line 564 of file geoms.cpp.
References trArea().
Referenced by xms::TrTinImpl::SwapEdge().
double xms::gmCross2D | ( | const double & | dx1, |
const double & | dy1, | ||
const double & | dx2, | ||
const double & | dy2 | ||
) |
Returns the cross product of two 2-d vectors.
dx1 | x component of vector 1 |
dy1 | y component of vector 1 |
dx2 | x component of vector 2 |
dy2 | y component of vector 2 |
Definition at line 579 of file geoms.cpp.
Referenced by gmAngleBetween2DVectors(), gmBisectingAngle(), and gmLinesCross().
Perform a cross product of Pt3d's.
[in] | a_vec1 | First vector to cross |
[in] | a_vec2 | Second vector to cross |
[out] | a_vec3 | Resulting cross product |
Definition at line 2034 of file geoms.cpp.
Referenced by gmIntersectTriangleAndLineSegment().
Perform a dot product of Pt3d's.
[in] | a_vec1 | First vector to dot |
[in] | a_vec2 | Second vector to dot |
Definition at line 2048 of file geoms.cpp.
Referenced by gmIntersectTriangleAndLineSegment().
Calculates the envelope of a vector of points.
[in] | a_pts | Array of points |
[out] | a_min | Min x,y,z of the points |
[out] | a_max | Max x,y,z of the points |
Definition at line 1606 of file geoms.cpp.
References XM_ENSURE_TRUE.
bool xms::gmEqualPointsXY | ( | double | x1, |
double | y1, | ||
double | x2, | ||
double | y2, | ||
double | tolerance | ||
) |
Returns true if the points are equal to within gmXyTol().
x1 | x of point 1. |
y1 | y of point 1. |
x2 | x of point 2. |
y2 | y of point 2. |
tolerance | tolerance. |
Definition at line 1308 of file geoms.cpp.
Referenced by xms::GmMultiPolyIntersectionSorterTerse::FixTValueAtDuplicateXy(), gmEqualPointsXY(), gmInsideOfLineWithTol(), gmInsideOrOnLineWithTol(), gmIntersectLineSegmentsWithTol(), gmOnLineWithTol(), and xms::GmMultiPolyIntersectorImpl::PointsOnSegment().
bool xms::gmEqualPointsXY | ( | double | x1, |
double | y1, | ||
double | x2, | ||
double | y2 | ||
) |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
x1 | x of point 1. |
y1 | y of point 1. |
x2 | x of point 2. |
y2 | y of point 2. |
Definition at line 1327 of file geoms.cpp.
References gmEqualPointsXY(), and gmXyTol().
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
a_pt1 | Point 1. |
a_pt2 | Point 2. |
tol | tolerance. |
Definition at line 1338 of file geoms.cpp.
References gmEqualPointsXY().
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
a_pt1 | Point 1. |
a_pt2 | Point 2. |
tol | tolerance. |
Definition at line 1349 of file geoms.cpp.
References gmEqualPointsXY().
bool xms::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.
x1 | x of point 1. |
y1 | y of point 1. |
z1 | z of point 1. |
x2 | x of point 2. |
y2 | y of point 2. |
z2 | z of point 2. |
tolerance | tolerance. |
Definition at line 1376 of file geoms.cpp.
Referenced by gmEqualPointsXYZ().
bool xms::gmEqualPointsXYZ | ( | double | x1, |
double | y1, | ||
double | z1, | ||
double | x2, | ||
double | y2, | ||
double | z2 | ||
) |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
x1 | x of point 1. |
y1 | y of point 1. |
z1 | z of point 1. |
x2 | x of point 2. |
y2 | y of point 2. |
z2 | z of point 2. |
Definition at line 1398 of file geoms.cpp.
References gmEqualPointsXYZ(), and gmXyTol().
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
pt1 | Point 1. |
pt2 | Point 2. |
tol | tolerance. |
Definition at line 1409 of file geoms.cpp.
References gmEqualPointsXYZ().
Get the 2D extents of a vector of points.
[in] | a_points | The points. |
[out] | a_min | Minimum point (xy) of bounding rectangle. |
[out] | a_max | Maximum point (xy) of bounding rectangle. |
Definition at line 1821 of file geoms.cpp.
References gmAddToExtents(), and XM_ENSURE_TRUE_VOID_NO_ASSERT.
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
[in] | a_points | Vector of points. |
[in,out] | a_min | Minimum value. |
[in,out] | a_max | Maximum value. |
Definition at line 1838 of file geoms.cpp.
References gmAddToExtents(), and XM_ENSURE_TRUE_VOID_NO_ASSERT.
Get the 3D extents of a vector of points.
[in] | a_points | The points. |
[out] | a_min | Minimum point (xyz) of bounding box. |
[out] | a_max | Maximum point (xyz) of bounding box. |
Definition at line 1856 of file geoms.cpp.
References gmAddToExtents(), and XM_ENSURE_TRUE_VOID_NO_ASSERT.
double xms::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.
[in] | a_pt1 | First point of segment |
[in] | a_pt2 | Second point of segment |
[in] | a_pt | Point used to find closest point on the segment |
[out] | a_newpt | The point on the segment (a_pt1, a_pt2) that is closest to a_pt. |
[in] | a_tol | Tolerance |
Definition at line 1714 of file geoms.cpp.
References gmPtDistanceAlongSegment().
Calculate convex hull using Monotone chain aka Andrew's algorithm.
[in] | a_pts | The input points. |
[out] | a_hull | The convex hull polygon. First point is NOT repeated. |
[in] | a_includeOn | Should points on the hull's line segments be included? |
Definition at line 2211 of file geoms.cpp.
References gmTurn().
Referenced by ConvexHull(), xms::XmUGrid::Impl::GetPlanViewPolygon3d(), and GeomsXmsngUnitTests::testConvexHull().
bool xms::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 oppositevertex.
[in] | a_vertex1 | : First point on the line segment |
[in] | a_vertex2 | : Second point on the line segment |
[in] | a_oppositevertex | : Point on one side of the line |
[in] | a_x | : X location of the point being tested |
[in] | a_y | : Y location of the point being tested |
[in] | a_tol | : Tolerance for degenerate cases |
Definition at line 1782 of file geoms.cpp.
References gmEqualPointsXY(), and sqr().
bool xms::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". ASSERTs in debug if "inpoint" is on the line (within tol).
[in] | p1 | Point 1. |
[in] | p2 | Point 2. |
[in] | inpoint | Point on same side of line as xy if true is returned. |
[in] | x | x of point. |
[in] | y | y of point. |
[in] | tol | Tolerance. |
[out] | dist | How far "outside" the point is. A negative distance indicates point is in by that distance. |
Definition at line 1450 of file geoms.cpp.
References gmEqualPointsXY(), Mfabs(), and sqr().
Referenced by gmPointInTriangleWithTol().
bool xms::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.
one1 | first point on segment 1 |
one2 | second point on segment 1 |
two1 | first point on segment 2 |
two2 | second point on segment 2 |
xi | x coord of intersection |
yi | y coord of intersection |
zi1 | z coord of intersection on segment 1 |
zi2 | z coord of intersection on segemnt 2 |
tol | tolerance for geometric comparison |
segment 1 = one1,one2 = one1 + lambda(one2 - one1) segment 2 = two1,two2 = two1 + mu (two2 - two1)
Definition at line 420 of file geoms.cpp.
References EQ_TOL(), gmColinearWithTol(), gmEqualPointsXY(), GT_TOL(), and LT_TOL().
Referenced by xms::TrBreaklineAdderImpl::FindIntersectingEdgeFromEdge(), xms::TrBreaklineAdderImpl::FindIntersectingEdgeFromPoint(), and xms::TrBreaklineAdderImpl::ProcessSegmentBySwapping().
int xms::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).
[in] | a_pt1 | First point defining a segment |
[in] | a_pt2 | Second point defining a segment |
[in] | a_t0 | First point defining a triangle |
[in] | a_t1 | Second point defining a triangle |
[in] | a_t2 | Third point defining a triangle |
[out] | a_IntersectPt | The point of intersection if the segment and triangle intersect |
Definition at line 2075 of file geoms.cpp.
References gmCross3D(), gmDot3D(), and XM_ZERO_TOL.
bool xms::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.
[in] | a_segment1Point1 | First point 3d of line segment 1 |
[in] | a_segment1Point2 | Second point 3d of line segment 1 |
[in] | a_segment2Point1 | First point 3d of line segment 2 |
[in] | a_segment2Point2 | Second point 3d of line segment 2 |
Definition at line 1013 of file geoms.cpp.
References gmCross2D().
Referenced by xms::XmUGrid::Impl::DoEdgesCrossWithPointChange(), DoLineSegmentsCross(), and GeomsXmsngUnitTests::testDoLineSegmentsCross().
bool xms::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)
one1 | first point on first segment |
one2 | second point on first segment |
two1 | first point on second segment |
two2 | second point on second segment |
Definition at line 966 of file geoms.cpp.
References GT_TOL(), LT_TOL(), and XM_ZERO_TOL.
double xms::gmMiddleZ | ( | const VecPt3d & | a_points | ) |
bool xms::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.
a_pt1 | first location defining segment |
a_pt2 | second location defining segment |
a_x | x coord |
a_y | y coord |
a_tol | tolerance for comparison |
Definition at line 767 of file geoms.cpp.
References gmOnLineWithTol(), Mmax(), and Mmin().
Referenced by xms::GmMultiPolyIntersectorImpl::PointsOnSegment().
bool xms::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't the same.
p1 | first location defining a line |
p2 | second location defining a line |
x | x coord of point to test |
y | y coord of point to test |
tol | tolerance for geometric comparison |
Definition at line 732 of file geoms.cpp.
References gmEqualPointsXY(), and sqr().
Referenced by gmColinearWithTol(), and gmOnLineAndBetweenEndpointsWithTol().
void xms::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 ordered counter clockwise based on the angle from the centroid where angle starts at point at startindex.
[in] | a_pts | Array of points |
[out] | a_ccwOrder | Array of indices to points |
[in] | a_startindex | Starting position in arrays |
Definition at line 1635 of file geoms.cpp.
Referenced by gmOrderPointsCounterclockwise().
void xms::gmOrderPointsCounterclockwise | ( | VecPt3d & | a_pts | ) |
Overload to gmOrderPointsCounterclockwise.
[out] | a_pts | Array of points |
Definition at line 1690 of file geoms.cpp.
References gmOrderPointsCounterclockwise().
int xms::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 vertices.
DON'T repeat the first point at the end of the polygon array
[in] | a_verts | The polygon in CW or CCW order (doesn't matter which). |
[in] | a_num | Number of verts in a_verts. |
[in] | a_x | X coordinate of point to test. |
[in] | a_y | Y coordinate of point to test. |
[in] | a_tol | Tolerance. |
Definition at line 1057 of file geoms.cpp.
Referenced by GmPointInPolyTester_gmPointInPolygon2D::CheckPoint(), gmPointInPolygon2D(), and GeomsXmsngIntermediateTests::test_gmPointInPolygon2D().
int xms::gmPointInPolygon2D | ( | const Pt3d * | a_verts, |
size_t | a_num, | ||
double | a_x, | ||
double | a_y | ||
) |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
[in] | a_verts | The polygon in CW or CCW order (doesn't matter which). |
[in] | a_num | Number of verts in a_verts. |
[in] | a_x | X coordinate of point to test. |
[in] | a_y | Y coordinate of point to test. |
Definition at line 1155 of file geoms.cpp.
References gmPointInPolygon2D(), and gmXyTol().
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
[in] | a_verts | The polygon in CW or CCW order (doesn't matter which). |
[in] | a_num | Number of verts in a_verts. |
[in] | a_pt | Point to test. |
Definition at line 1166 of file geoms.cpp.
References gmPointInPolygon2D(), and gmXyTol().
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
[in] | a_verts | The polygon in CW or CCW order (doesn't matter which). |
[in] | a_num | Number of verts in a_verts. |
[in] | a_pt | Point to test. |
Definition at line 1177 of file geoms.cpp.
References gmPointInPolygon2D(), and gmXyTol().
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
[in] | a_verts | The polygon in CW or CCW order (doesn't matter which). |
[in] | a_num | Number of verts in a_verts. |
[in] | a_pt | Point to test. |
Definition at line 1188 of file geoms.cpp.
References gmPointInPolygon2D(), and gmXyTol().
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
[in] | a_verts | The polygon in CW or CCW order (doesn't matter which). |
[in] | a_num | Number of verts in a_verts. |
[in] | a_pt | Point to test. |
Definition at line 1199 of file geoms.cpp.
References gmPointInPolygon2D(), and gmXyTol().
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
[in] | a_verts | The polygon in CW or CCW order (doesn't matter which). |
[in] | a_pt | Point to test. |
Definition at line 1209 of file geoms.cpp.
References gmPointInPolygon2D(), and gmXyTol().
Referenced by GeomsXmsngIntermediateTests::test_gmPointInPolygon2D().
template int xms::gmPointInPolygon2D< Pt2d > | ( | const Pt2d * | a_verts, |
const size_t | a_num, | ||
const double | a_x, | ||
const double | a_y, | ||
const double | a_tol | ||
) |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
[in] | a_verts | The polygon in CW or CCW order (doesn't matter which). |
[in] | a_num | Number of verts in a_verts. |
[in] | a_x | X coordinate of point to test. |
[in] | a_y | Y coordinate of point to test. |
[in] | a_tol | Tolerance. |
template int xms::gmPointInPolygon2D< Pt2i > | ( | const Pt2i * | a_verts, |
const size_t | a_num, | ||
const double | a_x, | ||
const double | a_y, | ||
const double | a_tol | ||
) |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
[in] | a_verts | The polygon in CW or CCW order (doesn't matter which). |
[in] | a_num | Number of verts in a_verts. |
[in] | a_x | X coordinate of point to test. |
[in] | a_y | Y coordinate of point to test. |
[in] | a_tol | Tolerance. |
template int xms::gmPointInPolygon2D< Pt3d > | ( | const Pt3d * | a_verts, |
const size_t | a_num, | ||
const double | a_x, | ||
const double | a_y, | ||
const double | a_tol | ||
) |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
[in] | a_verts | The polygon in CW or CCW order (doesn't matter which). |
[in] | a_num | Number of verts in a_verts. |
[in] | a_x | X coordinate of point to test. |
[in] | a_y | Y coordinate of point to test. |
[in] | a_tol | Tolerance. |
bool xms::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.
[in] | p1 | Point 1. |
[in] | p2 | Point 2. |
[in] | p3 | Point 3. |
[in] | x | x of point. |
[in] | y | y of point. |
[in] | tol | Tolerance. |
Definition at line 1423 of file geoms.cpp.
References gmInsideOrOnLineWithTol().
double xms::gmPolygonArea | ( | const Pt3d * | pts, |
size_t | npoints | ||
) |
Compute 2d planview projection of area of polygon.
pts | locations defining the polygon |
npoints | number of pts |
CCW = positive area, CW = negative. Don't repeat the last point.
Definition at line 1537 of file geoms.cpp.
Referenced by xms::XmUGrid::Impl::GetOrientationFromArea().
double xms::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.
[in] | a_pt1 | First point of segment |
[in] | a_pt2 | Second point of segment |
[in] | a_pt | Point used to find closest point on the segment |
[in] | a_tol | Tolerance |
Definition at line 1749 of file geoms.cpp.
Referenced by xms::GmMultiPolyIntersectorImpl::ComputeTValues(), and gmFindClosestPtOnSegment().
PtInOutOrOn_enum xms::gmPtInCircumcircle | ( | const Pt3d & | pt, |
Pt3d | circumcirclePts[3] | ||
) |
Is a given point inside a circumcircle defined by three points?
pt | The point. |
circumcirclePts | 3 points on the circumcircle. |
Definition at line 172 of file geoms.cpp.
References gmCircumcircleWithTol(), gmXyTol(), PT_ERROR, PT_IN, PT_ON, PT_OUT, and sqr().
Referenced by xms::TrTinImpl::PointIsInCircumcircle().
Determine if angle a_v1, a_v2, a_v3 is a left turn, a right turn, colinear 180 degrees, or colinear 0 degrees.
a_v1 | first of 3 locations defining an angle |
a_v2 | second of 3 locations defining an angle |
a_v3 | third of 3 locations defining an angle |
a_angtol | tolerance for comparing angle see notes |
Definition at line 863 of file geoms.cpp.
References TURN_COLINEAR_0, TURN_COLINEAR_180, TURN_LEFT, and TURN_RIGHT.
Referenced by gmGetConvexHull(), and xms::GmMultiPolyIntersectionSorterTerse::RemoveDuplicateEdges().
double xms::gmXyDistance | ( | double | x1, |
double | y1, | ||
double | x2, | ||
double | y2 | ||
) |
Compute the 2d distance between 2 points.
x1 | x coord of location 1 |
y1 | y coord of location 1 |
x2 | x coord of location 2 |
y2 | y coord of location 2 |
Definition at line 832 of file geoms.cpp.
References sqr().
Referenced by xms::GmMultiPolyIntersectorImpl::CalculateBuffer(), xms::TrTinImpl::CheckTriangle(), and gmComputeXyTol().
double xms::gmXyTol | ( | bool | a_set, |
double | a_value | ||
) |
Get or set (set first time) global xy tolerance for float operations.
a_set | True if setting. |
a_value | The value if setting. |
Definition at line 1279 of file geoms.cpp.
Referenced by xms::GmMultiPolyIntersectorImpl::ComputeTValues(), xms::GmMultiPolyIntersectionSorterTerse::FixTValueAtDuplicateXy(), gmEqualPointsXY(), gmEqualPointsXYZ(), gmPointInPolygon2D(), gmPtInCircumcircle(), and xms::GmMultiPolyIntersectorImpl::TraverseLineSegment().
double xms::gmZTol | ( | bool | a_set, |
double | a_value | ||
) |
void xms::iCartToBary | ( | const Pt3d & | a_pt1, |
const Pt3d & | a_pt2, | ||
const Pt3d & | a_pt3, | ||
BarycentricVals & | a_b | ||
) |
Converts from Cartesian coordinates to Barycentric coordinates (local triangle coords).
a_pt1 | Point 1. |
a_pt2 | Point 2. |
a_pt3 | Point 3. |
a_b | Barycentric values. |
Definition at line 188 of file GmTriSearch.cpp.
References xms::BarycentricVals::coef, xms::BarycentricVals::dir, gmBaryPrepare(), xms::BarycentricVals::orig, xms::Pt3< T >::x, and xms::Pt3< T >::y.
Referenced by xms::GmTriSearchImpl::GetTriBarycentricVals().
|
static |
Return the minimum distance from the point to a line on the ring.
[in] | a_ring | The ring. |
[in] | a_pt | A point. |
Definition at line 75 of file GmPolygon.cpp.
Referenced by xms::GmPolygonImpl::MinDistanceToBoundary().
void xms::iGetBarycentricCoords | ( | const Pt3d & | a_p, |
BarycentricVals & | a_b, | ||
Pt3d & | weights | ||
) |
Calculates weights at the 3 points of a triangle given the location "a_p" and the Barycentric coords of the triangle.
a_p | A point |
a_b | Barycentric values. |
weights | Weights. |
Definition at line 203 of file GmTriSearch.cpp.
References xms::BarycentricVals::coef, xms::BarycentricVals::dir, gmCartToBary(), xms::BarycentricVals::orig, xms::Pt3< T >::x, and xms::Pt3< T >::y.
Referenced by xms::GmTriSearchImpl::FindTriangle().
void xms::mxCopy4x4 | ( | double | copy[4][4], |
const double | matrix[4][4] | ||
) |
Copy a transformation matrix. Simply uses memcpy.
copy | matrix that will have result of copy |
matrix | 4x4 matrix to copy |
Definition at line 563 of file matrix.cpp.
Referenced by mxInvert4x4(), and mxMult4x4().
void xms::mxIdent4x4 | ( | double | matrix[4][4] | ) |
Sets a transformation matrix to identity.
matrix | matrix to change |
Definition at line 571 of file matrix.cpp.
Referenced by mxRotate4x4(), and mxTranslate4x4().
void xms::mxiLubksb3 | ( | const double | mat[3][3], |
const int | indx[3], | ||
double | b[3] | ||
) |
solve [mat]{x} = {b} by back substitution (mat = LU of mat)
mat | matrix[3][3] |
indx | |
b | vector of solution [mat]{x} = {b} |
Definition at line 231 of file matrix.cpp.
Referenced by mxSolve3x3().
bool xms::mxiLudcmp3 | ( | double | mat[3][3], |
int | indx[3], | ||
double * | d | ||
) |
Decomposes a 3 x 3 matrix into Upper & Lower trianglar(in place)
mat | matrix[3][3] |
indx | |
d | determinant of matrix |
Definition at line 157 of file matrix.cpp.
References xmlog::debug, and XM_LOG.
Referenced by mxSolve3x3().
int xms::mxInvert4x4 | ( | const double | matrix[4][4], |
double | inv[4][4] | ||
) |
Compute the inverse of a transformation matrix. Checks for special case of Orthogonal 4x4 and does easy inverse otherwise uses Gauss-Jordan Elimination with full pivoting and scaling with partial unraveling of loops for speed. See Numerical Recipes in C pages 36-37.
matrix | matrix[4][4] to invert |
inv | inverse of matrix |
Definition at line 422 of file matrix.cpp.
References xmlog::error, mxCopy4x4(), XM_FAILURE, XM_LOG, and XM_SUCCESS.
int xms::mxLUBcksub | ( | double ** | mat, |
int | n, | ||
const int * | indx, | ||
double * | b | ||
) |
solve [mat]{x} = {b} by back substitution (mat = LU of mat), from "Numerical Recipes in C" p 44.
mat | matrix[n][n] |
n | dimension of matrix mat |
indx | |
b | vector of solution [mat]{x} = {b} |
Definition at line 122 of file matrix.cpp.
References XM_SUCCESS.
Referenced by mxSolveNxN().
int xms::mxLUDecomp | ( | double ** | mat, |
int | n, | ||
int * | indx, | ||
double * | d | ||
) |
Decompose an n x n matrix into Upper & Lower trianglar (in place), from "Numerical Recipes in C" p 43.
mat | matrix[n][n] |
n | dimension of matrix mat |
indx | |
d | matrix determinant |
Definition at line 47 of file matrix.cpp.
References XM_FAILURE, and XM_SUCCESS.
Referenced by mxSolveNxN().
void xms::mxMult4x4 | ( | double | product[4][4], |
const double | matrix1[4][4], | ||
const double | matrix2[4][4] | ||
) |
Multiplies two transformation matrices.
product | result of matrix1*matrix2 |
matrix1 | [4][4] matrix multiplied by matrix2 |
matrix2 | [4][4] matrix multiplied by matrix1 |
Definition at line 585 of file matrix.cpp.
References mxCopy4x4().
void xms::mxPointMult | ( | Pt3d * | pt, |
const double | ctm[4][4] | ||
) |
Multiplies a point by a transformation matrix.
pt | point to multiply |
ctm | transform matrix[4][4]. |
Definition at line 549 of file matrix.cpp.
void xms::mxRotate4x4 | ( | double | xrot, |
double | yrot, | ||
double | zrot, | ||
double | rMatt[4][4] | ||
) |
return a rotation matrix with specified x, y, z rotations z - y - x order
xrot | x rotation amount |
yrot | y rotation amount |
zrot | z rotation amount |
rMatt | matrix to be built |
Definition at line 611 of file matrix.cpp.
References mxIdent4x4(), and XM_PI.
bool xms::mxSolve3x3 | ( | double | A[3][3], |
double | x[3], | ||
double | b[3] | ||
) |
Solves the matrix equation [A]*{x}={b} for {x} where [A] is 3x3 and {x} and {b} are 3x1 vectors.
A | matrix[3][3] |
x | vector[3] |
b | vector[3] of solution [A]{x} = {b} |
Definition at line 399 of file matrix.cpp.
References mxiLubksb3(), and mxiLudcmp3().
bool xms::mxSolveBandedEquations | ( | double ** | a, |
double * | x, | ||
double * | b, | ||
int | numeqs, | ||
int | halfbandwidth | ||
) |
Solves the matrix equation [a][x]=[b] using gauss elimination.
a | matrix[n][n] |
x | vector nx1 |
b | vector of solution [A]{x} = {b} |
numeqs | number of equations |
halfbandwidth |
Definition at line 325 of file matrix.cpp.
bool xms::mxSolveNxN | ( | double ** | A, |
double * | x, | ||
double * | b, | ||
int | n | ||
) |
Solves the matrix equation [A]*{x}={b} for {x} where [A] is NxN and {x} and {b} are Nx1 vectors. The matrix A is replaced with the LU decomposition of A.
A | matrix[n][n] |
x | vector nx1 |
b | vector of solution [A]{x} = {b} |
n | size of matrix A |
Definition at line 268 of file matrix.cpp.
References mxLUBcksub(), mxLUDecomp(), and XM_FAILURE.
void xms::mxTranslate4x4 | ( | const Pt3d & | a_translation, |
double | a_mx[4][4] | ||
) |
return a translation matrix with specified translation z - y - x order
a_translation | translation distances in x, y, and z. |
a_mx | matrix to be built |
Definition at line 638 of file matrix.cpp.
References mxIdent4x4().
std::shared_ptr< XmUGrid > xms::TEST_XmUBuild3DChevronUgrid | ( | ) |
Builds a UGrid with one 3D Chevron polyhedron.
Definition at line 3980 of file XmUGrid.cpp.
References xms::XmUGrid::New().
Referenced by XmUGridUnitTests::testGetCellPlanViewPolygon().
std::shared_ptr< XmUGrid > xms::TEST_XmUBuildHexahedronUgrid | ( | int | a_rows, |
int | a_cols, | ||
int | a_lays | ||
) |
Builds a UGrid of Quads at 1 spacing for rows & cols specified.
[in] | a_rows | number of rows in UGrid |
[in] | a_cols | number of columns in UGrid |
[in] | a_lays | number of layers in UGrid |
Definition at line 3817 of file XmUGrid.cpp.
Referenced by XmUGridUnitTests::testCell3dFunctionCaching(), XmUGridUnitTests::testGetCell3dFaceAdjacentCell(), XmUGridUnitTests::testGetCell3dFaceOrientationHexahedrons(), XmUGridUnitTests::testGetCellEdgeAdjacentCells(), XmUGridUnitTests::testGetExtents(), and XmUGridUnitTests::testIsCellValidWithPointChange().
std::shared_ptr< XmUGrid > xms::TEST_XmUBuildHexahedronUgrid | ( | int | a_rows, |
int | a_cols, | ||
int | a_lays, | ||
const xms::Pt3d & | a_origin | ||
) |
Builds a UGrid of Quads at 1 spacing for rows & cols specified.
[in] | a_rows | number of rows in UGrid |
[in] | a_cols | number of columns in UGrid |
[in] | a_lays | number of layers in UGrid |
[in] | a_origin | locatGetPointLocation UGrid begins (min x, min y, min z) |
Definition at line 3830 of file XmUGrid.cpp.
References xms::XmUGrid::New().
std::shared_ptr< XmUGrid > xms::TEST_XmUBuildPolyhedronUgrid | ( | int | a_rows, |
int | a_cols, | ||
int | a_lays | ||
) |
Builds a UGrid of Quads at 1 spacing for rows & cols specified.
[in] | a_rows | number of rows in UGrid |
[in] | a_cols | number of columns in UGrid |
[in] | a_lays | number of layers in UGrid |
Definition at line 3883 of file XmUGrid.cpp.
Referenced by XmUGridUnitTests::testGetCell3dFaceOrientationHexahedrons(), XmUGridUtilsTests::testLargeUGridBinarySpeed(), and XmUGridUnitTests::testLargeUGridLinkSpeed().
std::shared_ptr< XmUGrid > xms::TEST_XmUBuildPolyhedronUgrid | ( | int | a_rows, |
int | a_cols, | ||
int | a_lays, | ||
const xms::Pt3d & | a_origin | ||
) |
Builds a UGrid of Quads at 1 spacing for rows & cols specified.
[in] | a_rows | number of rows in UGrid |
[in] | a_cols | number of columns in UGrid |
[in] | a_lays | number of layers in UGrid |
[in] | a_origin | location wGetPointLocationd begins (min x, min y, min z) |
Definition at line 3896 of file XmUGrid.cpp.
References xms::XmUGrid::New().
std::shared_ptr< XmUGrid > xms::TEST_XmUBuildQuadUGrid | ( | int | a_rows, |
int | a_cols | ||
) |
Builds a UGrid of Quads at 1 spacing for rows & cols specified.
[in] | a_rows | number of rows in UGrid |
[in] | a_cols | number of columns in UGrid |
Definition at line 3769 of file XmUGrid.cpp.
std::shared_ptr< XmUGrid > xms::TEST_XmUBuildQuadUGrid | ( | int | a_rows, |
int | a_cols, | ||
const xms::Pt3d & | a_origin | ||
) |
Builds a UGrid of Quads at 1 spacing for rows & cols specified.
[in] | a_rows | number of rows in UGrid |
[in] | a_cols | number of columns in UGrid |
[in] | a_origin | location where the UGrid begins (min x, min y, min z) |
Definition at line 3781 of file XmUGrid.cpp.
References xms::XmUGrid::New().
std::shared_ptr< XmUGrid > xms::TEST_XmUGrid1Left90Tri | ( | ) |
Builds a 1 cell (left 90 degree triangle) 2D XmUGrid for testing.
Definition at line 3639 of file XmUGrid.cpp.
References xms::XmUGrid::New().
Referenced by XmUGridUnitTests::testOperators(), XmUGridUtilsTests::testReadBasicUGrid(), XmUGridUtilsTests::testReadReadyAtNextLine(), and XmUGridUtilsTests::testWriteBasicUGrid().
std::shared_ptr< XmUGrid > xms::TEST_XmUGrid2dLinear | ( | ) |
Builds an XmUGrid with supported 1D and 2D linear cells for testing.
[snip_test_2DShapes]
[snip_test_2DShapes]
Definition at line 3679 of file XmUGrid.cpp.
References xms::XmUGrid::New().
Referenced by XmUGridUnitTests::testGetCell3dFaceCount(), XmUGridUnitTests::testGetCell3dFacePointCount(), XmUGridUnitTests::testGetCell3dFacePoints(), XmUGridUnitTests::testGetCellAdjacentCells(), XmUGridUnitTests::testGetCellDimension(), XmUGridUnitTests::testGetCellEdge(), XmUGridUnitTests::testGetCellEdgeCount(), XmUGridUnitTests::testGetCellPlanViewPolygon(), XmUGridUnitTests::testGetCellPoints(), XmUGridUnitTests::testGetCellType(), XmUGridUnitTests::testGetExtents(), XmUGridUnitTests::testGetPointAdjacentCells(), XmUGridUtilsTests::testLinear2dWriteThenRead(), and XmUGridUtilsTests::testWriteLinear2dCells().
std::shared_ptr< XmUGrid > xms::TEST_XmUGrid3dLinear | ( | ) |
Builds an XmUGrid with supported 3D linear cells for testing.
[snip_test_3DShapes]
[snip_test_3DShapes]
Definition at line 3705 of file XmUGrid.cpp.
References xms::XmUGrid::New().
Referenced by XmUGridUnitTests::testCell3dFaceFunctions(), XmUGridUnitTests::testGetCell3dFaceCount(), XmUGridUnitTests::testGetCell3dFacePointCount(), XmUGridUnitTests::testGetCell3dFacePoints(), XmUGridUnitTests::testGetCellAdjacentCells(), XmUGridUnitTests::testGetCellDimension(), XmUGridUnitTests::testGetCellEdge(), XmUGridUnitTests::testGetCellEdgeCount(), XmUGridUnitTests::testGetCellPlanViewPolygon(), XmUGridUnitTests::testGetCellPoints(), XmUGridUnitTests::testGetCellType(), XmUGridUnitTests::testGetExtents(), XmUGridUnitTests::testGetPointAdjacentCells(), XmUGridUtilsTests::testLinear3dWriteThenRead(), XmUGridUtilsTests::testReadVersion1Dot0File(), XmUGridUtilsTests::testWriteLinear3dCells(), XmUGridUtilsTests::testWriteThenReadUGridBinary(), XmUGridUtilsTests::testWriteThenReadUGridFile(), and XmUGridUtilsTests::testWriteThenReadUGridFileToAscii().
std::shared_ptr< XmUGrid > xms::TEST_XmUGridHexagonalPolyhedron | ( | ) |
Builds a 1 cell hexagon with cell type polyhedron.
Definition at line 3740 of file XmUGrid.cpp.
References xms::XmUGrid::New().
Referenced by XmUGridUtilsTests::testReadPolyhedronUGrid(), and XmUGridUtilsTests::testWritePolyhedronUGrid().
std::shared_ptr< XmUGrid > xms::TEST_XmUGridSimpleQuad | ( | ) |
Builds a 2 cell (Quad) 2D XmUGrid for testing.
Definition at line 3662 of file XmUGrid.cpp.
References xms::XmUGrid::New().
Referenced by XmUGridUnitTests::testCellEdgeAdjacentCellFunctions(), XmUGridUnitTests::testCellEdges(), XmUGridUnitTests::testCellFunctions(), XmUGridUnitTests::testCellstreamFunctions(), XmUGridUnitTests::testEdgeAdjacentCells(), XmUGridUnitTests::testGetCellAdjacentCells(), XmUGridUnitTests::testGetCellEdge(), XmUGridUnitTests::testGetCellEdgeAdjacentCells(), XmUGridUnitTests::testGetCellPoints(), XmUGridUnitTests::testGetPointAdjacentPoints(), XmUGridUnitTests::testGetPointsAdjacentCells(), XmUGridUnitTests::testIsCellValidWithPointChange(), and XmUGridUnitTests::testPointFunctions().
Return the signed planar area of the triangle (CCW Positive).
This uses the formula: (-x2*y1 + x3*y1 + x1*y2 - x3*y2 - x1*y3 + x2*y3) / 2 reworked to ensure precision is not lost with large numbers by moving the triangle to the origin by subtracting the first vertex's location from all vertices.
[in] | a_pt1 | First point on triangle. |
[in] | a_pt2 | Second point on triangle. |
[in] | a_pt3 | Third point on triangle. |
Definition at line 46 of file triangles.cpp.
Referenced by gmCounterClockwiseTri(), xms::TrTriangulatorPoints::ReceiveTriangle(), xms::TrTinImpl::SwapEdge(), and xms::TrTinImpl::TriangleArea().
Build something like this:
/// (0, 0) /// 0-------------1-------------2 /// | | | /// | | | /// | 1 | 2 | /// | | | /// | | | /// 3-------------4-------------5 /// (20,-10) ///
[in] | rows | Number of rows |
[in] | cols | Number of columns |
[out] | pts | The points. |
[out] | polys | The polygons. |
Definition at line 132 of file triangles.cpp.
Referenced by trBuildGridTrianglePolys().
Create something like this:
/// (0, 0) /// 0-------------1-------------2 /// |\ 4 /|\ 8 /| /// | \ / | \ / | /// | \ / | \ / | /// | 1 3 3 | 5 4 7 | /// | / \ | / \ | /// | / \ | / \ | /// |/ 2 \|/ 6 \| /// 5-------------6-------------7 2 /// (20, -10) ///
[in] | rows | Number of rows |
[in] | cols | Number of columns |
[out] | a_points | The points. |
[out] | a_polys | The polygons. |
Definition at line 72 of file triangles.cpp.
References gmComputeCentroid(), and trBuildGridPolys().
Referenced by GmMultiPolyIntersectorUnitTests::testAlongEdgesInsideToInside(), GmMultiPolyIntersectorUnitTests::testAlongEdgesOutsideToOutside(), GmMultiPolyIntersectorUnitTests::testEdgeThroughOppositeVertexAtAngle(), GmMultiPolyIntersectorUnitTests::testEndAtEdgeFromAdjacent(), GmMultiPolyIntersectorUnitTests::testInsideToEdgeThenThroughAdjacent(), GmMultiPolyIntersectorUnitTests::testInsideToInside(), GmMultiPolyIntersector2IntermediateTests::testLargeNumPolys(), GmMultiPolyIntersector2IntermediateTests::testLargeNumPolysAndSegments(), GmMultiPolyIntersectorUnitTests::testOutsideToOutside(), GmMultiPolyIntersectorUnitTests::testStartAtEdgeThroughAdjacent(), GmMultiPolyIntersectorUnitTests::testTouchesEdge(), and GmMultiPolyIntersectorUnitTests::testTouchesVertex().
|
inline |
Faster than a % operation and we do this a lot.
i | The index. |
Definition at line 54 of file triangles.h.
Referenced by xms::TrTinImpl::CheckAndSwap(), xms::TrTinImpl::NextBoundaryPoint(), xms::TrTinImpl::SwapEdge(), and xms::TrTinImpl::TriangleAdjacentToEdge().
|
inline |
Faster than a % operation and we do this a lot.
i | The index. |
Definition at line 44 of file triangles.h.
Referenced by xms::TrTinImpl::AdjacentTriangle(), xms::TrTinImpl::CheckAndSwap(), xms::TrTinImpl::CheckTriangle(), xms::TrTinImpl::CommonEdgeIndex(), xms::TrBreaklineAdderImpl::FindIntersectingEdgeFromEdge(), xms::TrBreaklineAdderImpl::FindIntersectingEdgeFromPoint(), xms::TrTinImpl::GetBoundaryPoints(), xms::TrOuterTriangleDeleterImpl::MarkNeighbors(), xms::TrTinImpl::OptimizeTriangulation(), xms::TrTinImpl::PreviousBoundaryPoint(), xms::TrBreaklineAdderImpl::ProcessSegmentBySwapping(), xms::TrTinImpl::SwapEdge(), and xms::TrTinImpl::TriangleFromEdge().
Boost serialize function.
Renumbers items in a_vec as if deleting items in a_deleting.
Used to renumber m_tris when deleting points. Subtracts 1 from each point index that is greater than or equal to the point index that was deleted. Goes from greatest number being deleted to least so greater numbers get reduced more with each pass.
[in] | a_deleting | Set of numbered things in a_vec being deleted. |
[in,out] | a_vec | Vector of numbered things that gets renumbered. |
Definition at line 1376 of file TrTin.cpp.
Referenced by xms::TrTinImpl::DeletePoints().
bool xms::trTriangulateIt | ( | TrTriangulator & | a_Client | ) |
Do delaunay divide-and-conquer triangulation.
This algorithm can leave triangles that are VERY long and thin. The default ratio should be sufficient to get rid of these (on the exterior only). If you want to supply a bigger ratio to get rid of more triangles, you can, or you can do it on your own after you call this function.
[in,out] | a_Client | Class used to help do the triangulation. |
Definition at line 2656 of file triangulate.cpp.
References xmlog::error, and XM_LOG.
Referenced by TrTriangulatorPointsUnitTests::test1(), TrTriangulatorPointsUnitTests::test_bug12336(), and xms::TrTriangulator::Triangulate().
Test if two edges are the same ignoring direction.
a_edge1 | The first edge to compare against. |
a_edge2 | The second edge to compare against. |
Definition at line 158 of file XmEdge.cpp.
References xms::XmEdge::IsEquivalent().
Referenced by XmEdgeUnitTests::testIsEquivalent().
std::shared_ptr< XmUGrid > xms::XmReadUGridFromAsciiFile | ( | const std::string & | a_filePath | ) |
Read XmUGrid from an ASCII file.
[in] | a_filePath | filename to read including path, file name, and extension |
Definition at line 563 of file XmUGridUtils.cpp.
References XmReadUGridFromStream().
Referenced by XmUGridUtilsTests::testWriteThenReadUGridFileToAscii().
std::shared_ptr< XmUGrid > xms::XmReadUGridFromStream | ( | std::istream & | a_inStream | ) |
Read an XmUGrid from ASCII text from an input stream.
[in] | a_inStream | the stream to read |
Definition at line 574 of file XmUGridUtils.cpp.
References daReadLine(), xmlog::error, xms::XmUGrid::New(), xms::DaStreamReader::ReadVecInt(), xms::DaStreamReader::ReadVecPt3d(), and XM_LOG.
Referenced by GmMultiPolyIntersector2IntermediateTests::testBug12586(), XmUGridUtilsTests::testLargeUGridBinarySpeed(), XmUGridUtilsTests::testLinear2dWriteThenRead(), XmUGridUtilsTests::testLinear3dWriteThenRead(), XmUGridUtilsTests::testReadBasicUGrid(), XmUGridUtilsTests::testReadEmptyUGridAsciiFile(), XmUGridUtilsTests::testReadPolyhedronUGrid(), XmUGridUtilsTests::testReadReadyAtNextLine(), XmUGridUtilsTests::testReadVersion1Dot0File(), XmUGridUtilsTests::testWriteThenReadUGridBinary(), XmUGridUtilsTests::testWriteThenReadUGridFile(), and XmReadUGridFromAsciiFile().
void xms::XmWriteUGridToAsciiFile | ( | std::shared_ptr< XmUGrid > | a_ugrid, |
const std::string & | a_filePath | ||
) |
Write an XmUGrid to an ASCII file.
[in] | a_ugrid | The UGrid to write to the file |
[in] | a_filePath | filename to write including path, file name, and extension |
Definition at line 630 of file XmUGridUtils.cpp.
References XmWriteUGridToStream().
Referenced by XmUGridUtilsTests::testWriteThenReadUGridFileToAscii().
void xms::XmWriteUGridToStream | ( | std::shared_ptr< XmUGrid > | a_ugrid, |
std::ostream & | a_outStream | ||
) |
Save an XmUGrid ASCII text to output stream.
[in] | a_ugrid | the UGrid to save |
[in] | a_outStream | the stream to write |
Definition at line 640 of file XmUGridUtils.cpp.
Referenced by XmUGridUtilsTests::testLargeUGridBinarySpeed(), XmUGridUtilsTests::testLinear2dWriteThenRead(), XmUGridUtilsTests::testLinear3dWriteThenRead(), XmUGridUtilsTests::testWriteBasicUGrid(), XmUGridUtilsTests::testWriteEmptyUGrid(), XmUGridUtilsTests::testWriteLinear2dCells(), XmUGridUtilsTests::testWriteLinear3dCells(), XmUGridUtilsTests::testWritePolyhedronUGrid(), XmUGridUtilsTests::testWriteThenReadUGridFile(), and XmWriteUGridToAsciiFile().
void xms::XmWriteUGridToStream | ( | const XmUGrid & | a_ugrid, |
std::ostream & | a_outStream, | ||
bool | a_binary | ||
) |
Save an XmUGrid ASCII text to output stream.
[in] | a_ugrid | the UGrid to save |
[in] | a_outStream | the stream to write |
[in] | a_binary | should arrays be written in binary |
Definition at line 650 of file XmUGridUtils.cpp.