Multiple file access functions¶
The following group of functions should be the preferred method to access to the library. They allow to access to multiple ephemeris files at the same time, even by multiple threads.
When an error occurs, these functions execute error handlers according to the behavior defined by the function f90calceph_seterrorhandler()
.
Time notes¶
The functions f90calceph_compute()
, f90calceph_compute_unit()
, f90calceph_compute_order()
, f90calceph_orient_unit()
, ... only accept a date expressed in the same timescale as the ephemeris files, which can be retrieved using the function f90calceph_gettimescale()
. Ephemeris files are generally expressed using the timescale TDB.
If a date, expressed in the TT (Terrestrial Time) timescale, is supplied to them, these functions will return an erroneous position of the order of several tens of meters for the planets.
If a date, expressed in the Coordinated Universal Time (UTC), is supplied to them, these functions will return a very large erroneous position over several thousand kilometers for the planets.
Thread notes¶
If the standard I/O functions such as fread are not reentrant then the CALCEPH I/O functions using them will not be reentrant either.
It's safe for two threads to call the functions with the same handle of ephemeris object if and only if the function f90calceph_isthreadsafe()
returns a non-zero value. A previous call to the function f90calceph_prefetch()
is required for the function f90calceph_isthreadsafe()
to return a non-zero value.
It's safe for two threads to access simultaneously to the same ephemeris file with two different objects. In this case, each thread must open the same file.
Usage¶
The following examples, that can be found in the directory examples of the library sources, show the typical usage of this group of functions.
The example in Fortran 77/90/95 language is f77multiple.f
.
Functions¶
f90calceph_open¶
-
function
f90calceph_open
(eph, filename)¶ Parameters: - filename [CHARACTER(len=*), intent(in)] :: pathname of the file.
- eph [INTEGER(8), intent(out)] :: ephemeris descriptor
Return: f90calceph_open [INTEGER] :: 0 if an error occurs, otherwise non-zero value.
This function opens the file whose pathname is the string pointed to by filename, reads the two header blocks of this file and returns an ephemeris descriptor associated to it. This file must be compliant to the format specified by the 'original JPL binary' , 'INPOP 2.0 binary' or 'SPICE' ephemeris file. At the moment, supported SPICE files are the following :
- text Planetary Constants Kernel (KPL/PCK) files
- binary PCK (DAF/PCK) files.
- binary SPK (DAF/SPK) files containing segments of type 1, 2, 3, 5, 8, 9, 12, 13, 14, 17, 18, 20, 21, 102, 103 and 120.
- meta kernel (KPL/MK) files.
- frame kernel (KPL/FK) files. Only a basic support is provided.
Just after the call of f90calceph_open()
, the function f90calceph_prefetch()
should be called to accelerate future computations.
The function f90calceph_close()
must be called to free allocated memory by this function.
The following example opens the ephemeris file example1.dat
include 'f90calceph.h'
integer res
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
! ... computation ...
endif
call f90calceph_close(peph)
f90calceph_open_array¶
-
function
f90calceph_open_array
(eph, n, array_filename, len_filename)¶ Parameters: - eph [INTEGER(8), intent(out)] :: ephemeris descriptor
- n [INTEGER, intent(in)] :: number of files.
- array_filename [CHARACTER(len=*), dimension(*), intent(in)] :: array of pathname of the files
- len_filename [INTEGER, intent(in)] :: number of characters of each file's name.
Return: f90calceph_open_array [INTEGER] :: 0 if an error occurs, otherwise non-zero value.
This function opens n files whose pathnames are the string pointed to by array_filename, reads the header blocks of these files and returns an ephemeris descriptor associated to them.
These files must have the same type (e.g., all files are SPICE files or original JPL files). This file must be compliant to the format specified by the 'original JPL binary' , 'INPOP 2.0 or 3.0 binary' or 'SPICE' ephemeris file. At the moment, supported SPICE files are the following :
- text Planetary Constants Kernel (KPL/PCK) files
- binary PCK (DAF/PCK) files.
- binary SPK (DAF/SPK) files containing segments of type 1, 2, 3, 5, 8, 9, 12, 13, 14, 17, 18, 20, 21, 102, 103 and 120.
- meta kernel (KPL/MK) files.
- frame kernel (KPL/FK) files. Only a basic support is provided.
Just after the call of f90calceph_open_array()
, the function f90calceph_prefetch()
should be called to accelerate future computations.
The function f90calceph_close()
must be called to free allocated memory by this function.
The following example opens the ephemeris file example1.bsp and example1.tpc
integer*8 peph
integer res
character*80 filear(2)
filear(1) = "example1.bsp"
filear(2) = "example1.tpc"
res = f90calceph_open_array(peph, 2, filear, 80)
if (res.eq.1) then
res = f90calceph_prefetch(peph)
! ... computation ...
call calceph_close(peph)
endif
f90calceph_prefetch¶
-
function
f90calceph_prefetch
(eph)¶ Parameters: eph [INTEGER(8), intent(in)] :: ephemeris descriptor Return: f90calceph_prefetch [INTEGER] :: 0 if an error occurs, otherwise non-zero value.
This function prefetches to the main memory all files associated to the ephemeris descriptor eph.
This prefetching operation will accelerate the further computations performed with f90calceph_compute()
, f90calceph_compute_unit()
, f90calceph_compute_order()
, f90calceph_orient_unit()
, ... .
It requires that the file is smaller than the main memory. If multiple threads (e.g. threads of openMP or Posix Pthreads) prefetch the data for the same ephemeris file, the used memory will remain the same as if the prefetch operation was done by a single thread if and if the endianess of the file is the same as the computer and if the operating system, such as Linux, MacOS X other unix, supports the function mmap.
f90calceph_isthreadsafe¶
-
function
f90calceph_isthreadsafe
(eph)¶ Parameters: eph [INTEGER(8), intent(in)] :: ephemeris descriptor Return: f90calceph_isthreadsafe [INTEGER] :: returns 1 if multiple threads can access the same ephemeris ephemeris descriptor, otherwise 0.
This function returns 1 if multiple threads can access the same ephemeris ephemeris descriptor, otherwise 0.
A previous call to the function f90calceph_prefetch()
is required, and the library should be compiled with --enable-thread=yes on Unix-like operating system, for the function f90calceph_isthreadsafe()
to return a non-zero value. If the file is not encoded with the same endian as the current hardware, then function may return 0.
If this function returns 1, severals threads may use the same ephemeris descriptor for the computational functions f90calceph_compute()
, .... It allows to use the same object for parallel loops.
f90calceph_compute¶
-
function
f90calceph_compute
(eph, JD0, time, target, center, PV)¶ Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- JD0 [REAL(8), intent(in)] :: Integer part of the Julian date (TDB or TCB)
- time [REAL(8), intent(in)] :: Fraction part of the Julian date (TDB or TCB)
- target [INTEGER, intent(in)] :: The body or reference point whose coordinates are required (see the list, below).
- center [INTEGER, intent(in)] :: The origin of the coordinate system (see the list, below). If target is 14, 15, 16 or 17 (nutation, libration, TT-TDB or TCG-TCB), center must be 0.
- PV (6) [REAL(8), intent(out) :: Depending on the target value, an array to receive the cartesian position (x,y,z) and the velocity (xdot, ydot, zdot), or a time scale transformation value, or the angles of the librations of the Moon and their derivatives, or the nutation angles and their derivatives.
Return: f90calceph_compute [INTEGER] :: 0 if an error occurs, otherwise non-zero value.
This function reads, if needed, in the ephemeris file associated to eph and interpolates a single object, usually the position and velocity of one body (target) relative to another (center) for the time JD0+time and stores the results to PV. The ephemeris file associated to eph must have been previously opened with the function f90calceph_open()
.
The returned array PV has the following properties
- If the target is TT-TDB, only the first element of this array will get the result. The time scale transformation TT-TDB is expressed in seconds.
- If the target is TCG-TCB, only the first element of this array will get the result. The time scale transformation TCG-TCB is expressed in seconds.
- If the target is Librations, the array contains the angles of the librations of the Moon and their derivatives. The angles of the librations of the Moon are expressed in radians and their derivatives are expressed in radians per day.
- If the target is Nutations, the array contains the nutation angles and their derivatives. The nutation angles are expressed in radians and their derivatives are expressed in radians per day.
- Otherwise the returned values is the cartesian position (x,y,z), expressed in Astronomical Unit (au), and the velocity (xdot, ydot, zdot), expressed in Astronomical Unit per day (au/day).
The date (JD0, time) should be expressed in the same timescale as the ephemeris files, which can be retrieved using the function f90calceph_gettimescale()
. To get the best numerical precision for the interpolation, the time is splitted in two floating-point numbers. The argument JD0 should be an integer and time should be a fraction of the day. But you may call this function with time=0 and JD0, the desired time, if you don't take care about numerical precision.
Warning
If a date, expressed in the Coordinated Universal Time (UTC), is supplied to this function, a very large erroneous position will be returned.
The possible values for target and center are :
value | meaning |
---|---|
1 | Mercury Barycenter |
2 | Venus Barycenter |
3 | Earth |
4 | Mars Barycenter |
5 | Jupiter Barycenter |
6 | Saturn Barycenter |
7 | Uranus Barycenter |
8 | Neptune Barycenter |
9 | Pluto Barycenter |
10 | Moon |
11 | Sun |
12 | Solar Sytem barycenter |
13 | Earth-moon barycenter |
14 | Nutation angles |
15 | Librations |
16 | TT-TDB |
17 | TCG-TCB |
asteroid number + CALCEPH_ASTEROID | asteroid |
These accepted values by this function are the same as the value for the JPL function PLEPH, except for the values TT-TDB, TCG-TCB and asteroids.
For example, the value "CALCEPH_ASTEROID+4" for target or center specifies the asteroid Vesta.
The following example prints the heliocentric coordinates of Mars at time=2442457.5 and at 2442457.9
integer res
double precision jd0
double precision dt1, dt2
double precision PV(6)
jd0 = 2442457
dt1 = 0.5D0
dt2 = 0.9D0
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
res = f90calceph_compute(peph,jd0, dt1, 4, 11, PV)
write(*,*) PV
res = f90calceph_compute(peph,jd0, dt2, 4, 11, PV)
write(*,*) PV
call f90calceph_close(peph)
endif
f90calceph_compute_unit¶
-
function
f90calceph_compute_unit
(eph, JD0, time, target, center, unit, PV)¶ Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- JD0 [REAL(8), intent(in)] :: Integer part of the Julian date (TDB or TCB)
- time [REAL(8), intent(in)] :: Fraction part of the Julian date (TDB or TCB)
- target [INTEGER, intent(in)] :: The body or reference point whose coordinates are required. The numbering system depends on the parameter unit.
- center [INTEGER, intent(in)] :: The origin of the coordinate system. The numbering system depends on the parameter unit.
- unit [INTEGER, intent(in)] :: The units of PV.This integer is a sum of some unit constants (CALCEPH_UNIT_???) and/or the constant
CALCEPH_USE_NAIFID
.If the unit containsCALCEPH_USE_NAIFID
, the NAIF identification numbering system is used for the target and the center (NAIF identification numbers for the list).If the unit doesnot containCALCEPH_USE_NAIFID
, the old number system is used for the target and the center (see the list in the functionf90calceph_compute()
). - PV (6) [REAL(8), intent(out) :: Depending on the target value, an array to receive the cartesian position (x,y,z) and the velocity (xdot, ydot, zdot), or a time scale transformation value, or the angles of the librations of the Moon and their derivatives, or the nutation angles and their derivatives.
Return: f90calceph_compute_unit [INTEGER] :: 0 if an error occurs, otherwise non-zero value.
This function is similar to the function f90calceph_compute()
, except that the units of the output are specified.
This function reads, if needed, in the ephemeris file associated to eph and interpolates a single object, usually the position and velocity of one body (target) relative to another (center) for the time JD0+time and stores the results to PV. The ephemeris file associated to eph must have been previously opened with the function f90calceph_open()
.
The output values are expressed in the units specified by unit.
This function checks the units if invalid combinations of units are given to the function.
The date (JD0, time) should be expressed in the same timescale as the ephemeris files, which can be retrieved using the function f90calceph_gettimescale()
.
Warning
If a date, expressed in the Coordinated Universal Time (UTC), is supplied to this function, a very large erroneous position will be returned.
The returned array PV has the following properties
- If the target is the time scale transformation TT-TDB, only the first element of this array will get the result.
- If the target is the time scale transformation TCG-TCB, only the first element of this array will get the result.
- If the target is Librations, the array contains the angles of the librations of the Moon and their derivatives.
- If the target is Nutations, the array contains the nutation angles and their derivatives.
- Otherwise the returned value is the cartesian position (x,y,z) and the velocity (xdot, ydot, zdot).
The values stored in the array PV are expressed in the following units
- The position and velocity are expressed in Astronomical Unit (au) if unit contains
CALCEPH_UNIT_AU
.- The position and velocity are expressed in kilometers if unit contains
CALCEPH_UNIT_KM
.- The velocity, TT-TDB, TCG-TCB, the derivatives of the angles of the nutation, or the derivatives of the librations of the Moon or are expressed in days if unit contains
CALCEPH_UNIT_DAY
.- The velocity, TT-TDB, TCG-TCB, the derivatives of the angles of the nutation, or the derivatives of the librations of the Moon are expressed in seconds if unit contains
CALCEPH_UNIT_SEC
.- The angles of the librations of the Moon or the nutation angles are expressed in radians if unit contains
CALCEPH_UNIT_RAD
.
For example, to get the position and velocities expressed in kilometers and kilometers/seconds, the unit must be set to CALCEPH_UNIT_KM
+ CALCEPH_UNIT_SEC
.
The following example prints the heliocentric coordinates of Mars at time=2442457.5
integer*8 peph
integer res
double precision jd0
double precision dt1
double precision PV(6)
jd0 = 2442457
dt1 = 0.5D0
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
! the heliocentric coordinates of Mars in km and km/s
res = f90calceph_compute_unit(peph,jd0, dt1, 4, 11,
& CALCEPH_USE_NAIFID+CALCEPH_UNIT_KM+CALCEPH_UNIT_SEC,
& PV)
write(*,*) PV
! compute same quantity as the previous call using NAIF ID
res = f90calceph_compute_unit(peph,jd0, dt1,
& NAIFID_MARS_BARYCENTER, NAIFID_SUN,
& CALCEPH_USE_NAIFID+CALCEPH_UNIT_KM+CALCEPH_UNIT_SEC,
& PV)
write(*,*) PV
call f90calceph_close(peph)
endif
f90calceph_orient_unit¶
-
function
f90calceph_orient_unit
(eph, JD0, time, target, unit, PV)¶ Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- JD0 [REAL(8), intent(in)] :: Integer part of the Julian date (TDB or TCB)
- time [REAL(8), intent(in)] :: Fraction part of the Julian date (TDB or TCB)
- target [INTEGER, intent(in)] :: The body whose orientations are requested. The numbering system depends on the parameter unit.
- unit [INTEGER, intent(in)] :: The units of PV.This integer is a sum of some unit constants (CALCEPH_UNIT_???) and/or the constant
CALCEPH_USE_NAIFID
.If the unit containsCALCEPH_USE_NAIFID
, the NAIF identification numbering system is used for the target (NAIF identification numbers for the list).If the unit does not containCALCEPH_USE_NAIFID
, the old number system is used for the target (see the list in the functionf90calceph_compute()
). - PV (6) [REAL(8), intent(out) :: An array to receive the euler angles, or nutation angles, and their derivatives for the orientation of the body.
Return: f90calceph_orient_unit [INTEGER] :: 0 if an error occurs, otherwise non-zero value.
This function reads, if needed, in the ephemeris file associated to eph and interpolates the orientation of a single body (target) for the time JD0+time and stores the results to PV. The ephemeris file associated to eph must have been previously opened with the function f90calceph_open()
.
The output values are expressed in the units specified by unit.
The date (JD0, time) should be expressed in the same timescale as the ephemeris files, which can be retrieved using the function f90calceph_gettimescale()
.
This function checks the units if invalid combinations of units are given to the function.
The returned array PV has the following properties
- If unit contains
CALCEPH_OUTPUT_NUTATIONANGLES
, the array contains the nutation angles and their derivatives for the orientation of the body. At the present moment, only the nutation for the earth are supported in the original DE files.- If unit contains
CALCEPH_OUTPUT_EULERANGLES
, or doesnot containCALCEPH_OUTPUT_NUTATIONANGLES
, the array contains the euler angles and their derivatives for the orientation of the body.
The values stored in the array PV are expressed in the following units
- The derivatives of the angles are expressed in days if unit contains
CALCEPH_UNIT_DAY
.- The derivatives of the angles are expressed in seconds if unit contains
CALCEPH_UNIT_SEC
.- The angles and their derivatives are expressed in radians if unit contains
CALCEPH_UNIT_RAD
.
For example, to get the nutation angles of the Earth and their derivatives expressed in radian and radian/seconds using the NAIF identification numbering system, the target must be set to NAIFID_EARTH and the unit must be set to CALCEPH_OUTPUT_NUTATIONANGLES
+ CALCEPH_UNIT_RAD
+ CALCEPH_UNIT_SEC
.
The following example prints the angles of libration of the Moon at time=2442457.5
integer*8 peph
integer res
double precision jd0
double precision dt1
double precision PV(6)
jd0 = 2442457
dt1 = 0.5D0
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
res = f90calceph_orient_unit(peph,jd0, dt1, NAIFID_MOON,
& CALCEPH_USE_NAIFID+CALCEPH_UNIT_RAD+CALCEPH_UNIT_SEC,
& PV)
write(*,*) PV
call f90calceph_close(peph)
endif
f90calceph_rotangmom_unit¶
-
function
f90calceph_rotangmom_unit
(eph, JD0, time, target, unit, PV)¶ Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- JD0 [REAL(8), intent(in)] :: Integer part of the Julian date (TDB or TCB)
- time [REAL(8), intent(in)] :: Fraction part of the Julian date (TDB or TCB)
- target [INTEGER, intent(in)] :: The body whose orientations are requested. The numbering system depends on the parameter unit.
- unit [INTEGER, intent(in)] :: The units of PV.This integer is a sum of some unit constants (CALCEPH_UNIT_???) and/or the constant
CALCEPH_USE_NAIFID
.If the unit containsCALCEPH_USE_NAIFID
, the NAIF identification numbering system is used for the target (NAIF identification numbers for the list).If the unit does not containCALCEPH_USE_NAIFID
, the old number system is used for the target (see the list in the functionf90calceph_compute()
). - PV (6) [REAL(8), intent(out) :: An array to receive the angular momentum due to its rotation, divided by the product of the mass and of the square of the radius, and the derivatives, of the body.
Return: f90calceph_rotangmom_unit [INTEGER] :: 0 if an error occurs, otherwise non-zero value.
This function reads, if needed, in the ephemeris file associated to eph and interpolates the angular momentum vector due to the rotation of the body, divided by the product of the mass and of the square of the radius
, of a single body (target) for the time JD0+time and stores the results to PV. The ephemeris file associated to eph must have been previously opened with the function
f90calceph_open()
. The angular momentum , due to the rotation of the body, is defined as the product of the inertia matrix
by the angular velocity vector
. So the returned value is
The date (JD0, time) should be expressed in the same timescale as the ephemeris files, which can be retrieved using the function f90calceph_gettimescale()
.
The output values are expressed in the units specified by unit.
This function checks the units if invalid combinations of units are given to the function.
The values stored in the array PV are expressed in the following units
- The angular momentum and its derivative are expressed in days if unit contains
CALCEPH_UNIT_DAY
.- The angular momentum and its derivative are expressed in seconds if unit contains
CALCEPH_UNIT_SEC
.
The following example prints the angular momentum, due to its rotation, for the Earth at time=2451419.5
integer*8 peph
integer res
double precision jd0
double precision dt1
double precision G(6)
jd0 = 2451419
dt1 = 0.5D0
res = f90calceph_open(peph, "example2_rotangmom.dat")
if (res.eq.1) then
res = f90calceph_rotangmom_unit(peph,jd0, dt1, NAIFID_EARTH,
& CALCEPH_USE_NAIFID+CALCEPH_UNIT_SEC,
& G)
write(*,*) G
call f90calceph_close(peph)
endif
f90calceph_compute_order¶
-
function
f90calceph_compute_order
(eph, JD0, time, target, center, unit, order, PVAJ)¶ Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- JD0 [REAL(8), intent(in)] :: Integer part of the Julian date (TDB or TCB)
- time [REAL(8), intent(in)] :: Fraction part of the Julian date (TDB or TCB)
- target [INTEGER, intent(in)] :: The body or reference point whose coordinates are required. The numbering system depends on the parameter unit.
- center [INTEGER, intent(in)] :: The origin of the coordinate system. The numbering system depends on the parameter unit.
- unit [INTEGER, intent(in)] :: The units of PVAJ.This integer is a sum of some unit constants (CALCEPH_UNIT_???) and/or the constant
CALCEPH_USE_NAIFID
.If the unit containsCALCEPH_USE_NAIFID
, the NAIF identification numbering system is used for the target and the center (NAIF identification numbers for the list).If the unit doesnot containCALCEPH_USE_NAIFID
, the old number system is used for the target and the center (see the list in the functionf90calceph_compute()
). - order [INTEGER, intent(in)] ::
The order of derivatives
- = 0 , only the position is computed. The first three numbers of PVAJ are valid for the results.
- = 1 , only the position and velocity are computed. The first six numbers of PVAJ are valid for the results.
- = 2 , only the position, velocity and acceleration are computed. The first nine numbers of PVAJ are valid for the results.
- = 3 , the position, velocity and acceleration and jerk are computed. The first twelve numbers of PVAJ are valid for the results.
If order equals to 1, the behavior of
f90calceph_compute_order()
is the same asf90calceph_compute_unit()
. - PVAJ (12) [REAL(8), intent(out) :: Depending on the target value, an array to receive the cartesian position (x,y,z), the velocity (xdot, ydot, zdot), the acceleration and the jerk, or a time scale transformation value, or the angles of the librations of the Moon and their successive derivatives, or the nutation angles and their successive derivatives.
Return: f90calceph_compute_order [INTEGER] :: 0 if an error occurs, otherwise non-zero value.
This function is similar to the function f90calceph_compute_unit()
, except that the order of the computed derivatives is specified.
This function reads, if needed, in the ephemeris file associated to eph and interpolates a single object, usually the position and their derivatives of one body (target) relative to another (center) for the time JD0+time and stores the results to PVAJ. The ephemeris file associated to eph must have been previously opened with the function f90calceph_open()
.
The order of the derivatives are specified by order. The output values are expressed in the units specified by unit.
The returned array PVAJ has the following properties
- If the target is the time scale transformation TT-TDB, only the first elements of each component will get the result.
- If the target is the time scale transformation TCG-TCB, only the first elements of each component will get the result.
- If the target is Librations, the array contains the angles of the librations of the Moon and their successive derivatives.
- If the target is Nutations, the array contains the nutation angles and their successive derivatives.
- Otherwise the returned value is the cartesian position (x,y,z), the velocity (xdot, ydot, zdot), the jerk and the acceleration.
The returned array PVAJ must be large enough to store the results.
- PVAJ[1:3] contain the position (x,y,z) and is always valid.
- PVAJ[4:6] contain the velocity (dx/dt,dy/dt,dz/dt) and is only valid if order is greater or equal to 1.
- PVAJ[7:9] contain the acceleration (d^2x/dt^2,d^2y/dt^2,d^2z/dt^2) and is only valid if order is greater or equal to 2.
- PVAJ[10:12] contain the jerk (d^3x/dt^3,d^3y/dt^3,d^3z/dt^3) and is only valid if order is equal to 3.
The date (JD0, time) should be expressed in the same timescale as the ephemeris files, which can be retrieved using the function f90calceph_gettimescale()
.
Warning
If a date, expressed in the Coordinated Universal Time (UTC), is supplied to this function, a very large erroneous position will be returned.
The values stored in the array PVAJ are expressed in the following units
- The position, velocity, acceleration and jerk are expressed in Astronomical Unit (au) if unit contains
CALCEPH_UNIT_AU
.- The position, velocity, acceleration and jerk are expressed in kilometers if unit contains
CALCEPH_UNIT_KM
.- The velocity, acceleration, jerk, TT-TDB, TCG-TCB or the derivatives of the angles of the librations of the Moon are expressed in days if unit contains
CALCEPH_UNIT_DAY
.- The velocity, acceleration, jerk, TT-TDB, TCG-TCB or the derivatives of the angles of the librations of the Moon are expressed in seconds if unit contains
CALCEPH_UNIT_SEC
.- The angles of the librations of the Moon are expressed in radians if unit contains
CALCEPH_UNIT_RAD
.
For example, to get the positions, velocities, accelerations and jerks expressed in kilometers and kilometers/seconds, the unit must be set to CALCEPH_UNIT_KM
+ CALCEPH_UNIT_SEC
.
This function checks the units if invalid combinations of units are given to the function.
The following example prints the heliocentric coordinates of Mars at time=2442457.5
integer*8 peph
integer res
double precision jd0
double precision dt1
double precision P(3), PVAJ(12)
jd0 = 2442457
dt1 = 0.5D0
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
! compute only the heliocentric position of Mars in km
res = f90calceph_compute_order(peph, jd0, dt1,
& NAIFID_MARS_BARYCENTER,
& NAIFID_SUN,
& CALCEPH_USE_NAIFID+CALCEPH_UNIT_KM+CALCEPH_UNIT_SEC,
& 0, P);
write(*,*) P
! compute positions, velocities, accelerations and jerks of Mars in km and seconds
res = f90calceph_compute_order(peph, jd0, dt1,
& NAIFID_MARS_BARYCENTER,
& NAIFID_SUN,
& CALCEPH_USE_NAIFID+CALCEPH_UNIT_KM+CALCEPH_UNIT_SEC,
& 3, PVAJ);
write(*,*) PVAJ
call f90calceph_close(peph)
endif
f90calceph_orient_order¶
-
function
f90calceph_orient_order
(eph, JD0, time, target, unit, order, PVAJ)¶ Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- JD0 [REAL(8), intent(in)] :: Integer part of the Julian date (TDB or TCB)
- time [REAL(8), intent(in)] :: Fraction part of the Julian date (TDB or TCB)
- target [INTEGER, intent(in)] :: The body whose orientations are requested. The numbering system depends on the parameter unit.
- unit [INTEGER, intent(in)] :: The units of PV.This integer is a sum of some unit constants (CALCEPH_UNIT_???) and/or the constant
CALCEPH_USE_NAIFID
.If the unit containsCALCEPH_USE_NAIFID
, the NAIF identification numbering system is used for the target (NAIF identification numbers for the list).If the unit does not containCALCEPH_USE_NAIFID
, the old number system is used for the target (see the list in the functionf90calceph_compute()
). - order [INTEGER, intent(in)] ::
The order of derivatives.
- = 0 , only the angles is computed. The first three numbers of PVAJ are valid for the results.
- = 1 , only the angles and the first derivative are computed. The first six numbers of PVAJ are valid for the results.
- = 2 , only the angles and the first and second derivatives are computed. The first nine numbers of PVAJ are valid for the results.
- = 3 , the angles and the first, second and third derivatives are computed. The first twelve numbers of PVAJ are valid for the results.
If order equals to 1, the behavior of
f90calceph_orient_order()
is the same asf90calceph_orient_unit()
. - PVAJ (12) [REAL(8), intent(out) :: An array to receive the euler angles, or nutation angles, and their derivatives for the orientation of the body.
Return: f90calceph_orient_order [INTEGER] :: 0 if an error occurs, otherwise non-zero value.
This function is similar to the function f90calceph_orient_unit()
, except that the order of the computed derivatives is specified.
This function reads, if needed, in the ephemeris file associated to eph and interpolates the orientation of a single body (target) for the time JD0+time and stores the results to PVAJ.
The order of the derivatives are specified by order. The ephemeris file associated to eph must have been previously opened with the function f90calceph_open()
.
The output values are expressed in the units specified by unit.
This function checks the units if invalid combinations of units are given to the function.
The returned array PVAJ has the following properties
- If unit contains
CALCEPH_OUTPUT_NUTATIONANGLES
, the array contains the nutation angles and their successive derivatives for the orientation of the body. At the present moment, only the nutation for the earth are supported in the original DE files.- If unit contains
CALCEPH_OUTPUT_EULERANGLES
, or doesnot containCALCEPH_OUTPUT_NUTATIONANGLES
, the array contains the euler angles and their successive derivatives for the orientation of the body.
The returned array PVAJ must be large enough to store the results.
- PVAJ[1:3] contain the angles and is always valid.
- PVAJ[4:6] contain the first derivative and is only valid if order is greater or equal to 1.
- PVAJ[7:9] contain the second derivative and is only valid if order is greater or equal to 2.
- PVAJ[10:12] contain the third derivative and is only valid if order is equal to 3.
The values stored in the array PVAJ are expressed in the following units
- The derivatives of the angles are expressed in days if unit contains
CALCEPH_UNIT_DAY
.- The derivatives of the angles are expressed in seconds if unit contains
CALCEPH_UNIT_SEC
.- The angles and their derivatives are expressed in radians if unit contains
CALCEPH_UNIT_RAD
.
The date (JD0, time) should be expressed in the same timescale as the ephemeris files, which can be retrieved using the function f90calceph_gettimescale()
.
The following example prints only the angles of libration of the Moon at time=2442457.5
integer*8 peph
integer res
double precision jd0
double precision dt1
double precision P(3)
jd0 = 2442457
dt1 = 0.5D0
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
res = f90calceph_orient_order(peph,jd0, dt1, NAIFID_MOON,
& CALCEPH_USE_NAIFID+CALCEPH_UNIT_RAD+CALCEPH_UNIT_SEC,
& 0, P)
write(*,*) P
call f90calceph_close(peph)
endif
f90calceph_rotangmom_order¶
-
function
f90calceph_rotangmom_order
(eph, JD0, time, target, unit, order, PVAJ)¶ Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- JD0 [REAL(8), intent(in)] :: Integer part of the Julian date (TDB or TCB)
- time [REAL(8), intent(in)] :: Fraction part of the Julian date (TDB or TCB)
- target [INTEGER, intent(in)] :: The body whose orientations are requested. The numbering system depends on the parameter unit.
- unit [INTEGER, intent(in)] :: The units of PV.This integer is a sum of some unit constants (CALCEPH_UNIT_???) and/or the constant
CALCEPH_USE_NAIFID
.If the unit containsCALCEPH_USE_NAIFID
, the NAIF identification numbering system is used for the target (NAIF identification numbers for the list).If the unit does not containCALCEPH_USE_NAIFID
, the old number system is used for the target (see the list in the functionf90calceph_compute()
). - order [INTEGER, intent(in)] ::
The order of derivatives.
- = 0 , only the angular momentum is computed. The first three numbers of PVAJ are valid for the results.
- = 1 , only the angular momentum and the first derivative are computed. The first six numbers of PVAJ are valid for the results.
- = 2 , only the angular momentum and the first and second derivatives are computed. The first nine numbers of PVAJ are valid for the results.
- = 3 , the angular momentum and the first, second and third derivatives are computed. The first twelve numbers of PVAJ are valid for the results.
If order equals to 1, the behavior of
f90calceph_rotangmom_order()
is the same asf90calceph_rotangmom_unit()
. - PVAJ (12) [REAL(8), intent(out) :: An array to receive the angular momentum due to its rotation, divided by the product of the mass and of the square of the radius, and their different order of the derivatives, of the body.
Return: f90calceph_rotangmom_order [INTEGER] :: 0 if an error occurs, otherwise non-zero value.
This function is similar to the function f90calceph_orient_unit()
, except that the order of the computed derivatives is specified.
This function reads, if needed, in the ephemeris file associated to eph and interpolates the angular momentum vector due to the rotation of the body, divided by the product of the mass and of the square of the radius
, of a single body (target) for the time JD0+time and stores the results to PVAJ. The angular momentum
, due to the rotation of the body, is defined as the product of the inertia matrix
by the angular velocity vector
. So the returned value is
The order of the derivatives are specified by order. The ephemeris file associated to eph must have been previously opened with the function
f90calceph_open()
.
The output values are expressed in the units specified by unit.
This function checks the units if invalid combinations of units are given to the function.
The returned array PVAJ must be large enough to store the results.
- PVAJ[1:3] contain the angular momentum and is always valid.
- PVAJ[4:6] contain the first derivative and is only valid if order is greater or equal to 1.
- PVAJ[7:9] contain the second derivative and is only valid if order is greater or equal to 2.
- PVAJ[10:12] contain the third derivative and is only valid if order is equal to 3.
The values stored in the array PVAJ are expressed in the following units
- The angular momentum and its derivatives are expressed in days if unit contains
CALCEPH_UNIT_DAY
.- The angular momentum and its derivatives are expressed in seconds if unit contains
CALCEPH_UNIT_SEC
.
The date (JD0, time) should be expressed in the same timescale as the ephemeris files, which can be retrieved using the function f90calceph_gettimescale()
.
The following example prints only the angular momentum, due to its rotation, of the Earth at time=2451419.5
integer*8 peph
integer res
double precision jd0
double precision dt1
double precision G(3)
jd0 = 2451419
dt1 = 0.5D0
res = f90calceph_open(peph, "example2_rotangmom.dat")
if (res.eq.1) then
res = f90calceph_rotangmom_order(peph,jd0, dt1, NAIFID_EARTH,
& CALCEPH_USE_NAIFID+CALCEPH_UNIT_SEC,
& G)
write(*,*) G
call f90calceph_close(peph)
endif
f90calceph_getconstant¶
-
function
f90calceph_getconstant
(eph, name, value)¶ Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- name [CHARACTER(len=*), intent(in)] :: name of the constant
- value [REAL(8), intent(out)] :: first value of the constant
Return: f90calceph_getconstant [INTEGER] :: returns 0 if an error occurs, otherwise the number of values associated to the constant.
This function returns the value associated to the constant name in the header of the ephemeris file associated to eph. Only the first value is returned if multiple values are associated to a constant, such as a list of values.
This function is the same function as f90calceph_getconstantsd()
.
The following example prints the value of the astronomical unit stored in the ephemeris file
integer*8 peph
integer res
double precision AU
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
! print the value of AU
if (f90calceph_getconstant(peph, "AU", AU).eq.1) then
write (*,*) "AU=", AU
endif
call f90calceph_close(peph)
endif
f90calceph_getconstantsd¶
-
function
f90calceph_getconstantsd
(eph, name, value)¶ Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- name [CHARACTER(len=*), intent(in)] :: name of the constant
- value [REAL(8), intent(out)] :: first value of the constant
Return: f90calceph_getconstantsd [INTEGER] :: returns 0 if an error occurs, otherwise the number of values associated to the constant.
This function returns, as a floating-point number, the value associated to the constant name in the header of the ephemeris file associated to eph. Only the first value is returned if multiple values are associated to a constant, such as a list of values. The value must be a floating-point or integer number, otherwise an error is reported.
This function is the same function as f90calceph_getconstant()
.
The following example prints the value of the astronomical unit stored in the ephemeris file
integer*8 peph
integer res
double precision AU
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
! print the value of AU
if (f90calceph_getconstantsd(peph, "AU", AU).eq.1) then
write (*,*) "AU=", AU
endif
call f90calceph_close(peph)
endif
f90calceph_getconstantvd¶
-
function
f90calceph_getconstantvd
(eph, name, arrayvalue, nvalue)¶ Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- name [CHARACTER(len=*), intent(in)] :: name of the constant
- value [REAL(8), dimension(1:nvalue), intent(inout)] :: array of values for the constant
- nvalue [INTEGER, intent(in)] :: number of elements of the array
Return: f90calceph_getconstantvd [INTEGER] :: returns 0 if an error occurs, otherwise the number of values associated to the constant.
This function stores, to the array arrayvalue as floating-point numbers, the nvalue first values associated to the constant name in the header of the ephemeris file associated to eph. The integer value returned by the function is equal to the number of valid entries in the arrayvalue if nvalue is greater or equal to that integer value..
The required value nvalue to store all values can be determinated with the previous call to f90calceph_getconstantsd.
The values must be floating-point or integer numbers, otherwise an error is reported.
The following example prints the body radii of the earth stored in the ephemeris file
integer*8 peph
integer res, nvalue
double precision, allocatable :: radii
double precision svalue
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
! get the number of values
nvalue = calceph_getconstantsd(peph, "BODY399_RADII", svalue)
! fill the array
allocate(radii(1:nvalue))
res = calceph_getconstantvd(peph, "BODY399_RADII", radii, nvalue)
write(*,*) radii
call f90calceph_close(peph)
endif
f90calceph_getconstantss¶
-
function
f90calceph_getconstantss
(eph, name, value)¶ Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- name [CHARACTER(len=*), intent(in)] :: name of the constant
- value [CHARACTER(len=CALCEPH_MAX_CONSTANTNAME), intent(out)] :: first value of the constant
Return: f90calceph_getconstantss [INTEGER] :: returns 0 if an error occurs, otherwise the number of values associated to the constant.
This function returns, as a string of character, the value associated to the constant name in the header of the ephemeris file associated to eph. Only the first value is returned if multiple values are associated to a constant, such as a list of values. The value must be a string, otherwise an error is reported.
Trailing blanks are added to each value.
The following example prints the value of the unit stored in the ephemeris file
integer*8 peph
integer res
character(len=CALCEPH_MAX_CONSTANTVALUE) UNIT
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
! print the value of UNIT
if (f90calceph_getconstantss(peph, "UNIT", UNIT).eq.1) then
write (*,*) "UNIT=", trim(UNIT)
endif
call f90calceph_close(peph)
endif
f90calceph_getconstantvs¶
-
function
f90calceph_getconstantvs
(eph, name, arrayvalue, nvalue)¶ Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- name [CHARACTER(len=*), intent(in)] :: name of the constant
- value [CHARACTER(len=CALCEPH_MAX_CONSTANTNAME), dimension(1:nvalue), intent(inout)] :: array of values for the constant
- nvalue [INTEGER, intent(in)] :: number of elements of the array
Return: f90calceph_getconstantvs [INTEGER] :: returns 0 if an error occurs, otherwise the number of values associated to the constant.
This function stores, to the array arrayvalue as strings of characters, the nvalue first values associated to the constant name in the header of the ephemeris file associated to eph. The integer value returned by the function is equal to the number of valid entries in the arrayvalue if nvalue is greater or equal to that integer value.
The required value nvalue to store all values can be determinated with the previous call to f90calceph_getconstantss.
The values must be strings, otherwise an error is reported.
Trailing blanks are added to each value.
The following example prints the units of the mission stored in the ephemeris file
integer*8 peph
integer res, nvalue
character(len=CALCEPH_MAX_CONSTANTVALUE), allocatable :: mission_units
character(len=CALCEPH_MAX_CONSTANTVALUE) svalue
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
! get the number of values
nvalue = calceph_getconstantss(peph, "MISSION_UNITS", svalue)
! fill the array
allocate(mission_units(1:nvalue))
res = calceph_getconstantvs(peph, "MISSION_UNITS", mission_units, nvalue)
write(*,*) mission_units
call f90calceph_close(peph)
endif
f90calceph_getconstantcount¶
-
function
f90calceph_getconstantcount
(eph)¶ Parameters: eph [INTEGER(8), intent(in)] :: ephemeris descriptor Return: f90calceph_getconstantcount [INTEGER] :: number of constants. 0 if an error occurs, otherwise non-zero value.
This function returns the number of constants available in the header of the ephemeris file associated to eph.
The following example prints the number of available constants stored in the ephemeris file
integer*8 peph
integer res
integer n
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
n = f90calceph_getconstantcount(peph)
write (*,*) "number of constants", n
call f90calceph_close(peph)
endif
f90calceph_getconstantindex¶
-
function
f90calceph_getconstantindex
(eph, index, name, value)¶ Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- index [INTEGER, intent(in)] :: index of the constant, between 1 and
f90calceph_getconstantcount()
- name [CHARACTER(len=CALCEPH_MAX_CONSTANTNAME), intent(out)] :: name of the constant
- value [REAL(8), intent(out)] :: first value of the constant
Return: f90calceph_getconstantindex [INTEGER] :: returns 0 if an error occurs, otherwise the number of values associated to the constant.
This function returns the name and its value of the constant available at the specified index in the header of the ephemeris file associated to eph. The value of index must be between 1 and f90calceph_getconstantcount()
.
Only the first value is returned if multiple values are associated to a constant, such as a list of values. If the first value is not an floating-point number, such as a string, then the function returns 0 without raising an error.
Trailing blanks are added to the name of the constant.
The following example displays the name of the constants, stored in the ephemeris file, and their values
integer*8 peph
integer res
integer n
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
do j=1, f90calceph_getconstantcount(peph)
res = f90calceph_getconstantindex(peph,j,nameconstant,valueconstant)
write (*,*) nameconstant,"=",valueconstant
enddo
call f90calceph_close(peph)
endif
f90calceph_getfileversion¶
-
function
f90calceph_getfileversion
(eph, version)¶ Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- version [CHARACTER(len=CALCEPH_MAX_CONSTANTVALUE), intent(out)] :: version of the file
Return: f90calceph_getfileversion [INTEGER] :: returns 0 if the file version was not found, otherwise non-zero value.
This function returns the version of the ephemeris file, as a string. For example, the argument version will contain 'INPOP10B', 'EPM2017' or 'DE405', ... .
If the file is an original JPL binary planetary ephemeris, then the version of the file can always be determined. If the file is a spice kernel, the version of the file is retrieved from the constant INPOP_PCK_VERSION, EPM_PCK_VERSION, or PCK_VERSION.
The following example prints the version of the ephemeris file.
integer*8 peph
integer res
character(len=CALCEPH_MAX_CONSTANTVALUE) version
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
res = f90calceph_getfileversion(peph, version)
write (*,*) "The version of the file is ", version
call f90calceph_close(peph)
endif
f90calceph_gettimescale¶
-
function
f90calceph_gettimescale
(eph)¶ Parameters: eph [INTEGER(8), intent(in)] :: ephemeris descriptor Return: f90calceph_gettimescale [INTEGER] :: 0 if an error occurs, otherwise non-zero value.
- This function returns the timescale of the ephemeris file associated to eph :
- 1 if the quantities of all bodies are expressed in the TDB time scale.
- 2 if the quantities of all bodies are expressed in the TCB time scale.
The following example prints the time scale available in the ephemeris file
integer*8 peph
integer res
double precision AU
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
! print the time scale
timescale = calceph_gettimescale(peph)
write (*,*) "timescale=", timescale
call f90calceph_close(peph)
endif
f90calceph_gettimespan¶
-
function
f90calceph_gettimespan
(eph, firsttime, lasttime, continuous)¶ Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- firsttime [REAL(8), intent(out)] :: Julian date of the first time
- lasttime [REAL(8), intent(out)] :: Julian date of the last time
- continuous [INTEGER, intent(out)] :: information about the availability of the quantities over the time span
Return: f90calceph_gettimespan [INTEGER] :: 0 if an error occurs, otherwise non-zero value.
This function returns the first and last time available in the ephemeris file associated to eph. The Julian date for the first and last time are expressed in the time scale returned by f90calceph_gettimescale()
.
It returns the following value in the parameter continuous :
- 1 if the quantities of all bodies are available for any time between the first and last time.
- 2 if the quantities of some bodies are available on discontinuous time intervals between the first and last time.
- 3 if the quantities of each body are available on a continuous time interval between the first and last time, but not available for any time between the first and last time.
The following example prints the first and last time available in the ephemeris file
integer*8 peph
integer res
integer continuous
double precision firsttime, lasttime
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
if (f90calceph_gettimespan(peph, firsttime, lasttime, continuous).eq.1) then
write (*,*) firsttime, lasttime, countinuous
endif
call f90calceph_close(peph)
endif
f90calceph_getpositionrecordcount¶
-
function
f90calceph_getpositionrecordcount
(eph)¶ Parameters: eph [INTEGER(8), intent(in)] :: ephemeris descriptor Return: f90calceph_getpositionrecordcount [INTEGER] :: number of position's records. 0 if an error occurs, otherwise non-zero value.
This function returns the number of position's records available in the ephemeris file associated to eph. Usually, the number of records is equal to the number of bodies in the ephemeris file if the timespan is continuous. If the timespan is discontinuous for the target and center bodies, then each different timespan is counted as a different record. If the ephemeris file contain timescale transformations' records, such as TT-TDB or TCG-TCB, then these records are included in the returned value.
The following example prints the number of position's records available in the ephemeris file
integer*8 peph
integer res
integer n
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
n = f90calceph_getpositionrecordcount(peph)
write (*,*) "number of position's record", n
call f90calceph_close(peph)
endif
f90calceph_getpositionrecordindex¶
-
function
f90calceph_getpositionrecordindex
(eph, index, target, center, firsttime, lasttime, frame)¶ Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- index [INTEGER, intent(int)] :: index of the position's record, between 1 and
f90calceph_getpositionrecordcount()
- target [INTEGER, intent(out)] :: The target body
- center [INTEGER, intent(out)] :: The origin body
- firsttime [REAL(8), intent(out)] :: Julian date of the first time
- lasttime [REAL(8), intent(out)] :: Julian date of the last time
- frame [INTEGER, intent(out)] :: reference frame (see the list, below)
Return: f90calceph_getpositionrecordindex [INTEGER] :: 0 if an error occurs, otherwise non-zero value.
This function returns the target and origin bodies, the first and last time, and the reference frame available at the specified index for the position's records of the ephemeris file associated to eph.
The NAIF identification numbering system is used for the target and center integers (NAIF identification numbers for the list).
The Julian date for the first and last time are expressed in the time scale returned by f90calceph_gettimescale()
.
It returns the following value in the parameter frame :
value | Name |
---|---|
1 | ICRF |
The following example displays the position's records stored in the ephemeris file.
integer*8 peph
integer j, itarget, icenter, iframe
real*8 firsttime, lasttime
integer res
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
do j=1, f90calceph_getpositionrecordcount(peph)
res = f90calceph_getpositionrecordindex(peph,j,itarget, icenter, firsttime, lasttime, iframe)
write (*,*) itarget, icenter, firsttime, lasttime, iframe
enddo
call f90calceph_close(peph)
endif
f90calceph_getpositionrecordindex2¶
-
function
f90calceph_getpositionrecordindex2
(eph, index, target, center, firsttime, lasttime, frame, segid)¶ Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- index [INTEGER, intent(int)] :: index of the position's record, between 1 and
f90calceph_getpositionrecordcount()
- target [INTEGER, intent(out)] :: The target body
- center [INTEGER, intent(out)] :: The origin body
- firsttime [REAL(8), intent(out)] :: Julian date of the first time
- lasttime [REAL(8), intent(out)] :: Julian date of the last time
- frame [INTEGER, intent(out)] :: reference frame (see the list, below)
- segid [INTEGER, intent(out)] :: type of the segment.
Return: f90calceph_getpositionrecordindex2 [INTEGER] :: 0 if an error occurs, otherwise non-zero value.
This function returns the target and origin bodies, the first and last time, the reference frame, and the segment type available at the specified index for the position's records of the ephemeris file associated to eph.
The NAIF identification numbering system is used for the target and center integers (NAIF identification numbers for the list).
The Julian date for the first and last time are expressed in the time scale returned by f90calceph_gettimescale()
.
It returns the following value in the parameter frame :
value | Name |
---|---|
1 | ICRF |
It returns in the parameter segid one of the predefined constants CALCEPH_SEGTYPE_... (Constants).
The following example displays the position's records stored in the ephemeris file.
integer*8 peph
integer j, itarget, icenter, iframe, iseg
real*8 firsttime, lasttime
integer res
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
do j=1, f90calceph_getpositionrecordcount(peph)
res = f90calceph_getpositionrecordindex2(peph,j,itarget, icenter, firsttime, lasttime, iframe, iseg)
write (*,*) itarget, icenter, firsttime, lasttime, iframe, iseg
enddo
call f90calceph_close(peph)
endif
f90calceph_getorientrecordcount¶
-
function
f90calceph_getorientrecordcount
(eph)¶ Parameters: eph [INTEGER(8), intent(in)] :: ephemeris descriptor Return: f90calceph_getorientrecordcount [INTEGER] :: number of orientation's records. 0 if an error occurs, otherwise non-zero value.
This function returns the number of orientation's records available in the ephemeris file associated to eph. Usually, the number of records is equal to the number of bodies in the ephemeris file if the timespan is continuous. If the timespan is discontinuous for the target body, then each different timespan is counted as a different record.
The following example prints the number of orientation's records available in the ephemeris file
integer*8 peph
integer res
integer n
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
n = f90calceph_getorientrecordcount(peph)
write (*,*) "number of orientation's record", n
call f90calceph_close(peph)
endif
f90calceph_getorientrecordindex¶
-
function
f90calceph_getpositionrecordindex
(eph, index, target, firsttime, lasttime, frame) Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- index [INTEGER, intent(int)] :: index of the orientation's record, between 1 and
f90calceph_getorientrecordcount()
- target [INTEGER, intent(out)] :: The target body
- firsttime [REAL(8), intent(out)] :: Julian date of the first time
- lasttime [REAL(8), intent(out)] :: Julian date of the last time
- frame [INTEGER, intent(out)] :: reference frame (see the list, below)
Return: f90calceph_getorientrecordindex [INTEGER] :: 0 if an error occurs, otherwise non-zero value.
This function returns the target body, the first and last time, and the reference frame available at the specified index for the orientation's records of the ephemeris file associated to eph.
The NAIF identification numbering system is used for the target body (NAIF identification numbers for the list).
The Julian date for the first and last time are expressed in the time scale returned by f90calceph_gettimescale()
.
It returns the following value in the parameter frame :
value | Name |
---|---|
1 | ICRF |
The following example displays the orientation's records stored in the ephemeris file.
integer*8 peph
integer j, itarget, iframe
real*8 firsttime, lasttime
integer res
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
do j=1, f90calceph_getorientrecordcount(peph)
res = f90calceph_getorientrecordindex(peph,j,itarget, firsttime, lasttime, iframe)
write (*,*) itarget, firsttime, lasttime, iframe
enddo
call f90calceph_close(peph)
endif
f90calceph_getorientrecordindex2¶
-
function
f90calceph_getpositionrecordindex2
(eph, index, target, firsttime, lasttime, frame, segid) Parameters: - eph [INTEGER(8), intent(in)] :: ephemeris descriptor
- index [INTEGER, intent(int)] :: index of the orientation's record, between 1 and
f90calceph_getorientrecordcount()
- target [INTEGER, intent(out)] :: The target body
- firsttime [REAL(8), intent(out)] :: Julian date of the first time
- lasttime [REAL(8), intent(out)] :: Julian date of the last time
- frame [INTEGER, intent(out)] :: reference frame (see the list, below)
- segid [INTEGER, intent(out)] :: type of the segment.
Return: f90calceph_getorientrecordindex2 [INTEGER] :: 0 if an error occurs, otherwise non-zero value.
This function returns the target body, the first and last time, the reference frame and the segment type available at the specified index for the orientation's records of the ephemeris file associated to eph.
The NAIF identification numbering system is used for the target body (NAIF identification numbers for the list).
The Julian date for the first and last time are expressed in the time scale returned by f90calceph_gettimescale()
.
It returns the following value in the parameter frame :
value | Name |
---|---|
1 | ICRF |
It returns in the parameter segid one of the predefined constants CALCEPH_SEGTYPE_... (Constants).
The following example displays the orientation's records stored in the ephemeris file.
integer*8 peph
integer j, itarget, iframe, iseg
real*8 firsttime, lasttime
integer res
res = f90calceph_open(peph, "example1.dat")
if (res.eq.1) then
do j=1, f90calceph_getorientrecordcount(peph)
res = f90calceph_getorientrecordindex2(peph,j,itarget, firsttime, lasttime, iframe, iseg)
write (*,*) itarget, firsttime, lasttime, iframe, iseg
enddo
call f90calceph_close(peph)
endif