Low-level API

These functions are not meant to be used outside of the package. They are documented only to aid future developments of the package.

Ephemeris Provider

Ephemerides.get_dafFunction
daf(eph::EphemerisProvider)

Return the DAF files stored in the ephemeris provider.

source
daf(eph::EphemerisProvider, id::Int)

Return the DAF file in the ephemeris provider at index id.

source

DAF Routines

DAF Header

Ephemerides.DAFHeaderType
DAFHeader

The DAF header, or file record, is the first physical record in a DAF and stores general information about the content of the file.

Fields

  • ndInt32 number of double components in each array summary
  • niInt32 number of integer components in each array summary
  • fwdInt32 record number of initial summary record
  • bwdInt32 record number of final summary record
  • ffaInt32 first free address of the file
  • nameString internal name of the file
  • lendBool true if the file was generated in little endian

References

See Also

See also DAF and EphemerisProvider

source
Ephemerides.initial_recordFunction
initial_record(head::DAFHeader)

Return the record number of the initial summary record in the DAF

source
initial_record(daf::DAF)

Return the record number of the initial summary record in the DAF.

source
Ephemerides.final_recordFunction
final_record(head::DAFHeader)

Return the record number of the final summary record in the DAF

source
final_record(daf::DAF)

Return the record number of the final summary record in the DAF.

source
Ephemerides.free_addressFunction
free_address(head::DAFHeader)

Return the first free address in the file, i.e., the address at which the first element of the next array is to be added.

source
free_address(daf::DAF)

Return the first free address in the file, i.e., the address at which the first element of the next array is to be added.

source
Ephemerides.endianFunction
endian(head::DAFHeader)

Return true if the DAF file is in little-endian.

source
endian(daf::DAF)

Return true if the DAF is in little-endian.

source

DAF Descriptor

Ephemerides.DAFSegmentDescriptorType
DAFSegmentDescriptor

A container object to store both SPK and PCK descriptors information.

Fields

  • segtypeInt32 SPK/PCK segment type
  • tstartFloat64 initial segment type, in TDB seconds since J2000.0
  • tendFloat64 final segment type, in TDB seconds since J2000.0
  • tidInt32 target object NAIF ID
  • cidInt32 center object NAIF ID
  • axesidInt32 reference axes ID. Defaults to -1 for PCKs
  • iaaInt32 initial array address
  • faaInt32 final array address

References

source
Ephemerides.initial_timeFunction
initial_time(desc::DAFSegmentDescriptor)

Return the initial epoch of the interval for which ephemeris data are contained in the segment, in seconds since J2000.0

source
initial_time(link::SPKLink)

Return the initial epoch of the interval for which ephemeris data are contained in the segment associated to this link, in seconds since J2000.0

source
Ephemerides.final_timeFunction
final_time(desc::DAFSegmentDescriptor)

Return the final epoch of the interval for which ephemeris data are contained in the segment, in seconds since J2000.0

source
final_time(link::SPKLink)

Return the final epoch of the interval for which ephemeris data are contained in the segment associated to this link, in seconds since J2000.0

source
Ephemerides.centerFunction
center(desc::DAFSegmentDescriptor)

Return the NAIF integer code for the reference object or axes for SPK and PCK, respectively.

source
Ephemerides.targetFunction
target(desc::DAFSegmentDescriptor)

Return the NAIF integer code for the target object or axes for SPK and PCK, respectively.

source
Ephemerides.axesFunction
axes(desc::DAFSegmentDescriptor)

Return the NAIF integer code for the reference axes. It is valid only for SPK files and defaults to -1 for PCKs.

source
Ephemerides.get_segment_boundariesFunction
get_segment_boundaries(desclist::Vector{DAFSegmentDescriptor})

Parse all the segment descriptors of a given (center, target) pair and return a set of initial and final times, in TDB seconds since J2000, representing all the time sub-windows in which the ephemeris data for this pair is defined.

source

DAF

Ephemerides.DAFType
DAF

Container to hold the information of NAIF's Double precision Array File (DAF).

Fields

  • filepathString system filepath of the DAF
  • arrayVector{UInt8} binary content of the DAF
  • headerDAFHeader file record of the DAF
  • commentString text within the DAF comment area
  • ftypeInt file type, equals 1 for SPK and 2 for PCK
  • desc – DAF PCK/SPK segment descriptors
  • seglistSPKSegmentList list of the SPK/PCK segments within the DAF

References

See Also

See also DAFHeader, Ephemerides.SPKSegmentList and EphemerisProvider

source
Ephemerides.is_little_endianFunction
is_little_endian(array::Vector{UInt8})

Return true if the array corresponds to the string indicating a little-endian format.

source
Ephemerides.SPKLinkType
SPKLink

A link object to create a mapping between DAFSegmentDescriptor and its actual location within an EphemerisProvider object.

Fields

  • descDAFSegmentDescriptor for the segment associated to this link
  • fidInt index of the DAF containg the link data.
  • lidInt field number in the SPKSegmentList for this segment type.
  • eidInt index of the inner segment list that stores this SPK segment.
  • fctInt 1 or -1 depending on whether the (from, to) directions must be reversed.

See Also

See also SPKLinkTable, SPKSegmentList and add_spklinks!.

source
Ephemerides.list_idFunction
list_id(link::SPKLink)

Return the index of the list containing the segments of the given SPK/PCK type.

source
Ephemerides.reverse_linkFunction
reverse_link(link::SPKLink)

Reverse the sign, i.e. change the sign of the multiplicative factor, of the link.

source

SPK Segment List

Ephemerides.get_segmentFunction
get_segment(list::SPKSegmentList, lid::Int, eid::Int)

Return the segment contained in the lid list at index eid.

source

SPK Segment Types

Abstract SPK Types

SPK Type 1 and 21

Ephemerides.SPKSegmentHeader1Type
SPKSegmentHeader1 <: AbstractSPKHeader

Header instance for SPK segments of type 1 and 21.

Fields

  • nInt number of records in the segment
  • ndirsInt number of directory epochs
  • epochs – Storage for directory epochs or epochs (when ndirs = 0)
  • iaa - Int initial segment file address
  • etidInt initial address for the epoch table (after all the MDA records)
  • recsize - Int Number of double numbers stored in each MDA record
  • maxdim - Int MDA dimension (fixed to 15 for type 1)
source
Ephemerides.SPKSegmentCache1Type
SPKSegmentCache1 <: AbstractSPKCache

Cache instance for SPK segments of type 1 and 21. The fields contained within this cache are taken from the FORTRAN NAIF's SPICE implementation for type 1 SPK segments.

Fields

  • tl – Reference epoch of the difference line.
  • g – Stepsize function vector.
  • refpos – Reference position vector.
  • refvel – Reference velocity vector.
  • dt – Modified Divided Difference arrays, with size (maxdim, 3)
  • kqmax – Maximum integration order plus 1.
  • kq – Integration order array.
  • id – Index of the currently loaded logical record.
  • fc – Buffer for the MDA computations.
  • wc – Buffer for the MDA computations.
  • w – Buffer for the MDA computations.
  • vct – Buffer for the MDA computations.
source
Ephemerides.SPKSegmentType1Type
SPKSegmentType1 <: AbstractSPKSegment

Segment instance for SPK segments of type 1 and 21, which contain Modified Difference Arrays (MDA). This data type is normally used for spacecraft whose ephemerides are produced by JPL's principal trajectory integrator DPTRAJ.

Fields

  • head – Segment header
  • cache – Segment cache

References

source

SPK Type 2 and 3

Ephemerides.SPKSegmentHeader2Type
SPKSegmentHeader2 <: AbstractSPKHeader

Header instance for SPK segments of type 2 and 3.

Fields

  • tstartFloat64 initial epoch of the first record, in seconds since J2000
  • tlenFloat64 interval length covered by each record, in seconds
  • orderInt polynomial order
  • NInt number of coefficients in each window
  • nInt number of records in the segment
  • recsizeInt byte size of each logical record
  • ncompInt number of vector components
  • iaaInt initial segment file address
  • typeInt SPK segment type, either 2 or 2
source
Ephemerides.SPKSegmentCache2Type
SPKSegmentCache2 <: AbstractSPKCache

Cache instance for SPK segments of type 2 and 3.

Fields

  • A – Chebyshev's polynomial coefficients, with size (ncomp, order)
  • p – Stores the record mid point and radius and scale factor
  • buff – Stores the buffers for the Chebyshev polynomials
  • id – Index of the currently loaded logical record
source
Ephemerides.SPKSegmentType2Type
SPKSegmentType2 <: AbstractSPKSegment

Segment instance for SPK segments of type 2 and 3, which contain Chebyshev polynomial coefficients for the position and/or state of the body as function of time. This data type is normally used for planet barycenters, and for satellites whose ephemerides are integrated.

Fields

  • head – Segment header
  • cache – Segment cache

References

source

SPK Type 5

Ephemerides.SPKSegmentHeader5Type
SPKSegmentHeader5 <: AbstractSPKHeader

Header instance for SPK segments of type 5.

Fields

  • GMFloat64 Gravitational constant
  • nInt number of states
  • ndirsInt number of epoch directories
  • etidInt initial address for the epoch table (after all the state data)
  • epochs – Storage for directory epochs or epochs (when ndirs = 0)
  • iaa - Int initial segment file address
source
Ephemerides.SPKSegmentCache5Type
SPKSegmentCache5 <: AbstractSPKCache

Cache instance for SPK segments of type 5.

Fields

  • c1 – Twobody propagation cache for the left state.
  • c2 – Twobody propagation cache for the right state.
  • epochs – Epochs associated to the two states.
  • id – Index of the currently loaded logical record.
source

SPK Type 8 and 12

Ephemerides.SPKSegmentHeader8Type
SPKSegmentHeader8 <: AbstractSPKHeader

Header instance for SPK segments of type 8 and 12.

Fields

  • tstartFloat64 segment starting epoch, in TDB seconds since J2000
  • tlenFloat64 interval length, in seconds
  • orderInt interpolating polynomial degree
  • NInt group size
  • nInt number of states in the segment
  • iaa - Int initial segment file address
  • isevenBool true for even group size
  • typeInt SPK type (either 8 or 12)
source
Ephemerides.SPKSegmentCache8Type
SPKSegmentCache8 <: AbstractSPKCache

Cache instance for SPK segments of type 8 and 12.

Fields

  • states – Matrix storing the states of the interpolating points.
  • buff – Buffers to compute the interpolating polynomials.
  • id – Index of the currently loaded logical record.
source

SPK Type 9 and 13

Ephemerides.SPKSegmentHeader9Type
SPKSegmentHeader9 <: AbstractSPKHeader

Header instance for SPK segments of type 9 and 13.

Fields

  • nInt number of states in the segment
  • ndirsInt number of epoch directories
  • epochs – Storage for directory epochs or epochs (when ndirs = 0)
  • iaa - Int initial segment file address
  • etidInt initial address for the epoch table (after all the state data)
  • orderInt interpolating polynomial degree
  • NInt group size
  • isevenBool true for even group size
  • typeInt SPK type (either 9 or 13)
source
Ephemerides.SPKSegmentCache9Type
SPKSegmentCache9 <: AbstractSPKCache

Cache instance for SPK segments of type 9 and 13.

Fields

  • epochs – Epochs of the interpolating points.
  • states – Matrix storing the states of the interpolating points.
  • buff – Buffers to compute the interpolating polynomials.
  • id – Index of the currently loaded logical record.
source

SPK Type 14

Ephemerides.SPKSegmentHeader14Type
SPKSegmentHeader14 <: AbstractSPKHeader

Header instance for SPK segments of type 14.

Fields

  • orderInt interpolating polynomial degree
  • nInt number of packets in the segment
  • ndirsInt number of epoch directories
  • epochs – Storage for directory epochs or epochs (when ndirs = 0)
  • etidInt initial address for the epoch table (after all the state data)
  • ptidInt initial address for the packet table (after the constants)
  • pktsizeInt size of each data packet excluding the packet information area.
  • pktoffInt offset of the packet data from the packet start
  • ncompInt number of states coefficients (= 6 for SPK 14)
  • NInt number of polynomial coefficients
source
Note

SPK segments of type 14 have the same cache structure of SPK type 2 and 3.

SPK Type 15

Ephemerides.SPKSegmentHeader15Type
SPKSegmentHeader15 <: AbstractSPKHeader

Header instance for SPK segments of type 15.

Fields

  • epoch – Epoch of periapsis
  • tp – Trajectory pole, i.e., vector parallel to the angular momentum of the orbit
  • pv – Central body north pole unit vector
  • pa – Periapsis unit vector at epoch
  • p – Semi-latus rectum
  • ecc – Eccentricity
  • j2f – J2 processing flag
  • vj2 – J2 validation flag, true if the orbit shape is compliant with J2 pertubations.
  • GM – Central body gravitational constant (km³/s²)
  • J2 – Central body J2
  • R – Central body radius (km)
  • dmdt – Mean anomaly rate of change (rad/s)
  • kn – Gain factor for the regression of the nodes
  • kp – Gain factor for the precession of the pericenter
source
Note

The cache of SPK Type 15 segments is made of Ephemerides.TwoBodyUniversalCache objects.

SPK Type 17

Ephemerides.SPKSegmentHeader17Type
SPKSegmentHeader17 <: AbstractSPKHeader

Header instance for SPK segments of type 17.

Fields

  • epoch: epoch of periapsis (s)
  • sma: semi-major axis (km)
  • h: H term of the equinoctial elements
  • k: K term of the equinoctial elements
  • lon: mean longitude at epoch (rad)
  • p: P term of the equinoctial elements
  • q: Q term of the equinoctial elements
  • dlpdt: rate of longitude of the periapse (rad/s)
  • dmldt: mean longitude rate (mean motion rate), (rad/s)
  • dnodedt: longitude of the ascending node rate (rad/s)
  • ra: equatorial pole right ascension (rad)
  • de: equatorial pole declination (rad)
  • R: Rotation matrix from planetary equator to inertial reference frame
source
Note

SPK segments of type 17 do not require a cache structure.

SPK Type 18 and 19

Ephemerides.SPKSegmentHeader18Type
SPKSegmentHeader18 <: AbstractSPKHeader

Header instance for SPK segments of type 18.

Fields

  • nInt number of states in the segment
  • ndirsInt number of epoch directories
  • epochs – Storage for directory epochs or epochs (when ndirs = 0)
  • iaa - Int initial segment file address
  • etidInt initial address for the epoch table (after all the state data)
  • orderInt interpolating polynomial degree
  • NInt group size
  • subtypeInt type 18 subtype, either 0 (Hermite) or 1 (Lagrange)
  • packetsizeInt packet size for each point, either 12 (Hermite) or 6 (Lagrange)
source
Ephemerides.SPKSegmentCache18Type
SPKSegmentCache18 <: AbstractSPKCache

Cache instance for SPK segments of type 18.

Fields

  • p – Vector storing indexes of the first and last points as well as the window size.
  • epochs – Epochs of the interpolating points.
  • states – Matrix storing the states of the interpolating points.
  • buff – Buffers to compute the interpolating polynomials.
source
Ephemerides.SPKSegmentHeader19Type
SPKSegmentHeader19 <: AbstractSPKHeader

Header instance for SPK segments of type 19.

Fields

  • nInt number of states in the segment.
  • ndirsInt number of epoch directories.
  • times – Storage for interval directories or start times (when ndirs = 0).
  • iaa - Int initial segment file address.
  • etidInt byte address for the interval table (after all the minisegment data).
  • ptidInt byte for the pointer table.
  • usefirstBool boundary flag, true if the preceding segment should be used.
  • typeInt either type 18 or 19.
source
Ephemerides.SPKSegmentCache19Type
SPKSegmentCache19 <: AbstractSPKCache

Cache instance for SPK segments of type 19.

Fields

  • minihead – Header with the mini-segment properties.
  • minidata – Cache for the mini-segment.
  • id – Index of the currently loaded mini-segment.
source
Note

SPK segments of type 18 are the only ones that do not posses a dedicated SPK segment type structure, because they are treated as special cases (i.e., single minisegments) of the type 19.

SPK Type 20

Ephemerides.SPKSegmentHeader20Type
SPKSegmentHeader20 <: AbstractSPKHeader

Header instance for SPK segments of type 20.

Fields

  • dscaleFloat64 length conversion factor
  • tscaleFloat64 time conversion factor
  • tstartFloat64 initial epoch of the first record, in seconds since J2000
  • tlenFloat64 interval length covered by each record, in seconds
  • recsizeInt byte size of each logical record
  • orderInt polynomial order
  • NInt number of coefficients in each window
  • nInt number of records in the segment
  • iaaInt initial segment file address
source
Ephemerides.SPKSegmentCache20Type
SPKSegmentCache20 <: AbstractSPKCache

Cache instance for SPK segments of type 20.

Fields

  • id – Index of the currently loaded logical record
  • p – Stores the record position constants
  • A – Chebyshev's polynomial coefficients, with size (ncomp, order)
  • buff – Stores the buffers for the Chebyshev polynomials
source
Ephemerides.SPKSegmentType20Type
SPKSegmentType2 <: AbstractSPKSegment

Segment instance for SPK segments of type 20, which contain Chebyshev polynomial coefficients for the position and/or state of the body as function of time. This data type is normally used for planet barycenters, and for satellites whose ephemerides are integrated.

Fields

  • head – Segment header
  • cache – Segment cache

References

source

SPK Utility Functions

Ephemerides.normalise_timeFunction
normalise_time(cache::SPKSegmentCache2, time::Number)

Transform time in an interval between [-1, 1] for compliance with Chebyshev polynomials.

source
normalise_time(head::SPKSegmentHeader8, time::Number, index::Int)

Returned a normalised time that starts at 1 at the beginning of the interval.

source
normalise_time(head::SPKSegmentHeader20, time::Number, index::Int)
source
Ephemerides.find_logical_recordFunction
find_logical_record(daf::DAF, head::SPKSegmentHeader1, time::Number)
source
find_logical_record(head::SPKSegmentHeader2, time::Number)
source
find_logical_record(daf::DAF, head::SPKSegmentHeader5, time::Number)
source
find_logical_record(head::SPKSegmentHeader8, time::Number)
source
find_logical_record(daf::DAF, head::SPKSegmentHeader1, time::Number)
source
find_logical_record(daf::DAF, head::SPKSegmentHeader14, time::Number)
source
find_logical_record(daf::DAF, head::SPKSegmentHeader18, time::Number)
source
find_logical_record(head::SPKSegmentHeader20, time::Number)
source
Ephemerides.get_coefficients!Function
get_coefficients!(daf::DAF, head::SPKSegmentHeader1, cache::SPKSegmentCache1, index::Int)
source
get_coefficients!(daf::DAF, head::SPKSegmentHeader2, cache::SPKSegmentCache2, index::Int)
source
get_coefficients!(daf::DAF, head::SPKSegmentHeader5, cache::SPKSegmentCache5, index::Int)
source
get_coefficients!(daf::DAF, head::SPKSegmentHeader8, cache::SPKSegmentCache8, index::Int)
source
get_coefficients!(daf::DAF, head::SPKSegmentHeader9, cache::SPKSegmentCache9, index::Int)
source
get_coefficients!(daf::DAF, head::SPKSegmentHeader14, cache::SPKSegmentCache14, index::Int)
source
get_coefficients!(daf::DAF, head::SPKSegmentHeader18, cache::SPKSegmentCache18, first::Int, last::Int)
source
get_coefficients!(daf::DAF, head::SPKSegmentHeader20, cache::SPKSegmentCache20, index::Int)
source

Interpolating Functions

Caches

Ephemerides.get_bufferFunction
get_buffer(c::InterpCache, idx::int, x::Number)

Return the idx-th buffer from the corresponding DiffCache depending on the type of x.

source

Chebyshev Polynomials

Ephemerides.chebyshevFunction
chebyshev(cache::InterpCache, cₖ, t::Number, idx::Int, N::Int, ibuff::Int)

Evaluate a sum of Cheybyshev polynomials of the first kind at t using a recursive algorithm. It simultenously evalutes the 3 state components. idx is the index of the starting row (in 0-based notation) in the matrix of coefficients cₖ. ibuff is the index of the first free buffer.

Note

x is a re-work of the actual ascissa value that lies between [-1, 1]

See Also

See also ∂chebyshev, ∂²chebyshev and ∂³chebyshev

source
Ephemerides.∂chebyshevFunction
∂chebyshev(cache::InterpCache, cₖ, t::Number, idx::Int, N::Int, Δt, ibuff=1)

Evaluate a sum of Cheybyshev polynomials of the first kind and its derivative at t using a recursive algorithm. It simultenously evalutes the 3 state components. idx is the index of the starting row (in 0-based notation) in the matrix of coefficients cₖ. ibuff is the index of the first free buffer.

Note

x is a re-work of the actual ascissa value that lies between [-1, 1]

See Also

See also chebyshev, ∂²chebyshev and ∂³chebyshev

source
Ephemerides.∂²chebyshevFunction
∂²chebyshev(cache::InterpCache, cₖ, t::Number, idx::Int, N::Int, Δt, ibuff=1)

Evaluate a sum of Cheybyshev polynomials of the first kind and its two derivatives at t using a recursive algorithm. It simultenously evalutes the 3 state components. idx is the index of the starting row (in 0-based notation) in the matrix of coefficients cₖ. ibuff is the index of the first free buffer.

Note

x is a re-work of the actual ascissa value that lies between [-1, 1]

See Also

See also chebyshev, ∂chebyshev and ∂³chebyshev

source
Ephemerides.∂³chebyshevFunction
∂³chebyshev(cache::InterpCache, cₖ, t::Number, idx::Int, N::Int, Δt, ibuff=1)

Evaluate a sum of Cheybyshev polynomials of the first kind and its three derivatives at t using a recursive algorithm. It simultenously evalutes the 3 state components. idx is the index of the starting row (in 0-based notation) in the matrix of coefficients cₖ. ibuff is the index of the first free buffer.

Note

x is a re-work of the actual ascissa value that lies between [-1, 1]

See Also

See also chebyshev, ∂chebyshev and ∂²chebyshev

source
Ephemerides.∫chebyshevFunction
∫chebyshev(cache::InterpCache, cₖ, t::Number, N::Int, Δt, tlen, p₀)

Evaluate the integral of a sum of Cheybyshev polynomials of the first kind using a recursive algorithm. It simultenously evalutes the 3 state components. It assumes the Chebyshev polynomials up to degree N have already been computed and are stored in the buffer with index ibuff. tlen is the size of the record interval, Δt is the timescale factor, and p₀ is a vector containing the position coefficients at the midpoint (i.e., when the integral is evaluated at t = 0).

Note

x is a re-work of the actual ascissa value that lies between [-1, 1]

source
∫chebyshev(cache::InterpCache, cₖ, t::Number, N::Int, Δt, tlen, p₀)

Evaluate the integral of a sum of Cheybyshev polynomials of the first kind using a recursive algorithm. It simultenously evalutes the 3 state components. This function simultaneously computes both the Chebyshev polynomials as well as their integrals.

source

Lagrange Polynomials

Ephemerides.lagrangeFunction
lagrange(cache::InterpCache, states, x, idx::Int, N::Int)

Recursively evaluate a Lagrange polynomial at x by using Neville's algorithm. This function is valid only for equally-spaced polynomials. idx is the index of the desired state and N is the number of coefficients of the polynomial.

Note

x is a re-work of the actual ascissa value that starts at 1

See Also

See also ∂lagrange and ∂²lagrange.

source
lagrange(cache::InterpCache, states, epochs, x, idx::Int, N::Int)

Recursively evaluate a Lagrange polynomial at x by using Neville's algorithm. This function handles unequally-spaced polynomials, where the coefficients in states are interpolated at epochs.

source
Ephemerides.∂lagrangeFunction
∂lagrange(cache::InterpCache, states, x, idx::Int, N::Int, Δt::Number)

Recursively evaluate a Lagrange polynomial and its derivative at x by using Neville's algorithm. This function is valid only for equally-spaced polynomials. idx is the index of the desired state, N is the number of coefficients of the polynomial and Δtis the length of the polynomial interval

Note

x is a re-work of the actual ascissa value that starts at 1

See Also

See also lagrange and ∂²lagrange

source
∂lagrange(cache::InterpCache, states, epochs, x, idx::Int, N::Int)

Recursively evaluate a Lagrange polynomial and its derivative at x by using Neville's algorithm. This function handles unequally-spaced polynomials, where the coefficients in states are interpolated at epochs.

source
Ephemerides.∂²lagrangeFunction
∂²lagrange(cache::InterpCache, states, x, idx::Int, N::Int, Δt::Number)

Recursively evaluate a Lagrange polynomial and its two derivatives at x by using Neville's algorithm. This function is valid only for equally-spaced polynomials. idx is the index of the desired state, N is the number of coefficients of the polynomial and Δtis the length of the polynomial interval

Note

x is a re-work of the actual ascissa value that starts at 1

See Also

See also lagrange and ∂lagrange

source
∂²lagrange(cache::InterpCache, states, epochs, x, idx::Int, N::Int)

Recursively evaluate a Lagrange polynomial and its two derivatives at x by using Neville's algorithm. This function handles unequally-spaced polynomials, where the coefficients in states are interpolated at epochs.

See Also

See also lagrange and ∂lagrange

source

Hermite Polynomials

Ephemerides.hermiteFunction
hermite(cache::InterpCache, states, x, idx::Int, N::Int, Δt::Number)

Evalute a Hermite polynomial at x using a recursive algorithm. This function is valid only for equally-spaced polynomials. idx is the index of the desired state, N is the number of coefficients of the polynomial and Δt is the length of the polynomial interval.

Note

x is a re-work of the actual ascissa value that starts at 1

See Also

See also ∂hermite, ∂²hermite and ∂³hermite.

source
hermite(cache::InterpCache, states, epochs, x, idx::Int, N::Int)

Evalute a Hermite polynomial at x using a recursive algorithm. This function handles unequally-spaced polynomials, where the coefficients in states are interpolated at epochs.

source
Ephemerides.∂hermiteFunction
∂hermite(cache::InterpCache, states, x, idx::Int, N::Int, Δt::Number)

Evalute a Hermite polynomial and its derivative at x using a recursive algorithm. This function is valid only for equally-spaced polynomials. idx is the index of the desired state, N is the number of coefficients of the polynomial and Δt is the length of the polynomial interval.

Note

x is a re-work of the actual ascissa value that starts at 1

See Also

See also hermite, ∂²hermite and ∂³hermite.

source
∂hermite(cache::InterpCache, states, epochs, x, idx::Int, N::Int)

Evalute a Hermite polynomial and its derivative at x using a recursive algorithm. This function handles unequally-spaced polynomials, where the coefficients in states are interpolated at epochs.

source
Ephemerides.∂²hermiteFunction
∂²hermite(cache::InterpCache, states, x, idx::Int, N::Int, Δt::Number)

Evalute a Hermite polynomial and its two derivatives at x using a recursive algorithm. This function is valid only for equally-spaced polynomials. idx is the index of the desired state, N is the number of coefficients of the polynomial and Δt is the length of the polynomial interval.

Note

x is a re-work of the actual ascissa value that starts at 1

See Also

See also hermite, ∂hermite and ∂³hermite.

source
∂²hermite(cache::InterpCache, states, epochs, x, idx::Int, N::Int)

Evalute a Hermite polynomial and its two derivatives at x using a recursive algorithm. This function handles unequally-spaced polynomials, where the coefficients in states are interpolated at epochs.

source
Ephemerides.∂³hermiteFunction
∂³hermite(cache::InterpCache, states, x, idx::Int, N::Int, Δt::Number)

Evalute a Hermite polynomial and its three derivatives at x using a recursive algorithm. This function is valid only for equally-spaced polynomials. idx is the index of the desired state, N is the number of coefficients of the polynomial and Δt is the length of the polynomial interval.

Note

x is a re-work of the actual ascissa value that starts at 1

See Also

See also hermite, ∂hermite and ∂²hermite.

source
∂³hermite(cache::InterpCache, states, epochs, x, idx::Int, N::Int)

Evalute a Hermite polynomial and its three derivatives at x using a recursive algorithm. This function handles unequally-spaced polynomials, where the coefficients in states are interpolated at epochs.

source

Introspection

Ephemerides.initial_timesFunction
initial_times(record::AbstractEphemRecord)

Recover the initial times of each sub-window in which the ephemeris data of record is defined, expressed in TDB seconds since J2000

source
Ephemerides.final_timesFunction
final_times(record::AbstractEphemRecord)

Recover the final times of each sub-window in which the ephemeris data of record is defined, expressed in TDB seconds since J2000

source

TwoBody Routines

Ephemerides.TwoBodyUniversalCacheType
TwoBodyUniversalCache

A container to store precomputed quantities required for the two-body propagation based on universal variables. Only the quantities that depend on the initial state are computed.

source
Ephemerides.update_cache!Function
update_cache!(c::TwoBodyUniversalCache)

Update the precomputed values in the cache using the position and velocity in c.

source
Ephemerides.propagate_twobodyFunction
propagate_twobody(cache::TwoBodyUniversalCache, Δt::Number)

Propagate the state vector in cache of Δt using the universal variables formulation for Kepler's Equation and the Lagrange coefficients f and g.

Note

This routine is valid for any type of orbit and uses a bisection method to find the root of the universal variables Kepler's equation. The algorithm has been directly taken from the SPICE toolkit prob2b.f.

References

source
Ephemerides.stumpffFunction
stumpff(x::Number, p::AbstractVector)

Compute Stumpff's functions from C₀ up to C₃ at x.

Note

This routine uses the trigonometrical expressions of the functions when the absolute value of x is greater or equal to 1. If that is not the case, the C₂ and C₃ functions are computed from a truncated expression of the Maclaurin series at order 11, which guarantees a higher precision and avoid overflow errors when x is null.

References

source