Safe Haskell | Safe-Inferred |
---|---|
Language | Haskell2010 |
Types and functions for dealing with Kepler orbits.
Synopsis
- data Orbit a = Orbit {
- eccentricity :: !(Unitless a)
- periapsis :: !(Distance a)
- inclinationSpecifier :: !(InclinationSpecifier a)
- periapsisSpecifier :: !(PeriapsisSpecifier a)
- primaryGravitationalParameter :: !(Quantity ((:*) ((:^) Meter (Succ (Succ (Succ 'Zero)))) ((:^) Second (Pred (Pred 'Zero)))) a)
- data InclinationSpecifier a
- = Inclined {
- longitudeOfAscendingNode :: !(Angle a)
- inclination :: !(Angle a)
- | NonInclined
- = Inclined {
- data PeriapsisSpecifier a
- = Eccentric {
- argumentOfPeriapsis :: !(Angle a)
- | Circular
- = Eccentric {
- data Classification
- isValid :: (Ord a, Num a) => Orbit a -> Bool
- classify :: (Num a, Ord a) => Orbit a -> Classification
- normalizeOrbit :: (Floating a, Real a) => Orbit a -> Orbit a
- apoapsis :: (Fractional a, Ord a) => Orbit a -> Maybe (Distance a)
- meanMotion :: (Floating a, Ord a) => Orbit a -> Quantity ((:/) Radian Second) a
- period :: (Floating a, Ord a) => Orbit a -> Maybe (Time a)
- arealVelocity :: (Ord a, Floating a) => Orbit a -> Quantity ((:/) ((:^) Meter (Succ (Succ 'Zero))) Second) a
- semiMajorAxis :: (Fractional a, Ord a) => Orbit a -> Maybe (Distance a)
- semiMinorAxis :: (Floating a, Ord a) => Orbit a -> Distance a
- semiLatusRectum :: Num a => Orbit a -> Distance a
- hyperbolicApproachAngle :: (Floating a, Ord a) => Orbit a -> Maybe (Angle a)
- hyperbolicDepartureAngle :: (Floating a, Ord a) => Orbit a -> Maybe (Angle a)
- timeAtMeanAnomaly :: (Floating a, Ord a) => Orbit a -> Angle a -> Time a
- timeAtEccentricAnomaly :: (Floating a, Ord a) => Orbit a -> Angle a -> Maybe (Time a)
- timeAtHyperbolicAnomaly :: (Floating a, Ord a) => Orbit a -> AngleH a -> Maybe (Time a)
- timeAtTrueAnomaly :: (Real a, Floating a) => Orbit a -> Angle a -> Maybe (Time a)
- meanAnomalyAtTime :: (Floating a, Ord a) => Orbit a -> Time a -> Angle a
- meanAnomalyAtEccentricAnomaly :: (Floating a, Ord a) => Orbit a -> Angle a -> Maybe (Angle a)
- meanAnomalyAtHyperbolicAnomaly :: (Floating a, Ord a) => Orbit a -> AngleH a -> Maybe (Angle a)
- meanAnomalyAtTrueAnomaly :: (Real a, Floating a) => Orbit a -> Angle a -> Maybe (Angle a)
- eccentricAnomalyAtTime :: (Converge [a], Floating a, Real a) => Orbit a -> Time a -> Maybe (Angle a)
- eccentricAnomalyAtMeanAnomaly :: forall a. (Converge [a], Floating a, Real a) => Orbit a -> Angle a -> Maybe (Angle a)
- eccentricAnomalyAtMeanAnomalyFloat :: Orbit Float -> Angle Float -> Maybe (Angle Float)
- eccentricAnomalyAtTrueAnomaly :: (Floating a, Real a) => Orbit a -> Angle a -> Maybe (Angle a)
- hyperbolicAnomalyAtTime :: forall a. (Converge [a], RealFloat a) => Orbit a -> Time a -> Maybe (AngleH a)
- hyperbolicAnomalyAtMeanAnomaly :: forall a. (Converge [a], RealFloat a) => Orbit a -> Angle a -> Maybe (AngleH a)
- hyperbolicAnomalyAtMeanAnomalyDouble :: Orbit Double -> Angle Double -> Maybe (AngleH Double)
- hyperbolicAnomalyAtTrueAnomaly :: (Floating a, Ord a) => Orbit a -> Angle a -> Maybe (AngleH a)
- trueAnomalyAtTime :: forall a. (Converge [a], RealFloat a) => Orbit a -> Time a -> Angle a
- trueAnomalyAtMeanAnomaly :: (Converge [a], RealFloat a) => Orbit a -> Angle a -> Angle a
- trueAnomalyAtEccentricAnomaly :: RealFloat a => Orbit a -> Angle a -> Maybe (Angle a)
- trueAnomalyAtHyperbolicAnomaly :: (Ord a, Floating a) => Orbit a -> AngleH a -> Maybe (Angle a)
- specificAngularMomentum :: Floating a => Orbit a -> Quantity ((:*) ((:^) Meter (Succ (Succ 'Zero))) ((:^) Second (Pred 'Zero))) a
- specificOrbitalEnergy :: (Ord a, Floating a) => Orbit a -> Quantity ((:/) Joule ((:@) Kilo Gram)) a
- specificPotentialEnergyAtTrueAnomaly :: (Ord a, Floating a) => Orbit a -> Angle a -> Quantity ((:/) Joule ((:@) Kilo Gram)) a
- specificKineticEnergyAtTrueAnomaly :: (Ord a, Floating a) => Orbit a -> Angle a -> Quantity ((:/) Joule ((:@) Kilo Gram)) a
- speedAtTrueAnomaly :: (Ord a, Floating a) => Orbit a -> Angle a -> Speed a
- radiusAtTrueAnomaly :: (Ord a, Floating a) => Orbit a -> Angle a -> Distance a
- escapeVelocityAtDistance :: Floating a => Quantity ((:*) ((:^) Meter (Succ (Succ (Succ 'Zero)))) ((:^) Second (Pred (Pred 'Zero)))) a -> Distance a -> Speed a
- type Quantity u = MkQu_ULN u 'DefaultLCSU
- type Time = Quantity Second
- type Distance = Quantity Meter
- type Speed = Quantity ((:*) Meter ((:^) Second (Pred 'Zero)))
- type Mass = Quantity ((:@) Kilo Gram)
- type Angle = Quantity Radian
- type AngleH = Quantity RadianHyperbolic
- data RadianHyperbolic = RadianHyperbolic
- data PlaneAngleHyperbolic = PlaneAngleHyperbolic
- type Unitless = Quantity Number
- class Converge a
The Orbit data type and dependencies
Data type defining an orbit parameterized by the type used to represent values
Orbit | |
|
data InclinationSpecifier a Source #
Along with PeriapsisSpecifier
the InclinationSpecifier
describes
orbital elements extra to its geometry.
Inclined | The orbit does not lie exactly in the reference plane |
| |
NonInclined | The orbit lies in the reference plane |
Instances
Show a => Show (InclinationSpecifier a) Source # | |
Defined in Physics.Orbit showsPrec :: Int -> InclinationSpecifier a -> ShowS # show :: InclinationSpecifier a -> String # showList :: [InclinationSpecifier a] -> ShowS # | |
Eq a => Eq (InclinationSpecifier a) Source # | |
Defined in Physics.Orbit (==) :: InclinationSpecifier a -> InclinationSpecifier a -> Bool # (/=) :: InclinationSpecifier a -> InclinationSpecifier a -> Bool # |
data PeriapsisSpecifier a Source #
Along with InclinationSpecifier
the PeriapsisSpecifier
describes
orbital elements extra to its geometry.
Eccentric | The orbit is not circular |
| |
Circular | The orbit has an eccentricity of 0 so the
|
Instances
Show a => Show (PeriapsisSpecifier a) Source # | |
Defined in Physics.Orbit showsPrec :: Int -> PeriapsisSpecifier a -> ShowS # show :: PeriapsisSpecifier a -> String # showList :: [PeriapsisSpecifier a] -> ShowS # | |
Eq a => Eq (PeriapsisSpecifier a) Source # | |
Defined in Physics.Orbit (==) :: PeriapsisSpecifier a -> PeriapsisSpecifier a -> Bool # (/=) :: PeriapsisSpecifier a -> PeriapsisSpecifier a -> Bool # |
data Classification Source #
What form the orbit's geometry takes. This is dependant only on the
eccentricity
, e >= 0, of the orbit.
Elliptic | 0 <= e < 1 This includes circular orbits. |
Parabolic | e == 1 |
Hyperbolic | e > 1 |
Instances
Read Classification Source # | |
Defined in Physics.Orbit readsPrec :: Int -> ReadS Classification # readList :: ReadS [Classification] # | |
Show Classification Source # | |
Defined in Physics.Orbit showsPrec :: Int -> Classification -> ShowS # show :: Classification -> String # showList :: [Classification] -> ShowS # | |
Eq Classification Source # | |
Defined in Physics.Orbit (==) :: Classification -> Classification -> Bool # (/=) :: Classification -> Classification -> Bool # |
Functions for dealing with orbits
Utilities
isValid :: (Ord a, Num a) => Orbit a -> Bool Source #
Determines if the orbital elements are valid (e >= 0
etc...). The
behavior of all the other functions in this module is undefined when given
an invalid orbit.
normalizeOrbit :: (Floating a, Real a) => Orbit a -> Orbit a Source #
Return an equivalent orbit such that
- i ∈ [0..π)
- Ω ∈ [0..2π)
- ω ∈ [0..2π)
- inclinationSpecifier == NonInclined if i = 0
- periapsisSpecifier == Circular if e == 0 and ω == 0
Orbital elements
meanMotion :: (Floating a, Ord a) => Orbit a -> Quantity ((:/) Radian Second) a Source #
Calculate the mean motion, n, of an orbit
This is the rate of change of the mean anomaly with respect to time.
period :: (Floating a, Ord a) => Orbit a -> Maybe (Time a) Source #
Calculate the orbital period, p, of an elliptic orbit.
period
returns Nothing if given a parabolic or hyperbolic orbit.
arealVelocity :: (Ord a, Floating a) => Orbit a -> Quantity ((:/) ((:^) Meter (Succ (Succ 'Zero))) Second) a Source #
Calculate the areal velocity, A, of the orbit.
The areal velocity is the area swept out by the line between the orbiting body and the primary per second.
Geometry
semiMajorAxis :: (Fractional a, Ord a) => Orbit a -> Maybe (Distance a) Source #
semiMinorAxis :: (Floating a, Ord a) => Orbit a -> Distance a Source #
Calculate the semi-minor axis, b, of the Orbit
. Like semiMajorAxis
'semiMinorAxis' o
is negative when o
is a hyperbolic orbit. In the
case of a parabolic orbit semiMinorAxis
returns 0m
.
semiLatusRectum :: Num a => Orbit a -> Distance a Source #
Calculate the semiLatusRectum, l, of the Orbit
hyperbolicApproachAngle :: (Floating a, Ord a) => Orbit a -> Maybe (Angle a) Source #
Calculate the angle at which a body leaves the system when on a hyperbolic trajectory relative to the argument of periapsis. This is the limit of the true anomaly as time tends towards -infinity minus the argument of periapsis. The approach angle is in the closed range (-π..π/2).
This is the negation of the departure angle.
hyperbolicApproachAngle
returns Nothing when given a non-hyperbolic orbit
and -π when given a parabolic orbit.
hyperbolicDepartureAngle :: (Floating a, Ord a) => Orbit a -> Maybe (Angle a) Source #
Calculate the angle at which a body leaves the system when on an escape trajectory relative to the argument of periapsis. This is the limit of the true anomaly as time tends towards infinity minus the argument of periapsis. The departure angle is in the closed range (π/2..π).
This is the negation of the approach angle.
hyperbolicDepartureAngle
returns Nothing when given an elliptic orbit and
π when given a parabolic orbit.
Conversions
To time since periapse
timeAtMeanAnomaly :: (Floating a, Ord a) => Orbit a -> Angle a -> Time a Source #
Calculate the time since periapse, t, when the body has the given mean anomaly, M. M may be negative, indicating that the orbiting body has yet to reach periapse.
The sign of the time at mean anomaly M is the same as the sign of M.
The returned time is unbounded.
timeAtEccentricAnomaly :: (Floating a, Ord a) => Orbit a -> Angle a -> Maybe (Time a) Source #
Calculate the time since periapse, t, of an elliptic orbit when at eccentric anomaly E.
timeAtEccentricAnomaly
returns Nothing if given a parabolic or hyperbolic
orbit.
timeAtHyperbolicAnomaly :: (Floating a, Ord a) => Orbit a -> AngleH a -> Maybe (Time a) Source #
Calculate the time since periapse, t, of a hyperbolic orbit when at hyperbolic anomaly H.
Returns Nothing if given an elliptic or parabolic orbit.
timeAtTrueAnomaly :: (Real a, Floating a) => Orbit a -> Angle a -> Maybe (Time a) Source #
Calculate the time since periapse given the true anomaly, ν, of an orbiting body.
Returns Nothing
if the body never passed through the specified true
anomaly.
To mean anomaly
meanAnomalyAtTime :: (Floating a, Ord a) => Orbit a -> Time a -> Angle a Source #
Calculate the mean anomaly, M, at the given time since periapse, t. t may be negative, indicating that the orbiting body has yet to reach periapse.
The sign of the mean anomaly at time t is the same as the sign of t.
The returned mean anomaly is unbounded.
meanAnomalyAtEccentricAnomaly :: (Floating a, Ord a) => Orbit a -> Angle a -> Maybe (Angle a) Source #
Calculate the mean anomaly, M, of an elliptic orbit when at eccentric anomaly E
meanAnomalyAtEccentricAnomaly
returns Nothing if given a parabolic or
hyperbolic orbit.
The number of orbits represented by the anomalies is preserved;
i.e. M div
2π = E div
2π
meanAnomalyAtHyperbolicAnomaly :: (Floating a, Ord a) => Orbit a -> AngleH a -> Maybe (Angle a) Source #
Calculate the mean anomaly, M, of a hyperbolic orbit when at hyperbolic anomaly H
To eccentric anomaly
eccentricAnomalyAtTime :: (Converge [a], Floating a, Real a) => Orbit a -> Time a -> Maybe (Angle a) Source #
Calculate the eccentric anomaly, E, of an elliptic orbit at time t.
eccentricAnomalyAtTime
returns Nothing when given a parabolic or
hyperbolic orbit.
The number of orbits represented by the time is preserved;
i.e. t div
p = E div
2π
eccentricAnomalyAtMeanAnomaly :: forall a. (Converge [a], Floating a, Real a) => Orbit a -> Angle a -> Maybe (Angle a) Source #
Calculate the eccentric anomaly, E, of an elliptic orbit when at mean anomaly M. This function is considerably slower than most other conversion functions as it uses an iterative method as no closed form solution exists.
The number of orbits represented by the anomalies is preserved;
i.e. M div
2π = E div
2π
eccentricAnomalyAtMeanAnomaly
returns Nothing when given a parabolic or
hyperbolic orbit.
eccentricAnomalyAtMeanAnomalyFloat :: Orbit Float -> Angle Float -> Maybe (Angle Float) Source #
eccentricAnomalyAtMeanAnomaly
specialized to Float
.
This function is used to calculate the initial guess for
eccentricAnomalyAtMeanAnomaly
.
eccentricAnomalyAtTrueAnomaly :: (Floating a, Real a) => Orbit a -> Angle a -> Maybe (Angle a) Source #
To hyperbolic anomaly
hyperbolicAnomalyAtTime :: forall a. (Converge [a], RealFloat a) => Orbit a -> Time a -> Maybe (AngleH a) Source #
hyperbolicAnomalyAtMeanAnomaly :: forall a. (Converge [a], RealFloat a) => Orbit a -> Angle a -> Maybe (AngleH a) Source #
hyperbolicAnomalyAtMeanAnomalyDouble :: Orbit Double -> Angle Double -> Maybe (AngleH Double) Source #
Calculate the hyperbolic anomaly, H, at a given mean anomaly. Unline
eccentricAnomalyAtMeanAnomalyFloat
this uses double precision floats to
help avoid overflowing.
hyperbolicAnomalyAtTrueAnomaly :: (Floating a, Ord a) => Orbit a -> Angle a -> Maybe (AngleH a) Source #
To true anomaly
trueAnomalyAtTime :: forall a. (Converge [a], RealFloat a) => Orbit a -> Time a -> Angle a Source #
Calculate the true anomaly, ν, of a body at time since periapse, t.
trueAnomalyAtMeanAnomaly :: (Converge [a], RealFloat a) => Orbit a -> Angle a -> Angle a Source #
Calculate the true anomaly, ν, of an orbiting body when it has the given mean anomaly, _M.
trueAnomalyAtHyperbolicAnomaly :: (Ord a, Floating a) => Orbit a -> AngleH a -> Maybe (Angle a) Source #
Properties of orbits
specificAngularMomentum :: Floating a => Orbit a -> Quantity ((:*) ((:^) Meter (Succ (Succ 'Zero))) ((:^) Second (Pred 'Zero))) a Source #
Specific angular momentum, h, is the angular momentum per unit mass
specificOrbitalEnergy :: (Ord a, Floating a) => Orbit a -> Quantity ((:/) Joule ((:@) Kilo Gram)) a Source #
Specific orbital energy, ε, is the orbital energy per unit mass
specificPotentialEnergyAtTrueAnomaly :: (Ord a, Floating a) => Orbit a -> Angle a -> Quantity ((:/) Joule ((:@) Kilo Gram)) a Source #
Specific potential energy, εp, is the potential energy per unit mass at a particular true anomaly
specificKineticEnergyAtTrueAnomaly :: (Ord a, Floating a) => Orbit a -> Angle a -> Quantity ((:/) Joule ((:@) Kilo Gram)) a Source #
Specific kinetic energy, εk, is the kinetic energy per unit mass at a particular true anomaly
speedAtTrueAnomaly :: (Ord a, Floating a) => Orbit a -> Angle a -> Speed a Source #
What is the speed, v, of a body at a particular true anomaly
radiusAtTrueAnomaly :: (Ord a, Floating a) => Orbit a -> Angle a -> Distance a Source #
The distance, r, from the primary body to the orbiting body at a particular true anomaly.
Other utilities
escapeVelocityAtDistance :: Floating a => Quantity ((:*) ((:^) Meter (Succ (Succ (Succ 'Zero)))) ((:^) Second (Pred (Pred 'Zero)))) a -> Distance a -> Speed a Source #
The escape velocity for a primary with specified gravitational parameter at a particular distance.
Unit synonyms
type Quantity u = MkQu_ULN u 'DefaultLCSU Source #
type Speed = Quantity ((:*) Meter ((:^) Second (Pred 'Zero))) Source #
A measure in meters per second.
type AngleH = Quantity RadianHyperbolic Source #
A measure in radians (hyperbolic)
data RadianHyperbolic Source #
Instances
Show RadianHyperbolic Source # | |
Defined in Physics.Orbit.Metrology showsPrec :: Int -> RadianHyperbolic -> ShowS # show :: RadianHyperbolic -> String # showList :: [RadianHyperbolic] -> ShowS # | |
Unit RadianHyperbolic Source # | |
Defined in Physics.Orbit.Metrology type BaseUnit RadianHyperbolic # type DimOfUnit RadianHyperbolic # type UnitFactorsOf RadianHyperbolic :: [Factor Type] # | |
type BaseUnit RadianHyperbolic Source # | |
Defined in Physics.Orbit.Metrology | |
type DimOfUnit RadianHyperbolic Source # | |
Defined in Physics.Orbit.Metrology | |
type UnitFactorsOf RadianHyperbolic Source # | |
Defined in Physics.Orbit.Metrology |
data PlaneAngleHyperbolic Source #
Instances
Dimension PlaneAngleHyperbolic Source # | |
Defined in Physics.Orbit.Metrology type DimFactorsOf PlaneAngleHyperbolic :: [Factor Type] # | |
type DimFactorsOf PlaneAngleHyperbolic Source # | |
Defined in Physics.Orbit.Metrology | |
type DefaultUnitOfDim PlaneAngleHyperbolic Source # | |
Defined in Physics.Orbit.Metrology |
Reexported from CReal
If a type is an instance of Converge then it represents a stream of values which are increasingly accurate approximations of a desired value
Instances
Converge [CReal n] | The overlapping instance for It's important to note when the error function reaches zero this function
behaves like Find where log x = π using Newton's method
|
Eq a => Converge [a] | Every list of equatable values is an instance of |