SatcomCoordinates.jl
Documentation for SatcomCoordinates.jl
.
SatcomCoordinates.AER Type
AER{CRS <: AbstractCRS, T}
Represents a position in the Azimuth-Elevation-Range (AER) coordinate system, which is a local coordinate system centered at a point on or above the surface of an Ellipsoid. The Elevation and Azimuth angles are always defined w.r.t. the ENU CRS with the same origin. More specifically:
Note
This is actually a type alias for SphericalCRS{ENU{CRS, T}, AzEl{ENU{CRS, T}}}
, but is provided for convenience.
Properties/Units/Aliases
az
=> u"°" => (azimuth,)el
=> u"°" => (elevation,)r
=> u"m" => (distance, range)
Constructor
AER(args...) = SphericalCRS(AzEl(ENU(args...)))
The constructor simply forwards all the arguments to the ENU
constructor and then wraps the result in a SphericalCRS
with the AzEl
pointing CRS.
See also: ENU
, NED
, AzEl
, SphericalCRS
SatcomCoordinates.Abstract2DPointingCRS Type
Abstract2DPointingCRS{CRS <: AbstractCRS} <: AbstractPointingCRS{CRS}
Abstract type representing any pointing type defined over a 3D Cartesian CRS that is defined by two coordinates (e.g. theta/phi, azimuth/elevation, etc.).
See also: AbstractPointingCRS
, UV
, ThetaPhi
, AzOverEl
, ElOverAz
, AzEl
SatcomCoordinates.AbstractCRS Type
AbstractCRS
Abstract type representing an arbitrary coordinate reference system. This is intended to be used as type parameter for all coordinate types to allow consistency check between operation among coordinates and in some cases simplify conversion of a given coordinate to different CRSs.
All concrete subtypes of AbstractCRS
are expected to either by directly a 3-dimensional cartesian CRS (i.e. CRS whose units are all lengths) or be based on a 3-dimensional cartesian CRS.
The reference 3-dimensional cartesian CRS should always be accessible from any custom subtype of AbstractCRS
via the cartesiancrs
function.
By default, a CRS is considered cartesian if all its properties have units of length. Specific customization of whether a custom CRS is considered or not cartesian can be done by adding a custom method to the iscartesiancrs
function (which must return either true
or false
).
No-Argument Constructor
Each concrete subtype of AbstractCRS
is expected to define a no-argument constructor without any type parameter (i.e. CRS()
) which returns the default instance of said CRS.
Extended Help
A subtype of an AbstractCRS
shall satisfy the following conditions:
It can be cartesian or not cartesian (and should have a valid method for
iscartesiancrs
representing this)It can be a root crs, or a derived crs. In the latter case, it should have the wrapped CRS type as it's first type parameter, and contain the specific instance of the wrapped CRS as a field within it.
Derived CRSs should only depend on a single wrapped CRS. Complex derivations shall be implemented by nesting CRSs rather than by having multiple different wrapped CRSs as fields.
Derived CRSs should be wrapping a Cartesian CRS. Any kind of -non-cartesian CRS shall only be the last level of nesting.
Derived CRSs shall have a valid method for the
raw_linkedcrs_transform
function, which should return a transformation to express a coordinate in the derived CRS to the equivalent ones expressed in the wrapped CRS.
Additionally, for each CRS whose type parameters can be completely identified by their fields, only one inner constructor without explicit type parameters should be defined. .
sourceSatcomCoordinates.AbstractCRSTransform Type
AbstractCRSTransform
Abstract type representing a coordinate transform between two CRSs with numbertype T
.
SatcomCoordinates.AbstractLinkedCRS Type
AbstractLinkedCRS{CRS <: AbstractCRS} <: AbstractCRS
Abstract type representing a CRS that is linked to another CRS.
Each instance of a concrete subtype of AbstractLinkedCRS
should have at least one field whose type subtypes AbstractCRS
.
The first field represents an AbstractCRS
is automatically returned by the [linkedcrs](@ref)(crs::AbstractLinkedCRS)
function.
SatcomCoordinates.AbstractPointingCRS Type
AbstractPointingCRS{CRS <: AbstractCRS} <: AbstractCRS
Abstract type representing any pointing type defined over a 3D Cartesian CRS. Different pointing types identify different ways of identifying a position on the unitary sphere (over the specified CRS)
Although these are not strictly speaking CRSs themselves, they are considered a subtype of AbstractCRS
to better fit within the package interface.
See also: AbstractCRS
, AngularPointingCRS
SatcomCoordinates.AbstractRawCRSTransform Type
AbstractRawCRSTransform <: Transform
Abstract type representing a raw transform between two CRSs. A raw transform is the lower level CRS transformation that operate directly on the raw coordinates of the two CRSs (i.e. the one returned by tuplecoords(coords::Coordinate{CRS})
)
A raw transform is used inernally by the user-facing AbstractCRSTransform
and is expected to receive as input a NTuple{N, T <: AbstractFloat}
and produce as output a NTuple{M, T <:AbstractFloat}
, where N
and M
are the number of coordinates (i.e. the output of ncoords(crs)
) of the input and output CRSs respectively.
See extended help for more details on necessary methods for new concrete subtypes of AbstractRawCRSTransform
.
Extended Help
Each new concrete subtype of AbstractRawCRSTransform
(indicated with type name RAW_T
in the rest of this section) is expected to implement a method for the following two methods:
ncoords_out(::Type{<:RAW_T})
: Returns the number of coordinates in the output NTuple (i.e. theM
above)ncoords_in(::Type{<:RAW_T})
: Returns the number of coordinates in the input NTuple (i.e. theN
above)
For transforms where N === M
, it is sufficient to just implement ncoords(::Type{<:RAW_T})
as that is the default fallback for both ncoords_out
and ncoords_in
.
Additionally, new Raw transforms must implement the appropriate interface functions from TransformsBase.jl
which as a minimum requires a valid method for:
TransformsBase.apply(t::RAW_T, tup::NTuple{N, <:AbstractFloat}) where N
: Remember that this function must provide 2 separate outputs, the first is the actual transformed coordinate asNTuple{M, <:AbstractFloat}
and the second is thecache
output of theTransformsBase
interface. For all expected subtypes ofAbstractRawCRSTransform
, this second output should simply benothing
.
SatcomCoordinates.AbstractSatcomCoordinate Type
AbstractSatcomCoordinate{CRS <: AbstractCRS, T, N}
General abstract type identifying a coordinate with N
dimension and defined with respect to a specific Coordinate Reference System CRS
. The parameter T
represent the underlying numeric type (e.g. machine precision) of the coordinate components.
Note
The term coordinate is used here in a loose sense, identifying both position in space as well as pointing directions
SatcomCoordinates.AffineCartesian Type
AffineCartesian{CRSₗ, CRS, T} <: AbstractLinkedCRS{CRS}
A Cartesian CRS that is linked to another Cartesian CRS through a user-defined affine transform.
Fields
linked::CRSₗ <: AbstractCRS
: The linked CRSbase::CRS <: AbstractCRS
: The base CRStransform::T <: RawAffineTransform
: The affine transform that links the base CRS to the linked CRS
The two CRSs used to define the AffineLinkedCRS instance must be Cartesian CRSs (i.e. they must have iscartesiancrs(crs) == true
).
For linking a non-cartesian base CRS, simply use the AffineCartesian
as the base for further nesting. As an example, if one desires to use an LLA
CRS as the base for an AffineCartesian
CRS, simply create an AffineCartesian
CRS with the underlying ECEF
CRS and use that for wrapping the LLA
CRS.
SatcomCoordinates.AnyN Type
struct AnyN end
Singletone structure just used to match any integer number
sourceSatcomCoordinates.AzEl Type
AzEl{CRS <: AbstractCRS} <: Abstract2DPointingCRS{CRS}
Object specifying a pointing direction in "Elevation/Azimuth" coordinates, defined following the convention used for Azimuth-Elevation-Range (AER
) coordinates used by MATLAB and by this package.
This represents the azimuth/elevation definition often used for describing pointing from a user terminal on ground towards a satellite. They are very similar to the ElOverAz
definition but with a rotation of the underlying CRS such that elevation is 90° in the direction of the +Z axis and azimuth is computed on the XY plane from the +Y axis towards the +X axis.
Assuming u
, v
, and w
to be direction cosines of the pointing versor ̂p
, their relation with the El
and Az
angles is:
u = cos(El) * sin(Az)
v = cos(El) * cos(Az)
w = sin(El)
Properties/Units/Aliases
az
=> u"°" => (azimuth,)el
=> u"°" => (elevation,)
See also: ThetaPhi
, UV
, ElOverAz
, AzOverEl
SatcomCoordinates.AzOverEl Type
AzOverEl{CRS <: AbstractCRS} <: Abstract2DPointingCRS{CRS}
Object specifying a pointing direction in "Azimuth over Elevation" coordinates, which specify the elevation and azimuth angles that needs to be fed to an azimuth-over-elevation positioner for pointing to a target towards the pointing direction ̂p.
Following the convention used in most Antenna-related literature, the elevation and azimuth are 0° in the direction of the +Z axis of the reference frame.
Assuming u
, v
, and w
to be direction cosines of the pointing versor ̂p
, their relation with the Az
and El
angles is:
u = -sin(Az)
v = cos(Az) * sin(El)
w = cos(Az) * cos(El)
Note
The equations above are used to represent the "Azimuth over Elevation" coordinates in GRASP. Some textbooks, however, use the opposite convention, meaning that the same equations (with a possible flip in the az/u sign) are used to describe an "Elevation over Azimuth" coordinate system. This is for example the case in the book "Theory and Practice of Modern Antenna Range Measurements" by Clive Parini et al.
Properties/Units/Aliases
az
=> u"°" => (azimuth,)el
=> u"°" => (elevation,)
Note
The fields of AzOverEl
objects can also be accessed via getproperty
using the azimuth
and elevation
aliases.
SatcomCoordinates.Cartesian Type
Cartesian <: AbstractCRS
Generic Cartesian CRS, for use in cases that do not require any specific identification of a CRS/Position
sourceSatcomCoordinates.Coordinate Type
Coordinate{CRS <: AbstractCRS, T, N} <: AbstractSatcomCoordinate{CRS, T, N}
A generic coordinate type, which wraps a CRS and a tuple of coordinates.
The interpretation of the coordinate depends on the specific CRS it is defined on (e.g. a Coordinate over a PointingCRS can be considered a Pointing direction)
The type parameter CRS
is the CRS type, T
is the type of the coordinates, and N
is the number of coordinates.
The type parameter CRS
is a subtype of AbstractCRS
, and the type parameter T
is a subtype of AbstractFloat
.
Constructors
Coordinate(crs::AbstractCRS, coords...)
(crs::AbstractCRS)(coords...)
Coordinates can be constructed directly with the Coordinate
constructor by providing a CRS instance as first argument and the specific numerical coordinates (either as tuple/svector or as separate values) as remaining argument[s]
Alternatively (and usually easier) an instance of a specific CRS can be directly used with the numerical coordinates to construct a Coordinate
instance
Examples
julia> using SatcomCoordinates
julia> cartcrs = Cartesian();
julia> cartcrs(1, 2, 3)
Coordinate{Cartesian}:
x = 1.0 m
y = 2.0 m
z = 3.0 m
julia> cartcrs((1,2,3)) == Coordinate(cartcrs, 1, 2, 3)
true
Extended Help
The coordinates over the specific CRS are stored internally as a plain NTuple{N, T}
(where T <: AbstractFloat
is currently enforced during construction) which holds the raw coordinates (i.e. without unit and transformed to corresponding SI unit (for angles we transform to radians). This means that for example if a specific CRS has the unit of u"km"
, the raw values stored internally are the numerical value in meters)
The constructor basically goes over the following 3 steps:
Eventually strip units and normalize to unitless SI floating point values using the
preprocess_input_coords
function.- This function should not need to be customized for custom CRSs and only operates on the specific units associated to each properties specified with the
@define_properties
macro.
- This function should not need to be customized for custom CRSs and only operates on the specific units associated to each properties specified with the
Perform additional checks and manipulation on the raw coordinates based on the specific CRS instance provided to the constructor, using the
process_unitless_coords
function.- Not all custom CRSs require special processing and the default method for the
process_unitless_coords
function simply returns the input tuple as is. An example of CRSs that need to do some further processing are the concrete subtypes ofAbstractPointingCRS
which wraps the angles to always have them within the range[-π, π]
- Not all custom CRSs require special processing and the default method for the
SatcomCoordinates.ECEF Type
ECEF{ID} <: AbstractCRS
This represents an Ellipsoid-Centered Ellipsoid-Fixed (ECEF) CRS which is uniquely identified by its only field id::ID
.
This is a generalization of the conventianal ECEF acronym which is only referred to earth-centered coordinates.
The id
field can be any object that uniquely identifies a specific instance of ECEF CRS
When not specified, the default ID for ECEF CRSs is an instance of the singleton type EarthDefault
.
Properties/Units/Aliases
x
=> u"m" => (,): X coordinate in the ECEF framey
=> u"m" => (,): Y coordinate in the ECEF framez
=> u"m" => (,): Z coordinate in the ECEF frame
Constructor
ECEF(crs)
ECEF(x, y, z)
See also: frameid
, ellipsoidparams
, EarthDefault
, ECI
, LLA
Extended Help
Conversion between ECEF and LLA
To allow proper conversion betwen ECEF and LLA, a valid method for:
SatcomCoordinates.ellipsoidparams(id::ID)
must be defined. See the docstrings of ellipsoidparams
for more details.
Additionally, the conversion is implement inside an extension and relies on the package SatelliteToolboxTransformations
to be loaded to perform the actual conversion.
Conversion between ECEF and ECI
For conversion between coordinates in ECEF and ECI CRSs, the Rotation matrix between the two frames at a specific instane must be provided.
If that is already available (if it was already computed by other means), this can be directly passed to the
change_crs
function as a keyword argumentR_eci_to_ecef
.Alternatively, for custom ECI/ECEF frames, a custom function to compute this rotation matrix based on the IDs of the specific ECI and ECEF frames needs defined as an additional method to the
SatcomCoordinates.eci_to_ecef_rotation
function.
This package only defines a single method of this function where both ECI and ECEF frames are based on the EarthDefault
frame. Check the docstring of SatComCoordinates.eci_to_ecef_rotation
for more details.
SatcomCoordinates.ECI Type
ECI{ID} <: AbstractCRS
This represents an Ellipsoid-Centered Inertial (ECI) CRS which is uniquely identified by its only field id::ID
.
This is a generalization of the conventianal ECI acronym which is only referred to earth-centered coordinates.
The id
field can be any object that uniquely identifies a specific instance of ECI CRS
When not specified, the default ID for ECEF CRSs is an instance of the singleton type EarthDefault
.
Conversion between ECEF and ECI
For conversion between coordinates in ECEF and ECI CRSs, the Rotation matrix between the two frames at a specific instane must be provided.
If that is already available (if it was already computed by other means), this can be directly passed to the
change_crs
function as a keyword argumentR_eci_to_ecef
.Alternatively, for custom ECI/ECEF frames, a custom function to compute this rotation matrix based on the IDs of the specific ECI and ECEF frames needs defined as an additional method to the
SatcomCoordinates.eci_to_ecef_rotation
function.
This package only defines a single method of this function where both ECI and ECEF frames are based on the EarthDefault
frame. Check the docstring of SatComCoordinates.eci_to_ecef_rotation
for more details.
Properties/Units/Aliases
x
=> u"m" => (,): X coordinate in the ECI framey
=> u"m" => (,): Y coordinate in the ECI framez
=> u"m" => (,): Z coordinate in the ECI frame
See also: frameid
, ellipsoidparams
, EarthDefault
, ECEF
, eci_to_ecef_rotation
SatcomCoordinates.ENU Type
ENU{CRS <: AbstractCRS, T} <: AbstractTopocentricCRS{CRS, T}
Represents a position in the East-North-Up (ENU) coordinate system, which is a local coordinate system centered at a point on or above the surface of an Ellipsoid.
This type stores inside of it it's origin (both in ECEF and LLA) as well as the rotation matrix to convert between a point from ENU to it's linked ECEF CRS.
Properties/Units/Aliases
e
=> u"m" => (east, x)n
=> u"m" => (north, y)u
=> u"m" => (up, z)
Constructor
ENU(origin::Coordinate)
Constructs a ENU
CRS by providing it's origin as Coordinate expressed either in an ECEF or LLA CRS.
Convenience methods
Each custom CRS which satisfies the topocentric trait (for which istopocentriccrs(crs) == true
) must implement the following methods:
ecef_origin(crs::CRS)
: This function shall return the coordinate (in the linked ECEF CRS) of the origin of the providedcrs
.lla_origin(crs::CRS)
: This function shall return the coordinate (in the linked LLA CRS) of the origin of the providedcrs
.
Examples
enu_crs = ENU(LLA(0,0,1200km)) # ENU centered around a point at 1200km
enu_coord = enu_crs(0, 0, 100km) # A point 0m East, 0m North, 100km Up from the origin of the CRS
change_crs(LLA(), enu_coord) == LLA(0, 0, 1300km) # The same point expressed in LLA CRS
See also: istopocentriccrs
, ecef_origin
, lla_origin
, NED
SatcomCoordinates.EarthDefault Type
EarthDefault
This represents the default frame ID representing the Earth for both the ECEF and ECI CRSs
It assumes the WGS84 ellipsoid and assumes that the frames from conversion between ECEF and ECI are the following:
ECEF: ITRF
ECI: GCRS
See also: ellipsoidparams
, ECEF
, ECI
SatcomCoordinates.ElOverAz Type
ElOverAz{CRS <: AbstractCRS} <: Abstract2DPointingCRS{CRS}
Object specifying a pointing direction in "Elevation over Azimuth" coordinates, which specify the azimuth and elevation angles that needs to be fed to an elevation-over-azimuth positioner for pointing to a target towards the pointing direction ̂p.
Following the convention used in most Antenna-related literature, the elevation and azimuth are 0° in the direction of the +Z axis of the reference frame.
Assuming u
, v
, and w
to be direction cosines of the pointing versor ̂p
, their relation with the El
and Az
angles is:
u = -sin(Az) * cos(El)
v = sin(El)
w = cos(Az) * cos(El)
Note
The equations above are used to represent the "Elevation over Azimuth" coordinates in GRASP. Some textbooks, however, use the opposite convention, meaning that the same equations (with sometimes an additional flip in the az/u sign) are used to describe an "Azimuth over Elevation" coordinate system. This is for example the case in the book "Theory and Practice of Modern Antenna Range Measurements" by Clive Parini et al.
Properties/Units/Aliases
az
=> u"°" => (azimuth,)el
=> u"°" => (elevation,)
See also: AzOverEl
, ThetaPhi
, UV
SatcomCoordinates.FieldOrCoordinate Type
abstract type FieldOrCoordinate end
Abstract type representing either a field or a coordinate.
This is wrap all possible types within this package that are not CRSs themselves. It is used to contraint custom methods we have for Base
functions.
SatcomCoordinates.LLA Type
LLA{CRS<:AbstractCRS} <: AbstractLinkedCRS{CRS}
A coordinate system that represents a point in latitude, longitude and altitude over a specific ellipsoid. It is derived from an ECEF CRS and wraps it.
Properties/Units/Aliases
lat
=> u"°" => (latitude, lat, l): Latitudelon
=> u"°" => (longitude, lon, l): Longitudealt
=> u"m" => (altitude, alt, h, height): Altitude
Constructor
LLA(wrapped_crs::AbstractCRS)
The constructor expects a CRS as input, which must be an ECEF CRS.
See also: ECEF
, frameid
, ellipsoidparams
, EarthDefault
Extended Help
Conversion between ECEF and LLA
To allow proper conversion betwen ECEF and LLA, a valid method for:
SatcomCoordinates.ellipsoidparams(id::ID)
must be defined. See the docstrings of ellipsoidparams
for more details.
Additionally, the conversion is implement inside an extension and relies on the package SatelliteToolboxTransformations
to be loaded to perform the actual conversion.
SatcomCoordinates.NED Type
NED{CRS <: AbstractCRS, T} <: AbstractTopocentricCRS{CRS, T}
Represents a position in the North-East-Down (NED) coordinate system, which is a local coordinate system centered at a point on or above the surface of an Ellipsoid.
This type stores inside of it it's origin (both in ECEF and LLA) as well as the rotation matrix to convert between a point from NED to it's linked ECEF CRS.
Properties/Units/Aliases
n
=> u"m" => (north, x)e
=> u"m" => (east, y)d
=> u"m" => (down, z)
Constructor
NED(origin::Coordinate)
Constructs a NED
CRS by providing it's origin as Coordinate expressed either in an ECEF or LLA CRS.
Convenience methods
Each custom CRS which satisfies the topocentric trait (for which istopocentriccrs(crs) == true
) must implement the following methods:
ecef_origin(crs::CRS)
: This function shall return the coordinate (in the linked ECEF CRS) of the origin of the providedcrs
.lla_origin(crs::CRS)
: This function shall return the coordinate (in the linked LLA CRS) of the origin of the providedcrs
.
See also: istopocentriccrs
, ecef_origin
, lla_origin
, ENU
SatcomCoordinates.PointingTransform Type
abstract type PointingTransform <: AbstractRawCRSTransform end
Abstract type for all transformations that used to convert between a 2D pointing CRS to the DirectionCosines one (and vice-versa).
These are raw transforms and are also used to go to map to the
sourceSatcomCoordinates.Raw Type
Raw{C <: FieldOrCoordinate} <: FieldOrCoordinate
Structure that is only used to wrap a FieldOrCoordinate
object and allow to access its raw values (i.e. normalized and without units) directly via Base.getproperty
.
Example
using SatcomCoordinates
# Define our custom Cartesian CRS
struct CustomKM <: AbstractCRS end
# We specify that this has x,y,z properties which have default unit of km (i.e. plain numbers are interpreted as km). The raw coordinates (i.e. how they are stored internally in the coordinate instances) are actually floating point values represented in meters (as that is the SI unit for length)
@define_properties CustomKM [
x => u"km"
y => u"km"
z => u"km"
]
# Create a position in the custom CRS at x = 1km, y = 2km, z = 3km
pkm = Position(CartesianKM(), 1,2,3)
# Normal property access
pkm.x === 1.0u"km" # true
# If you want to directly access the raw unitless data (e.g. inside of hot loops), you can wrap the coordinate in `Raw`
(; y) = Raw(pkm) # This uses the julia desctructuring synthax that relies on `Base.getproperty`
y === 2000.0 # true
SatcomCoordinates.RawAffineTransform Type
RawAffineTransform{R, T} <: Transform
A transform that is defined by a rotation and a translation that is intended to operate on SVectors or NTuples of size N
Arguments
rotation::Union{Identity, RotMatrix{N, P}}
: The rotation of the transform.translation::Union{Identity, NTuple{N, P}, SVector{N, P}}
: The translation of the transform.
SatcomCoordinates.RawComposedTransform Type
RawComposedTransform{T1 <: Transform, T2 <: Transform} <: AbstractRawCRSTransform
A transform that is the composition of two raw transforms.
This type is automatically created when composing two raw (but non Identity) transforms that are not both affine (i.e. that do not satisfy isaffinetransform(t) == true
).
This is for example the case when composing an affine transform with a transformation between LLA and ECEF or between Spherical and Cartesian coordinates.
Fields
t1::T1
: The first transform.t2::T2
: The second transform.
When applied to an object, the operation order is t2(t1(obj))
SatcomCoordinates.SphericalCRS Type
SphericalCRS{CRS <: AbstractPointingCRS, PT <: Abstract2DPointingCRS{CRS}} <: AbstractLinkedCRS{CRS}
A generic spherical CRS, which wraps a Cartesian CRS as root and has properties that are specified based on the pointing CRS PT
Constructor
SphericalCRS(pointing_crs::Abstract2DPointingCRS)
The only constructor expects directly a 2D pointing CRS as input and will automatically extracts its base Cartesian CRS.
Example
# Generate a Spherical CRS based on ThetaPhi over a plain Cartesian CRS
SphericalCRS() # This is equivalent to SphericalCRS(ThetaPhi())
# Generate one where the pointing direction is `AzEl`
SphericalCRS(AzEl())
# Make the so-called `AER` CRS which is AzEl over an ENU frame
aer_crs = SphericalCRS(AzEl(ENU(LLA(0,0,1200km))))
SatcomCoordinates.ThetaPhi Type
ThetaPhi{CRS <: AbstractCRS} <: Abstract2DPointingCRS{CRS}
An object specifying a pointing direction in ThetaPhi coordinates over the Cartesian CRS CRS
, defined as the θ and φ in the (ISO/Physics definition) spherical coordinates representation
Assuming u
, v
, and w
to be direction cosines of the pointing versor ̂p
, their relation with the θ
and φ
angles is:
u = sin(θ) * cos(φ)
v = sin(θ) * sin(φ)
w = cos(θ)
Properties/Units/Aliases
θ
=> u"°" => (theta, t): The polar angleφ
=> u"°" => (phi, p, ϕ): The azimuthal angle
While the field name use the greek letters, the specific fields of an arbitrary ThetaPhi
object tp
can be accessed with alternative symbols:
tp.θ
,tp.theta
andtp.t
can be used to access theθ
fieldtp.φ
,tp.ϕ
,tp.phi
andtp.p
can be used to access theφ
field
See also: DirectionCosines
, UV
SatcomCoordinates.UV Type
UV{CRS <: AbstractCRS} <: Abstract2DPointingCRS{CRS}
Specify a pointing direction in UV coordinates over the cartesian CRS CRS
(which must be a Cartesian CRS). U,V Coordinates are equivalent to the direction cosines with respect to the X
and Y
axis of the reference frame CRS
. They can also be related to the spherical coordinates (ISO/Physics) spherical coordinates representation by the following equations:
u = sin(θ) * cos(φ)
v = sin(θ) * sin(φ)
Valid Directions
UV coordinates can only be used to represent pointing direction in the half-hemisphere containing the cartesian +Z axis.
Constructor Tolerance
The inputs values must also satisfy u^2 + v^2 <= 1 + SatcomCoordinates.UV_CONSTRUCTOR_TOLERANCE[]
or an error will be thrown. The SatcomCoordinates.UV_CONSTRUCTOR_TOLERANCE
is a Ref{Float64}
which defaults to 1e-5 (In case u^2 + v^2 > 1
the inputs are normalized to ensure u^2 + v^2 = 1
).
Properties/Units/Aliases
u
=> NoUnitsv
=> NoUnits
See also: AbstractPointingCRS
, ThetaPhi
SatcomCoordinates.add_unit Method
add_unit(userunit::Unitful.Units, refunit::Unitful.Units, val::Real)
Take a real value, add the refunit
to it and then convert it to userunit
using Unitful.uconvert
.
This is used internally in the coords
function to transform the raw coordinates stored in a Coordinate
type into unitful values according to the units of the specific CRS.
SatcomCoordinates.cartesiancrs Method
cartesiancrs(C::Type{<:AbstractCRS})
Returns the first Cartesian CRS type associated to the provided CRS type C
.
A Cartesian CRS is a CRS which is defined in a 3D cartesian space.
Any new custom Cartesian CRS type should define a method for this function that returns the Cartesian CRS type when given the custom CRS type as input.
See also: cartesiancrs
, rootcrs
, linkedcrs
SatcomCoordinates.compose Method
compose(t1, t2, ts...)
Compose the transforms t1
, t2
, ts...
as if they were applied in argument order (e.g. first t1
then t2
then ts...
).
SatcomCoordinates.ecef_origin Method
ecef_origin(crs::AbstractCRS)
ecef_origin(obj::FieldOrCoordinate)
Return the coordinate (in ECEF CRS) of the origin of the provided obj
or crs
.
Custom topocentric CRSs must implement a specific method for this function.
For all other CRSs, this function will simply traverse the chaing of linked CRSs starting from crs
until it either finds a valid topocentric CRS or it reaches the end of the chain (In which case it will simply throw an error).
See also lla_origin
SatcomCoordinates.ecefcrs Method
ecefcrs
CRS Trait function to represent CRSs which are Ellipsoid-Centered-Ellipsoid-Fixed (ECEF).
sourceSatcomCoordinates.eci_to_ecef_rotation Function
eci_to_ecef_rotation(eci_id, ecef_id; kwargs...)
This function is expected to return the rotation matrix that transforms ECI coordinates into ECEF coordinates.
Returned Matrix
The returned matrix is automatically converted to Rotations.RotMatrix3
so the ouptut of this should either by already a RotMatrix3
or something that can be implicitly converted to it..
This function is called when trying to transform between ECI and ECEF coordinates.
Any translation between ECI and ECEF coordinates will extract the rotation matrix by calling this function with the id
of the ECI and ECEF frames respectively and then simply apply it to the raw tuple coordinates (directly from ECI to ECEF, inverting the matrix before for ECEF to ECI)
All the auxiliary data required to compute the rotation matrix shall be provided via keyword arguments which are simply propagated from the change_crs
call (and via the lower_level transform_tuplecoords
one).
SatcomCoordinates.ecicrs Method
ecicrs
CRS Trait function to represent CRSs which are Ellipsoid-Centered-Inertial (ECI).
sourceSatcomCoordinates.ellipsoidparams Method
ellipsoidparams(frame_id)
Function that shall have a valid method for all valid id
(of either ECEF or ECI CRS) instances and shall return a NamedTuple with the following fields representing the useful ellipsoid parameters:
a
: Semimajor axisf
: Flatteningb
: Semiminor axise²
: First eccentricity squaredel²
: Second eccentricity squared
Example
Calling this function on the default Earth frames provides the output for the WGS84 ellipsoid:
julia> using SatcomCoordinates: ellipsoidparams, EarthDefault
julia> ellipsoidparams(EarthDefault())
(a = 6.378137e6, f = 0.0033528106647474805, b = 6.356752314245179e6, e² = 0.0066943799901413165, el² = 0.006739496742276434)
SatcomCoordinates.frameid Method
frameid(crs::AbstractCRS)
Extracts the frame identifier associated to the CRS crs
. This is currently only used to extract the id
field from the ECI and ECEF CRSs for further processing (e.g. extracting the ellipsoid parameters associated to the frame id)
See also: ellipsoidparams
, EarthDefault
SatcomCoordinates.getcrstransform Method
getcrstransform(traitfunc::Function, crsᵢ::AbstractCRS)
getcrstransform(traitfunc::Function, coord::FieldOrCoordinate)
Return the transformation taking as input a Coordinate in the input crs crsᵢ
(eventually extracted from the coordinate coord
) and returning as output a Coordinate in the output CRS identified as the first CRS within the nested CRS tree of crsᵢ
which satisfies the provided CRS trait traitfunc
.
If a CRS satisfying traitfunc
can not be found within the CRS tree of crsᵢ
, an error is thrown.
The returned transformation is an instance of CRSTransform
.
Note
This function simply calls getcrstransform_raw
and wraps the result in a CRSTransform
object.
Example
using SatcomCoordinates
# We first create a NED CRS at a specific location above Earth
enu_crs = ENU(LLA(0,0,1200km))
# We then create a Spherical CRS (AzEl) that is linked to the ENU CRS. This is a double nested CRS as it's itself based on a ENU which is based on an ECEF CRS.
aer_crs = SphericalCRS(AzEl(enu_crs))
getcrs(linkedcrs, aer_crs) == enu_crs # The linked CRS is the one immediately below the provided CRS, which is the ENU CRS
getcrs(rootcrs, aer_crs) == ECEF() # The root CRS is the one at the bottom of the nested CRS, which is the ECEF CRS
# Extract the transformation towards the linked CRS (NED)
tlinked = getcrstransform(linkedcrs, aer_crs)
# 90° el, 10m distance is (0,0,10) in ENU
tlinked(aer_crs(0,90,10)) ≈ enu_crs(0,0,10) # The transformation is applied to the coordinate
# We now extract the transformation towards the root CRS (ECEF)
troot = getcrstransform(rootcrs, aer_crs)
# Check that the origin is correctly transformed
troot(aer_crs(0,90,0)) ≈ ecef_origin(aer_crs)
See also: CRSTransform
, getcrstransform_raw
SatcomCoordinates.getcrstransform_raw Method
getcrstransform_raw(traitfunc::Function, crsᵢ::AbstractCRS)
getcrstransform_raw(traitfunc::Function, coord::FieldOrCoordinate)
Return the raw transformation that translates raw coordinates (i.e. tuples of numbers) from the input crs crsᵢ
(eventually extracted from the coordinate coord
) to the first CRS within the nested CRS tree of crsᵢ
which satisfies the provided trait traitfunc
.
If the CRS crsᵢ
is already satisfying the trait traitfunc
, this function returns the identity transformation.
If a CRS satisfying traitfunc
can not be found within the CRS tree of crsᵢ
, an error is thrown.
See also: getcrstransform
SatcomCoordinates.hascrstrait Method
hascrstrait(traitfunc::Function, crs::AbstractCRS)
This is an helper function which should be used in the definition of the default method for a trait function. For example, the istopocentriccrs
trait function's default method is actually defined as follows:
istopocentriccrs(obj) = hascrstrait(istopocentriccrs, obj)
hascrstrait(::typeof(istopocentriccrs), ::Type{<:AbstractCRS}) = false
This is mainly done to simplify the handling of complex linked CRSs (like AffineCartesian) which should actually forward all the trait functions to their base CRS.
Custom CRSs that simply need to add a method to an existing CRS trait function do not need to use this function, and this is only needed to be extended by a downstream user/developer when defining completely new CRS traits.
The three basic CRS traits below are not defined using this function as they rely on core properties of all CRSs for their default implementation:
sourceSatcomCoordinates.is_same_crs Method
is_same_crs(crs1::AbstractCRS, crs2::AbstractCRS)
Returns true
or false
to indicate whether the two provided CRSs are the same.
Defaults to true for CRSs of the same type and false otherwise.
This is done in place of simply doing crs1 == crs2
because it's slightly faster when the check can be done at compile time (on the types only).
Custom CRSs which may be different despite having the same type (e.g. the Topocentric CRSs) should override this function accordingly.
sourceSatcomCoordinates.isaffinetransform Method
isaffinetransform(t::Transform)
isaffinetransform(t::Type{<:Transform})
Function that returns true if the transform is an affine transform, meaning that it is a composition of a rotation and a translation.
Any transform which is considered affine should have a valid method for the following two functions:
SatcomCoordinates.raw_rotation(t::Transform)
: Returns the rotation matrix associated tot
as aRotMatrix
SatcomCoordinates.raw_translation(t::Transform)
: Returns the translation vector associated tot
as aSVector
SatcomCoordinates.israwtransform Method
israwtransform(t::Transform)
israwtransform(t::Type{<:Transform})
Function that returns true if the transform is considered raw, meaning that it directy operates and returns raw vectors or NTuples without units and is not specifically tied to input and output CRSs.
Raw transforms are supposed to be internally used as fields to build concrete subtypes of AbstractCRSTransform
.
Expected Input/Output
Raw transforms are always assumed to get a NTuple
of floating numbers as input and produce a NTuple
of floating numbers as output.
See CRSTransform
for an example of the only concrete subtype implemented in this package.
SatcomCoordinates.linkedcrs Function
linkedcrs
CRS Trait function identifying that a CRS is linked to another CRS.
For this specific trait, the getcrs
and getcrstype
functions will return the instance or type of the linked CRS.
SatcomCoordinates.lla_origin Method
lla_origin(crs::AbstractCRS)
lla_origin(obj::FieldOrCoordinate)
Return the coordinate (in LLA CRS) of the origin of the provided obj
or crs
.
Custom topocentric CRSs must implement a specific method for this function.
For all other CRSs, this function will simply traverse the chaing of linked CRSs starting from crs
until it either finds a valid topocentric CRS or it reaches the end of the chain (In which case it will simply throw an error).
See also ecef_origin
, AbstractTopocentricCRS
SatcomCoordinates.llacrs Method
llacrs
CRS Trait function to represent CRSs which store coordinates in Latitude/Longitude/Altitude (LLA).
sourceSatcomCoordinates.preprocess_input_coords Method
preprocess_input_coords(CRS::Type{<:AbstractCRS}, T::Type{<:AbstractFloat}, coords::Point{N, Any}) where N
This function take input coordinates and process them by eventually removing units they come with (ensuring consistency with the units expected from a CRS) and converting them to the specific machine precision specified by the type parameter T
.
SatcomCoordinates.process_unitless_coords Method
process_unitless_coords(::Type{C}, crs::AbstractCRS, coords::NTuple{<:Any, T}) where {C <: FieldOrCoordinate, T}
This function is the last step in the pipeline for creating a coordinate as part of the create_coordinate
function.
It takes as input the NTuple
already stripped of eventual units (and normalized to the reference/preferred unit before being stripped) and should just do final checks on the inputs and return an eventualy modified NTuple (e.g. wrap angles in radians if beyond [-π, π]), which will be then used to create the coordinate instance via the constructor_without_checks
function.
The default method for this function simply returns the input NTuple
as is. Custom CRSs that require specific checks or modification of the unitless inputs should add a method to this function.
SatcomCoordinates.raw_isapprox Method
raw_isapprox(C, crs1, crs2, coords1, coords2; kwargs...)
Low-level function that need to be extended for customizing the behavior of Base.isapprox
for two instances c1
and c2
of the same base subtype of AbstractSatcomCoordinate
Arguments
C::Type{<:AbstractSatcomCoordinate}
: The base type of the coordinate to be compared. E.g.Coordinate
crs1::AbstractCRS
: The CRS of the first coordinatecrs2::AbstractCRS
: The CRS of the second coordinatecoords1::NTuple{N, <:AbstractFloat}
: The raw coordinates of the first coordinate (i.e. the output oftuplecoords(c1)
)coords2::NTuple{N, <:AbstractFloat}
: The raw coordinates of the second coordinate (i.e. the output oftuplecoords(c2)
)kwargs...
: Additional keyword arguments to be passed to theisapprox
function
This function has a default implementation for all comparisons between coordinates in the same CRS (i.e. is_same_crs(crs1, crs2) == true
) which simply transforms the raw coordiantes in SVectors and calls isapprox
on them.
SatcomCoordinates.raw_linkedcrs_transform Function
raw_linkedcrs_transform(crs::AbstractCRS)
Returns the raw transform that goes from the provided crs
to its linked one.
This function should return a raw transformation (i.e. a transformation operating directly the output of tuplecoords(coordinate)
rather than on the coordinate itself).
A raw transformation shall expects a NTuple{N, <:AbstractFloat} as input (where N
is the number of dimensions of the CRS) and return a NTuple{N, <:AbstractFloat} as output.
SatcomCoordinates.remove_unit Method
remove_unit(userunit::Unitful.Units, refunit::Unitful.Units, val::Number)
Take a number (with or without unit), convert or interpret it (depending on whether it has or not a unit) to the userunit
, convert again to refunit
and then strip the unit.
This is used internally in the constructor of Coordinate
s to convert user inputs into the raw coordinate for the specific CRS.
SatcomCoordinates.rootcrs Method
rootcrs(C::Type{<:AbstractCRS})
Returns the Root CRS type associated to the provided CRS type C
A root CRS is a CRS which is not derived or linked to any other CRS.
All root CRSs should also be cartesian CRSs
See also: rootcrs
, cartesiancrs
SatcomCoordinates.topocrs Method
topocrs
CRS Trait function to represent CRSs which are topocentric.
sourceSatcomCoordinates.transform_tuplecoords Method
transform_tuplecoords(crsₒ::AbstractCRS, crsᵢ::AbstractCRS, tup::NTuple{N, <:AbstractFloat}; kwargs...)
Low-level function that is used to convert to an output CRS (CRSₒ
) the raw coordinates tup
(coming as output of tuplecords
) of a coordinate in an input CRS (CRSᵢ
).
This function must return a ntuple with M
coordinates (where M == ncoords(CRSₒ)
).
Custom CRSs should implement a specific method of this function to enable conversion with other CRSs via the change_crs
user facing function.
SatcomCoordinates.transform_tuplecoords Method
transform_tuplecoords(ecef_crs::ECEF, eci_crs::ECI, tup::NTuple{3,T}; R_eci_to_ecef = NotProvided(), kwargs...) where T<:AbstractFloat
transform_tuplecoords(eci_crs::ECI, ecef_crs::ECEF, tup::NTuple{3,T}; R_eci_to_ecef = NotProvided(), kwargs...) where T<:AbstractFloat
Transforms a tuple of coordinates between an ECEF and ECI CRS. If the rotation matrix to go from the specific ECI CRS to the ECEF CRS is explicitly provided via the R_eci_to_ecef
keyword argument, it will be used to directly transform the coordinates. Otherwise, the eci_to_ecef_rotation
function will be used to compute it using the frameid
of both ECI and ECEF CRSs.
SatcomCoordinates.upreferred Method
upreferred(CRS::Type, unit::Unitful.Units)
Allow to customize the preferred unit on different types of CRSs. Defaults to Unitful.upreferred(unit)
.
Note
This is not the same function as Unitful.upreferred
but just shares the same name as they basically have the same end goal. It is nonetheless redefined internally to SatComCoordinates
to avoid polluting the methods of Unitful.upreferred
.
SatcomCoordinates.@define_properties Macro
@define_properties CRS [
propname1 => unit1
propname2 => unit2[ => (aliases2...)]
...
]
This macro simplifies the definition of properties for a custom CRS
type and association between each property and its user-facing unit as well as optional property aliases (supported by Base.getproperty
on coordinates, i.e. instances of AbstractSatcomCoordinate
).
The arguments expected by the macro are the following:
1st argument: The desired custom CRS type to extend
2nd argument: A vector of pairs (one per property of coordinates of the custom CRS)
The pairs are of the form propname => unit
(or propname => unit => (aliases...)
in case aliases are desired for specific properties) where:
propname
: is the name of the property/coordinate of the custom CRSunit
: is the unit of the property/coordinate. This MUST be aUnit
(and not aQuantity
) fromUnitful.jl
. In case of properties without a unit (e.g. the coords for UV pointing),NoUnits
must be used.aliases
: [OPTIONAL] A tuple of names the corresponding property can be accessed to viaBase.getproperty
on coordinates defined over the custom CRS.
Example
julia> using SatcomCoordinates
julia> struct CustomCartesian <: AbstractCRS end;
julia> SatcomCoordinates.@define_properties CustomCartesian [
x => u"km" # For some reason, we want x to be shown and parsed in km by default
y => u"m" => (y2,) # And for some other reason, we want to be able to access y also with the name y2
z => u"m"
]
julia> p = Position(CustomCartesian(), 1,2,3);
julia> p.x
1.0 km
julia> (; y2) = p;
julia> y2
2.0 m
When to use this macro
By default, CRSs (i.e. subtypes of AbstractCRS
) are cartesian CRSs and have the following default properties:
x
y
z
which are all associated to the u"m"
unit and no custom aliases. For all other CRSs (or if one wants to customize either the units or add aliases for custom CRSs) this macro must be used on the custom CRS type properly use the other types and function of SatcomCoordinates.jl
.
See the extended help section below for more details and advanced usage.
See also: AbstractSatcomCoordinate
, AbstractCRS
, Position
Extended Help
Generated Code
This macro automatically adds custom methods for the provided CRS (and following the provided properties and associated units) to the following functions which are internal to SatcomCoordinates.jl
:
SatcomCoordinates.units
(CRS::Type{<:AbstractCRS})SatcomCoordinates.resolve_property
(CRS::Type{<:AbstractCRS}, propname::Symbol)
SatcomCoordinates.units
The SatcomCoordinates.units
function operates on a CRS type and must returns a NamedTuple
with the properties as keys and the associated units as values.
In the default case for custom CRSs, the function method is the following:
SatcomCoordinates.units(CRS::Type{<:AbstractCRS}) = (; x = u"m", y = u"m", z = u"m")
And all additional methods for custom CRSs are expected to be of the same form (i.e. returning a NamedTuple
with the properties (the baseline ones, not the aliases) as keys and the associated units as values).
SatcomCoordinates.resolve_property
The SatcomCoordinates.resolve_property
function operates on a CRS type and a property name and must return the name of the property to be used for the getproperty
call on any coordinate defined over the extended CRS.
Let's consider as an example the the standard (ISO) spherical CRS, which can be mocked up as follows (this package actually defines and export a more generic SphericalCRS
which includes the ISO one, the one below is just a simplified mockup):
using SatcomCoordinates
struct ISOSpherical <: AbstractCRS end
SatcomCoordinates.@define_properties CustomCRS [
θ => u"°" => (theta, t)
φ => u"°" => (phi, p, ϕ)
r => u"m" => (distance, range)
]
In the above code, we want to be able to access e.g. the polar angle θ
also with the alternative aliases theta
or t
.
The SatcomCoordinates.resolve_property
function generated by the macro in the snippet above is the following:
Base.@constprop :aggressive function SatcomCoordinates.resolve_property(CRS::Type{<:ISOSpherical}, propname::Symbol)
if propname in (:θ, :theta, :t)
:θ
elseif propname in (:φ, :phi, :p, :ϕ)
:φ
elseif propname in (:r, :distance, :range)
:r
else
:__could_not_resolve_property__
end
end
The Base.@constprop :aggressive
macro is used to encourage constant propagation of symbols when accessed via getproperty
resulting in no cost access to the same property via multiple aliases.
Advanced Usage
The macro actually supports also a more advance signature for cases of CRSs that wrap a parent CRS and simply want to reuse the same properties of the parent CRS.
This is for example the case for the SphericalCRS
type, which is a wrapper around another CRS specifying the pointing type (e.g. ThetaPhi
or AzEl
). The SphericalCRS
needs to inherit the first 2 properties from it's wrapped pointing CRS and just add the third one (r
).
It is not possible to statically define all the CRS properties in this case, and so the following synthax is used:
@define_properties SphericalCRS [
pointingcrs(_)...
r => u"m" => (distance, range)
]
In this alternative synthax, the macro looks for any occurrence of the ...
at the end of an expression
This expression (without the ...
, so pointingcrs(_)
above) need to return the wrapped CRS (i.e. the pointing CRS in this example) type.
For convenience, the expression can contain the _
placeholder to represent the CRS type being extended.
In the specific case above, pointingcrs
is a function of SatComCoordinates
which returns the underlying pointing CRS type when called with a SphericalCRS
type as input.
This special synthax currently only supports a single ...
within a @define_properties
call, and creates the following generated code (in the case of the example above):
function SatcomCoordinates.units(CRS::Type{<:SphericalCRS})
(; SatcomCoordinates.units(SatcomCoordinates.pointingcrs(CRS))..., r = u"m")
end
Base.@constprop :aggressive function SatcomCoordinates.resolve_property(CRS::Type{<:SphericalCRS}, propname::Symbol)
if propname in (:r, :distance, :range)
:r
else
SatcomCoordinates.resolve_property(SatcomCoordinates.pointingcrs(CRS), propname)
end
end