Skip to content

Commit

Permalink
Merge pull request #52 from aplavin/unitful
Browse files Browse the repository at this point in the history
support Unitful quantities and units
  • Loading branch information
mileslucas authored Aug 1, 2023
2 parents b2e06a4 + 6b8de18 commit 7b7bac2
Show file tree
Hide file tree
Showing 4 changed files with 47 additions and 1 deletion.
6 changes: 5 additions & 1 deletion Project.toml
Original file line number Diff line number Diff line change
@@ -1,7 +1,7 @@
name = "SkyCoords"
uuid = "fc659fc5-75a3-5475-a2ea-3da92c065361"
authors = "Kyle Barbary, Mosé Giordano, and contributors"
version = "1.3.0"
version = "1.4.0"

[deps]
AstroAngles = "5c4adb95-c1fc-4c53-b4ea-2a94080c53d2"
Expand All @@ -12,17 +12,21 @@ StaticArrays = "90137ffa-7385-5640-81b9-e52037218182"

[weakdeps]
Accessors = "7d9f7c33-5ae7-4f3b-8dc6-eff91059b697"
Unitful = "1986cc42-f94f-5a68-af5c-568840ba703d"

[extensions]
AccessorsExt = "Accessors"
UnitfulExt = "Unitful"

[compat]
Accessors = "0.1"
AstroAngles = "0.1"
ConstructionBase = "1"
Rotations = "1"
StaticArrays = "0.8, 0.9, 1"
Unitful = "1.15"
julia = "1.6"

[extras]
Accessors = "7d9f7c33-5ae7-4f3b-8dc6-eff91059b697"
Unitful = "1986cc42-f94f-5a68-af5c-568840ba703d"
18 changes: 18 additions & 0 deletions ext/UnitfulExt.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,18 @@
module UnitfulExt

using Unitful
using SkyCoords

_COORDTYPES_LATLON = Union{ICRSCoords, GalCoords, FK5Coords, EclipticCoords}

(::Type{T})(lon::Quantity, lat) where {T<:_COORDTYPES_LATLON} = T(ustrip(u"rad", lon), lat)
(::Type{T})(lon, lat::Quantity) where {T<:_COORDTYPES_LATLON} = T(lon, ustrip(u"rad", lat))
(::Type{T})(lon::Quantity, lat::Quantity) where {T<:_COORDTYPES_LATLON} = T(ustrip(u"rad", lon), ustrip(u"rad", lat))

SkyCoords.lon(u::Unitful.Units, c) = SkyCoords.lon(c) * u"rad" |> u
SkyCoords.lat(u::Unitful.Units, c) = SkyCoords.lat(c) * u"rad" |> u

SkyCoords.separation(u::Unitful.Units, c1, c2) = SkyCoords.separation(c1, c2) * u"rad" |> u
SkyCoords.position_angle(u::Unitful.Units, c1, c2) = SkyCoords.position_angle(c1, c2) * u"rad" |> u

end
1 change: 1 addition & 0 deletions test/Project.toml
Original file line number Diff line number Diff line change
Expand Up @@ -9,3 +9,4 @@ PythonCall = "6099a3de-0909-46bc-b1f4-468b9a2dfc0d"
StableRNGs = "860ef19b-820b-49d6-a774-d7a799459cd3"
Statistics = "10745b16-79ce-11e8-11f9-7d13ad32a3b2"
Test = "8dfed614-e22c-5e08-85e1-65c5234f0b40"
Unitful = "1986cc42-f94f-5a68-af5c-568840ba703d"
23 changes: 23 additions & 0 deletions test/runtests.jl
Original file line number Diff line number Diff line change
@@ -1,5 +1,6 @@
using AstroAngles
using Accessors
using Unitful
using ConstructionBase: setproperties
using DelimitedFiles
using LinearAlgebra: normalize
Expand Down Expand Up @@ -189,6 +190,28 @@ VERSION > v"1.9-DEV" && @testset "Accessors" begin
end
end

VERSION > v"1.9-DEV" && @testset "Unitful" begin
@test ICRSCoords(1u"rad", 0.5) === ICRSCoords(1, 0.5)
@test GalCoords(1u"rad", 0.5u"rad") === GalCoords(1, 0.5)
@test FK5Coords{2000}(1u"°", 0.5) === FK5Coords{2000}(deg2rad(1), 0.5)
@test EclipticCoords{2000}(1u"°", 0.5u"°") === EclipticCoords{2000}(deg2rad(1), deg2rad(0.5))

@test SkyCoords.lat(u"rad", ICRSCoords(1, 0.5)) === 0.5u"rad"
@test SkyCoords.lon(u"°", ICRSCoords(1, 0.5)) === rad2deg(1)u"°"

@test separation(u"rad", ICRSCoords(1, 0.5), ICRSCoords(1, -0.2)) === 0.7u"rad"
@test separation(u"°", ICRSCoords(1, 0.5), ICRSCoords(1, -0.2)) === rad2deg(0.7)u"°"

@test position_angle(u"rad", ICRSCoords(1, 0.5), ICRSCoords(1, -0.2)) === Float64(π)*u"rad"
@test position_angle(u"°", ICRSCoords(1, 0.5), ICRSCoords(1, -0.2)) === 180.0u"°"

# offset() works without special Unitful support
@test offset(ICRSCoords(1, 0.5), 0.1u"rad", 2) === offset(ICRSCoords(1, 0.5), 0.1, 2)
@test offset(ICRSCoords(1, 0.5), 0.1u"rad", 2u"rad") === offset(ICRSCoords(1, 0.5), 0.1, 2)
@test offset(ICRSCoords(1, 0.5), 0.1, 100u"°") === offset(ICRSCoords(1, 0.5), 0.1, deg2rad(100))
@test offset(ICRSCoords(1, 0.5), 0.1u"°", 100u"°") === offset(ICRSCoords(1, 0.5), deg2rad(0.1), deg2rad(100))
end

@testset "equality" begin
@testset for T in [ICRSCoords, GalCoords, FK5Coords{2000}, EclipticCoords{2000}]
c1 = T(1., 2.)
Expand Down

0 comments on commit 7b7bac2

Please sign in to comment.