Common Argument Definitions

Throughout the library, numerous input and output arguments retain the same meaning and format. This section provides the definitions for those arguments. The particular definitions given here are for the functions that perform a particular calculation for multiple points (ntime, up to NTIME_MAX, defined in IRBEM_GET_NTIME_MAX()). For a few functions, there is also a pitch-angle dimension given by Nipa, which is up to NALPHA_MAX = 25. For functions that can compute for only a single point or angle the NTIME_MAX argument can be ignored.

ntime: long integer number of time in arrays (max allowed is NTIME_MAX)

nipa: long integer number of pitch angles in arrays (max allowed is NALPHA_MAX)

kext: long integer to select external magnetic field

0   = no external field
1   = Mead & Fairfield [1975] (uses 0≤Kp≤9 - Valid for rGEO≤17. Re)
2   = Tsyganenko short [1987] (uses 0≤Kp≤9 - Valid for rGEO≤30. Re)
3   = Tsyganenko long [1987] (uses 0≤Kp≤9 - Valid for rGEO≤70. Re)
4   = Tsyganenko [1989c] (uses 0≤Kp≤9 - Valid for rGEO≤70. Re)
5   = Olson & Pfitzer quiet [1977] (default - Valid for rGEO≤15. Re)
6   = Olson & Pfitzer dynamic [1988] (uses 5.≤dens≤50., 300.≤velo≤500., -100.≤Dst≤20. - Valid for rGEO≤60. Re)
7   = Tsyganenko [1996] (uses -100.≤Dst (nT)≤20., 0.5≤Pdyn (nPa)≤10., |ByIMF| (nT)≤10., |BzIMF| (nT)≤10. - Valid for rGEO≤40. Re)
8   = Ostapenko & Maltsev [1997] (uses dst,Pdyn,BzIMF, Kp)
9   = Tsyganenko [2001] (uses -50.≤Dst (nT)≤20., 0.5≤Pdyn (nPa)≤5., |ByIMF| (nT)≤5., |BzIMF| (nT)≤5., 0.≤G1≤10., 0.≤G2≤10. - Valid for xGSM≥-15. Re)
10 =Tsyganenko [2001] storm  (uses Dst, Pdyn, ByIMF, BzIMF, G2, G3 - there is no upper or lower limit for those inputs - Valid for xGSM≥-15. Re)
11 =Tsyganenko [2004] storm  (uses Dst, Pdyn, ByIMF, BzIMF, W1, W2, W3, W4, W5, W6 - there is no upper or lower limit for those inputs - Valid for xGSM≥-15. Re)
12 =Alexeev [2000] - also known as Paraboloid model - Submitted to ISO  (uses Dens, velo, Dst, BzIMF, AL)


Notes:
  1. when the magnetic field model inputs are outside the allowed range bad data values are returned.
  2. When solar wind inputs are required they must be taken in the vicinity of the day side magnetopause and not at L1.


options:
array(5) of long integer to set some control options on computed values

  1. options(1st element):  0 - don't compute L* or Φ ;  1 - compute L*; 2- compute Φ
  2. options(2nd element): 0 - initialize IGRF field once per year (year.5);  n - n is the  frequency (in days) starting on January 1st of each year (i.e. if options(2nd element)=15 then IGRF will be updated on the following days of the year: 1, 15, 30, 45 ...)
  3. options(3rd element): resolution to compute L* (0 to 9) where 0 is the recomended value to ensure a good ratio precision/computation time (i.e. an error of ~2% at L=6). The higher the value the better will be the precision, the longer will be the computing time. Generally there is not much improvement for values larger than 4. Note that this parameter defines the integration step (θ) along the field line such as dθ=(π)/(720*[options(3rd element)+1])
  4. options(4th element): resolution to compute L* (0 to 9). The higher the value the better will be the precision, the longer will be the computing time. It is recommended to use 0 (usually sufficient) unless L* is not computed on a LEO orbit. For LEO orbit higher values are recommended. Note that this parameter defines the integration step (φ) along the drift shell such as dφ=(2π)/(25*[options(4th element)+1])
  5. options(5th element): allows to select an internal magnetic field model (default is set to IGRF)

sysaxes, sysaxesIN, and sysaxesOUT: long integer to define which coordinate system is provided in

0: GDZ (alti, lati, East longi - km,deg.,deg)
1: GEO (cartesian) - Re
2: GSM (cartesian) - Re
3: GSE (cartesian) - Re
4: SM (cartesian) - Re
5: GEI (cartesian) - Re
6: MAG (cartesian) - Re
7: SPH (geo in spherical) - (radial distance, lati, East longi - Re, deg., deg.)
8: RLL  (radial distance, lati, East longi - Re, deg., deg. - prefered to 7)

iyear, iyr or iyearsat: array(NTIME_MAX) of long integer year when measurements are given. /font>
idoy: array(NTIME_MAX) of long integer doy when measurements are given
UT or secs: array(NTIME_MAX) of double, UT in seconds
Note: In matlab, year/doy/UT arguments are replaced by the single matlabd argument, which is a Matlab date numbers (a double precision day counter created with datenum). They can be generated, if needed, using the helper function [iyear,idoy,UT] = onera_desp_lib_matlabd2yds(matlabd).

x1: array(NTIME_MAX) of double, first coordinate according to sysaxes. If sysaxes is 0 then altitude has to be in km otherwise use dimensionless variables (in Re)

x2: array(NTIME_MAX) of double, second coordinate according to sysaxes. If sysaxes is 0 then latitude has to be in degrees otherwise use dimensionless variables (in Re)

x3: array(NTIME_MAX) of double, third coordinate according to sysaxes. If sysaxes is 0 then longitude has to be in degrees otherwise use dimensionless variables (in Re).

alpha: array(NALPHA_MAX) double, pitch-angle at input location (in degres).

R0: double, radius, in Re, of the minimum allowed radial distance along the drift orbit (use R0 < 1) for the drift loss cone.

maginput: array (25,NTIME_MAX) of double to specify magnetic fields inputs such as:

  1. maginput(1st element,*) =Kp: value of Kp as in OMNI2 files but has to be double instead of integer type. (NOTE, consistent with OMNI2, this is Kp*10, and it is in the range 0 to 90)
  2. maginput(2nd element,*) =Dst: Dst index (nT)
  3. maginput(3rd element,*) =dens: Solar Wind density (cm-3)
  4. maginput(4th element,*) =velo: Solar Wind velocity (km/s)
  5. maginput(5th element,*) =Pdyn: Solar Wind dynamic pressure (nPa)
  6. maginput(6th element,*) =ByIMF: GSM y component of IMF mag. field (nT)
  7. maginput(7th element,*) =BzIMF: GSM z component of IMF mag. field (nT)
  8. maginput(8th element,*) =G1:  G1=< Vsw*(Bperp/40)2/(1+Bperp/40)*sin3(theta/2) > where the <> mean an average over the previous 1 hour, Vsw is the solar wind speed, Bperp is the transverse IMF component (GSM) and theta it's clock angle.
  9. maginput(9th element,*) =G2: G2=< a*Vsw*Bs > where the <> mean an average over the previous 1 hour, Vsw is the solar wind speed, Bs=|IMF Bz| when IMF Bz < 0 and Bs=0 when IMF Bz > 0, a=0.005.
  10. maginput(10th element,*) =G3:  G3=< Vsw*Dsw*Bs /2000.>
    where the <> mean an average over the previous 1 hour, Vsw is the solar wind speed, Dsw is the solar wind density, Bs=|IMF Bz| when IMF Bz < 0 and Bs=0 when IMF Bz > 0.
  11. maginput(11th element,*) =W1 see definition in (JGR-A, v.110(A3), 2005.) (PDF 1.2MB)
  12. maginput(12th element,*) =W2 see definition in (JGR-A, v.110(A3), 2005.) (PDF 1.2MB)
  13. maginput(13th element,*) =W3 see definition in (JGR-A, v.110(A3), 2005.) (PDF 1.2MB)
  14. maginput(14th element,*) =W4 see definition in (JGR-A, v.110(A3), 2005.) (PDF 1.2MB)
  15. maginput(15th element,*) =W5 see definition in (JGR-A, v.110(A3), 2005.) (PDF 1.2MB)
  16. maginput(16th element,*) =W6 see definition in (JGR-A, v.110(A3), 2005.) (PDF 1.2MB)
  17. maginput(17th element,*) =AL the auroral index
  18. maginput(18th element,*) to maginput(25th element,*): for future use

Lm: L McIlwain (array(NTIME_MAX,NALPHA_MAX) of double)

IMPORTANT: When using a non dipole field Lm is pitch-angle dependent because of the McIlwain L definition which has some limitation but this dependency has nothing to do with shell-splitting!

Lstar or Φ : L Roederer  or Φ=2π*Bo*/Lstar [nT Re2] (array(NTIME_MAX,NALPHA_MAX) of double)
Blocal: magnitude of magnetic field at point (array(NTIME_MAX,NALPHA_MAX) of double) - [nT]
Bmir or Bmirror: magnitude of magnetic field at mirror point (array(NTIME_MAX,NALPHA_MAX) of double) - [nT]
Bmin: magnitude of magnetic field at equator (array(NTIME_MAX,NALPHA_MAX) of double) - [nT]
XJ: I, related to second adiabatic invariant (array(NTIME_MAX,NALPHA_MAX) of double) - [Re]
MLT: magnetic local time in hours, (array(NTIME_MAX,NALPHA_MAX) of double) - [hour]

POSIT: (note: size varies. Changes from what's given here will be provided in the function definitions). xGEO,yGEO and zGEO along the drift shell, (array(3,1000,48) of double)
1st element: x, y, z GEO coord, 2nd element: points along field line, 3rd element: number of field lines. Note, sometimes the long dimension is 3000 rather than 1000; see specific functions for clarification.
Nposit: (note: size varies. Changes from what's given here will be provided in the function definitions) long integer array (48) providing the number of points along the field line for each field line traced in 2nd element of POSIT array (max 1000 for some functions, 3000 for others).

Note: Posit is organized as follow: array(xyz,along bounce,drift index). Nposit tells how many of the 1000 (or 3000) points available are used to described one bounce.




Comments on adiabatic invariants:
J2=2*p*I/(2*π) where p is the particle momemtum - [kg m2 s-2 Re]
J3=q/(2*π) *Φ - [C nT Re2]


Below is provided an chart explaining the logic for the coding of Lm and Lstar from the routine:

 


MAKE_LSTAR

DESCRIPTION:

This function allows one to compute magnetic coordinate at any s/c position, i.e. L, L*, Blocal/Bmirror, Bequator. A set of internal/external field can be selected.

INPUTS:

ntime, kext, options, sysaxes, iyear, idoy, UT, x1, x2, x3, maginput

see Common Argument Definitions

IMPORTANT: all inputs must be present. For those which are not used a dummy value can be provided.
 

OUTPUTS:

Lm, Lstar or Φ, Blocal, Bmin, XJ, MLT

see Common Argument Definitions


CALLING SEQUENCE FROM MATLAB:

[Lm,Lstar,Blocal,Bmin,J,MLT] = onera_desp_lib_make_lstar(kext,options,sysaxes,matlabd,x1,x2,x3,maginput)

CALLING SEQUENCE FROM IDL:

result = call_external(lib_name, 'make_lstar_', ntime,kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3, maginput, lm,lstar,blocal,bmin,xj,mlt, /f_value)
 

CALLING SEQUENCE FROM FORTRAN:

call make_lstar1(ntime,kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3, maginput, lm,lstar,blocal,bmin,xj,mlt)


LANDI2LSTAR

DESCRIPTION:

This function allows one to compute magnetic coordinate at any s/c position, i.e. L, L*, Blocal/Bmirror, Bequator. A set of internal/external field can be selected. Right now this routine differs from MAKE_LSTAR because L* is deduced from Lm, I and day of year empirically. This is set only for IGRF+Olson-Pfitzer quiet field model (kext can only be 5). The errors in L* values are less than 2%.

INPUTS:

ntime, kext, options, sysaxes, iyear, idoy, UT, x1, x2, x3, maginput

see Common Argument Definitions

IMPORTANT: all inputs must be present. For those which are not used a dummy value can be provided.

OUTPUTS:

Lm, Lstar or Φ, Blocal, Bmin, XJ, MLT

see Common Argument Definitions

CALLING SEQUENCE FROM MATLAB:

[Lm,Lstar,Blocal,Bmin,J,MLT] = onera_desp_lib_landi2lstar(kext,options,sysaxes,matlabd,x1,x2,x3,maginput)

CALLING SEQUENCE FROM IDL:

result = call_external(lib_name, 'landi2lstar_', ntime,kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3, maginput, lm,lstar,blocal,bmin,xj,mlt, /f_value)
 

CALLING SEQUENCE FROM FORTRAN:

call landi2lstar1(ntime,kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3, maginput, lm,lstar,blocal,bmin,xj,mlt)


EMPIRICALLSTAR

DESCRIPTION:

This function allows one to compute L* empirically being given Lm, I and day of year. This is set only for IGRF+Olson-Pfitzer quiet field model (kext can only be 5) so Lm and I provided as input must be computed with this field. The errors in L* values are less than 2%.

INPUTS:

ntime, kext, options, iyearsat, idoy, maginput, Lm, XJ

see Common Argument Definitions

IMPORTANT: all inputs must be present. For those which are not used a dummy value can be provided.

OUTPUTS:

Lstar or Φ

see Common Argument Definitions

CALLING SEQUENCE FROM MATLAB:

Lstar = onera_desp_lib_empiricallstar(kext,options,matlabd,maginput,Lm,J)

CALLING SEQUENCE FROM IDL:

result = call_external(lib_name, 'empiricalLstar_', ntime,kext,options,iyearsat,idoy, maginput, lm,xj,lstar, /f_value)
 

CALLING SEQUENCE FROM FORTRAN:

call EmpiricalLstar1(ntime,kext,options,iyearsat,idoy,maginput, lm,xj,lstar)



MAKE_LSTAR_SHELL_SPLITTING

DESCRIPTION:

This function allows one to compute the magnetic coordinates at any s/c position and pitch-angle, i.e. L, L*, Bmirror, Bequator. A set of internal/external field can be selected.

INPUTS:

ntime, Nipa, kext, options, sysaxes, iyear, idoy, UT, x1, x2, x3, alpha, maginput

see Common Argument Definitions

IMPORTANT: all inputs must be present. For those which are not used a dummy value can be provided.
 

OUTPUTS:

Lm, Lstar or Φ, Blocal, Bmin, XJ, MLT

see Common Argument Definitions

CALLING SEQUENCE FROM MATLAB:

[Lm,Lstar,Bmirror,Bmin,J,MLT] = onera_desp_lib_make_lstar_shell_splitting(kext,options,sysaxes,matlabd,x1,x2,x3,alpha,maginput)

CALLING SEQUENCE FROM IDL:

result = call_external(lib_name, 'make_lstar_shell_splitting_', ntime,Npa,kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3,alpha,maginput,lm,lstar,blocal,bmin,xj,mlt, /f_value)
 

CALLING SEQUENCE FROM FORTRAN:

call make_lstar_shell_splitting1(ntime,Npa,kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3, alpha,maginput,lm,lstar,blocal,bmin,xj,mlt)



LANDI2LSTAR_SHELL_SPLITTING

DESCRIPTION:

This function allows one to compute the magnetic coordinates at any s/c position and pitch-angle, i.e. L, L*, Bmirror, Bequator. A set of internal/external field can be selected. Right now this routine differs from MAKE_LSTAR_SHELL_SPLITTING because L* is deduced from Lm, I and day of year empirically. This is set only for IGRF+Olson-Pfitzer quiet field model (kext can only be 5). The errors in L* values are less than 2%.

INPUTS:

ntime, Nipa, kext, options, sysaxes, iyear, idoy, UT, x1, x2, x3, alpha, maginput

see Common Argument Definitions

IMPORTANT: all inputs must be present. For those which are not used a dummy value can be provided.
 

OUTPUTS:

Lm, Lstar or Φ, Bmirror, Bmin, XJ, MLT

see Common Argument Definitions

CALLING SEQUENCE FROM MATLAB:

[Lm,Lstar,Bmirror,Bmin,J,MLT] = onera_desp_lib_landi2lstar_shell_splitting(kext,options,sysaxes,matlabd,x1,x2,x3,alpha,maginput)

CALLING SEQUENCE FROM IDL:

result = call_external(lib_name, 'landi2lstar_shell_splitting_', ntime,Npa,kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3,alpha,maginput,lm,lstar,blocal,bmin,xj,mlt, /f_value)
 

CALLING SEQUENCE FROM FORTRAN:

call landi2lstar_shell_splitting1(ntime,Npa,kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3, alpha,maginput,lm,lstar,blocal,bmin,xj,mlt)



DRIFT_SHELL

DESCRIPTION:

This function traces a full drift shell for particles that have their mirror point at the input location.  The output is a full array of positions of the drift shell, usefull for plotting and visualisation (for just the points on the drift-bounce orbit, use DRIFT_BOUNCE_ORBIT). A set of internal/external field can be selected.

INPUTS:

INPUTS:

kext, options, sysaxes, iyear, idoy, UT, x1, x2, x3, maginput

see Common Argument Definitions

IMPORTANT: all inputs must be present. For those which are not used a dummy value can be provided.
 

OUTPUTS:

Lm, Lstar or Φ, Blocal, Bmin, XJ, POSIT, Nposit
Blocal is (1000,48), Lm, Lstar/Φ, Bmin, and XJ are a scalars

see Common Argument Definitions

CALLING SEQUENCE FROM MATLAB:

[Lm,Lstar,Blocal,Bmin,J,POSIT] = onera_desp_lib_drift_shell(kext,options,sysaxes,matlabd,x1,x2,x3,maginput)

CALLING SEQUENCE FROM IDL:

result = call_external(lib_name, 'drift_shell_',  kext,options,sysaxes,iyear,idoy,ut,  x1,x2,x3, maginput, lm,lstar,blocal,bmin,xj,posit,Nposit,  /f_value)
 

CALLING SEQUENCE FROM FORTRAN:

call drift_shell1(kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3, maginput, lm,lstar,blocal,bmin,xj,posit,Nposit)



DRIFT_BOUNCE_ORBIT

DESCRIPTION:

This function traces a full drift-bounce orbit for particles with a specified pitch angle at the input location.  The output is a full array of positions of the drift-bounce orbit, usefull for computing drift-bounce averages. Key differences from DRIFT_SHELL: (1) only positions between mirror points are returned, (2) 25 rather than 48 azimuths are returned, (3) Lstar accuracy respects options(3) and options(4), (4) a new parameter R0 is required which specifies the minimum radial distance allowed along the drift path (usually R0=1, but use R0<1 in the drift loss cone), (5) hmin and hmin_lon outputs provide the altitude and longitude (GDZ) of the minimum altitude point along the orbit (among those traced, not just those returned). A set of internal/external field can be selected.

INPUTS:

kext, options, sysaxes, iyear, idoy, UT, x1, x2, x3, alpha, maginput, R0

see Common Argument Definitions

IMPORTANT: all inputs must be present. For those which are not used a dummy value can be provided.
 

OUTPUTS:

Lm, Lstar, Blocal, Bmin, Bmirror, XJ, POSIT, Nposit, hmin, hmin_lon
Note: POSIT is size (3,1000,25) and Nposit is size (25)
Blocal is (1000,25), Lm, Lstar/Φ, Bmin, Bmirror, XJ, hmin, and hmin_lon are a scalars

see Common Argument Definitions

hmin: double, minimum altitude, GDZ, km, along the drift orbit
hmin_lon: double, longitude, GDZ, degrees, at hmin.

CALLING SEQUENCE FROM MATLAB:

[Lm,Lstar,Blocal,Bmin,Bmir,J,POSIT,hmin,hmin_lon] = onera_desp_lib_drift_bounce_orbit(kext,options,sysaxes,matlabd,x1,x2,x3,alpha,maginput,R0)

CALLING SEQUENCE FROM IDL:

result = call_external(lib_name, 'drift_bounce_orbit2_',  kext,options,sysaxes,iyear,idoy,ut,  x1,x2,x3,alpha, maginput,R0, lm,lstar,blocal,bmin,bmir,xj,posit,Nposit,hmin,hmin_lon,  /f_value)
 

CALLING SEQUENCE FROM FORTRAN:

call drift_bounce_orbit2_1(kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3,alpha, maginput,R0 lm,lstar,blocal,bmin,bmir,xj,posit,Nposit,hmin,hmin_lon)



FIND_MIRROR_POINT

DESCRIPTION:

This function finds the magnitude and location of the mirror point along a field line traced from any given location and local pitch-angle for a set of internal/external field to be selected. 

INPUTS:

kext, options, sysaxes, iyear, idoy, UT, x1, x2, x3, alpha, maginput

see Common Argument Definitions

IMPORTANT: all inputs must be present. For those which are not used a dummy value can be provided.
 

OUTPUTS:

Blocal, Bmir, POSIT
Blocal and Bmir are scalars. POSIT is size (3) and provides the GEO coordinaets of the mirror point in Re. see Common Argument Definitions

CALLING SEQUENCE FROM MATLAB:

[Blocal,Bmirror,xGEO] = onera_desp_lib_find_mirror_point(kext,options,sysaxes,matlabd,x1,x2,x3,alpha,maginput)

CALLING SEQUENCE FROM IDL:

result = call_external(lib_name,  'find_mirror_point_', kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3,alpha,maginput, blocal,bmir,posit, /f_value)
 

CALLING SEQUENCE FROM FORTRAN:

call find_mirror_point1(kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3,alpha, maginput,blocal,bmir,posit)



FIND_MAGEQUATOR

DESCRIPTION:

This function finds the GEO coordinates of the magnetic equator along the field line starting from input location for a set of internal/external field to be selected. 

INPUTS:

kext, options, sysaxes, iyear, idoy, UT, x1, x2, x3, maginput

see Common Argument Definitions

IMPORTANT: all inputs must be present. For those which are not used a dummy value can be provided.
 

OUTPUTS:

Bmin: magnitude of magnetic field at equator (double)

POSIT: xGEO,yGEO and zGEO of the magnetic equator location, (array(3) of double)

CALLING SEQUENCE FROM MATLAB:

[Bmin,xGEO] = onera_desp_lib_find_magequator(kext,options,sysaxes,matlabd,x1,x2,x3,maginput)

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'find_magequator_',  kext,options,sysaxes,iyear,idoy,ut,  x1,x2,x3, maginput, bmin,posit,  /f_value)
 

CALLING SEQUENCE from FORTRAN:

call find_magequator1(kext,options,sysaxes,iyear,idoy,ut,  x1,x2,x3, maginput,bmin,posit)



FIND_FOOT_POINT

DESCRIPTION:

This function finds the of the field line crossing a specified altitude in a specified hemisphere for a set of internal/external field to be selected. 

INPUTS:

kext, options, sysaxes, iyear, idoy, UT, x1, x2, x3, stop_alt, hemi_flag, maginput

See Common Argument Definitions

stop_alt: double, desired altitude of field-line crossing, km.

hemi_flag: long integer, hemisphere flag.

0   = same magnetic hemisphere as starting point
+1   = northern magnetic hemisphere
-1   = southern magnetic hemisphere
+2   = opposite magnetic hemisphere as starting point

IMPORTANT: all inputs must be present. For those which are not used a dummy value can be provided.
 

OUTPUTS:

XFOOT: location of foot point, GDZ, (array(3) of double)
BFOOT: magnetic field vector at foot point, GEO, nT, (array(3) of double)
BFOOTMAG: magnetic field magnitude at foot point, GEO, nT, (array(3) of double)
 

CALLING SEQUENCE FROM MATLAB:

[Xfoot,Bfoot,BfootMag] = onera_desp_lib_find_foot_point(kext,options,sysaxes,matlabd,x1,x2,x3,stop_alt,hemi_flag,maginput)

CALLING SEQUENCE FROM IDL:

result = call_external(lib_name,  'find_foot_point_', kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3,stop_alt,hemi_flag,maginput, xfoot,bfoot,bfootmag, /f_value)
 

CALLING SEQUENCE FROM FORTRAN:

call find_foot_point1(kext,options,sysaxes,iyearsat,idoy,UT,xIN1,xIN2,xIN3,stop_alt,hemi_flag,maginput,XFOOT,BFOOT,BFOOTMAG)



GET_FIELD

DESCRIPTION:

This function is deprecated. Use GET_FIELD_MULTI instead.

GET_FIELD_MULTI

DESCRIPTION:

This function computes the GEO vector of the magnetic field at input location for a set of internal/external magnetic field to be selected. 
 

INPUTS:

ntime, kext, options, sysaxes, iyear, idoy, UT, x1, x2, x3, maginput

See Common Argument Definitions
IMPORTANT: all inputs must be present. For those which are not used a dummy value can be provided.

OUTPUTS:

Bgeo: BxGEO,ByGEO and BzGEO of the magnetic field (nT), (array(3,NTIME_MAX) of double)
Bl: magnitude of magnetic field (array(NTIME_MAX) of double) (nT)
 

CALLING SEQUENCE FROM MATLAB:

[Bgeo,B] = onera_desp_lib_get_field(kext,options,sysaxes,matlabd,x1,x2,x3,maginput)

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'get_field_multi_idl_',ntime,kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3, maginput,Bgeo, Bl,  /f_value)
 

CALLING SEQUENCE from FORTRAN:

call GET_FIELD_MULTI(ntime,kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3, maginput,Bgeo,Bl)



GET_BDERIVS

DESCRIPTION:

This function computes the magnetic field and its 1st-order derivatives at each input location for a set of internal/external magnetic field to be selected. 
 

INPUTS:

ntime, kext, options, sysaxes, dX, iyear, idoy, UT, x1, x2, x3, maginput

See Common Argument Definitions

dX: dX is the step size, in RE for the numerical derivatives (recommend 1E-3), double scalar.

IMPORTANT: all inputs must be present. For those which are not used a dummy value can be provided.

OUTPUTS:

Bgeo: BxGEO,ByGEO and BzGEO of the magnetic field (nT), (array(3,NTIME_MAX) of double)
Bmag: magnitude of magnetic field (array(NTIME_MAX) of double) (nT)
gradBmag: gradient (GEO) of magnitude of magnetic field (array(3,NTIME_MAX) of double) (nT)
diffB: derivatives of magnetic field vector (array(3,3,NTIME_MAX) of double) (nT)
diffB(i,j,t) = dB_i/dx_j at t'th location. GEO coordinates.

 

CALLING SEQUENCE FROM MATLAB:

[Bgeo,B,gradBmag,diffB] = onera_desp_lib_get_bderivs(kext,options,sysaxes,matlabd,x1,x2,x3,maginput,...)

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'get_bderivs_idl',ntime,kext,options,sysaxes,dX,iyear,idoy,ut, x1,x2,x3, maginput,Bgeo,Bmag,gradBmag,diffB,  /f_value)
 

CALLING SEQUENCE from FORTRAN:

call GET_BDERIVS(ntime,kext,options,sysaxes,dX,iyear,idoy,ut, x1,x2,x3, maginput,Bgeo,Bmag,gradBmag,diffB)



COMPUTE_GRAD_CURV_CURL

DESCRIPTION:

This function computes gradient and curvature force factors and div/curl of B from the outputs of GET_BDERIVS 
 

INPUTS:

ntime, Bgeo, Bmag, gradBmag, diffB defined as for GET_BDERIVS

OUTPUTS:

grad_par: gradient of Bmag along B nT/RE, (array(NTIME_MAX) of double)
grad_perp: gradient of Bmag perpendicular to B nT/RE (array(3,NTIME_MAX) of double)
grad_drift: (bhat x grad_perp)/B, 1/RE (part of gradient drift velocity) (array(3,NTIME_MAX) of double)
curvature: (bhat dot grad)bhat, 1/RE (part of curvature force) (array(3,NTIME_MAX) of double)
Rcurv: 1/|curvature| RE (radius of curvature) (array(3,NTIME_MAX) of double)
curv_drift: (bhat x curvature), 1/RE (part of curvature drift) (array(NTIME_MAX) of double)
curlB: curl of B (nT/RE) (part of electrostatic current term) (array(3,NTIME_MAX) of double)
divB: divergence of B (nT/RE) (should be zero!) (array(NTIME_MAX) of double)

 

CALLING SEQUENCE FROM MATLAB:

[grad_par,grad_perp,grad_drift,curvature,Rcurv,curv_drift,curlB,divB] = onera_desp_lib_compute_grad_curv_curl(Bgeo,B,gradBmag,diffB)

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'compute_grad_curv_idl',ntime,Bgeo,Bmag,gradBmag,diffB, grad_par,grad_perp,grad_drift,curvature,Rcurv,curv_drift,curlB,divB,   /f_value)
 

CALLING SEQUENCE from FORTRAN:

call COMPUTE_GRAD_CURV_CURL(ntime,Bgeo,Bmag,gradBmag,diffB, grad_par,grad_perp,grad_drift,curvature,Rcurv,curv_drift,curlB,divB)



TRACE_FIELD_LINE2

DESCRIPTION:

This function traces a full field line which crosses the input position.  The output is a full array of positions of the field line, usefull for plotting and visualisation for a set of internal/external fields to be selected. A new option (R0) for TRACE_FIELD_LINE2 allows user to specify the radius (RE) of the reference surface between which the line is traced (R0=1 in TRACE_FIELD_LINE)  

INPUTS:

kext, options, sysaxes, iyear, idoy, UT, x1, x2, x3, maginput, R0

See Common Argument Definitions

R0: double, specifies radius of reference surface between which field line is traced.


IMPORTANT: all inputs must be present. For those which are not used a dummy value can be provided.
 

OUTPUTS:

Lm, Blocal, Bmin, XJ, POSIT, Nposit

See Common Argument Definitions
POSIT is (3,3000), Nposit is a scalar.  


CALLING SEQUENCE FROM MATLAB:

[Lm,Blocal,Bmin,J,POSIT] = onera_desp_lib_trace_field_line(kext,options,sysaxes,matlabd,x1,x2,x3,maginput,R0)

CALLING SEQUENCE FROM IDL:

result = call_external(lib_name, 'trace_field_line2_',  kext,options,sysaxes,iyear,idoy,ut,  x1,x2,x3,maginput,R0,lm,blocal,bmin,xj,posit,Nposit,  /f_value)
 

CALLING SEQUENCE FROM FORTRAN:

call trace_field_line2_1(kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3, maginput,R0, lm,blocal,bmin,xj,posit,Nposit)



TRACE_FIELD_LINE_TOWARD_EARTH

DESCRIPTION:

This function traces a field line from the input position to the Earth surface.  The output is a full array of positions of the field line, usefull for plotting and visualisation for a set of internal/external fields to be selected. 

INPUTS:

kext, options, sysaxes, iyear, idoy, UT, x1, x2, x3, maginput, ds

See Common Argument Definitions

ds: double, Integration step along field line (in Re).


IMPORTANT: all inputs must be present. For those which are not used a dummy value can be provided.
 

OUTPUTS:

POSIT, Nposit

See Common Argument Definitions
POSIT is (3,3000), Nposit is a scalar.

CALLING SEQUENCE FROM MATLAB:

POSIT = onera_desp_lib_trace_field_line_towards_earth(kext,options,sysaxes,matlabd,x1,x2,x3,maginput,ds)

CALLING SEQUENCE FROM IDL:

result = call_external(lib_name, 'trace_field_line_towards_earth_',  kext,options,sysaxes,iyear,idoy,ut,  x1,x2,x3, maginput, ds,posit,Nposit,  /f_value)
 

CALLING SEQUENCE FROM FORTRAN:

call trace_field_line_towards_earth1(kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3, maginput, ds,posit,Nposit)


GET_MLT

DESCRIPTION:

Routine to get Magnetic Local Time (MLT) from a Cartesian GEO position and date
 

INPUTS:

iyr : long integer year
idoy : long integer day of year
secs : UT in seconds - double
xGEO : 3D array (double) of cartesian position in GEO (Re)
 

OUTPUTS:

MLT: Magnetic Local Time (hours) - double

 

CALLING SEQUENCE FROM MATLAB:

MLT = onera_desp_lib_get_mlt(matlabd,xGEO)

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'get_mlt_', iyr,idoy,secs,xGEO,MLT, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call get_mlt1(iyr,idoy,secs,xGEO,MLT)


GET_HEMI

DESCRIPTION:

This routine is deprecated. Use GET_HEMI_MULTI instead.

GET_HEMI_MULTI

DESCRIPTION:

This function computes in which magnetic hemisphere is the input location for a set of internal/externalmagnetic field to select. 

 

INPUTS:

ntime, kext, options, sysaxes, iyear, idoy, UT, x1, x2, x3, maginput

See Common Argument Definitions

IMPORTANT: all inputs must be present. For those which are not used a dummy value can be provided.

OUTPUTS:

xHEMI: +1 means Northern magnetic hemisphere, -1 means Southern magnetic hemisphere, 0 means cannot be defined (magnetic field not valid) - long 
 

CALLING SEQUENCE FROM MATLAB:

[xHEMI] = onera_desp_lib_get_hemi(kext,options,sysaxes,matlabd,x1,x2,x3,maginput)

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'get_hemi_multi_idl_',ntime,kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3, maginput,xHEMI,  /f_value)
 

CALLING SEQUENCE from FORTRAN:

call GET_HEMI_MULTI(ntime,kext,options,sysaxes,iyear,idoy,ut, x1,x2,x3, maginput,xHEMI)




LSTAR_PHI

DESCRIPTION:

This function allows one to convert L* to Φ or vice versa.

INPUTS:

ntime, whichinv, options, iyear, idoy

See Common Argument Definitions
NOTE: options does not support user-supplied internal field model.

whichinv: long integer to select which conversion to perform 

1   = Lstar to Φ 
2   = Φ to Lstar

IMPORTANT: all inputs must be present. For those which are not used a dummy value can be provided.
 

INPUTS or OUTPUTS (according to whichinv):

Lstar, Phi See Common Argument Definitions

Phi: Φ=2π*Bo/Lstar [nT Re2] (array(GET_IRBEM_NTIME_MAX()) of double) 

CALLING SEQUENCE FROM MATLAB:

out = onera_desp_lib_lstar_phi(which,options,matlabd,in)

CALLING SEQUENCE FROM IDL:

result = call_external(lib_name, lstar_phi_', ntime,whichinv,options,iyear,idoy,lstar,phi, /f_value)
 

CALLING SEQUENCE FROM FORTRAN:

call lstar_phi1(ntime,whichinv,options,iyear,idoy,lstar,phi)


COORD_TRANS

DESCRIPTION:

Deprecated. Use COORD_TRANS_VEC

COORD_TRANS_VEC

DESCRIPTION:

Generic coordinates transformation from one Earth or Heliospheric coordinate system to another  Earth or Heliospheric coordinates system. Handle up to GET_IRBEM_NTIME_MAX() times at once.
 

INPUTS:

ntime, sysaxesIN, sysaxesOUT, iyr, idoy, secs, xIN

see Common Argument Definitions

sysaxesIN: long integer indicating input coordinates system 
sysaxesOUT : long integer indicating output coordinates system
xIN : array (3,NTIME_MAX) of double with position in input coordinates system.

IMPORTANT: all inputs must be present. For those which are not used a dummy value can be provided.

Input

GDZ (0)

GEO (1)

GSM (2)

GSE (3)

SM(4)

GEI(5)

MAG (6)

SPH(7)

RLL(8)

HEE (9)

HAE(10)

HEEQ (11)

Output

GDZ (0)

GEO (1)

GSM (2)

GSE (3)

SM (4)

GEI (5)

MAG (6)

SPH (7)

RLL (8)

HEE (9)

HAE (10)

HEEQ (11)

Table: Available coordinates transformation:
 ⇒Transformation implemented
 ⇒Will not crash but not very usefull 
 ⇒Not implemented

OUTPUTS:

xOUT : array (3,NTIME_MAX) of double with position in output coordinates system. (systems defined by sysaxesOUT)
 

CALLING SEQUENCE FROM MATLAB:

Y=onera_desp_lib_coord_trans(X,rotation,matlabd) 

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'coord_trans_vec_', ntime,sysaxesIN,sysaxesOUT,iyr,idoy,secs,xIN,xOUT, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call coord_trans_vec1(ntime,sysaxesIN,sysaxesOUT,iyr,idoy,secs,xIN,xOUT)


GEO2GSM

DESCRIPTION:

Routine to transform Cartesian GEO to cartesian GSM coordinates
 

INPUTS:

iyr : long integer year
idoy : long integer day of year
secs : UT in seconds - double
xGEO : 3D array (double) of cartesian position in GEO (Re)
 

OUTPUTS:

psi : angle for GSM coordinate - double
xGSM :
3D array (double) of cartesian position in GSM (Re)

 

CALLING SEQUENCE FROM MATLAB:

xGSM = onera_desp_lib_rotate(xGEO,'geo2gsm',matlabd);
[xGSM,psi] = onera_desp_lib_rotate(xGEO,'geo2gsm',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'geo2gsm_', iyr,idoy,secs,psi,xGEO,xGSM, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call geo2gsm1(iyr,idoy,secs,psi,xGEO,xGSM)




GSM2GEO

DESCRIPTION:

Routine to transform Cartesian GSM to cartesian GEO coordinates

INPUTS:

iyr : long integer year
idoy : long integer day of year
secs : UT in seconds - double
xGSM : 3D array (double) of cartesian position in GSM (Re)
 

OUTPUTS:

psi : angle for GSM coordinate - double
xGEO : 3D array (double) of cartesian position in GEO (Re)
 

CALLING SEQUENCE FROM MATLAB:

xGEO = onera_desp_lib_rotate(xGSM,'gsm2geo',matlabd);
[xGEO,psi] = onera_desp_lib_rotate(xGSM,'gsm2geo',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'gsm2geo_', iyr,idoy,secs,psi,xGSM,xGEO, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call gsm2geo1(iyr,idoy,secs,psi,xGSM,xGEO)



GEO2GSE

DESCRIPTION:

Routine to transform Cartesian GEO to cartesian GSM coordinates

INPUTS:

iyr : long integer year
idoy : long integer day of year
secs : UT in seconds - double
xGEO : 3D array (double) of cartesian position in GEO (Re)
 

OUTPUTS:

xGSE : 3D array (double) of cartesian position in GSE (Re)
 

CALLING SEQUENCE FROM MATLAB:

xGSE = onera_desp_lib_rotate(xGEO,'geo2gse',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'geo2gse_', iyr,idoy,secs,xGEO,xGSE, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call geo2gse1(iyr,idoy,secs,xGEO,xGSE)



GSE2GEO

DESCRIPTION:

Routine to transform Cartesian GSE to cartesian GEO coordinates
 

INPUTS:

iyr : long integer year
idoy : long integer day of year
secs : UT in seconds - double
xGSE : 3D array (double) of cartesian position in GSE (Re)
 

OUTPUTS:

xGEO : 3D array (double) of cartesian position in GEO (Re)
 

CALLING SEQUENCE FROM MATLAB:

xGEO = onera_desp_lib_rotate(xGSE,'gse2geo',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'gse2geo_', iyr,idoy,secs,xGSE,xGEO, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call gse2geo1(iyr,idoy,secs,xGSE,xGEO)



GDZ2GEO

DESCRIPTION:

Routine to transform GEODEZIC coordinates to cartesian GEO coordinates
 

INPUTS:

lati : latitude (degres) (double)
longi : East longitude (degres) (double)
alti : altitude (km) (double)
 

OUTPUTS:

xx : xGEO (Re) (double)
yy : yGEO (Re) (double)
zz : zGEO (Re) (double)
 

CALLING SEQUENCE FROM MATLAB:

******NOTE ORDER OF GDZ COORDINATES********
xGEO = onera_desp_lib_rotate([alti(:), lati(:), longi(:)],'gdz2geo');
xx = xGEO(:,1); yy = xGEO(:,2); zz = xGEO(:,3);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'gdz2geo_', lati,longi,alti,xx,yy,zz, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call gdz_geo(lati,longi,alti,xx,yy,zz)



GEO2GDZ

DESCRIPTION:

Routine to transform cartesian GEO coordinates to GEODEZIC coordinates
 

INPUTS:

xx : xGEO (Re) (double)
yy : yGEO (Re) (double)
zz : zGEO (Re) (double)
 

OUTPUTS:

lati : latitude (degres) (double)
longi : East longitude (degres) (double)
alti : altitude (km) (double)
 

CALLING SEQUENCE FROM MATLAB:

******NOTE ORDER OF GDZ COORDINATES********
xGDZ = onera_desp_lib_rotate([xx(:) yy(:) zz(:)],'geo2gdz');
alti = xGDZ(:,1); lati = xGDZ(:,2); longi = xGDZ(:,3);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'geo2gdz_', xx,yy,zz,lati,longi,alti, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call geo_gdz(xx,yy,zz,lati,longi,alti)



GEO2GEI

DESCRIPTION:

Routine to transform Cartesian GEO to cartesian GEI coordinates
 

INPUTS:

iyr : long integer year
idoy : long integer day of year
secs : UT in seconds -double
xGEO : 3D array (double) of cartesian position in GEO (Re)
 

OUTPUTS:

xGEI: 3D array (double) of cartesian position in GEI (Re)
 

CALLING SEQUENCE FROM MATLAB:

xGEI = onera_desp_lib_rotate(xGEO,'geo2gei',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'geo2gei_', iyr,idoy,secs,xGEO,xGEI, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call geo2gei1(iyr,idoy,secs,xGEO,xGEI)



GEI2GEO

DESCRIPTION:

Routine to transform Cartesian GEI to cartesian GEO coordinates
 

INPUTS:

iyr : long integer year
idoy : long integer day of year
secs : UT in seconds - double
xGEI : 3D array (double) of cartesian position in GEI (Re)
 

OUTPUTS:

xGEO: 3D array (double) of cartesian position in GEO (Re)
 

CALLING SEQUENCE FROM MATLAB:

xGEO = onera_desp_lib_rotate(xGEI,'gei2geo',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'gei2geo_', iyr,idoy,secs,xGEI,xGEO, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call gei2geo1(iyr,idoy,secs,xGEI,xGEO)



GEO2SM

DESCRIPTION:

Routine to transform Cartesian GEO to cartesian SM coordinates
 

INPUTS:

iyr : long integer year
idoy : long integer day of year
secs : UT in seconds - double
xGEO : 3D array (double) of cartesian position in GEO (Re)
 

OUTPUTS:

xSM: 3D array (double) of cartesian position in SM (Re)
 

CALLING SEQUENCE FROM MATLAB:

xSM = onera_desp_lib_rotate(xGEO,'geo2sm',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'geo2sm_', iyr,idoy,secs,xGEO,xSM, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call geo2sm1(iyr,idoy,secs,xGEO,xSM)



SM2GEO

DESCRIPTION:

Routine to transform Cartesian SM to cartesian GEO coordinates
 

INPUTS:

iyr : long integer year
idoy : long integer day of year
secs : UT in seconds - double
xSM : 3D array (double) of cartesian position in SM (Re)
 

OUTPUTS:

xGEO: 3D array (double) of cartesian position in GEO (Re)
 

CALLING SEQUENCE FROM MATLAB:

xGEO = onera_desp_lib_rotate(xSM,'sm2geo',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'sm2geo_', iyr,idoy,secs,xSM,xGEO, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call sm2geo1(iyr,idoy,secs,xSM,xGEO)



GSM2SM

DESCRIPTION:

Routine to transform Cartesian GSM to cartesian SM coordinates
 

INPUTS:

iyr : long integer year
idoy : long integer day of year
secs : UT in seconds - double
xGSM : 3D array (double) of cartesian position in GSM (Re)
 

OUTPUTS:

xSM: 3D array (double) of cartesian position in SM (Re)
 

CALLING SEQUENCE FROM MATLAB:

xSM = onera_desp_lib_rotate(xGSM,'gsm2sm',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'gsm2sm_', iyr,idoy,secs,xGSM,xSM, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call gsm2sm1(iyr,idoy,secs,xGSM,xSM)



SM2GSM

DESCRIPTION:

Routine to transform Cartesian SM to cartesian GSM coordinates
 

INPUTS:

iyr : long integer year
idoy : long integer day of year
secs : UT in seconds - double
xSM : 3D array (double) of cartesian position in SM (Re)
 

OUTPUTS:

xGSM: 3D array (double) of cartesian position in GSM (Re)
 

CALLING SEQUENCE FROM MATLAB:

xGSM = onera_desp_lib_rotate(xSM,'sm2gsm',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'sm2gsm_', iyr,idoy,secs,xSM,xGSM, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call sm2gsm1(iyr,idoy,secs,xSM,xGSM)



GEO2MAG

DESCRIPTION:

Routine to transform Cartesian GEO to cartesian MAG coordinates
 

INPUTS:

iyr : long integer year
xGEO : 3D array (double) of cartesian position in GEO (Re)
 

OUTPUTS:

xMAG: 3D array (double) of cartesian position in MAG (Re)
 

CALLING SEQUENCE FROM MATLAB:

xMAG = onera_desp_lib_rotate(xGEO,'geo2mag',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'geo2mag_', iyr,xGEO,xMAG, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call geo2mag1(iyr,xGEO,xMAG)



MAG2GEO

DESCRIPTION:

Routine to transform Cartesian MAG to cartesian GEO coordinates
 

INPUTS:

iyr : long integer year
xMAG : 3D array (double) of cartesian position in MAG (Re)
 

OUTPUTS:

xGEO: 3D array (double) of cartesian position in GEO (Re)
 

CALLING SEQUENCE FROM MATLAB:

xGEO = onera_desp_lib_rotate(xMAG,'mag2geo',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'mag2geo_', iyr,xMAG,xGEO, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call mag2geo1(iyr,xMAG,xGEO)



SPH2CAR

DESCRIPTION:

Routine to transform spherical coordinates to cartesian
 

INPUTS:

r : radial distance (double)
lati: (degree)
longi: (degree)
 

OUTPUTS:

x: 3D array (double) of cartesian position (same unit as r)
 

CALLING SEQUENCE FROM MATLAB:

xCAR = onera_desp_lib_rotate([r(:), lati(:), longi(:)],'sph2car');

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'sph2car_', r,lati,longi,x, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call SPH_CAR(r,lati,longi,x)



CAR2SPH

DESCRIPTION:

Routine to transform cartesian coordinates to spherical
 

INPUTS:

x : 3D array (double) of cartesian position
 

OUTPUTS:

r : radial distance (double) - same unit as x
lati : (degree)
longi : (degree)
 

CALLING SEQUENCE FROM MATLAB:

xSPH = onera_desp_lib_rotate(xCAR,'car2sph');
r = xSPH(:,1); lati = xSPH(:,2); longi = xSPH(:,3);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'car2sph_', x,r,lati,longi, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call CAR_SPH(x,r,lati,longi)



RLL2GDZ

DESCRIPTION:

Routine to transform radial distance, latitude, longitude to altitude, latitude, longitude
 

INPUTS:

r : radial distance (double) - Re
lati : (double- degree)
longi : (double- degree)
 

OUTPUTS:

alti : altitude (double) - km

 

CALLING SEQUENCE FROM MATLAB:

******NOTE ORDER OF GDZ COORDINATES********
xGDZ = onera_desp_lib_rotate([r(:), lati(:), longi(:)],'rll2gdz');
alti = xGDZ(:,1); lati = xGDZ(:,2); longi = xGDZ(:,3);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'rll2gdz_', r,lati,longi,alti, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call RLL_GDZ(r,lati,longi,alti)



GSE2HEE

DESCRIPTION:

Routine to transform geocentric coordinates GSE to heliospheric coordinates HEE (Heliocentric Earth Ecliptic also sometime known as Heliospheric Solar Ecliptic (HSE)
 

INPUTS:

iyr : year (long) 
idoy: day of year - January 1st is idoy=1 (long) 
UT : Universal Time (double- seconds)
xGSE : 3D array (double) of cartesian position in GSE (Re)

OUTPUTS:

xHEE : 3D array (double) of cartesian position in HEE (AU)

 

CALLING SEQUENCE FROM MATLAB:

xHEE = onera_desp_lib_rotate(xGSE,'gse2hee',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'gse2hee_', iyr,idoy,UT,xGSE,xHEE, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call gse2hee1(iyr,idoy,UT,xGSE,xHEE)



HEE2GSE

DESCRIPTION:

Routine to transform heliospheric coordinates HEE (Heliocentric Earth Ecliptic also sometime known as Heliospheric Solar Ecliptic (HSE) to geocentric coordinates GSE
 

INPUTS:

iyr : year (long) 
idoy: day of year - January 1st is idoy=1 (long) 
UT : Universal Time (double- seconds)
xHEE : 3D array (double) of cartesian position in HEE (AU)

OUTPUTS:

xGSE : 3D array (double) of cartesian position in GSE (Re)

 

CALLING SEQUENCE FROM MATLAB:

xGSE = onera_desp_lib_rotate(xHEE,'hee2gse',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'hee2gse_', iyr,idoy,UT,xHEE,xGSE, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call hee2gse1(iyr,idoy,UT,xHEE,xGSE)



HEE2HAE

DESCRIPTION:

Routine to transform heliospheric coordinates HEE (Heliocentric Earth Ecliptic also sometime known as Heliospheric Solar Ecliptic (HSE) to heliospheric coordinates HAE (Heliocentric Aries Ecliptic also sometime known as Heliospheric Solar Ecliptic (HSE)
 

INPUTS:

iyr : year (long) 
idoy: day of year - January 1st is idoy=1 (long) 
UT : Universal Time (double- seconds)
xHEE : 3D array (double) of cartesian position in HEE (AU)

OUTPUTS:

xHAE : 3D array (double) of cartesian position in HAE (AU)

 

CALLING SEQUENCE FROM MATLAB:

xHAE = onera_desp_lib_rotate(xHEE,'hee2hae',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'hee2hae_', iyr,idoy,UT,xHEE,xHAE, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call hee2hae1(iyr,idoy,UT,xHEE,xHAE)



HAE2HEE

DESCRIPTION:

Routine to transform heliospheric coordinates HAE (Heliocentric Aries Ecliptic also sometime known as Heliospheric Solar Ecliptic (HSE) to heliospheric coordinates HEE (Heliocentric Earth Ecliptic also sometime known as Heliospheric Solar Ecliptic (HSE)
 

INPUTS:

iyr : year (long) 
idoy: day of year - January 1st is idoy=1 (long) 
UT : Universal Time (double- seconds)
xHAE : 3D array (double) of cartesian position in HAE (AU)

OUTPUTS:

xHEE : 3D array (double) of cartesian position in HEE (AU)

 

CALLING SEQUENCE FROM MATLAB:

xHEE = onera_desp_lib_rotate(xHAE,'hae2hee',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'hae2hee_', iyr,idoy,UT,xHAE,xHEE, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call hae2hee1(iyr,idoy,UT,xHAE,xHEE)



HAE2HEEQ

DESCRIPTION:

Routine to transform heliospheric coordinates HAE (Heliocentric Aries Ecliptic also sometime known as Heliospheric Solar Ecliptic (HSE) to heliospheric coordinates HEEQ (Heliocentric Earth Equatorial also sometime known as Heliospheric Solar (HS)
 

INPUTS:

iyr : year (long) 
idoy: day of year - January 1st is idoy=1 (long) 
UT : Universal Time (double- seconds)
xHAE : 3D array (double) of cartesian position in HAE (AU)

OUTPUTS:

xHEEQ : 3D array (double) of cartesian position in HEEQ (AU)

 

CALLING SEQUENCE FROM MATLAB:

xHEEQ = onera_desp_lib_rotate(xHAE,'hae2heeq',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'hae2heeq_', iyr,idoy,UT,xHAE,xHEEQ, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call hae2heeq1(iyr,idoy,UT,xHAE,xHEEQ)



HEEQ2HAE

DESCRIPTION:

Routine to transform heliospheric coordinates HEEQ (Heliocentric Earth Equatorial also sometime known as Heliospheric Solar (HS) to heliospheric coordinates HAE (Heliocentric Aries Ecliptic also sometime known as Heliospheric Solar Ecliptic (HSE)
 

INPUTS:

iyr : year (long) 
idoy: day of year - January 1st is idoy=1 (long) 
UT : Universal Time (double- seconds)
xHEEQ : 3D array (double) of cartesian position in HEEQ (AU)

OUTPUTS:

xHAE : 3D array (double) of cartesian position in HAE (AU)

 

CALLING SEQUENCE FROM MATLAB:

xHAE = onera_desp_lib_rotate(xHEEQ,'heeq2hae',matlabd);

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'heeq2hae_', iyr,idoy,UT,xHEEQ,xHAE, /f_value)
 

CALLING SEQUENCE from FORTRAN:

call heeq2hae1(iyr,idoy,UT,xHEEQ,xHAE)



JULDAY

DESCRIPTION:

Calculate the Julian Day Number for a given month, day, and year.
 

INPUTS:

MONTH: Number of the desired month (1 = January, ..., 12 = December). - long

DAY: Number of day of the month. - long

YEAR: Number of the desired year.Year parameters must be valid values from the civil calendar. Years B.C.E. are represented as negative integers. Years in the common era are represented as positive integers. In particular, note that there is no year 0 in the civil calendar. 1 B.C.E. (-1) is followed by 1 C.E. (1). - long
 

OUTPUTS:

JULDAY returns the Julian Day Number (which begins at noon) of the specified calendar date. - long
 

CALLING SEQUENCE FROM MATLAB:

Not available, see MATLAB build-in routines

CALLING SEQUENCE from IDL:

Not available, see IDL build-in routines
 

CALLING SEQUENCE from FORTRAN:

Result = JULDAY(Year,Month, Day)



CALDAT

DESCRIPTION:

Return the calendar date given julian day. This is the inverse of the function JULDAY.
 

INPUTS:

JULIAN contains the Julian Day Number (which begins at noon) of the specified calendar date. It should be a long integer.
 

OUTPUTS:

MONTH: Number of the desired month (1 = January, ..., 12 = December). - long

DAY: Number of day of the month. - long

YEAR: Number of the desired year. - long
 

CALLING SEQUENCE FROM MATLAB:

Not available, see MATLAB build-in routines

CALLING SEQUENCE from IDL:

Not available, see IDL build-in routines
 

CALLING SEQUENCE from FORTRAN:

CALL CALDAT(julian, year,month, day)



GET_DOY

DESCRIPTION:

Calculate the day of year for a given month, day, and year.
 

INPUTS:

MONTH: Number of the desired month (1 = January, ..., 12 = December).

DAY: Number of day of the month.

YEAR: Number of the desired year.Year parameters must be valid values from the civil calendar. Years B.C.E. are represented as negative integers. Years in the common era are represented as positive integers. In particular, note that there is no year 0 in the civil calendar. 1 B.C.E. (-1) is followed by 1 C.E. (1).
 

OUTPUTS:

GET_DOY returns the day of year of the specified calendar date.
 

CALLING SEQUENCE FROM MATLAB:

Not available, see MATLAB build-in routines

CALLING SEQUENCE from IDL:

Not available, see IDL build-in routines
 

CALLING SEQUENCE from FORTRAN:

Result = GET_DOY(Year,Month, Day)



DECY2DATE_AND_TIME

DESRIPTION:

Calculate the date and time (yea,month,day of month, day of year, hour, minute and second and Universal Time)
 

INPUTS:

DEC_Y : Decimal year where yyyy.0d0 is January 1st at 00:00:00 - double
 

OUTPUTS:

YEAR: Number year.Year parameters must be valid values from the civil calendar. Years B.C.E. are represented as negative integers. Years in the common era are represented as positive integers. In particular, note that there is no year 0 in the civil calendar. 1 B.C.E. (-1) is followed by 1 C.E. (1). - long

MONTH: Number month (1 = January, ..., 12 = December). - long

DAY: Number of day of the month. - long

DOY: Number of day of year (DOY=1 is for January 1st) - long

HOUR, MINUTE and SECOND: Universal time in the day - long

UT: Univeral time in seconds - double
 

CALLING SEQUENCE FROM MATLAB:

Not available, see MATLAB build-in routines

CALLING SEQUENCE from IDL:

Not available, see IDL build-in routines

CALLING SEQUENCE from FORTRAN:

CALL DECY2DATE_AND_TIME(Dec_y,Year,Month, Day, doy, hour,minute,second,UT)


 


DATE_AND_TIME2DECY

DESRIPTION:

Calculate the decimal year from date and time (yea,month,day of month, hour, minute and second)
 

INPUTS:

YEAR: Number year.Year parameters must be valid values from the civil calendar. Years B.C.E. are represented as negative integers. Years in the common era are represented as positive integers. In particular, note that there is no year 0 in the civil calendar. 1 B.C.E. (-1) is followed by 1 C.E. (1). - long

MONTH: Number month (1 = January, ..., 12 = December). - long

DAY: Number of day of the month. - long

HOUR, MINUTE and SECOND: Universal time in the day - long


OUTPUTS:

DEC_Y : Decimal year where yyyy.0d0 is January 1st at 00:00:00 - double

CALLING SEQUENCE FROM MATLAB:

Not available, see MATLAB build-in routines

CALLING SEQUENCE from IDL:

Not available, see IDL build-in routines

CALLING SEQUENCE from FORTRAN:

CALL DATE_AND_TIME2DECY(Year,Month,Day,hour,minute,second,Dec_y)


 

DOY_AND_UT2DATE_AND_TIME

DESRIPTION:

Calculate month, day, year from year and day of year
Calculate time (hour, minute, second) from UT

 

INPUTS:

YEAR: Number year.Year parameters must be valid values from the civil calendar. Years B.C.E. are represented as negative integers. Years in the common era are represented as positive integers. In particular, note that there is no year 0 in the civil calendar. 1 B.C.E. (-1) is followed by 1 C.E. (1). - long

DOY: Number of day of year (DOY=1 is for January 1st) - long

UT: Universal time in seconds - double


OUTPUTS:

Month : Number month (1 = January, ..., 12 = December). - long

Day : Number of day of the month - long

HOUR, MINUTE and SECOND: Universal time in the day - long

CALLING SEQUENCE FROM MATLAB:

Not available, see MATLAB build-in routines

CALLING SEQUENCE from IDL:

Not available, see IDL build-in routines

CALLING SEQUENCE from FORTRAN:

CALL DOY_AND_UT2DATE_AND_TIME(Year,Doy,UT,Month, Day, hour,minute,second)





FLY_IN_NASA_AEAP

DESCRIPTION:

This function allows one to fly any spacecraft in NASA AE8 min/max and AP8 min/max models.
The output can be differential flux or flux within an energy range or integral flux.

INPUTS:

ntime: long integer number of time in arrays (max allowed is GET_IRBEM_NTIME_MAX())

sysaxes: long integer to define which coordinate system is provided in

0: GDZ (alti, lati, East longi - km,deg.,deg)
1: GEO (cartesian) - Re
2: GSM (cartesian) - Re
3: GSE (cartesian) - Re
4: SM (cartesian) - Re
5: GEI (cartesian) - Re
6: MAG (cartesian) - Re
7: SPH (geo in spherical) - (radial distance, lati, East longi - Re, deg., deg.)
8: RLL  (radial distance, lati, East longi - Re, deg., deg. - prefered to 7)

whichm: long integer to select in which NASA model to fly

whichmModel
1=AE8 MIN
2=AE8 MAX
3=AP8 MIN
4=AP8 MAX
-1=AE8 MIN - ESA Interpolation
-2=AE8 MAX - ESA Interpolation
-3=AP8 MIN - ESA Interpolation
-4=AP8 MAX - ESA Interpolation

The ESA Interpolation scheme provides for better flux interpolations at low altitudes.
See: E. J. Daly, et al., "Problems with models of the radiation belts",
IEEE Trans. Nucl. Sci, Vol 43, No. 2, Apr. 1996. DOI 10.1109/23.490889

whatf: long integer to select what flux to compute

1=differential flux (MeV-1 cm-2 s-1) at energy(1,*)
2=flux within an ernergy range (MeV-1 cm-2 s-1) - energy(1,*) to energy(2,*)
3=integral flux (cm-2 s-1) at energy(1,*)

Nene: long integer number of energies in array energy (max allowed is 25)

energy: array(2,25) of double. If whatf=1 or 3 then energy(2,*) is not considered.
       Note: if energy is in MeV then differential flux will be per MeV.

iyear: array(GET_IRBEM_NTIME_MAX()) of long integer year when positions are given. Can be a dummy value for the following sysaxes values: 0,1,6,7,8.

idoy: array(GET_IRBEM_NTIME_MAX()) of long integer doy when positions are given. Can be a dummy value for the following sysaxes values: 0,1,6,7,8. 

UT: array(GET_IRBEM_NTIME_MAX()) of double, UT in seconds when positions are given. Can be a dummy value for the following sysaxes values: 0,1,6,7,8 

x1: array(GET_IRBEM_NTIME_MAX()) of double, first coordinate according to sysaxes. If sysaxes is 0 then altitude has to be in km otherwise use dimensionless variables (in Re)

x2: array(GET_IRBEM_NTIME_MAX()) of double, second coordinate according to sysaxes. If sysaxes is 0 then latitude has to be in degrees otherwise use dimensionless variables (in Re)

x3: array(GET_IRBEM_NTIME_MAX()) of double, third coordinate according to sysaxes. If sysaxes is 0 then East longitude has to be in degrees otherwise use dimensionless variables (in Re).


OUTPUTS:

flux: flux according to selection of whatf for all times and energies (array(GET_IRBEM_NTIME_MAX(),25) of double)

CALLING SEQUENCE FROM MATLAB:

Flux = onera_desp_lib_fly_in_nasa_aeap(sysaxes,whichm,energy,matlabd,x1,x2,x3)

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'fly_in_nasa_aeap_', ntime,sysaxes,whichm,whatf,Nene,energy,iyear,idoy, UT,x1,x2,x3,flux, /f_value)

CALLING SEQUENCE from FORTRAN:

CALL fly_in_nasa_aeap1(ntime,sysaxes,whichm,whatf,Nene,energy,iyear,idoy, UT,x1,x2,x3,flux)



GET_AE8_AP8_FLUX

DESCRIPTION:

This function allows one to compute NASA AE8 min/max and AP8 min/max flux for any B/Bo, L position.
The output can be differential flux or flux within an energy range or integral flux.

INPUTS:

ntime: long integer number of points in arrays (max allowed is GET_IRBEM_NTIME_MAX())

whichm: long integer to select in which NASA model to fly

whichmModel
1=AE8 MIN
2=AE8 MAX
3=AP8 MIN
4=AP8 MAX
-1=AE8 MIN - ESA Interpolation
-2=AE8 MAX - ESA Interpolation
-3=AP8 MIN - ESA Interpolation
-4=AP8 MAX - ESA Interpolation

The ESA Interpolation scheme provides for better flux interpolations at low altitudes.
See: E. J. Daly, et al., "Problems with models of the radiation belts",
IEEE Trans. Nucl. Sci, Vol 43, No. 2, Apr. 1996. DOI 10.1109/23.490889

whatf: long integer to select what flux to compute

1=differential flux (MeV-1 cm-2 s-1) at energy(1,*)
2=flux within an ernergy range (MeV-1 cm-2 s-1) - energy(1,*) to energy(2,*)
3=integral flux (cm-2 s-1) at energy(1,*)

Nene: long integer number of energies in array energy (max allowed is 25)

energy: array(2,25) of double. If whatf=1 or 3 then energy(2,*) is not considered.
       Note: if energy is in MeV then differential flux will be per MeV.

BBo: array(GET_IRBEM_NTIME_MAX()) of double. Provide B/Bequator for all position where to compute the fluxes. Note that the Jensen and Cain 1960 magnetic field model must be used for any call to AE8min, AE8max, AP8min and GSFC 1266 (extended to year 1970) for any call to AP8max.

L: array(GET_IRBEM_NTIME_MAX()) of double. Provide McIlwain L for all position where to compute the fluxes.  Note that the Jensen and Cain 1960 magnetic field model must be used for any call to AE8min, AE8max, AP8min and GSFC 1266 (extended to year 1970) for any call to AP8max.


OUTPUTS:

flux: flux according to selection of whatf for all times and energies (array(GET_IRBEM_NTIME_MAX(),25) of double)

CALLING SEQUENCE FROM MATLAB:

Flux = onera_desp_lib_get_ae8_ap8_flux(whichm,energy,BBo,L)

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'get_ae8_ap8_flux_idl_', ntime,whichm,whatf,Nene,energy,BBo,L,,flux, /f_value)

CALLING SEQUENCE from FORTRAN:

CALL get_ae8_ap8_flux (ntime,whichm,whatf,Nene,energy,BBo,L,flux )

 


FLY_IN_AFRL_CRRES

DESCRIPTION:

This function allows one to fly any spacecraft in AFRL CRRESPRO and CRRESELE models.
The output can be differential flux or flux within an energy range or integral flux.

Caution: integral flux for electron are from E to  5.75 MeV and for proton are from E to 81.3 MeV.

To run CRRES models you have to set the full path for the source code directory (where CRRES data files are located).

INPUTS:

ntime: long integer number of time in arrays (max allowed is GET_IRBEM_NTIME_MAX())

sysaxes: long integer to define which coordinate system is provided in

0: GDZ (alti, lati, East longi - km,deg.,deg)
1: GEO (cartesian) - Re
2: GSM (cartesian) - Re
3: GSE (cartesian) - Re
4: SM (cartesian) - Re
5: GEI (cartesian) - Re
6: MAG (cartesian) - Re
7: SPH (geo in spherical) - (radial distance, lati, East longi - Re, deg., deg.)
8: RLL  (radial distance, lati, East longi - Re, deg., deg. - prefered to 7)

whichm: long integer to select in which AFRL CRRES model to fly

1=CRRESPRO QUIET
2=CRRESPRO ACTIVE
3=CRRESELE AVERAGE
4=CRRESELE WORST CASE
5=CRRESELE Ap15

whatf: long integer to select what flux to compute

1=differential flux (MeV-1 cm-2 s-1) at energy(1,*)
2=flux within an ernergy range (MeV-1 cm-2 s-1) - energy(1,*) to energy(2,*)
3=integral flux (cm-2 s-1) at energy(1,*)

Nene: long integer number of energies in array energy (max allowed is 25)

energy: array(2,25) of double (MeV). If whatf=1 or 3 then energy(2,*) is not considered.

iyear: array(GET_IRBEM_NTIME_MAX()) of long integer year when positions are given. Can be a dummy value for the following sysaxes values: 0,1,6,7,8.

idoy: array(GET_IRBEM_NTIME_MAX()) of long integer doy when positions are given. Can be a dummy value for the following sysaxes values: 0,1,6,7,8. 

UT: array(GET_IRBEM_NTIME_MAX()) of double, UT in seconds when positions are given. Can be a dummy value for the following sysaxes values: 0,1,6,7,8. 

x1: array(GET_IRBEM_NTIME_MAX()) of double, first coordinate according to sysaxes. If sysaxes is 0 then altitude has to be in km otherwise use dimensionless variables (in Re)

x2: array(GET_IRBEM_NTIME_MAX()) of double, second coordinate according to sysaxes. If sysaxes is 0 then latitude has to be in degrees otherwise use dimensionless variables (in Re)

x3: array(GET_IRBEM_NTIME_MAX()) of double, third coordinate according to sysaxes. If sysaxes is 0 then longitude has to be in degrees otherwise use dimensionless variables (in Re).

Ap15: array(GET_IRBEM_NTIME_MAX()) of double, preceding 15-day running average of Ap index assuming a one day delay. Array can be set to 0 if whichm non equal to 5.

path: byte array where the number of elements is the length of the string to be converted to byte, provides the path to locate the files which describe CRRES models. Note that the path should end by a "/" or "\" depending on your operating system. Before submitting the path it must me converted to byte array where each element is the ascii code for the corresponding character in the path. 

path_len: long integer, provide here the length of the path string


OUTPUTS:

flux: flux according to selection of whatf for all times and energies (array(GET_IRBEM_NTIME_MAX(),25) of double)

CALLING SEQUENCE FROM MATLAB:

Flux = onera_desp_lib_fly_in_afrl_crres(sysaxes,whichm,energy,matlabd,x1,x2,x3,Ap15,crres_path)

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'fly_in_afrl_crres_', ntime,sysaxes,whichm,whatf,Nene,energy,iyear,idoy, UT,x1,x2,x3,Ap15,flux,path,path_len, /f_value)

CALLING SEQUENCE from FORTRAN:

CALL fly_in_afrl_crres1(ntime,sysaxes,whichm,whatf,Nene,energy,iyear,idoy, UT,x1,x2,x3,Ap15,flux,path,path_len)


GET_CRRES_FLUX

DESCRIPTION:

This function allows one to compute AFRL CRRESPRO and CRRESELE flux for any B/Bo, L position.
The output can be differential flux or flux within an energy range or integral flux.

Caution: integral flux for electron are from E to  5.75 MeV and for proton are from E to 81.3 MeV.

To run CRRES models you have to set the full path for the source code directory (where CRRES data files are located).

INPUTS:

ntime: long integer number of points in arrays (max allowed is GET_IRBEM_NTIME_MAX())

whichm: long integer to select in which AFRL CRRES model to fly

1=CRRESPRO QUIET
2=CRRESPRO ACTIVE
3=CRRESELE AVERAGE
4=CRRESELE WORST CASE
5=CRRESELE Ap15

whatf: long integer to select what flux to compute

1=differential flux (MeV-1 cm-2 s-1) at energy(1,*)
2=flux within an ernergy range (MeV-1 cm-2 s-1) - energy(1,*) to energy(2,*)
3=integral flux (cm-2 s-1) at energy(1,*)

Nene: long integer number of energies in array energy (max allowed is 25)

energy: array(2,25) of double (MeV). If whatf=1 or 3 then energy(2,*) is not considered.

BBo: array(GET_IRBEM_NTIME_MAX()) of double. Provide B/Bequator for all position where to compute the fluxes. Note that the IGRF1985 magnetic field model must be used.

L: array(GET_IRBEM_NTIME_MAX()) of double. Provide McIlwain L for all position where to compute the fluxes.  Note that the IGRF1985 magnetic field model must be used.

Ap15: array(GET_IRBEM_NTIME_MAX()) of double, preceding 15-day running average of Ap index assuming a one day delay. Array can be set to 0 if whichm non equal to 5.

path: byte array where the number of elements is the length of the string to be converted to byte, provides the path to locate the files which describe CRRES models. Note that the path should end by a "/" or "\" depending on your operating system. Before submitting the path it must me converted to byte array where each element is the ascii code for the corresponding character in the path. 

path_len: long integer, provide here the length of the path string


OUTPUTS:

flux: flux according to selection of whatf for all times and energies (array(GET_IRBEM_NTIME_MAX(),25) of double)

CALLING SEQUENCE FROM MATLAB:

Flux = onera_desp_lib_get_crres_flux(whichm,energy,BBo,L,Ap15,crres_path)

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'get_crres_flux_idl_', ntime,whichm,whatf,Nene,energy,BBo,L,Ap15,flux,path,path_len, /f_value)

CALLING SEQUENCE from FORTRAN:

CALL get_crres_flux(ntime,whichm,whatf,Nene,energy,BBo,L,Ap15,flux,path,path_len)


FLY_IN_IGE

DESCRIPTION:

This function allows one to fly any geostationary spacecraft in IGE (International Geostationary Electron) models. The use of the model is limited to geostationary altitude as it is based on LANL-GEO  bird series (1976 to 2006) and JAXA-DRTS spacecraft (added in IGE2006). Three versions of the model are provided:

1: known as POLE-V1, it covers electron energies from 30 keV-1.3 MeV (issued in 2003)
Published in Boscher D., S. Bourdarie, R. Friedel and D. Belian, A model for the geostationary electron environment : POLE, IEEE Trans. Nuc. Sci., 50 (6), 2278-2283, Dec. 2003.

2: known as POLE-V2, it covers electron energies from 30 keV-5.2 MeV (issued in 2005)
Published in Sicard-Piet A., S. Bourdarie, Boscher D. and R. Friedel, A model for the geostationary electron environment : POLE from 30 keV to 5.2 MeV, IEEE Trans. Nuc. Sci., 53 (4), 1844-1850, Aug. 2006.

3: known as IGE-2006, it covers electron energies from 0.9 keV-5.2 MeV (issued in 2006)
Submitted in Sicard-Piet A., S. Bourdarie, Boscher D., R. Friedel M. Thomsen, T. Goka, H. Matsumoto, H. Koshiishi, A new international geostationary model: IGE-2006 from 1 keV to 5.2 MeV, JGR-Space Weather, 2007.

Note that POLE stands for Particle-ONERA-LANL-Electrons

The output can be differential flux or flux within an energy range or integral flux.

Caution: integral flux for electron are from E to Emax of the given selected model. So one should be very carrefull when looking at integral fluxes with POLE-V1 .

INPUTS:

launch_year: year of spacecraft start of life in space - long integer

mission_duration: duration of the mission (number of years) - long integer

whichm: long integer to select in which IGE model to fly

1=POLE-V1
2=POLE-V2
3=IGE-2006

whatf: long integer to select what flux to compute

1=differential flux (MeV-1 cm-2 s-1 sr-1) at energy(1,*)
2=flux within an ernergy range (MeV-1 cm-2 s-1 sr-1) - energy(1,*) to energy(2,*)
3=integral flux (cm-2 s-1 sr-1) at energy(1,*)

Nene: long integer number of energies in array energy (max allowed is 50). If Nene =0 then a default energy grid is being used. The number of default channels is then returned in Nene. 

energy: array(2,50) of double (MeV). If whatf=1 or 3 then energy(2,*) is not considered. If Nene=0 then the default energies are returned (in energy(1,1:Nene) in fortran or in energy(0,0:Nene-1)).

OUTPUTS:

Lower_flux: Lower flux according to selection of whatf for all energies averaged over entire mission duration - This has to be considered as a lower envelop to bound expected flux at GEO for any solar cycle (array(50) of double)

Mean_flux: Mean flux according to selection of whatf for all energies averaged over entire mission duration - This spectrum is an averaged expected flux at GEO for any solar cycle , no margins are included at this point(array(50) of double)

Upper_flux: Upper flux according to selection of whatf for all energies averaged over entire mission duration - This has to be considered as an upper envelop for expected flux at GEO for any solar cycle, this spectrum can be used for any conservative approach as margins are included at this point. Note that the margins are enrgy depended and can be assesed by looking at Upper_flux/Mean_flux (array(50) of double)

Note: all flux are expressed in MeV-1 cm-2 s-1 sr-1 for differential flux or in cm-2 s-1 sr-1 for integrated flux. To derive omnidirectional flux at GEO one should multiply these flux values by 4π.


 

CALLING SEQUENCE FROM MATLAB:

[Lower_flux,Mean_flux,Upper_flux] = onera_desp_lib_fly_in_ige(launch_year,mission_duration,whichm,whatf,energy)

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'fly_in_ige_', launch_year,mission_duration,whichm,whatf,Nene,energy,Lower_flux,Mean_flux,Upper_flux, /f_value)

CALLING SEQUENCE from FORTRAN:

CALL fly_in_ige1(launch_year,mission_duration,whichm,whatf,Nene,energy,Lower_flux,Mean_flux,Upper_flux)



FLY_IN_MEO_GNSS

DESCRIPTION:

This function allows one to fly any MEO GNSS type spacecraft in MEO ONERA models. The use of the model is limited to GPS altitude (~20000 km - 55° inclination) as it is based on LANL-GPS  bird series (1990 to 2006). Two versions of the model are provided:

1: known as MEO-V1, it covers electron energies from 280 keV-1.12 MeV (issued in 2006). Note that in this model there is no solar cycle variations which is to say that the model should be used for long term studies on the order of a solar cycle duration (i.e. 11 years).
Published in Sicard-Piet A., S. Bourdarie,D. Boscher, R. Friedel and T. Cayton, Solar cycle electron environment at GNSS like altitudes, IAC-06-D5.2.04,2006..

2: known as MEO-V2, it covers electron energies from 280 keV-2.24 MeV (issued in 2007)
Not published yet.


The output can be differential flux or flux within an energy range or integral flux.


INPUTS:

launch_year: year of spacecraft start of life in space - long integer

mission_duration: duration of the mission (number of years) - long integer

whichm: long integer to select in which MEO model to fly

1=MEO-V1
2=MEO-V2

whatf: long integer to select what flux to compute

1=differential flux (MeV-1 cm-2 s-1 sr-1) at energy(1,*)
2=flux within an ernergy range (MeV-1 cm-2 s-1 sr-1) - energy(1,*) to energy(2,*)
3=integral flux (cm-2 s-1 sr-1) at energy(1,*)

Nene: long integer number of energies in array energy (max allowed is 50). If Nene =0 then a default energy grid is being used. The number of default channels is then returned in Nene. 

energy: array(2,50) of double (MeV). If whatf=1 or 3 then energy(2,*) is not considered. If Nene=0 then the default energies are returned (in energy(1,1:Nene) in fortran or in energy(0,0:Nene-1)).

OUTPUTS:

Lower_flux: Lower flux according to selection of whatf for all energies averaged over entire mission duration - This has to be considered as a lower envelop to bound expected flux at MEO-GNSS for any solar cycle (array(50) of double)

Mean_flux: Mean flux according to selection of whatf for all energies averaged over entire mission duration - This spectrum is an averaged expected flux at MEO-GNSS for any solar cycle , no margins are included at this point(array(50) of double)

Upper_flux: Upper flux according to selection of whatf for all energies averaged over entire mission duration - This has to be considered as an upper envelop for expected flux at MEO-GNSS for any solar cycle, this spectrum can be used for any conservative approach as margins are included at this point. Note that the margins are enrgy depended and can be assesed by looking at Upper_flux/Mean_flux (array(50) of double)

Note: all flux are expressed in MeV-1 cm-2 s-1 sr-1 for differential flux or in cm-2 s-1 sr-1 for integrated flux. To derive omnidirectional flux at MEO-GNSS one should multiply these flux values by 4π.


 

CALLING SEQUENCE FROM MATLAB:

[Lower_flux,Mean_flux,Upper_flux] = onera_desp_lib_fly_in_meo_gnss(launch_year,mission_duration,whichm,energy)

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'fly_in_meo_gnss_', launch_year,mission_duration,whichm,whatf,Nene,energy,Lower_flux,Mean_flux,Upper_flux, /f_value)

CALLING SEQUENCE from FORTRAN:

CALL fly_in_meo_gnss1(launch_year,mission_duration,whichm,whatf,Nene,energy,Lower_flux,Mean_flux,Upper_flux)



MSIS86

DESCRIPTION:

The Mass-Spectrometer-Incoherent-Scatter-1986 (MSIS-86) neutral atmosphere model describes the neutral temperature and the densities of He, O, N2, O2, Ar, H, and N. The MSIS model is based on the extensive data compilation and analysis work of A. E. Hedin and his collaborators [A. E. Hedin et al., J. Geophys. Res. 82, 2139-2156, 1977; A. E. Hedin, J. Geophys. Res. 88, 10170- 10188, 1983; A. E. Hedin, J. Geophys. Res. 92, 4649, 1987]. MSIS-86 constitutes the upper part of the COSPAR International Reference Atmosphere (CIRA-86).

Data sources for the present model include temperature and density measurements from several rockets, satellites (OGO-6, San Marco 3, Aeros-A, AE-C, AE-D, AE-E, ESRO 4 and DE-2) and incoherent scatter radars (Millstone Hill, St. Santin, Arecibo, Jicamarca, and Malvern). Since the MSIS-83 model, terms were added or changed to better represent seasonal variations in the polar regions under both quiet and magnetically disturbed conditions and local time variations in the magnetic activity effect. In addition a new species, atomic nitrogen, was added to the list of species covered by the model.

INPUTS:

Ntime: long integer number of time in arrays (max allowed is GET_IRBEM_NTIME_MAX())

WhichAp: long integer to select the kind of Ap input:

1=Only daily Ap magnetic index is provided in Ap array (i.e. 1st element)
2=All fields in Ap array are provided (i.e. the 7th elements)

DOY: array(GET_IRBEM_NTIME_MAX()) of long integer, Number of day of year (DOY=1 is for January 1st)

UT: array(GET_IRBEM_NTIME_MAX()) of double, UT in seconds when positions are given..

Alt: array(GET_IRBEM_NTIME_MAX()) of double, Altitude in km (greater than 85 km).

Lat: array(GET_IRBEM_NTIME_MAX()) of double, Geodetic latitude (Deg.)

Long: array(GET_IRBEM_NTIME_MAX()) of double, Geodetic longitude (Deg.)

F107A: array(GET_IRBEM_NTIME_MAX()) of double, 3 month average of F10.7 flux.

F107: array(GET_IRBEM_NTIME_MAX()) of double, daily F10.7 flux for previous day.

Ap: array(7,GET_IRBEM_NTIME_MAX()) of double where:

1st element=Daily Ap
2nd element=3 hours Ap index for current time
3rd element=3 hours Ap index for 3 hours before current time
4th element=3 hours Ap index for 6 hours before current time
5th element=3 hours Ap index for 9 hours before current time
6th element=Average of eight 3 hours Ap indices from 12 to 33 hours prior to current time
7th element=Average of eight 3 hours Ap indices from 36 to 59 hours prior to current time

OUTPUTS:

Dens: array(8,GET_IRBEM_NTIME_MAX()) of double where:

Dens(1st element,*)=He number density (cm-3)
Dens(2nd element,*)=O number density (cm-3)
Dens(3rd element,*)=N2 number density (cm-3)
Dens(4th element,*)=O2 number density (cm-3)
Dens(5th element,*)=AR number density (cm-3)
Dens(6th element,*)=Total mass density (g/cm-3)
Dens(7th element,*)=H number density (cm-3)
Dens(8th element,*)=N number density (cm-3)

Temp: array(2,GET_IRBEM_NTIME_MAX()) of double where:

Temp(1st element,*)=Exospheric temperature (K)
Temp(2nd element,*)=Temperature at alt (K) 

 

CALLING SEQUENCE FROM MATLAB:

out = onera_desp_lib_msis('msis86',date,X,sysaxes,F107A,F107,Ap)

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'msis86_idl_', ntime,whichAp,DOY,UT,Alt,Lat,Lon,F107A,F107,Ap,Dens,Temp, /f_value)

CALLING SEQUENCE from FORTRAN:

CALL msis86(ntime,whichAp,DOY,UT,Alt,Lat,Lon,F107A,F107,Ap,Dens,Temp)



MSISE90

DESCRIPTION:

The MSISE model describes the neutral temperature and densities in Earth's atmosphere from ground to thermospheric heights. Below 72.5 km the model is primarily based on the MAP Handbook (Labitzke et al., 1985) tabulation of zonal average temperature and pressure by Barnett and Corney, which was also used for the CIRA-86. Below 20 km these data were supplemented with averages from the National Meteorological Center (NMC). In addition, pitot tube, falling sphere, and grenade sounder rocket measurements from 1947 to 1972 were taken into consideration. Above 72.5 km MSISE-90 is essentially a revised MSIS-86 model taking into account data derived from space shuttle flights and newer incoherent scatter results. For someone interested only in the thermosphere (above 120 km), the author recommends the MSIS-86 model. MSISE is also not the model of preference for specialized tropospheric work. It is rather for studies that reach across several atmospheric boundaries.

INPUTS:

Ntime: long integer number of time in arrays (max allowed is GET_IRBEM_NTIME_MAX())

WhichAp: long integer to select the kind of Ap input:

1=Only daily Ap magnetic index is provided in Ap array (i.e. 1st element)
2=All fields in Ap array are provided (i.e. the 7th elements)

DOY: array(GET_IRBEM_NTIME_MAX()) of long integer, Number of day of year (DOY=1 is for January 1st)

UT: array(GET_IRBEM_NTIME_MAX()) of double, UT in seconds when positions are given..

Alt: array(GET_IRBEM_NTIME_MAX()) of double, Altitude (km).

Lat: array(GET_IRBEM_NTIME_MAX()) of double, Geodetic latitude (Deg.)

Long: array(GET_IRBEM_NTIME_MAX()) of double, Geodetic longitude (Deg.)

F107A: array(GET_IRBEM_NTIME_MAX()) of double, 3 month average of F10.7 flux.

F107: array(GET_IRBEM_NTIME_MAX()) of double, daily F10.7 flux for previous day.

Ap: array(7,GET_IRBEM_NTIME_MAX()) of double where:

1st element=Daily Ap
2nd element=3 hours Ap index for current time
3rd element=3 hours Ap index for 3 hours before current time
4th element=3 hours Ap index for 6 hours before current time
5th element=3 hours Ap index for 9 hours before current time
6th element=Average of eight 3 hours Ap indices from 12 to 33 hours prior to current time
7th element=Average of eight 3 hours Ap indices from 36 to 59 hours prior to current time

OUTPUTS:

Dens: array(8,GET_IRBEM_NTIME_MAX()) of double where:

Dens(1st element,*)=He number density (cm-3)
Dens(2nd element,*)=O number density (cm-3)
Dens(3rd element,*)=N2 number density (cm-3)
Dens(4th element,*)=O2 number density (cm-3)
Dens(5th element,*)=AR number density (cm-3)
Dens(6th element,*)=Total mass density (g/cm-3)
Dens(7th element,*)=H number density (cm-3)
Dens(8th element,*)=N number density (cm-3)

Temp: array(2,GET_IRBEM_NTIME_MAX()) of double where:

Temp(1st element,*)=Exospheric temperature (K)
Temp(2nd element,*)=Temperature at alt (K) 

 

CALLING SEQUENCE FROM MATLAB:

out = onera_desp_lib_msis('msise90',date,X,sysaxes,F107A,F107,Ap)

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'msise90_idl_', ntime,whichAp,DOY,UT,Alt,Lat,Lon,F107A,F107,Ap,Dens,Temp, /f_value)

CALLING SEQUENCE from FORTRAN:

CALL msise90(ntime,whichAp,DOY,UT,Alt,Lat,Lon,F107A,F107,Ap,Dens,Temp)



NRLMSIS00

DESCRIPTION:

The NRLMSIS-00 empirical atmosphere model was developed by Mike Picone, Alan Hedin, and Doug Drob based on the MSISE90 model. The main differences to MSISE90 are noted in the comments at the top of the computer code. They involve (1) the extensive use of drag and accelerometer data on total mass density, (2) the addition of a component to the total mass density that accounts for possibly significant contributions of O+ and hot oxygen at altitudes above 500 km, and (3) the inclusion of the SMM UV occultation data on [O2]. The MSISE90 model describes the neutral temperature and densities in Earth's atmosphere from ground to thermospheric heights. Below 72.5 km the model is primarily based on the MAP Handbook (Labitzke et al., 1985) tabulation of zonal average temperature and pressure by Barnett and Corney, which was also used for the CIRA-86. Below 20 km these data were supplemented with averages from the National Meteorological Center (NMC). In addition, pitot tube, falling sphere, and grenade sounder rocket measurements from 1947 to 1972 were taken into consideration. Above 72.5 km MSISE-90 is essentially a revised MSIS-86 model taking into account data derived from space shuttle flights and newer incoherent scatter results. For someone interested only in the thermosphere (above 120 km), the author recommends the MSIS-86 model. MSISE is also not the model of preference for specialized tropospheric work. It is rather for studies that reach across several atmospheric boundaries.

INPUTS:

Ntime: long integer number of time in arrays (max allowed is GET_IRBEM_NTIME_MAX())

WhichAp: long integer to select the kind of Ap input:

1=Only daily Ap magnetic index is provided in Ap array (i.e. 1st element)
2=All fields in Ap array are provided (i.e. the 7th elements)

DOY: array(GET_IRBEM_NTIME_MAX()) of long integer, Number of day of year (DOY=1 is for January 1st)

UT: array(GET_IRBEM_NTIME_MAX()) of double, UT in seconds when positions are given..

Alt: array(GET_IRBEM_NTIME_MAX()) of double, Altitude in km (greater than 85 km).

Lat: array(GET_IRBEM_NTIME_MAX()) of double, Geodetic latitude (Deg.)

Long: array(GET_IRBEM_NTIME_MAX()) of double, Geodetic longitude (Deg.)

F107A: array(GET_IRBEM_NTIME_MAX()) of double, 3 month average of F10.7 flux.

F107: array(GET_IRBEM_NTIME_MAX()) of double, daily F10.7 flux for previous day.

Ap: array(7,GET_IRBEM_NTIME_MAX()) of double where:

1st element=Daily Ap
2nd element=3 hours Ap index for current time
3rd element=3 hours Ap index for 3 hours before current time
4th element=3 hours Ap index for 6 hours before current time
5th element=3 hours Ap index for 9 hours before current time
6th element=Average of eight 3 hours Ap indices from 12 to 33 hours prior to current time
7th element=Average of eight 3 hours Ap indices from 36 to 59 hours prior to current time

OUTPUTS:

Dens: array(9,GET_IRBEM_NTIME_MAX()) of double where:

Dens(1st element,*)=He number density (cm-3)
Dens(2nd element,*)=O number density (cm-3)
Dens(3rd element,*)=N2 number density (cm-3)
Dens(4th element,*)=O2 number density (cm-3)
Dens(5th element,*)=AR number density (cm-3)
Dens(6th element,*)=Total mass density (g/cm-3)
Dens(7th element,*)=H number density (cm-3)
Dens(8th element,*)=N number density (cm-3)
Dens(9th element,*)=Anomalous oxygen number density (cm-3)

Temp: array(2,GET_IRBEM_NTIME_MAX()) of double where:

Temp(1st element,*)=Exospheric temperature (K)
Temp(2nd element,*)=Temperature at alt (K) 

 

CALLING SEQUENCE FROM MATLAB:

out = onera_desp_lib_msis('nrlmsise00',date,X,sysaxes,F107A,F107,Ap) 

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'nrlmsise00_idl_', ntime,whichAp,DOY,UT,Alt,Lat,Lon,F107A,F107,Ap,Dens,Temp, /f_value)

CALLING SEQUENCE from FORTRAN:

CALL nrlmsise00(ntime,whichAp,DOY,UT,Alt,Lat,Lon,F107A,F107,Ap,Dens,Temp)



SGP4_TLE

DESCRIPTION:

SGP4 stands for Simplified General Perturbation -4 and consists in an orbit propagator. This function allows one to propagate NORAD two lines elements (TLE sets can be found at http://celestrak.com/). More about SGP4 can be found at: http://celestrak.com/NORAD/documentation/spacetrk.pdf and at http://www.centerforspace.com/downloads/files/pubs/AIAA 2006-6753 Revisiting Spacetrack Report 3.pdf

Important notice: those who are running the library on Unix or linux system have to convert TLE files from NORAD from DOS to UNIX using the command : dos2unix file.tle newfile.tle . Also be aware that some TLE files contains random errors, like elements from another satellite, repeated elements, ... there are no specific filters implemented in the SGP4 function.

For those who are familiar with SGP4 code, be aware for one difference between orginal SGP4 code and the one provided here: input start and stop time and time step are not in minutes but in seconds. This was chosen for convenience.

INPUTS:

runtype: long integer to select in which mode SGP4_TLE will run 

0: defines start and stop time to propagate each TLE according to input file 
1: propagates each TLE according to user start and stop time 

startsfe: double giving start time in seconds from date provided in each TLE. This number can be negative. Note that this value is not use if runtype is equal to 0.

stopsfe: double giving stop time in seconds from date provided in each TLE. This number can be negative. Note that this value is not use if runtype is equal to 0.

deltasec: double giving step time in seconds to propagateTLE.

InFileByte: byte array where the number of elements is the length of the string to be converted to byte, provides the path and name to locate the input TLE to be propagated. Before submitting the path it must me converted to byte array where each element is the ascii code for the corresponding character in the path. 

strlenIn: long integer providing the length of InFileByte string

OutFileByte: byte array where the number of elements is the length of the string to be converted to byte, provides the path and name to locate the output file. Before submitting the path it must me converted to byte array where each element is the ascii code for the corresponding character in the path. 

strlenOut: long integer providing here the length of OutFileByte string 
 

OUTPUTS:

None: They are provided in OutFileByte. This file is composed of 6 columns: 

1: date (dd/mm/yyyy) 
2: time (hh:mm:ss)
3: decimal year
4: altitude (km)
5: latitude (deg)
6: East longitude (deg) 

CALLING SEQUENCE FROM MATLAB:

The Matlab wrapper handles calculation of runtype strlenIn, OutfileByte, strlenOut. OutfileByte and strlenOut will be set appropriately for a temporary file (onera_desp_lib_sgp4_tle.tmp or onera_desp_lib_sgp4_tle.tmp.####). The wrapper reads the library routine's output from the temporary file and passes it back as a Matlab structure. The temporary file is deleted, so the Matlab structure is the only result of a successful call to sgp4_tle. In order to avoid a memory overflow for very long TLE runs, it is possible to leave the output in a text file (as in the FORTRAN call): the syntax is then onera_desp_lib_sgp4_tle(InFile,OutFile)

pos = onera_desp_lib_sgp4_tle(InFileByte) %implies runtype=0

pos = onera_desp_lib_sgp4_tle(startsfe,stopsfe,deltasec,InFileByte) %implies runtype=1

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'sgp4_tle_', runtype,startsfe,stopsfe,deltasec,InFileByte,strlenIn,OutFileByte,strlenOut, /f_value)

CALLING SEQUENCE from FORTRAN:

CALL sgp4_tle1(runtype,startsfe,stopsfe,deltasec,InFileByte,strlenIn,OutFileByte,strlenOut)

 


SGP4_ELE

DESCRIPTION:

SGP4 stands for Simplified General Perturbation -4 and consists in an orbit propagator. This function allows one to produce orbit coordinates from different types of orbital elements... More about SGP4 can be found at: http://celestrak.com/NORAD/documentation/spacetrk.pdf and at http://www.centerforspace.com/downloads/files/pubs/AIAA 2006-6753 Revisiting Spacetrack Report 3.pdf

For those who are familiar with SGP4 code, be aware for one difference between orginal SGP4 code and the one provided here: input start and stop time and time step are not in minutes but in seconds. This was chosen for convenience.

INPUTS:

sysaxesOUT: long integer to define which coordinate system is provided out

0: GDZ (alti, lati, East longi - km,deg.,deg)
1: GEO (cartesian) - Re
2: GSM (cartesian) - Re
3: GSE (cartesian) - Re
4: SM (cartesian) - Re
5: GEI (cartesian) - Re
6: MAG (cartesian) - Re
7: SPH (geo in spherical) - (radial distance, lati, East longi - Re, deg., deg.)
8: RLL  (radial distance, lati, East longi - Re, deg., deg. - prefered to 7)

year: long integer giving start year of start date

month: long integer giving start month of start date

day: long integer giving start day of month of start date.

hour: long integer giving start hour of start time.

minute: long integer giving start minute of start time.

sec: double giving start seconds of start time.

e1: double see detailed definition according to options array.

e2: double see detailed definition according to options array.

e2: double see detailed definition according to options array.

e4: double see detailed definition according to options array.

e5: double see detailed definition according to options array.

e6: double see detailed definition according to options array.

options: array(5) of long integer to define which type of element is provided in
Case options(1st element) of:

1: ONERA type elements: 
e1 = Inclination - deg
e2 = geocentric altitude of perigee - km
e3 = geocentric altitude of apogee - km
e4 = longitude of the ascending node - deg
Case options(2nd element) of:
1: e5 = argument of perigee - deg
2: e5 = longitude of perigee - deg
Case options(3rd element) of:
1: e6 = time since perigee passage - sec
2: e6 = true anomaly at epoch - deg
3: e6 = argument of latitude at epoch - deg
4: e6 = true longitude at epoch - deg
5: e6 = mean anomaly at epoch - deg

2: Classical type elements: 
e1 = semimajor axis - Re
e2 = eccentricity
e3 = inclination - deg
e4 = longitude of ascending node - deg
Case options(2nd element) of:
1: e5 = argument of perigee - deg
2: e5 = longitude of perigee - deg
Case options(3rd element) of:
1: e6 = time since perigee passage - sec
2: e6 = true anomaly at epoch - deg
3: e6 = argument of latitude at epoch - deg
4: e6 = true longitude at epoch - deg
5: e6 = mean anomaly at epoch - deg

3: RV type elements: 
e1 = xGEI - km
e2 = yGEI - km
e3 = zGEI - km
e4 = VxGEI - km/s
e5 = VyGEI - km/s
e6 = VzGEI - km/s

4: SOLAR type elements: 
e1 = inclination - deg
e2 = geocentric altitude of perigee - km
e3 = geocentric altitude of apogee - km
e4 = local time of apogee - hours
e5 = local time of maximum inclination - hours
Case options(3rd element) of:
1: e6 = time since perigee passage - sec
2: e6 = true anomaly at epoch - deg
3: e6 = argument of latitude at epoch - deg
4: e6 = true longitude at epoch - deg
5: e6 = mean anomaly at epoch - deg

5: MEAN type elements: 
e1 = mean motion - rev/day
e2 = eccentricity
e3 = inclination - deg
e4 = longitude of the ascending node - deg
Case options(2nd element) of:
1: e5 = argument of perigee - deg
2: e5 = longitude of perigee - deg
Case options(3rd element) of:
1: e6 = time since perigee passage - sec
2: e6 = true anomaly at epoch - deg
3: e6 = argument of latitude at epoch - deg
4: e6 = true longitude at epoch - deg
5: e6 = mean anomaly at epoch - deg

startsfe: double giving start time in seconds from date provided before. This number can be negative.

stopsfe: double giving stop time in seconds from date provided before. This number can be negative.

deltasec: double giving step time in seconds to produce orbit outputs.

OUTPUTS:

OUTyear: array of long integer(GET_IRBEM_NTIME_MAX()) giving output year for each orbital locations

OUTdoy: array of long integer(GET_IRBEM_NTIME_MAX()) giving output day of year for each orbital locations

UT: array of double(GET_IRBEM_NTIME_MAX()) giving output UT (seconds) for each orbital locations

X1: array of double(GET_IRBEM_NTIME_MAX()) giving first coordinate of orbit according to sysaxes

X2: array of double(GET_IRBEM_NTIME_MAX()) giving second coordinate of orbit according to sysaxes

X3: array of double(GET_IRBEM_NTIME_MAX()) giving third coordinate of orbit according to sysaxes

CALLING SEQUENCE FROM MATLAB:

pos = onera_desp_lib_sgp4_ele([e1,e2,e3,e4,e5,e6],startdate,enddate,deltasec,sysaxesOUT) 

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'sgp4_ele_', sysaxesOUT,year,month,day,hour,minute,sec, e1,e2,e3,e4,e5,e6,options,startsfe,stopsfe,deltasec,OUTyear,OUTdoy,UT,X1,X2,X3, /f_value)

CALLING SEQUENCE from FORTRAN:

CALL sgp4_ele1(sysaxesOUT,year,month,day,hour,minute,sec, e1,e2,e3,e4,e5,e6,options,startsfe,stopsfe,deltasec,OUTyear,OUTdoy,UT,X1,X2,X3)


 

RV2COE

DESCRIPTION:

This function finds the classical orbital elements given the Geocentric Equatorial Position and Velocity vectors. It comes from SGP4 distribution. SGP4 stands for Simplified General Perturbation -4 and consists in an orbit propagator. More about SGP4 can be found at: http://celestrak.com/NORAD/documentation/spacetrk.pdf and at http://www.centerforspace.com/downloads/files/pubs/AIAA 2006-6753 Revisiting Spacetrack Report 3.pdf

INPUTS:

R: array(3) of double. Position vector in GEI (km) 

V: array(3) of double. Velocity vector in GEI (km/s).

OUTPUTS:

P: double. SemiLatus rectum (km) 

A: double. Semimajor axis (km) 

Ecc: double. Eccentricity 

Incl: double. Inclination (rad) 

Omega: double. Longitude of ascending node (rad)

Argp: double. Argument of perigee (rad) 

Nu: double. True anomaly (rad) 

M: double. Mean anomaly (rad) 

ArgLat: double. Argument of latitude (rad)

LamTrue: double. True longitude (rad) 

LonPer: double. Longitude of Periapsis (rad) 

CALLING SEQUENCE FROM MATLAB:

elements=onera_desp_lib_rv2coe(R,V) 

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'rv2coe_idl_', R,V,P,A,Ecc,Incl,Omega,Argp,Nu,M,ArgLat,LamTrue,LonPer, /f_value)

CALLING SEQUENCE from FORTRAN:

CALL rv2coe(R,V,P,A,Ecc,Incl,Omega,Argp,Nu,M,ArgLat,LamTrue,LonPer)


 



IRBEM_FORTRAN_VERSION

DESCRIPTION:

Provides the repository version number of the fortran source code.
 

INPUTS:

NONE

OUTPUTS:

VERSN Repository version number (long integer).
 

CALLING SEQUENCE FROM MATLAB:

VERSN = onera_desp_lib_fortran_version

CALLING SEQUENCE FROM IDL:

result = call_external(lib_name, 'irbem_fortran_version_',versn, /f_value)

CALLING SEQUENCE from FORTRAN:

call IRBEM_FORTRAN_VERSION1(VERSN)



IRBEM_FORTRAN_RELEASE

DESCRIPTION:

Provides the repository release tag of the fortran source code.
 

INPUTS:

NONE

OUTPUTS:

RLS Repository release tag (character(80)).
 

CALLING SEQUENCE FROM MATLAB:

RLS = onera_desp_lib_fortran_release

CALLING SEQUENCE FROM IDL:

result = call_external(lib_name, 'irbem_fortran_release_', rls, /f_value)
Tip: In IDL, string arguments should be defined as byte arrays, i.e.:    rls = BYTARR(80)

CALLING SEQUENCE from FORTRAN:

call IRBEM_FORTRAN_RELEASE1(RLS)



GET_IRBEM_NTIME_MAX

DESCRIPTION:

Returns size of time dimension in preallocated fortran arrays.
 

INPUTS:

NONE

OUTPUTS:

NTIME_MAXTime dimension size (long integer).
 

CALLING SEQUENCE FROM MATLAB:

ntime_max = onera_desp_lib_ntime_max

CALLING SEQUENCE FROM IDL:

result = call_external(lib_name, 'get_irbem_ntime_max_',ntime_max, /f_value)

CALLING SEQUENCE from FORTRAN:

call GET_IRBEM_NTIME_MAX1(ntime_max1)
 

SHIELDOSE2

DESCRIPTION:

SHIELDOSE2 [Seltzer, 1994] is a computer code for space-shielding radiation dose calculations. It determines the absorbed dose as a function of depth in aluminium shielding material of spacecraft, given the electron and proton fluences encountered in orbit. The code makes use of precalculated, mono-energetic depth-dose data for an isotropic, broad-beam fluence of radiation incident on uniform aluminium plane media. Such data are particularly suitable for routine dose predictions in situations where the geometrical and compositional complexities of the spacecraft are not known. Furthermore, the restriction to these rather simple geometries has allowed for the development of accurate electron and electron-bremsstrahlung data sets based on detailed transport calculations rather than on more approximate methods.
SHIELDOSE2 calculates, for arbitrary proton and electron incident spectra, the dose absorbed in small volumes of different detector materials for the following aluminium shield geometries:

  1. in a semi-infinite plane medium, as a function of depth; irradiation is from one side only (the assumed infinite backing effectively insures this).
  2. at the transmission surface of a plane slab, as a function of slab thickness; irradiation is from one side only.
  3. at the centre of a solid sphere, as a function of sphere radius; irradiation is from all directions.

INPUTS:

IDET: long integer to set the index of the detector type:

  1. Al detector
  2. Graphite detector
  3. Si detector
  4. Air detector
  5. Bone detector
  6. Calcium Fluoride detector
  7. Gallium Arsenide detector
  8. Lithium Fluoride detector
  9. Silicon Dioxide detector
  10. Tissue detector
  11. Water detector
INUC: long integer to indicate the nuclear interaction to account for:
  1. No nuclear attenuation for protons in Al
  2. Nuclear attenuation, local charged-secondary energy deposition
  3. Nuclear attenuation, local charged-secondary energy deposition, and approx exponential distribution of neutron dose
IMAX: long integer to set the number of shielding depth (max allowed=71). It is recommended to set this number close to maximum allowed value to compute accurate doses in semi-infinite aluminium medium and at center of aluminium spheres.

IUNT: long integer to set the shielding depth unit
  1. Mils
  2. g/cm2
  3. mm
Zin: double array(71) of thickness in unit of IUNT. Below is provided a reasonable thickness array in g/cm2 with IMAX=70
data Zin /1.000E-06,2.000E-06,5.000E-06,1.000E-05,2.000E-05,5.000E-05,1.000E-04,2.000E-04,&
5.000E-04,1.000E-03,1.000E-01,2.000E-01,5.000E-01,7.000E-01,1.000E+00,1.250E+00,&
1.500E+00,1.750E+00,2.000E+00,2.500E+00,3.000E+00,3.500E+00,4.000E+00,4.500E+00,&
5.000E+00,6.000E+00,7.000E+00,8.000E+00,9.000E+00,1.000E+01,1.100E+01,1.200E+01,&
1.300E+01,1.400E+01,1.500E+01,1.600E+01,1.700E+01,1.800E+01,1.900E+01,2.000E+01,&
2.100E+01,2.200E+01,2.300E+01,2.400E+01,2.500E+01,2.600E+01,2.700E+01,2.800E+01,&
2.900E+01,3.000E+01,3.100E+01,3.200E+01,3.300E+01,3.400E+01,3.500E+01,3.600E+01,&
3.700E+01,3.800E+01,3.900E+01,4.000E+01,4.100E+01,4.200E+01,4.300E+01,4.400E+01,&
4.500E+01,4.600E+01,4.700E+01,4.800E+01,4.900E+01,5.000E+01,0.000E+00/


EminS: min energy of solar protons spectrum (double) [MeV]

EmaxS: max energy of solar protons spectrum (double) [MeV]

EminP: min energy of trapped protons spectrum (double) [MeV]

EmaxP: max energy of trapped protons spectrum (double) [MeV]

NPTSP: number of spectrum points which divides proton spectra for integration. A value of 1001 is recommended (long integer).

EminE: min energy of trapped electrons spectrum (double) [MeV]

EmaxE: max energy of trapped electrons spectrum (double) [MeV]

NPTSE: Number of spectrum points which divides electron spectra for integration. A value of 1001 is recommended (long integer).

JSMAX: Number of points in falling spectrum of solar protons, max allowed=301 (long integer).

JPMAX: Number of points in falling spectrum of trapped protons, max allowed=301 (long integer).

JEMAX: Number of points in falling spectrum of trapped electrons, max allowed=301 (long integer).

EUNIT: Conversion factor from /energy to /MeV; e.g. EUNIT = 1000 if flux is /keV (double).

DURATN: Mission duration in multiples of unit time [s] (double).

ESin: Energy array(301) of solar proton spectrum [MeV] (double).

SFLUXin: Solar flare flux array(301) for solar protons [in /energy/cm2] (double).

EPin: Energy array(301) of trapped proton spectrum [MeV] (double).

PFLUXin: Incident omnidirectional flux array(301) for trapped protons [in /energy/cm2/s] (double).

EEin: Energy array(301) of trapped electron spectrum [MeV] (double).

EFLUXin: Incident omnidirectional flux array(301) for trapped electrons [in /energy/cm2/s] (double).

OUTPUTS:

SolDose: Dose profile array (71,3) for solar protons [rads] (double)

  1. Dose in semi-infinite aluminium medium
  2. Dose at transmission surface of finite aluminium slab shields
  3. 1/2 Dose at center of aluminium spheres

ProtDose: Dose profile array (71,3) for trapped protons [rads] (double)

  1. Dose in semi-infinite aluminium medium
  2. Dose at transmission surface of finite aluminium slab shields
  3. 1/2 Dose at center of aluminium spheres

ElecDose: Dose profile array (71,3) for trapped electrons [rads] (double)

  1. Dose in semi-infinite aluminium medium
  2. Dose at transmission surface of finite aluminium slab shields
  3. 1/2 Dose at center of aluminium spheres

BremDose: Dose profile array (71,3) for Bremsstrahlung [rads] (double)

  1. Dose in semi-infinite aluminium medium
  2. Dose at transmission surface of finite aluminium slab shields
  3. 1/2 Dose at center of aluminium spheres

TotDose: Total dose profile array (71,3) [rads] (double)

  1. Dose in semi-infinite aluminium medium
  2. Dose at transmission surface of finite aluminium slab shields
  3. 1/2 Dose at center of aluminium spheres

CALLING SEQUENCE FROM MATLAB:

[ProtDose,ElecDose,BremDose,SolDose,TotDose] = onera_desp_lib_shieldose2(ProtSpect,ElecSpect,SolSpect,Target,...)  

CALLING SEQUENCE from IDL:

result = call_external(lib_name, 'shieldose2idl_', IDET,INUC,IMAX,IUNT,Zin,EMINS,EMAXS,EMINP,EMAXP,NPTSP,EMINE,EMAXE,NPTSE,JSMAX,JPMAX,JEMAX,EUNIT,DURATN,ESin,SFLUXin,EPin,PFLUXin,EEin,EFLUXin,SolDose,ProtDose,ElecDose,BremDose,TotDose,/f_value)

CALLING SEQUENCE from FORTRAN:

CALL shieldose2(IDET,INUC,IMAX,IUNT,Zin,EMINS,EMAXS,EMINP,EMAXP,NPTSP,EMINE,EMAXE,NPTSE,JSMAX,JPMAX,JEMAX,EUNIT,DURATN,ESin,SFLUXin,EPin,PFLUXin,EEin,EFLUXin,SolDose,ProtDose,ElecDose,BremDose,TotDose)