API
Types
Initial Conditions
NbodyGradient.CartesianIC — TypeCartesianIC{T<:AbstractFloat} <: InitialConditions{T}Initial conditions, specified by the Cartesian coordinates and masses of each body.
Fields
x::Matrix{T}: Positions of each body [dimension, body].v::Matrix{T}: Velocities of each body [dimension, body].m::Vector{T}: masses of each body.nbody::Int64: Number of bodies in system.t0::T: Initial time.
NbodyGradient.Elements — TypeElements{T<:AbstractFloat} <: AbstractInitialConditionsOrbital elements of a binary, and mass of a 'outer' body. See Tutorials for units and conventions.
Fields
m::T: Mass of outer body.P::T: Period [Days].t0::T: Initial time of transit [Days].ecosω::T: Eccentricity vector x-component (eccentricity times cosine of the argument of periastron)esinω::T: Eccentricity vector y-component (eccentricity times sine of the argument of periastron)I::T: Inclination, as measured from sky-plane [Radians].Ω::T: Longitude of ascending node, as measured from +x-axis [Radians].a::T: Orbital semi-major axis [AU].e::T: Eccentricity.ω::T: Argument of periastron [Radians].tp::T: Time of periastron passage [Days].
NbodyGradient.Elements — MethodElements(m,P,t0,ecosω,esinω,I,Ω)Main Elements constructor. May use keyword arguments, see Tutorials.
NbodyGradient.ElementsIC — TypeElementsIC{T<:AbstractFloat} <: InitialConditions{T}Initial conditions, specified by a hierarchy vector and orbital elements.
Fields
elements::Matrix{T}: Masses and orbital elements.ϵ::Matrix{T}: Matrix of Jacobi coordinatesamat::Matrix{T}: 'A' matrix from Hamers and Portegies Zwart 2016.nbody::Int64: Number of bodies.m::Vector{T}: Masses of bodies.t0::T: Initial time [Days].
NbodyGradient.ElementsIC — MethodElementsIC(t0,H,elems)Collects Elements and produces an ElementsIC struct.
Arguments
t0::T: Initial time [Days].H: Hierarchy specification.elems: The orbital elements and masses of the system.
There are a number of way to specify the initial conditions. Below we've described the arguments ElementsIC takes. Any combination of H and elems may be used. For a concrete example see Tutorials.
Elements
elems...: A sequence ofElements{T}. Elements should be passed in the order they appear in the hierarchy (left to right).elems::Vector{Elements}: A vector ofElements. As above, Elements should be in order.elems::Matrix{T}: An matrix containing the masses and orbital elements.elems::String: Name of a file containing the masses and orbital elements.
Each method is simply populating the ElementsIC.elements field, which is a Matrix{T}.
Hierarchy
- Number of bodies: 
H::Int64: The system will be given by a 'fully-nested' Keplerian. 
H = 4 corresponds to:
3        ____|____
        |         |
2    ___|___      d
    |       |
1 __|__     c
 |     |
 a     b- Hierarchy Vector: 
H::Vector{Int64}: The first elements is the number of bodies. Each subsequent is the number of binaries on a level of the hierarchy. 
H = [4,2,1]. Two binaries on level 1, one on level 2.
2    ____|____
    |         |
1 __|__     __|__
 |     |   |     |
 a     b   c     d- Full Hierarchy Matrix: 
H::Matrix{<:Real}: Provide the hierarchy matrix, directly. 
H = [-1 1 0 0; 0 0 -1 1; -1 -1 1 1; -1 -1 -1 -1]. Produces the same system as H = [4,2,1].
State
NbodyGradient.State — TypeState{T<:AbstractFloat} <: AbstractStateCurrent state of simulation.
Fields (relevant to the user)
x::Matrix{T}: Positions of each body [dimension, body].v::Matrix{T}: Velocities of each body [dimension, body].t::Vector{T}: Current time of simulation.m::Vector{T}: Masses of each body.jac_step::Matrix{T}: Current Jacobian.dqdt::Vector{T}: Derivative with respect to time.
NbodyGradient.State — MethodState(ic)Constructor for State type.
Arguments
ic::InitialConditions{T}: Initial conditions for the system.
Integrator
NbodyGradient.Integrator — TypeIntegrator{T<:AbstractFloat}Integrator. Used as a functor to integrate a State.
Fields
scheme::Function: The integration scheme to use.h::T: Step size.t0::T: Initial time.tmax::T: Duration of simulation.
TransitTiming
NbodyGradient.TransitParameters — TypeTransitParameters{T<:AbstractFloat} <: TransitOutput{T}Transit times, impact parameters, sky-velocities, and derivatives.
(User-facing) Fields
ttbv::Matrix{T}: The transit times, impact parameter, and sky-velocity of each body.dtbvdq0::Array{T,5}: Derivatives of the transit times, impact parameters, and sky-velocities with respect to the initial Cartesian coordinates and masses.dtbvdelements::Array{T,5}: Derivatives of the transit times, impact parameters, and sky-velocities with respect to the initial orbital elements and masses.
NbodyGradient.TransitParameters — MethodTransitParameters(tmax, ic; ti)Constructor for TransitParameters type.
Arguments
tmax::T: Expected total elapsed integration time. (Allocates arrays accordingly)ic::ElementsIC{T}: Initial conditions for the system
Optional
ti::Int64=1: Index of the body with respect to which transits are measured. (Default is the central body)
NbodyGradient.TransitTiming — TypeTransitTiming{T<:AbstractFloat} <: TransitOutput{T}Transit times and derivatives.
(User-facing) Fields
tt::Matrix{T}: The transit times of each body.dtdq0::Array{T,4}: Derivatives of the transit times with respect to the initial Cartesian coordinates and masses.dtdelements::Array{T,4}: Derivatives of the transit times with respect to the initial orbital elements and masses.
NbodyGradient.TransitTiming — MethodTransitTiming(tmax, ic; ti)Constructor for TransitTiming type.
Arguments
tmax::T: Expected total elapsed integration time. (Allocates arrays accordingly)ic::ElementsIC{T}: Initial conditions for the system
Optional
ti::Int64=1: Index of the body with respect to which transits are measured. (Default is the central body)