1.1.0
CCC
 
GeometryHelper Class Reference

#include <GeometryHelper.hpp>

Inheritance diagram for GeometryHelper:
GH

Static Public Member Functions

static void rotationGetSinCos (double angle, double *cosa, double *sina)
 
static void rotationMatrixIdentityInPlace (int ndim, VectorDouble &rot)
 
static void rotation2DMatrixInPlace (double angle, VectorDouble &rot)
 
static void rotation3DMatrixInPlace (double alpha, double beta, double gamma, VectorDouble &rot)
 
static void rotationMatrixInPlace (int ndim, const VectorDouble &angles, VectorDouble &rot)
 
static VectorDouble rotationMatrix (int ndim, const VectorDouble &angles)
 
static void rotationGetAnglesInPlace (int ndim, const double *rot, double *angles)
 
static void rotationGetAnglesInPlace (const VectorDouble &rot, VectorDouble &angles)
 
static void rotationCopy (int ndim, const double *rotin, double *rotout)
 
static bool rotationIsIdentity (int ndim, double *rot, double eps=EPSILON10)
 
static MatrixSquareGeneral EulerToRotation (const VectorDouble &angles, const ERotation &convrot=ERotation::fromKey("SXYZ"))
 
static void rotationGetRandomDirection (double ct, double st, double *a, double *codir)
 
static void rotationGetDirection2D (const VectorDouble &angles, VectorDouble &codir)
 
static void rotationGetDirectionDefault (int ndim, VectorDouble &codir)
 
static void rotationGetAnglesFromCodirInPlace (const VectorDouble &codir, VectorDouble &angles)
 
static VectorDouble rotationGetAngles (const VectorDouble &codir, bool flagResize=false)
 
static VectorDouble rotationToEuler (const MatrixSquareGeneral &mat, const ERotation &convrot=ERotation::fromKey("SXYZ"), double eps=EPSILON10)
 
static void mergeBoxes (VectorDouble &mini1, VectorDouble &maxi1, VectorDouble &mini2, VectorDouble &maxi2)
 
static double distancePointToSegment (double x0, double y0, double x1, double y1, double x2, double y2, double *xd, double *yd, int *nint)
 
static bool segmentIntersect (double xd1, double yd1, double xe1, double ye1, double xd2, double yd2, double xe2, double ye2, double *xint, double *yint)
 
static double geodeticAngularDistance (double long1, double lat1, double long2, double lat2, double radius=1.)
 
static void geodeticAngles (double long1, double lat1, double long2, double lat2, double long3, double lat3, double *a, double *b, double *c, double *A, double *B, double *C)
 
static double geodeticTrianglePerimeter (double long1, double lat1, double long2, double lat2, double long3, double lat3)
 
static double geodeticTriangleSurface (double long1, double lat1, double long2, double lat2, double long3, double lat3)
 
static bool isInSphericalTriangle (double *coor, double surface, double *pts1, double *pts2, double *pts3, double *wgts, double eps=EPSILON6)
 
static bool isInSphericalTriangleOptimized (const double *coor, double *ptsa, double *ptsb, double *ptsc, double *wgts, double eps=EPSILON6)
 
static VectorVectorDouble convertLongLat (const VectorDouble &longitude, const VectorDouble &latitude, double dilate=1., double radius_arg=1.)
 
static void convertCart2Sph (double x, double y, double z, double *rlong, double *rlat, double radius_arg=1.)
 
static void convertSph2Cart (double rlong, double rlat, double *x, double *y, double *z, double radius_arg=1.)
 
static MatrixSquareGeneral gradXYToRotmat (double dzoverdx, double dzoverdy)
 
static MatrixRectangulargetDirectionsInR3 (const MatrixRectangular *U)
 
static MatrixRectangulargetDirectionsInRn (const MatrixRectangular *U)
 
static double formatAngle (double anglein, double basis=360.)
 
static VectorDouble formatAngles (const VectorDouble &anglesin, double basis=360.)
 
static VectorDouble rayTriangleIntersect (const VectorDouble &dir, const VectorDouble &v0, const VectorDouble &v1, const VectorDouble &v2)
 
static VectorVectorDouble sphBarCoord (const VectorVectorDouble &sphPts, const MatrixRectangular &apices, const MatrixInt &meshes)
 
static double getCosineAngularTolerance (double tolang)
 

Member Function Documentation

void GeometryHelper::convertCart2Sph ( double  x,
double  y,
double  z,
double *  rlong,
double *  rlat,
double  radius_arg = 1. 
)
static

Convert the cartesian coordinates into spherical coordinates

Parameters
[in]xFirst cartesian coordinate
[in]ySecond cartesian coordinate
[in]zThird cartesian coordinate
[in]radius_argRadius of the sphere (Earth if TEST)
[out]rlongLongitude (in degrees)
[out]rlatLatitude (in degrees)
VectorVectorDouble GeometryHelper::convertLongLat ( const VectorDouble longitude,
const VectorDouble latitude,
double  dilate = 1.,
double  radius_arg = 1. 
)
static

Returns the Vector of Sample coordinates in 3-D from Longitude-Latitude

Parameters
longitudeArray of longitude values
latitudeArray of latitude values
dilateDilation applied to radius
radius_argRadius (if note defined, taken from variety definition)
Returns
void GeometryHelper::convertSph2Cart ( double  rlong,
double  rlat,
double *  x,
double *  y,
double *  z,
double  radius_arg = 1. 
)
static

Convert the spherical coordinates into cartesian coordinates

Parameters
[in]rlongLongitude (in degrees)
[in]rlatLatitude (in degrees)
[in]radius_argradius of the sphere (Earth if TEST)
[out]xFirst cartesian coordinate
[out]ySecond cartesian coordinate
[out]zThird cartesian coordinate
double GeometryHelper::distancePointToSegment ( double  x0,
double  y0,
double  x1,
double  y1,
double  x2,
double  y2,
double *  xd,
double *  yd,
int *  nint 
)
static

Find the shortest distance between the point (x0,y0) and the segment with the two end points (x1,y1) and (x2,y2)

Returns
Minimum algebraic distance (positive or negative)
Parameters
[in]x0,y0Coordinates of the target point
[in]x1,y1Coordinate of the first end-point of the segment
[in]x2,y2Coordinate of the second end-point of the segment
[out]xd,ydCoordinates of the closest point
[out]nint=1 if the projection belongs to the segment =0 if it is set to one of the segment vertices
MatrixSquareGeneral GeometryHelper::EulerToRotation ( const VectorDouble angles,
const ERotation &  convrot = ERotation::fromKey("SXYZ") 
)
static

Returns the Rotation matrix, starting from the Euler angles

Parameters
anglesOrdered list of Euler angles
convrotRotation convention
Returns
Remarks
The code is coming from the following reference (BSD license)
https://github.com/matthew-brett/transforms3d/blob/master/transforms3d/euler.py
double GeometryHelper::formatAngle ( double  anglein,
double  basis = 360. 
)
static

Format an Angle to be lying in [0, basis]

Parameters
angleinInput angle value
basisBasis (should be 360 by default; could be 180 for variogram or covariance)
Returns
VectorDouble GeometryHelper::formatAngles ( const VectorDouble anglesin,
double  basis = 360. 
)
static
void GeometryHelper::geodeticAngles ( double  long1,
double  lat1,
double  long2,
double  lat2,
double  long3,
double  lat3,
double *  a,
double *  b,
double *  c,
double *  A,
double *  B,
double *  C 
)
static

Calculate all geodetic angles from a spherical triangle

Parameters
[in]long1Longitude of the first point (in degrees)
[in]lat1Latitude of the first point (in degrees)
[in]long2Longitude of the second point (in degrees)
[in]lat2Latitude of the second point (in degrees)
[in]long3Longitude of the third point (in degrees)
[in]lat3Latitude of the third point (in degrees)
[out]aAngle (P2,O,P3)
[out]bAngle (P3,O,P1)
[out]cAngle (P1,O,P2)
[out]AAngle (P2,P1,P3)
[out]BAngle (P3,P2,P1)
[out]CAngle (P1,P3,P2)
double GeometryHelper::geodeticAngularDistance ( double  long1,
double  lat1,
double  long2,
double  lat2,
double  radius = 1. 
)
static

Calculate the geodetic angular distance between two points on the sphere

Returns
Angular distance
Parameters
[in]long1Longitude of the first point (in degrees)
[in]lat1Latitude of the first point (in degrees)
[in]long2Longitude of the second point (in degrees)
[in]lat2Latitude of the second point (in degrees)
[in]radiusRadius of the sphere
double GeometryHelper::geodeticTrianglePerimeter ( double  long1,
double  lat1,
double  long2,
double  lat2,
double  long3,
double  lat3 
)
static

Calculate the perimeter of the spherical triangle

Returns
The Perimeter of the spherical triangle
Parameters
[in]long1Longitude of the first point (in degrees)
[in]lat1Latitude of the first point (in degrees)
[in]long2Longitude of the second point (in degrees)
[in]lat2Latitude of the second point (in degrees)
[in]long3Longitude of the third point (in degrees)
[in]lat3Latitude of the third point (in degrees)
double GeometryHelper::geodeticTriangleSurface ( double  long1,
double  lat1,
double  long2,
double  lat2,
double  long3,
double  lat3 
)
static

Calculate the surface of the spherical triangle

Returns
The Surface of the spherical triangle (with unit radius)
Parameters
[in]long1Longitude of the first point (in degrees)
[in]lat1Latitude of the first point (in degrees)
[in]long2Longitude of the second point (in degrees)
[in]lat2Latitude of the second point (in degrees)
[in]long3Longitude of the third point (in degrees)
[in]lat3Latitude of the third point (in degrees)
double GeometryHelper::getCosineAngularTolerance ( double  tolang)
static

Returns the cosine of the angular tolerance

Parameters
[in]tolangAngular tolerance
MatrixRectangular * GeometryHelper::getDirectionsInR3 ( const MatrixRectangular U)
static
MatrixRectangular * GeometryHelper::getDirectionsInRn ( const MatrixRectangular U)
static

Function to compute directions in Rd

Parameters
Ua matrix [n, d] of uniform values between [0,1]. 'n' is the number of direction and 'd' is the dimension of the space
Returns
A vector of a matrix [n, d] of the coordinates of the 'n' directions in the Euclidean space Rd.
MatrixSquareGeneral GeometryHelper::gradXYToRotmat ( double  dzoverdx,
double  dzoverdy 
)
static

Calculate the rotation matrix starting from the partial derivatives along X and Y of a tilted plane

Parameters
dzoverdxPartial derivative along X
dzoverdyPartial derivative along Y
bool GeometryHelper::isInSphericalTriangle ( double *  coor,
double  surface,
double *  pts1,
double *  pts2,
double *  pts3,
double *  wgts,
double  eps = EPSILON6 
)
static

Is a point inside a spherical triangle

Returns
1 if the point belongs to the spherical triangle; 0 otherwise
Parameters
[in]coorCoordinates of the target point (long,lat)
[in]surfaceSurface of the spherical triangle
[in]pts1Coordinates of the first point of the triangle
[in]pts2Coordinates of the second point of the triangle
[in]pts3Coordinates of the third point of the triangle
[in]epsTolerance
[out]wgtsArray of weights
bool GeometryHelper::isInSphericalTriangleOptimized ( const double *  coor,
double *  ptsa,
double *  ptsb,
double *  ptsc,
double *  wgts,
double  eps = EPSILON6 
)
static

Is a point inside a spherical triangle

Returns
True if the point belongs to the spherical triangle; False otherwise
Parameters
[in]coorCoordinates of the target point (long,lat)
[in]ptsaCoordinates of the first point of the triangle
[in]ptsbCoordinates of the second point of the triangle
[in]ptscCoordinates of the third point of the triangle
[in]epsTolerance
[out]wgtsArray of weights
void GeometryHelper::mergeBoxes ( VectorDouble mini1,
VectorDouble maxi1,
VectorDouble mini2,
VectorDouble maxi2 
)
static

Merge the extensions of the boxes (parallel to main axes)

Parameters
[in]mini1Input array containing the minimum along each axis
[in]maxi1Input array containing the maximum along each axis
[in]mini2Output array containing the minimum along each axis
[in]maxi2Output array containing the maximum along each axis
VectorDouble GeometryHelper::rayTriangleIntersect ( const VectorDouble dir,
const VectorDouble v0,
const VectorDouble v1,
const VectorDouble v2 
)
static
void GeometryHelper::rotation2DMatrixInPlace ( double  angle,
VectorDouble rot 
)
static

Calculates the 2-D rotation matrix

Parameters
[in]angleRotation angle (in degrees)
[out]rotRotation matrix (Dimension = 4)
void GeometryHelper::rotation3DMatrixInPlace ( double  alpha,
double  beta,
double  gamma,
VectorDouble rot 
)
static

Calculates the 3-D rotation matrix

Parameters
[in]alphaangle (in degrees) / oz
[in]betaangle (in degrees) / oy'
[in]gammaangle (in degrees) / ox''
[out]rotdirect rotation matrix (Dimension = 9)
void GeometryHelper::rotationCopy ( int  ndim,
const double *  rotin,
double *  rotout 
)
static

Copy a rotation matrix

Parameters
[in]ndimSpace dimension
[in]rotinInput rotation matrix
[out]rotoutOutput rotation matrix (already allocated)
VectorDouble GeometryHelper::rotationGetAngles ( const VectorDouble codir,
bool  flagResize = false 
)
static

From the vector of direction coefficients (codir) returns the vector of angles

Parameters
codirInput vector giving the direction coefficients
flagResizeWhen TRUE (and if in 2-D) the returned vector is resized to 1
Returns
void GeometryHelper::rotationGetAnglesFromCodirInPlace ( const VectorDouble codir,
VectorDouble angles 
)
static

Calculates the rotation angle from the direction coefficient

Parameters
[in]codirDirection vector (Dimension = ndim)
[out]anglesRotation angles (Dimension = ndim)
void GeometryHelper::rotationGetAnglesInPlace ( int  ndim,
const double *  rot,
double *  angles 
)
static

Calculates the rotation angles from the rotation matrix

Parameters
[in]ndimSpace dimension
[in]rotRotation matrix (Dimension = 9)
[out]anglesRotation angles (Dimension = ndim)
void GeometryHelper::rotationGetAnglesInPlace ( const VectorDouble rot,
VectorDouble angles 
)
static
void GeometryHelper::rotationGetDirection2D ( const VectorDouble angles,
VectorDouble codir 
)
static

Convert angles to a set of Directions

Parameters
[in]anglesVector giving the angles characteristics (in degrees) As this provides rotation in 2D, 'angles' is dimensioned to 'ndir' (one angle requires a single value)
[out]codirVector of the direction (Dim: ndir * ndim)
Remarks
'ndir' is given by the dimension of 'angles', 'ndim' is given by 'codir'
void GeometryHelper::rotationGetDirectionDefault ( int  ndim,
VectorDouble codir 
)
static

Create a Direction (used as default)

Parameters
[in]ndimNumber of space dimensions
[out]codirVector of the direction (Dim: ndim)
void GeometryHelper::rotationGetRandomDirection ( double  ct,
double  st,
double *  a,
double *  codir 
)
static

Rotation of a Direction in 3-D

Parameters
[in]ct,stCosine and Sine of the rotation angle
[in]aRandom direction
[in,out]codirDirection to be rotated
void GeometryHelper::rotationGetSinCos ( double  angle,
double *  cosa,
double *  sina 
)
static

Calculates the trigonometric features

Parameters
[in]angleinput angle (in degrees)
[out]cosacosine function
[out]sinasine function
bool GeometryHelper::rotationIsIdentity ( int  ndim,
double *  rot,
double  eps = EPSILON10 
)
static

Starting from a rotation matrix, check it is different from the Identity

Returns
False if a rotation is defined; True if it is an Identity
Parameters
[in]ndimSpace dimension
[in]rotRotation matrix
[in]epsTolerance
VectorDouble GeometryHelper::rotationMatrix ( int  ndim,
const VectorDouble angles 
)
static

Calculates the rotation matrix. Returns the rotation matrix as a VectorDouble

Parameters
[in]ndimSpace dimension
[in]anglesArray of angles
void GeometryHelper::rotationMatrixIdentityInPlace ( int  ndim,
VectorDouble rot 
)
static

Initialize a rotation matrix

Parameters
[in]ndimSpace dimension
[out]rotRotation matrix
void GeometryHelper::rotationMatrixInPlace ( int  ndim,
const VectorDouble angles,
VectorDouble rot 
)
static

Calculates the rotation matrix

Parameters
[in]ndimSpace dimension
[in]anglesArray of angles
[out]rotdirect rotation matrix (dimensionned to ndim*ndim)
VectorDouble GeometryHelper::rotationToEuler ( const MatrixSquareGeneral M,
const ERotation &  convrot = ERotation::fromKey(                                          "SXYZ"),
double  eps = EPSILON10 
)
static

Returns the Euler angles, starting from a rotation matrix

Parameters
MInput matrix
convrotRotation convention
epsTolerance
Returns
Remarks
The code is coming from the following reference (BSD license)
https://github.com/matthew-brett/transforms3d/blob/master/transforms3d/euler.py
bool GeometryHelper::segmentIntersect ( double  xd1,
double  yd1,
double  xe1,
double  ye1,
double  xd2,
double  yd2,
double  xe2,
double  ye2,
double *  xint,
double *  yint 
)
static

Check if two 2-D segments intersect

Returns
True if there is an intersection; False otherwise
Parameters
[in]xd1,yd1Starting point for the first segment
[in]xe1,ye1Ending point for the first segment
[in]xd2,yd2Starting point for the second segment
[in]xe2,ye2Ending point for the second segment
[out]xint,yintCoordinates of the intersection
VectorVectorDouble GeometryHelper::sphBarCoord ( const VectorVectorDouble sphPts,
const MatrixRectangular apices,
const MatrixInt meshes 
)
static

Calculate the Barycenter coordinates of points with in the Spherical Meshing

Parameters
sphPtsCoordinates of target samples (dim = np * 3)
apicesCoordinates of the apices of the Meshing
meshesMesh information of the Meshing
Returns
A vector of barycenters (dimension = 4 * np) where
  • 1: rank of the triangle
  • 2, 3, 4: barycentric coordinates

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