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_daf — Functiondaf(eph::EphemerisProvider)Return the DAF files stored in the ephemeris provider.
daf(eph::EphemerisProvider, id::Int)Return the DAF file in the ephemeris provider at index id.
Ephemerides.spk_links — Functionspk_links(eph::EphemerisProvider)Return the [SPKLinkTable] for the SPK segments.
Ephemerides.pck_links — Functionpck_links(eph::EphemerisProvider)Return the SPKLinkTable for the PCK segments.
DAF Routines
DAF Header
Ephemerides.DAFHeader — TypeDAFHeaderThe DAF header, or file record, is the first physical record in a DAF and stores general information about the content of the file.
Fields
nd–Int32number of double components in each array summaryni–Int32number of integer components in each array summaryfwd–Int32record number of initial summary recordbwd–Int32record number of final summary recordffa–Int32first free address of the filename–Stringinternal name of the filelend–Booltrue if the file was generated in little endian
References
See Also
See also DAF and EphemerisProvider
Ephemerides.initial_record — Functioninitial_record(head::DAFHeader)Return the record number of the initial summary record in the DAF
initial_record(daf::DAF)Return the record number of the initial summary record in the DAF.
Ephemerides.final_record — Functionfinal_record(head::DAFHeader)Return the record number of the final summary record in the DAF
final_record(daf::DAF)Return the record number of the final summary record in the DAF.
Ephemerides.free_address — Functionfree_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.
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.
Ephemerides.endian — Functionendian(head::DAFHeader)Return true if the DAF file is in little-endian.
endian(daf::DAF)Return true if the DAF is in little-endian.
Ephemerides.filename — Functionfilename(head::DAFHeader)Return the internal description of the DAF.
Ephemerides.summary_size — Functionsummary_size(head::DAFHeader)Compute the size of a single summary record of a DAF file, in bytes.
DAF Descriptor
Ephemerides.DAFSegmentDescriptor — TypeDAFSegmentDescriptorA container object to store both SPK and PCK descriptors information.
Fields
segtype–Int32SPK/PCK segment typetstart–Float64initial segment type, in TDB seconds since J2000.0tend–Float64final segment type, in TDB seconds since J2000.0tid–Int32target object NAIF IDcid–Int32center object NAIF IDaxesid–Int32reference axes ID. Defaults to -1 for PCKsiaa–Int32initial array addressfaa–Int32final array address
References
Ephemerides.segment_type — Functionsegment_type(desc::DAFSegmentDescriptor)Return the SPK/PCK segment type.
Ephemerides.initial_time — Functioninitial_time(desc::DAFSegmentDescriptor)Return the initial epoch of the interval for which ephemeris data are contained in the segment, in seconds since J2000.0
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
Ephemerides.final_time — Functionfinal_time(desc::DAFSegmentDescriptor)Return the final epoch of the interval for which ephemeris data are contained in the segment, in seconds since J2000.0
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
Ephemerides.center — Functioncenter(desc::DAFSegmentDescriptor)Return the NAIF integer code for the reference object or axes for SPK and PCK, respectively.
Ephemerides.target — Functiontarget(desc::DAFSegmentDescriptor)Return the NAIF integer code for the target object or axes for SPK and PCK, respectively.
Ephemerides.axes — Functionaxes(desc::DAFSegmentDescriptor)Return the NAIF integer code for the reference axes. It is valid only for SPK files and defaults to -1 for PCKs.
Ephemerides.initial_address — Functioninitial_address(desc::DAFSegmentDescriptor)Return the initial address of the segment array in the DAF.
Ephemerides.final_address — Functionfinal_address(desc::DAFSegmentDescriptor)Return the final address of the segment array in teh DAF.
Ephemerides.get_segment_boundaries — Functionget_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.
Ephemerides.parse_spk_segment_descriptor — Functionparse_spk_segment_descriptor(summary::Vector{UInt8}, lend::Bool)Create a DAFSegmentDescriptor object by parsing a binary SPK segment descriptor.
References
Ephemerides.parse_pck_segment_descriptor — Functionparse_pck_segment_descriptor(summary::Vector{UInt8}, lend::Bool)Create a DAFSegmentDescriptor object by parsing a binary PCK segment descriptor. A default value of -1 is used to fill the reference frame field. The target and center fields are used for the actual target and center axes.
References
DAF
Ephemerides.DAF — TypeDAFContainer to hold the information of NAIF's Double precision Array File (DAF).
Fields
filepath–Stringsystem filepath of the DAFarray–Vector{UInt8}binary content of the DAFheader–DAFHeaderfile record of the DAFcomment–Stringtext within the DAF comment areaftype–Intfile type, equals 1 for SPK and 2 for PCKdesc– DAF PCK/SPK segment descriptorsseglist–SPKSegmentListlist of the SPK/PCK segments within the DAF
References
See Also
See also DAFHeader, Ephemerides.SPKSegmentList and EphemerisProvider
Ephemerides.DAF_RECORD_LENGTH — ConstantEphemerides.FTPSTR — ConstantEphemerides.comment — Functionget_comment(daf::DAF)Return the comment written in the DAF comment section.
Ephemerides.header — Functionheader(spk::AbstractSPKSegment)Return the segment header.
get_header(daf::DAF)Return the DAFHeader header of the DAF.
Ephemerides.array — Functionget_array(daf::DAF)Return the byte content of the DAF file.
Ephemerides.descriptors — Functionget_descriptors(daf::DAF)Return the SPK/PCK segment descriptors contained in the DAF.
Ephemerides.segment_list — Functionget_segment_list(daf::DAF)Return the Ephemerides.SPKSegmentList list of segments stored in the DAF.
Ephemerides.filepath — Functionfilepath(daf::DAF)Return the system path of the DAF.
Ephemerides.is_spk — Functionis_spk(daf::DAF)Return true if the DAF stores SPK data.
Ephemerides.is_pck — Functionis_pck(daf::DAF)Return true if the DAF stores PCK data.
Ephemerides.parse_daf_comment — Functionparse_daf_comment(array::Vector{UInt8}, header::DAFHeader)Retrieve the comment section of a binary DAF.
References
See Also
See also DAF, DAFHeader and parse_daf_summaries
Ephemerides.parse_daf_summaries — Functionparse_daf_summaries(array::Vector{UInt8}, head::DAFHeader)Parse the DAF binary content and retrieve all the summary records.
References
See Also
See also DAF, DAFHeader and parse_daf_comment
Ephemerides.initialise_segments! — Functioninitialise_segments!(daf::DAF)Fill the Ephemerides.SPKSegmentList by initialising the SPK/PCK segments associated to all the descriptors stores within the DAF.
See Also
See also DAF and create_spk_segment
Ephemerides.create_spk_segment — Functioncreate_spk_segment(daf::DAF, desc::DAFSegmentDescriptor)Initialise an SPK segment according to the segment type defined in the DAFSegmentDescriptor desc.
Ephemerides.get_record — Functionget_record(array::Vector{UInt8}, index::Integer)Retrieve a whole DAF record at position index.
Ephemerides.is_little_endian — Functionis_little_endian(array::Vector{UInt8})Return true if the array corresponds to the string indicating a little-endian format.
SPK Links
Ephemerides.SPKLink — TypeSPKLinkA link object to create a mapping between DAFSegmentDescriptor and its actual location within an EphemerisProvider object.
Fields
desc–DAFSegmentDescriptorfor the segment associated to this linkfid–Intindex of the DAF containg the link data.lid–Intfield number in theSPKSegmentListfor this segment type.eid–Intindex of the inner segment list that stores this SPK segment.fct–Int1 or -1 depending on whether the (from, to) directions must be reversed.
See Also
See also SPKLinkTable, SPKSegmentList and add_spklinks!.
Ephemerides.SPKLinkTable — TypeSPKLinkTableDictionary object providing all the SPKLink available between a set of (from, to) objects
Ephemerides.descriptor — Functiondescriptor(link::SPKLink)Return the SPK/PCK segment descriptor associated to this link.
Ephemerides.file_id — Functionfile_id(link::SPKLink)Return the DAF file index.
Ephemerides.list_id — Functionlist_id(link::SPKLink)Return the index of the list containing the segments of the given SPK/PCK type.
Ephemerides.element_id — Functionelement_id(link::SPKLink)Return the segment index in the inner SPK/PCK segment list.
Ephemerides.factor — Functionfactor(link::SPKLink)Return the direction multiplicative factor.
Ephemerides.reverse_link — Functionreverse_link(link::SPKLink)Reverse the sign, i.e. change the sign of the multiplicative factor, of the link.
Ephemerides.create_linktables — Functioncreate_linktables(dafs::Vector{DAF})Create the SPK and PCK SPKLinkTable for all the segments stores in the input DAFs.
Ephemerides.add_spklinks! — Functionadd_spklinks!(table::SPKLinkTable, daf::DAF, fid::Int)Insert in the input SPKLinkTable all the SPK or PCK links associated to the segment descriptors of the input DAF.
SPK Segment List
Ephemerides.SPKSegmentList — TypeSPKSegmentListA container object to efficiently store all the different SPK segments that are contained within a single DAF file.
SPKSegmentList()Initialises an empty SPKSegmentList object.
See also
See also Ephemerides.add_segment!
Ephemerides.add_segment! — Functionadd_segment!(list::SPKSegmentList, spk::AbstractSPKSegment)Add the SPK segment to the proper vector within the given Ephemerides.SPKSegmentList list
Ephemerides.get_segment — Functionget_segment(list::SPKSegmentList, lid::Int, eid::Int)Return the segment contained in the lid list at index eid.
Ephemerides.SPK_SEGMENTLIST_MAPPING — ConstantSPK_SEGMENT_MAPPINGA dictionary mapping SPK segment types to the field index of the SPKSegmentList.
Ephemerides.TCB_SEGMENTS — ConstantTCB_SEGMENTSList of the SPK/PCK segment types for which the time argument is expressed in the TCB scale.
Ephemerides.TDB_SEGMENTS — ConstantTDB_SEGMENTSList of the SPK/PCK segment types for which the time argument is expressed in the TDB scale.
SPK Segment Types
Abstract SPK Types
Ephemerides.AbstractSPKHeader — TypeAbstractSPKHeaderAbstract type for all SPK segment type headers.
Ephemerides.AbstractSPKCache — TypeAbstractSPKCacheAbstract type for all SPK segment type caches.
Ephemerides.AbstractSPKSegment — TypeAbstractSPKSegmentAbstract type for all SPK segment types.
Ephemerides.cache — Functioncache(spk::AbstractSPKSegment)Return the segment cache data.
Ephemerides.spk_field — Functionspk_field(spk::AbstractSPKSegment)Return the field number in the Ephemerides.SPKSegmentList associated to the given SPK segment type.
SPK Type 1 and 21
Ephemerides.SPKSegmentHeader1 — TypeSPKSegmentHeader1 <: AbstractSPKHeaderHeader instance for SPK segments of type 1 and 21.
Fields
n–Intnumber of records in the segmentndirs–Intnumber of directory epochsepochs– Storage for directory epochs or epochs (when ndirs = 0)iaa-Intinitial segment file addressetid–Intinitial address for the epoch table (after all the MDA records)recsize-IntNumber of double numbers stored in each MDA recordmaxdim-IntMDA dimension (fixed to 15 for type 1)
Ephemerides.SPKSegmentCache1 — TypeSPKSegmentCache1 <: AbstractSPKCacheCache 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.
Ephemerides.SPKSegmentType1 — TypeSPKSegmentType1 <: AbstractSPKSegmentSegment 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 headercache– Segment cache
References
Ephemerides.compute_mda_position — Functioncompute_mda_position(cache::SPKSegmentCache1, Δ::Number)Ephemerides.compute_mda_velocity — Functioncompute_mda_velocity(cache::SPKSegmentCache1, Δ::Number)SPK Type 2 and 3
Ephemerides.SPKSegmentHeader2 — TypeSPKSegmentHeader2 <: AbstractSPKHeaderHeader instance for SPK segments of type 2 and 3.
Fields
tstart–Float64initial epoch of the first record, in seconds since J2000tlen–Float64interval length covered by each record, in secondsorder–Intpolynomial orderN–Intnumber of coefficients in each windown–Intnumber of records in the segmentrecsize–Intbyte size of each logical recordncomp–Intnumber of vector componentsiaa–Intinitial segment file addresstype–IntSPK segment type, either 2 or 2
Ephemerides.SPKSegmentCache2 — TypeSPKSegmentCache2 <: AbstractSPKCacheCache 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 factorbuff– Stores the buffers for the Chebyshev polynomialsid– Index of the currently loaded logical record
Ephemerides.SPKSegmentType2 — TypeSPKSegmentType2 <: AbstractSPKSegmentSegment 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 headercache– Segment cache
References
SPK Type 5
Ephemerides.SPKSegmentHeader5 — TypeSPKSegmentHeader5 <: AbstractSPKHeaderHeader instance for SPK segments of type 5.
Fields
GM–Float64Gravitational constantn–Intnumber of statesndirs–Intnumber of epoch directoriesetid–Intinitial address for the epoch table (after all the state data)epochs– Storage for directory epochs or epochs (when ndirs = 0)iaa-Intinitial segment file address
Ephemerides.SPKSegmentCache5 — TypeSPKSegmentCache5 <: AbstractSPKCacheCache 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.
Ephemerides.SPKSegmentType5 — TypeSPKSegmentType5 <: AbstractSPKSegmentSegment instance for SPK segments of type 5.
Fields
head– Segment headercache– Segment cache
References
SPK Type 8 and 12
Ephemerides.SPKSegmentHeader8 — TypeSPKSegmentHeader8 <: AbstractSPKHeaderHeader instance for SPK segments of type 8 and 12.
Fields
tstart–Float64segment starting epoch, in TDB seconds since J2000tlen–Float64interval length, in secondsorder–Intinterpolating polynomial degreeN–Intgroup sizen–Intnumber of states in the segmentiaa-Intinitial segment file addressiseven–Booltrue for even group sizetype–IntSPK type (either 8 or 12)
Ephemerides.SPKSegmentCache8 — TypeSPKSegmentCache8 <: AbstractSPKCacheCache 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.
Ephemerides.SPKSegmentType8 — TypeSPKSegmentType8 <: AbstractSPKSegmentSegment instance for SPK segments of type 8 and 12.
Fields
head– Segment headercache– Segment cache
References
SPK Type 9 and 13
Ephemerides.SPKSegmentHeader9 — TypeSPKSegmentHeader9 <: AbstractSPKHeaderHeader instance for SPK segments of type 9 and 13.
Fields
n–Intnumber of states in the segmentndirs–Intnumber of epoch directoriesepochs– Storage for directory epochs or epochs (when ndirs = 0)iaa-Intinitial segment file addressetid–Intinitial address for the epoch table (after all the state data)order–Intinterpolating polynomial degreeN–Intgroup sizeiseven–Booltrue for even group sizetype–IntSPK type (either 9 or 13)
Ephemerides.SPKSegmentCache9 — TypeSPKSegmentCache9 <: AbstractSPKCacheCache 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.
Ephemerides.SPKSegmentType9 — TypeSPKSegmentType9 <: AbstractSPKSegmentSegment instance for SPK segments of type 9 and 13.
Fields
head– Segment headercache– Segment cache
References
SPK Type 14
Ephemerides.SPKSegmentHeader14 — TypeSPKSegmentHeader14 <: AbstractSPKHeaderHeader instance for SPK segments of type 14.
Fields
order–Intinterpolating polynomial degreen–Intnumber of packets in the segmentndirs–Intnumber of epoch directoriesepochs– Storage for directory epochs or epochs (when ndirs = 0)etid–Intinitial address for the epoch table (after all the state data)ptid–Intinitial address for the packet table (after the constants)pktsize–Intsize of each data packet excluding the packet information area.pktoff–Intoffset of the packet data from the packet startncomp–Intnumber of states coefficients (= 6 for SPK 14)N–Intnumber of polynomial coefficients
Ephemerides.SPKSegmentType14 — TypeSPKSegmentType14 <: AbstractSPKSegmentSegment instance for SPK segments of type 14.
Fields
head– Segment headercache– Segment cache
References
SPK segments of type 14 have the same cache structure of SPK type 2 and 3.
SPK Type 15
Ephemerides.SPKSegmentHeader15 — TypeSPKSegmentHeader15 <: AbstractSPKHeaderHeader instance for SPK segments of type 15.
Fields
epoch– Epoch of periapsistp– Trajectory pole, i.e., vector parallel to the angular momentum of the orbitpv– Central body north pole unit vectorpa– Periapsis unit vector at epochp– Semi-latus rectumecc– Eccentricityj2f– J2 processing flagvj2– J2 validation flag, true if the orbit shape is compliant with J2 pertubations.GM– Central body gravitational constant (km³/s²)J2– Central body J2R– Central body radius (km)dmdt– Mean anomaly rate of change (rad/s)kn– Gain factor for the regression of the nodeskp– Gain factor for the precession of the pericenter
Ephemerides.SPKSegmentType15 — TypeSPKSegmentType15 <: AbstractSPKSegmentSegment instance for SPK segments of type 15.
Fields
head– Segment headercache– Segment cache
References
The cache of SPK Type 15 segments is made of Ephemerides.TwoBodyUniversalCache objects.
SPK Type 17
Ephemerides.SPKSegmentHeader17 — TypeSPKSegmentHeader17 <: AbstractSPKHeaderHeader instance for SPK segments of type 17.
Fields
epoch: epoch of periapsis (s)sma: semi-major axis (km)h: H term of the equinoctial elementsk: K term of the equinoctial elementslon: mean longitude at epoch (rad)p: P term of the equinoctial elementsq: Q term of the equinoctial elementsdlpdt: 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
Ephemerides.SPKSegmentType17 — TypeSPKSegmentType17 <: AbstractSPKSegmentSegment instance for SPK segments of type 17.
Fields
head– Segment header
SPK segments of type 17 do not require a cache because they do not extract any additional coefficients at runtime.
References
SPK segments of type 17 do not require a cache structure.
SPK Type 18 and 19
Ephemerides.SPKSegmentHeader18 — TypeSPKSegmentHeader18 <: AbstractSPKHeaderHeader instance for SPK segments of type 18.
Fields
n–Intnumber of states in the segmentndirs–Intnumber of epoch directoriesepochs– Storage for directory epochs or epochs (when ndirs = 0)iaa-Intinitial segment file addressetid–Intinitial address for the epoch table (after all the state data)order–Intinterpolating polynomial degreeN–Intgroup sizesubtype–Inttype 18 subtype, either 0 (Hermite) or 1 (Lagrange)packetsize–Intpacket size for each point, either 12 (Hermite) or 6 (Lagrange)
Ephemerides.SPKSegmentCache18 — TypeSPKSegmentCache18 <: AbstractSPKCacheCache 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.
Ephemerides.reset_indexes! — Functionreset_indexes!(cache::SPKSegmentCache18)Reset the cache indexes to force the coefficients reload.
Ephemerides.update_header! — Functionupdate_header!(head::SPKSegmentHeader18, daf::DAF, iaa, faa, type)Update the header of type 18 segments.
Ephemerides.SPKSegmentHeader19 — TypeSPKSegmentHeader19 <: AbstractSPKHeaderHeader instance for SPK segments of type 19.
Fields
n–Intnumber of states in the segment.ndirs–Intnumber of epoch directories.times– Storage for interval directories or start times (when ndirs = 0).iaa-Intinitial segment file address.etid–Intbyte address for the interval table (after all the minisegment data).ptid–Intbyte for the pointer table.usefirst–Boolboundary flag, true if the preceding segment should be used.type–Inteither type 18 or 19.
Ephemerides.SPKSegmentCache19 — TypeSPKSegmentCache19 <: AbstractSPKCacheCache 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.
Ephemerides.SPKSegmentType19 — TypeSPKSegmentType19 <: AbstractSPKSegmentSegment instance for SPK segments of type 18 and 19. Type 18 segments are treated as special cases of a type 19 with a single mini-segment.
Fields
head– Segment headercache– Segment cache
References
Ephemerides.find_minisegment — Functionfind_minirecord(daf::DAF, head::SPKSegmentHeader19, time::Number)Ephemerides.load_minisegment! — Functionload_minisegment!(daf::DAF, head::SPKSegmentHeader19, cache::SPKSegmentCache19, index::Int)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.SPKSegmentHeader20 — TypeSPKSegmentHeader20 <: AbstractSPKHeaderHeader instance for SPK segments of type 20.
Fields
dscale–Float64length conversion factortscale–Float64time conversion factortstart–Float64initial epoch of the first record, in seconds since J2000tlen–Float64interval length covered by each record, in secondsrecsize–Intbyte size of each logical recordorder–Intpolynomial orderN–Intnumber of coefficients in each windown–Intnumber of records in the segmentiaa–Intinitial segment file address
Ephemerides.SPKSegmentCache20 — TypeSPKSegmentCache20 <: AbstractSPKCacheCache instance for SPK segments of type 20.
Fields
id– Index of the currently loaded logical recordp– Stores the record position constantsA– Chebyshev's polynomial coefficients, with size (ncomp, order)buff– Stores the buffers for the Chebyshev polynomials
Ephemerides.SPKSegmentType20 — TypeSPKSegmentType2 <: AbstractSPKSegmentSegment 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 headercache– Segment cache
References
SPK Utility Functions
Ephemerides.normalise_time — Functionnormalise_time(cache::SPKSegmentCache2, time::Number)Transform time in an interval between [-1, 1] for compliance with Chebyshev polynomials.
normalise_time(head::SPKSegmentHeader8, time::Number, index::Int)Returned a normalised time that starts at 1 at the beginning of the interval.
normalise_time(head::SPKSegmentHeader20, time::Number, index::Int)Ephemerides.find_logical_record — Functionfind_logical_record(daf::DAF, head::SPKSegmentHeader1, time::Number)find_logical_record(head::SPKSegmentHeader2, time::Number)find_logical_record(daf::DAF, head::SPKSegmentHeader5, time::Number)find_logical_record(head::SPKSegmentHeader8, time::Number)find_logical_record(daf::DAF, head::SPKSegmentHeader1, time::Number)find_logical_record(daf::DAF, head::SPKSegmentHeader14, time::Number)find_logical_record(daf::DAF, head::SPKSegmentHeader18, time::Number)find_logical_record(head::SPKSegmentHeader20, time::Number)Ephemerides.get_coefficients! — Functionget_coefficients!(daf::DAF, head::SPKSegmentHeader1, cache::SPKSegmentCache1, index::Int)get_coefficients!(daf::DAF, head::SPKSegmentHeader2, cache::SPKSegmentCache2, index::Int)get_coefficients!(daf::DAF, head::SPKSegmentHeader5, cache::SPKSegmentCache5, index::Int)get_coefficients!(daf::DAF, head::SPKSegmentHeader8, cache::SPKSegmentCache8, index::Int)get_coefficients!(daf::DAF, head::SPKSegmentHeader9, cache::SPKSegmentCache9, index::Int)get_coefficients!(daf::DAF, head::SPKSegmentHeader14, cache::SPKSegmentCache14, index::Int)get_coefficients!(daf::DAF, head::SPKSegmentHeader18, cache::SPKSegmentCache18, first::Int, last::Int)get_coefficients!(daf::DAF, head::SPKSegmentHeader20, cache::SPKSegmentCache20, index::Int)Interpolating Functions
Caches
Ephemerides.InterpCache — TypeInterpCache{T}Ephemerides.get_buffer — Functionget_buffer(c::InterpCache, idx::int, x::Number)Return the idx-th buffer from the corresponding DiffCache depending on the type of x.
Chebyshev Polynomials
Ephemerides.chebyshev — Functionchebyshev(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.
x is a re-work of the actual ascissa value that lies between [-1, 1]
See Also
See also ∂chebyshev, ∂²chebyshev and ∂³chebyshev
Ephemerides.∂chebyshev — Function∂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.
x is a re-work of the actual ascissa value that lies between [-1, 1]
See Also
See also chebyshev, ∂²chebyshev and ∂³chebyshev
Ephemerides.∂²chebyshev — Function∂²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.
x is a re-work of the actual ascissa value that lies between [-1, 1]
See Also
See also chebyshev, ∂chebyshev and ∂³chebyshev
Ephemerides.∂³chebyshev — Function∂³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.
x is a re-work of the actual ascissa value that lies between [-1, 1]
See Also
See also chebyshev, ∂chebyshev and ∂²chebyshev
Ephemerides.∫chebyshev — Function∫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).
x is a re-work of the actual ascissa value that lies between [-1, 1]
∫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.
Lagrange Polynomials
Ephemerides.lagrange — Functionlagrange(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.
x is a re-work of the actual ascissa value that starts at 1
See Also
See also ∂lagrange and ∂²lagrange.
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.
Ephemerides.∂lagrange — Function∂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
x is a re-work of the actual ascissa value that starts at 1
See Also
See also lagrange and ∂²lagrange
∂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.
Ephemerides.∂²lagrange — Function∂²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
x is a re-work of the actual ascissa value that starts at 1
See Also
∂²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
Hermite Polynomials
Ephemerides.hermite — Functionhermite(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.
x is a re-work of the actual ascissa value that starts at 1
See Also
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.
Ephemerides.∂hermite — Function∂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.
x is a re-work of the actual ascissa value that starts at 1
See Also
∂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.
Ephemerides.∂²hermite — Function∂²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.
x is a re-work of the actual ascissa value that starts at 1
See Also
∂²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.
Ephemerides.∂³hermite — Function∂³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.
x is a re-work of the actual ascissa value that starts at 1
See Also
∂³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.
Introspection
Ephemerides.AbstractEphemRecord — TypeAbstractEphemRecordAbstract type for ephemeris segment records.
Ephemerides.initial_times — Functioninitial_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
Ephemerides.final_times — Functionfinal_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
Ephemerides.analyse_timespan — Functionanalyse_timespan(records)Analyse a set of AbstractEphemRecord, returning the minimum and maximum covered times, in TDB seconds since J2000, together with a continuity parameter.
References
- CALCEPH C++ library
See Also
See also ephem_spk_timespan and ephem_pck_timespan.
TwoBody Routines
Ephemerides.TwoBodyUniversalCache — TypeTwoBodyUniversalCacheA 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.
Ephemerides.update_cache! — Functionupdate_cache!(c::TwoBodyUniversalCache)Update the precomputed values in the cache using the position and velocity in c.
Ephemerides.propagate_twobody — Functionpropagate_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.
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
Ephemerides.stumpff — Functionstumpff(x::Number, p::AbstractVector)Compute Stumpff's functions from C₀ up to C₃ at x.
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