The Tiny HTM Library
Data Structures | Functions
include/tinyhtm/geometry.h File Reference

Minimalistic functions and types for spherical geometry. More...

#include "common.h"
+ Include dependency graph for geometry.h:
+ This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Data Structures

struct  htm_v3
 Cartesian coordinates for a point in R3. More...
struct  htm_sc
 Spherical coordinates (in degrees) for a point in S2. More...
struct  htm_s2ellipse
 An ellipse on the sphere, defined as the points p such that angsep(p, f1) + angsep(p, f2) <= 2a, where a is the semi-major axis angle of the ellipse and f1/f2 are its focii. More...
struct  htm_s2cpoly
 A convex polygon on the sphere. More...

Functions

HTM_INLINE enum htm_errcode htm_v3_init (struct htm_v3 *out, double x, double y, double z)
 Stores the given cartesian coordinates in out.
HTM_INLINE enum htm_errcode htm_sc_init (struct htm_sc *out, double lon_deg, double lat_deg)
 Stores the given spherical coordinates in out; coordinates are assumed to be in degrees.
HTM_INLINE void htm_v3_add (struct htm_v3 *out, const struct htm_v3 *v1, const struct htm_v3 *v2)
 Stores the vector sum v1 + v2 in out.
HTM_INLINE void htm_v3_sub (struct htm_v3 *out, const struct htm_v3 *v1, const struct htm_v3 *v2)
 Stores the vector difference v1 - v2 in out.
HTM_INLINE void htm_v3_neg (struct htm_v3 *out, const struct htm_v3 *v)
 Stores the vector v * -1 in out.
HTM_INLINE void htm_v3_mul (struct htm_v3 *out, const struct htm_v3 *v, double s)
 Stores the vector-scalar product v * s in out.
HTM_INLINE void htm_v3_div (struct htm_v3 *out, const struct htm_v3 *v, double s)
 Stores the vector-scalar quotient v / s in out.
HTM_INLINE double htm_v3_dot (const struct htm_v3 *v1, const struct htm_v3 *v2)
 Returns the dot product of the 3-vectors v1 and v2.
HTM_INLINE void htm_v3_cwise_mul (struct htm_v3 *out, const struct htm_v3 *v1, const struct htm_v3 *v2)
 Stores the component-wise product of the 3-vectors v1 and v2 in out.
HTM_INLINE double htm_v3_norm2 (const struct htm_v3 *v)
 Returns the squared norm of the 3-vector v, which must not be NULL.
HTM_INLINE double htm_v3_norm (const struct htm_v3 *v)
 Returns the norm of the 3-vector v, which must not be NULL.
HTM_INLINE void htm_v3_normalize (struct htm_v3 *out, const struct htm_v3 *v)
 Stores a normalized copy of v in out.
HTM_INLINE void htm_v3_rcross (struct htm_v3 *out, const struct htm_v3 *v1, const struct htm_v3 *v2)
 Stores twice the vector cross product of v1 and v2 in out.
HTM_INLINE void htm_v3_cross (struct htm_v3 *out, const struct htm_v3 *v1, const struct htm_v3 *v2)
 Stores the vector cross product of v1 and v2 in out.
enum htm_errcode htm_v3_ne (struct htm_v3 *north, struct htm_v3 *east, const struct htm_v3 *v)
 Computes the N,E basis unit vectors at a position v on the unit sphere.
enum htm_errcode htm_v3_tanrot (double *angle, const struct htm_v3 *v1, const struct htm_v3 *v2, double r)
 Computes the angle by which the normal of the plane defined by the origin, v1 and v2 should be rotated about the axis v1 - v2, such that the resulting plane is tangent to two small circles of radius r centered at v1 and v2.
enum htm_errcode htm_v3_rot (struct htm_v3 *out, const struct htm_v3 *v, const struct htm_v3 *k, double angle_deg)
 Rotates vector v around the axis k by the specified number of degrees.
enum htm_errcode htm_v3_centroid (struct htm_v3 *cen, const struct htm_v3 *points, size_t n)
 Computes the normalized sum of the given list of vectors.
enum htm_errcode htm_sc_tov3 (struct htm_v3 *out, const struct htm_sc *p)
 Converts the spherical coordinate pair p to a unit 3-vector and stores the results in out.
enum htm_errcode htm_v3_tosc (struct htm_sc *out, const struct htm_v3 *v)
 Converts the 3-vector v to spherical coordinates and stores the results in out.
double htm_sc_dist2 (const struct htm_sc *p1, const struct htm_sc *p2)
 Returns the square of the distance between the unit vectors corresponding to points p1 and p2.
double htm_sc_angsep (const struct htm_sc *p1, const struct htm_sc *p2)
 Returns the angular separation (in degrees) between the points p1 and p2.
HTM_INLINE double htm_v3_dist2 (const struct htm_v3 *v1, const struct htm_v3 *v2)
 Returns the square of the distance betwen vectors v1 and v2.
double htm_v3_angsepu (const struct htm_v3 *v1, const struct htm_v3 *v2)
 Returns the angular separation (in degrees) between unit vectors v1 and v2.
double htm_v3_angsep (const struct htm_v3 *v1, const struct htm_v3 *v2)
 Returns the angular separation (in degrees) between vectors v1 and v2, which need not have unit norm.
double htm_v3_edgedist2 (const struct htm_v3 *v, const struct htm_v3 *v1, const struct htm_v3 *v2, const struct htm_v3 *e)
 Returns the minimum square distance between v, and points on the edge from v1 to v2 (where e is a vector parallel to the cross product of v1 and v2).
enum htm_errcode htm_s2ellipse_init (struct htm_s2ellipse *e, const struct htm_v3 *f1, const struct htm_v3 *f2, double a)
 Initializes a spherical ellipse from the given focii (which must be unit vectors) and semi-major axis angle a (in degrees).
enum htm_errcode htm_s2ellipse_init2 (struct htm_s2ellipse *e, const struct htm_v3 *cen, double a, double b, double angle)
 Initializes a spherical ellipse from the given center (which must be a unit vector), axis angles a and b, and orientation (in degrees).
HTM_INLINE int htm_s2ellipse_cv3 (const struct htm_s2ellipse *e, const struct htm_v3 *v)
 Returns 1 if the spherical ellipse e contains vector v, and 0 otherwise.
HTM_INLINE struct htm_v3htm_s2cpoly_verts (struct htm_s2cpoly *poly)
 Returns a pointer to the vertices of poly.
HTM_INLINE struct htm_v3htm_s2cpoly_edges (struct htm_s2cpoly *poly)
 Returns a pointer to the edges of poly.
struct htm_s2cpolyhtm_s2cpoly_init (const struct htm_v3 *verts, size_t n, enum htm_errcode *err)
 Creates a polygon from a list of at least 3 vertices.
struct htm_s2cpolyhtm_s2cpoly_box (const struct htm_v3 *cen, double width, double height, double angle, enum htm_errcode *err)
 Creates a polygon corresponding to the box with the given width, height, center and orientation.
struct htm_s2cpolyhtm_s2cpoly_ngon (const struct htm_v3 *cen, double r, size_t n, enum htm_errcode *err)
 Creates an n-gon inscribed in the circle of the given radius and center.
struct htm_s2cpolyhtm_s2cpoly_line (const struct htm_v3 *v1, const struct htm_v3 *v2, double r, enum htm_errcode *err)
 Creates a polygon corresponding to the quadrilateral completely enclosing two small circles with the given centers and a radius of r degrees.
int htm_s2cpoly_cv3 (const struct htm_s2cpoly *poly, const struct htm_v3 *v)
 Returns 1 if the spherical convex polygon poly contains vector v, and 0 otherwise.
double htm_s2cpoly_area (const struct htm_s2cpoly *poly)
 Returns the area (in steradians) enclosed by the given spherical convex polygon.
struct htm_s2cpolyhtm_s2cpoly_clone (const struct htm_s2cpoly *poly)
 Returns a copy of the given polygon, or NULL if poly is NULL or the required memory allocation failed.
enum htm_errcode htm_s2cpoly_pad (struct htm_s2cpoly *poly, double r)
 Pads poly by the angle r (degrees).
int htm_v3_hemispherical (const struct htm_v3 *points, size_t n, enum htm_errcode *err)
 Returns 1 if the given set of points (which need not consist of unit vectors) is hemispherical and 0 otherwise.
int htm_v3_convex (const struct htm_v3 *points, size_t n, enum htm_errcode *err)
 Tests whether an ordered list of points form a spherical convex polygon.
struct htm_s2cpolyhtm_s2cpoly_hull (const struct htm_v3 *points, size_t n, enum htm_errcode *err)
 Creates a polygon corresponding to the convex hull of the given point set.

Detailed Description

Minimalistic functions and types for spherical geometry.

Authors:
Serge Monkewitz

Definition in file geometry.h.


Function Documentation

double htm_s2cpoly_area ( const struct htm_s2cpoly poly)

Returns the area (in steradians) enclosed by the given spherical convex polygon.

If poly is NULL, 0.0 is returned.

Definition at line 606 of file geometry.c.

References htm_v3_dot(), htm_v3_norm(), htm_v3_rcross(), htm_s2cpoly::n, and htm_s2cpoly::ve.

struct htm_s2cpoly* htm_s2cpoly_box ( const struct htm_v3 cen,
double  width,
double  height,
double  angle,
enum htm_errcode err 
) [read]

Creates a polygon corresponding to the box with the given width, height, center and orientation.

The orientation is specified as a rotation angle in degrees around the box center. An angle of zero degrees will result in a box aligned to the north (height) and east (width) axes at the box center.

To release resources for the polygon, call free() on the returned pointer.

Returns:
  • a newly allocated polygon on success.
  • NULL if an error occurs or the inputs are invalid. In this case, *err is additionally set to indicate the reason for the failure.

Definition at line 429 of file geometry.c.

References HTM_EANG, HTM_ENULLPTR, HTM_EZERONORM, htm_s2cpoly_init(), htm_v3_ne(), htm_v3_norm2(), htm_v3_normalize(), htm_v3_rcross(), and htm_v3_rot().

struct htm_s2cpoly* htm_s2cpoly_clone ( const struct htm_s2cpoly poly) [read]

Returns a copy of the given polygon, or NULL if poly is NULL or the required memory allocation failed.

To release resources for the polygon, call free() on the returned pointer.

Definition at line 628 of file geometry.c.

References htm_s2cpoly::n.

int htm_s2cpoly_cv3 ( const struct htm_s2cpoly poly,
const struct htm_v3 v 
)

Returns 1 if the spherical convex polygon poly contains vector v, and 0 otherwise.

Arguments must not be NULL pointers.

Definition at line 593 of file geometry.c.

References htm_v3_dot(), htm_s2cpoly::n, and htm_s2cpoly::ve.

Referenced by htm_tree_s2cpoly_count(), and htm_tree_s2cpoly_scan().

struct htm_s2cpoly* htm_s2cpoly_hull ( const struct htm_v3 points,
size_t  n,
enum htm_errcode err 
) [read]

Creates a polygon corresponding to the convex hull of the given point set.

Points must be specified as unit vectors.

To release resources for the polygon, call free() on the returned pointer.

Returns:
  • a newly allocated polygon on success.
  • NULL if an error occurs or the inputs are invalid. In this case, *err is additionally set to indicate the reason for the failure.

Definition at line 1389 of file geometry.c.

References HTM_EDEGEN, HTM_EHEMIS, HTM_ELEN, HTM_ENOMEM, HTM_ENULLPTR, HTM_OK, htm_v3_add(), htm_v3_centroid(), htm_v3_dist2(), htm_v3_div(), htm_v3_dot(), htm_v3_hemispherical(), htm_v3_norm(), htm_v3_norm2(), htm_v3_rcross(), htm_s2cpoly::n, v, htm_s2cpoly::ve, htm_s2cpoly::vsum, htm_v3::x, htm_v3::y, and htm_v3::z.

struct htm_s2cpoly* htm_s2cpoly_init ( const struct htm_v3 verts,
size_t  n,
enum htm_errcode err 
) [read]

Creates a polygon from a list of at least 3 vertices.

Vertices can be in clockwise or counter-clockwise order, but:

  • vertices must be unit vectors
  • vertices must to be hemispherical
  • vertices must define edges that do not intersect except at vertices
  • vertices must define edges forming a convex polygon

To release resources for the polygon, call free() on the returned pointer.

Returns:
  • a newly allocated polygon on success.
  • NULL if an error occurs or the inputs are invalid. In this case, *err is additionally set to indicate the reason for the failure.

Definition at line 369 of file geometry.c.

References HTM_ELEN, HTM_ENOMEM, HTM_ENULLPTR, HTM_OK, htm_v3_add(), htm_v3_dot(), htm_v3_neg(), htm_v3_rcross(), htm_s2cpoly::n, htm_s2cpoly::ve, and htm_s2cpoly::vsum.

Referenced by htm_s2cpoly_box(), htm_s2cpoly_line(), and htm_s2cpoly_ngon().

struct htm_s2cpoly* htm_s2cpoly_line ( const struct htm_v3 v1,
const struct htm_v3 v2,
double  r,
enum htm_errcode err 
) [read]

Creates a polygon corresponding to the quadrilateral completely enclosing two small circles with the given centers and a radius of r degrees.

The quadrilateral is oriented in the direction v2 - v1.

To release resources for the polygon, call free() on the returned pointer.

Returns:
  • a newly allocated polygon on success.
  • NULL if an error occurs or the inputs are invalid. In this case, *err is additionally set to indicate the reason for the failure.

Definition at line 553 of file geometry.c.

References HTM_OK, htm_s2cpoly_init(), htm_v3_normalize(), htm_v3_rcross(), htm_v3_rot(), htm_v3_sub(), and htm_v3_tanrot().

struct htm_s2cpoly* htm_s2cpoly_ngon ( const struct htm_v3 cen,
double  r,
size_t  n,
enum htm_errcode err 
) [read]

Creates an n-gon inscribed in the circle of the given radius and center.

To release resources for the polygon, call free() on the returned pointer.

Returns:
  • a newly allocated polygon on success.
  • NULL if an error occurs or the inputs are invalid. In this case, *err is additionally set to indicate the reason for the failure.

Definition at line 489 of file geometry.c.

References HTM_EANG, HTM_ELEN, HTM_ENOMEM, HTM_ENULLPTR, HTM_EZERONORM, HTM_RAD_PER_DEG, htm_s2cpoly_init(), htm_v3_ne(), htm_v3_norm2(), htm_v3_normalize(), htm_v3::x, htm_v3::y, and htm_v3::z.

HTM_INLINE int htm_s2ellipse_cv3 ( const struct htm_s2ellipse e,
const struct htm_v3 v 
)

Returns 1 if the spherical ellipse e contains vector v, and 0 otherwise.

Arguments must not be NULL pointers.

Definition at line 398 of file geometry.h.

Referenced by htm_tree_s2ellipse_count(), and htm_tree_s2ellipse_scan().

enum htm_errcode htm_s2ellipse_init2 ( struct htm_s2ellipse e,
const struct htm_v3 cen,
double  a,
double  b,
double  angle 
)

Initializes a spherical ellipse from the given center (which must be a unit vector), axis angles a and b, and orientation (in degrees).

The orientation is the position angle (east of north) of the first axis with respect to the north pole.

Definition at line 298 of file geometry.c.

References htm_s2ellipse::a, htm_s2ellipse::cen, HTM_EANG, HTM_ENULLPTR, HTM_OK, HTM_RAD_PER_DEG, htm_v3_add(), htm_v3_mul(), htm_v3_ne(), htm_v3_sub(), htm_v3::x, htm_s2ellipse::xx, htm_s2ellipse::xy, htm_s2ellipse::xz, htm_v3::y, htm_s2ellipse::yy, htm_s2ellipse::yz, htm_v3::z, and htm_s2ellipse::zz.

double htm_sc_angsep ( const struct htm_sc p1,
const struct htm_sc p2 
)

Returns the angular separation (in degrees) between the points p1 and p2.

Arguments must not not be NULL pointers, but may alias.

Definition at line 197 of file geometry.c.

References HTM_DEG_PER_RAD, htm_sc_dist2(), and htm_v3::x.

double htm_sc_dist2 ( const struct htm_sc p1,
const struct htm_sc p2 
)

Returns the square of the distance between the unit vectors corresponding to points p1 and p2.

Arguments must not be NULL pointers, but may alias.

Definition at line 183 of file geometry.c.

References HTM_RAD_PER_DEG, htm_sc::lat, htm_sc::lon, htm_v3::x, htm_v3::y, and htm_v3::z.

Referenced by htm_sc_angsep().

double htm_v3_angsep ( const struct htm_v3 v1,
const struct htm_v3 v2 
)

Returns the angular separation (in degrees) between vectors v1 and v2, which need not have unit norm.

Arguments must not be NULL pointers, but may alias.

Definition at line 216 of file geometry.c.

References HTM_DEG_PER_RAD, htm_v3_cross(), htm_v3_dot(), and htm_v3_norm().

Referenced by htm_tri_init(), htm_v3_convex(), and htm_v3_tanrot().

double htm_v3_angsepu ( const struct htm_v3 v1,
const struct htm_v3 v2 
)

Returns the angular separation (in degrees) between unit vectors v1 and v2.

Arguments must not be NULL pointers, but may alias.

Definition at line 206 of file geometry.c.

References HTM_DEG_PER_RAD, htm_v3_dist2(), and htm_v3::x.

Referenced by htm_s2ellipse_init().

int htm_v3_convex ( const struct htm_v3 points,
size_t  n,
enum htm_errcode err 
)

Tests whether an ordered list of points form a spherical convex polygon.

Returns:
  • 1: points form a spherical convex polygon and are in counter-clockwise order.
  • 0: points do not form a spherical convex polygon, or an error occurred. In case of error, *err is set to indicate the reason for the failure.
  • -1: points form a spherical convex polygon and are in clockwise order.

Definition at line 1189 of file geometry.c.

References HTM_EDEGEN, HTM_ELEN, HTM_OK, htm_v3_angsep(), htm_v3_centroid(), htm_v3_div(), htm_v3_dot(), htm_v3_hemispherical(), htm_v3_norm2(), and htm_v3_rcross().

HTM_INLINE double htm_v3_dist2 ( const struct htm_v3 v1,
const struct htm_v3 v2 
)

Returns the square of the distance betwen vectors v1 and v2.

Arguments must not be NULL pointers, but may alias.

Definition at line 323 of file geometry.h.

Referenced by htm_s2cpoly_hull(), htm_tree_s2circle_count(), htm_tree_s2circle_scan(), htm_v3_angsepu(), and htm_v3_edgedist2().

double htm_v3_edgedist2 ( const struct htm_v3 v,
const struct htm_v3 v1,
const struct htm_v3 v2,
const struct htm_v3 e 
)

Returns the minimum square distance between v, and points on the edge from v1 to v2 (where e is a vector parallel to the cross product of v1 and v2).

The vectors v, v1, and v2 are assumed to be normalized, but e need not have unit norm.

Definition at line 231 of file geometry.c.

References htm_v3_cross(), htm_v3_dist2(), htm_v3_dot(), htm_v3_norm2(), htm_v3::x, and htm_v3::y.

int htm_v3_hemispherical ( const struct htm_v3 points,
size_t  n,
enum htm_errcode err 
)

Returns 1 if the given set of points (which need not consist of unit vectors) is hemispherical and 0 otherwise.

If an error occurs, 0 is returned and *err is set to indicate the reason for the failure.

Definition at line 1061 of file geometry.c.

References HTM_ELEN, HTM_ENOMEM, HTM_ENULLPTR, and HTM_OK.

Referenced by htm_s2cpoly_hull(), htm_s2cpoly_pad(), and htm_v3_convex().