pygplates.GreatCircleArc
- class pygplates.GreatCircleArc
Bases:
instance
A great-circle arc on the surface of the unit globe.
Great circle arcs are equality (
==
,!=
) comparable (but not hashable - cannot be used as a key in adict
).A GreatCircleArc can also be pickled.
Changed in version 0.42: Added pickle support.
- __init__(start_point, end_point)
Create a great circle arc from two points.
- Parameters:
start_point (
PointOnSphere
orLatLonPoint
or tuple (latitude,longitude), in degrees, or tuple (x,y,z)) – the start point of the arc.end_point (
PointOnSphere
orLatLonPoint
or tuple (latitude,longitude), in degrees, or tuple (x,y,z)) – the end point of the arc.
- Raises:
IndeterminateResultError if points are antipodal (opposite each other)
An arc is specified by a start-point and an end-point:If these two points are not antipodal, a unique great-circle arc (with angle-span less than PI radians) will be determined between them. If they are antipodal then IndeterminateResultError will be raised. Note that an error is not raised if the two points are coincident.great_circle_arc = pygplates.GreatCircleArc(start_point, end_point)
Methods
__init__
(start_point, end_point)Create a great circle arc from two points.
get_arc_direction
(...)Return the direction along the arc at a point on the arc.
Returns the arc length of this great circle arc (in radians).
get_arc_point
(...)Return a point on this arc.
Return the arc's end point geometry.
Return the unit vector normal direction of the great circle this arc lies on.
Return the rotation axis of the arc as a 3D vector.
Return the (latitude, longitude) equivalent of
get_rotation_axis()
.Return the arc's start point geometry.
Return whether this great circle arc is of zero length.
to_tessellated
(tessellate_radians)Returns a list of
points
tessellated from this great circle arc such that adjacent points are separated by no more than tessellate_radians on the globe.to_uniform_points
(point_spacing_radians, ...)Returns a sequence of points uniformly spaced along this great circle arc.
- get_arc_direction(normalised_distance_from_start_point)
Return the direction along the arc at a point on the arc.
- Parameters:
normalised_distance_from_start_point (float) – distance from start point where zero is the start point, one is the end point and between zero and one are points along the arc
- Returns:
the unit-length 3D vector
- Return type:
- Raises:
ValueError if arc normalised_distance_from_start_point is not in the range [0,1]
- Raises:
IndeterminateGreatCircleArcDirectionError if arc is zero length
The returned direction is tangential to the Earth’s surface and is aligned with the direction of the great circle arc (in the direction going from the start point towards the end point). This direction is perpendicular to the great circle normal direction (see
get_great_circle_normal()
).The direction at the midpoint of an arc:
if not arc.is_zero_length(): arc_midpoint_direction = arc.get_arc_direction(0.5)
If normalised_distance_from_start_point is zero then the direction at start point is returned. If normalised_distance_from_start_point is one then the direction at end point is returned. Values of normalised_distance_from_start_point between zero and one return directions at points on the arc. If normalised_distance_from_start_point is outside the range from zero to one then then ValueError is raised.
- get_arc_length()
Returns the arc length of this great circle arc (in radians).
- Return type:
float
To convert to distance, multiply the result by the Earth radius (see
Earth
).
- get_arc_point(normalised_distance_from_start_point)
Return a point on this arc.
- Parameters:
normalised_distance_from_start_point (float) – distance from start point where zero is the start point, one is the end point and between zero and one are points along the arc
- Return type:
- Raises:
ValueError if arc normalised_distance_from_start_point is not in the range [0,1]
The midpoint of an arc:
arc_midpoint = arc.get_arc_point(0.5)
If normalised_distance_from_start_point is zero then the start point is returned. If normalised_distance_from_start_point is one then the end point is returned. Values of normalised_distance_from_start_point between zero and one return points on the arc. If normalised_distance_from_start_point is outside the range from zero to one then then ValueError is raised.
- get_end_point()
Return the arc’s end point geometry.
- Return type:
- get_great_circle_normal()
Return the unit vector normal direction of the great circle this arc lies on.
- Returns:
the unit-length 3D vector
- Return type:
- Raises:
IndeterminateGreatCircleArcNormalError if arc is zero length
if not arc.is_zero_length(): normal = arc.get_great_circle_normal()
Note
This returns the same (x, y, z) result as
get_rotation_axis()
, but in the form of aVector3D
instead of an (x, y, z) tuple.Note
The normal to the great circle can be considered to be the tangential direction (to the Earth’s surface) at any point along the great circle arc that is most pointing away from (perpendicular to) the direction of the arc (from start point to end point - see
get_arc_direction()
).The normal vector is the same direction as the
cross product
of the start point and the end point. In fact it is equivalent topygplates.Vector3D.cross(arc.start_point().to_xyz(), arc.end_point().to_xyz()).to_normalised()
.If the arc start and end points are the same (if
is_zero_length()
isTrue
) then IndeterminateGreatCircleArcNormalError is raised.See also
- get_rotation_axis()
Return the rotation axis of the arc as a 3D vector.
- Returns:
the unit-length 3D vector (x,y,z)
- Return type:
the tuple (float, float, float)
- Raises:
IndeterminateArcRotationAxisError if arc is zero length
if not arc.is_zero_length(): axis_x, axis_y, axis_z = arc.get_rotation_axis()
Note
This returns the same (x, y, z) result as
get_great_circle_normal()
, but in the form of an (x, y, z) tuple instead of aVector3D
.The rotation axis is the unit-length 3D vector (x,y,z) returned in the tuple.
The rotation axis direction is such that it rotates the start point towards the end point along the arc (assuming a right-handed coordinate system).
If the arc start and end points are the same (if
is_zero_length()
isTrue
) then IndeterminateArcRotationAxisError is raised.See also
- get_rotation_axis_lat_lon()
Return the (latitude, longitude) equivalent of
get_rotation_axis()
.- Returns:
the axis as (latitude, longitude)
- Return type:
the tuple (float, float)
- Raises:
IndeterminateArcRotationAxisError if arc is zero length
if not arc.is_zero_length(): axis_lat, axis_lon = arc.get_rotation_axis_lat_lon()
The rotation axis is the (latitude, longitude) returned in the tuple.
The rotation axis direction is such that it rotates the start point towards the end point along the arc (assuming a right-handed coordinate system).
If the arc start and end points are the same (if
is_zero_length()
isTrue
) then IndeterminateArcRotationAxisError is raised.
- get_start_point()
Return the arc’s start point geometry.
- Return type:
- is_zero_length()
Return whether this great circle arc is of zero length.
- Return type:
bool
If this arc is of zero length, it will not have a determinate rotation axis and a call to
get_rotation_axis()
will raise an error.
- to_tessellated(tessellate_radians)
Returns a list of
points
tessellated from this great circle arc such that adjacent points are separated by no more than tessellate_radians on the globe.- Parameters:
tessellate_radians (float) – maximum tessellation angle (in radians)
- Return type:
list
points
- Raises:
ValueError if tessellate_radians is negative or zero
Note
If this great circle arc subtends an angle less than tessellate_radians then only its
start point
andend point
are returned. For example, this applies to azero length
arc. Otherwise tessellated points within this arc are also returned.Tessellate a great circle arc to 2 degrees:
tessellation_points = great_circle_arc.to_tessellated(math.radians(2))
Note
Since a GreatCircleArc is immutable it cannot be modified. Which is why a tessellated list of PointOnSphere is returned.
See also
- to_uniform_points(point_spacing_radians[, first_point_spacing_radians=0.0][, return_segment_interpolations=False])
Returns a sequence of points uniformly spaced along this great circle arc.
- Parameters:
point_spacing_radians (float) – spacing between points (in radians)
first_point_spacing_radians (float) – Spacing of first uniform point from this arc’s start point (in radians). By default the first uniform point coincides with this arc’s start point. Ideally this is non-negative (but, for example, if it’s slightly negative then the first uniform point will be slightly off this arc near its start point but still on its great circle).
return_segment_interpolations (bool) – whether to also return information about the polyline segment that each uniform point is on - default is
False
- Returns:
list of points, or (if return_segment_interpolations is
True
) a 2-tuple containing a list of points and a list of segment interpolations (where each uniform point is located, on this great circle arc, in the range [0,1])- Return type:
list of
PointOnSphere
, or tuple (list ofPointOnSphere
, list of float) if return_segment_interpolations isTrue
- Raises:
ValueError if point_spacing_radians is negative or zero
Note
The distance (along the arc) between the last uniform point and the arc’s end point can be less than point_spacing_radians (since the length of the arc minus first_point_spacing_radians might not be an integer multiple of point_spacing_radians).
Note
If first_point_spacing_radians is greater than thearc's length
then no uniform points will be generated.And if the arc iszero length
and first_point_spacing_radians is zero then a single uniform point will be generated.Create points uniformly spaced by 1 degree along a great circle arc starting 0.5 degrees from its start point:
uniform_points = arc.to_uniform_points( math.radians(1), first_point_spacing_radians = math.radians(0.5))
Next, we extend the above example by associating an arc direction, tangential to the globe, at each uniform point:
uniform_points, uniform_point_segment_interpolations = arc.to_uniform_points( math.radians(1), first_point_spacing_radians = math.radians(0.5)) return_segment_interpolations = True) # We end up with a list of 3D direction vectors (with a list length equal to the number of uniform points). uniform_point_arc_directions = [arc.get_arc_direction(segment_interpolation) for segment_interpolation in uniform_point_segment_interpolations]
See also
Added in version 0.47.