spacepy.irbempy¶
Wrapper for the fortran library irbem_lib
Reference: https://github.com/PRBEM/IRBEM
Boscher, S. Bourdarie, T. P. O’Brien, T. Guild, IRBEM library V4.3, 2004-2008
Most functions in this module use an options list to define the models used and the settings that define the quality level of the result. The options list is a 5-element list and is defined as follows.
Options¶
1st element: 0 - don’t compute L* or phi ; 1 - compute L*; 2- compute phi
- 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 …)
- 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 (theta) along the field line such as dtheta=(2pi)/(720*[options(3rd element)+1])
- 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 (phi) along the drift shell such as dphi=(2pi)/(25*[options(4th element)+1])
- 5th element: allows to select an internal magnetic field model (default is set to IGRF)
0 = IGRF
1 = Eccentric tilted dipole
2 = Jensen&Cain 1960
3 = GSFC 12/66 updated to 1970
4 = User-defined model (Default: Centred dipole + uniform [Dungey open model] )
5 = Centred dipole
The routines also require specification of the external magnetic field model. The default is the Tsyganenko 2001 storm-time model. The external model is always specified using the extMag keyword and the following options exist.
extMag¶
‘0’ = No external field model
‘MEAD’ = Mead and Fairfield
‘T87SHORT’ = Tsyganenko 1987 short (inner magnetosphere)
‘T87LONG’ = Tsyganenko 1987 long (valid in extended tail region)
‘T89’ = Tsyganenko 1989
‘OPQUIET’ = Olsen-Pfitzer static model for quiet conditions
‘OPDYN’ = Olsen-Pfitzer static model for active conditions
‘T96’ = Tsyganenko 1996
‘OSTA’ = Ostapenko and Maltsev
‘T01QUIET’ = Tsyganenko 2001 model for quiet conditions
‘T01STORM’ = Tsyganenko 2001 model for active conditions
‘T05’ = Tsyganenko and Sitnov 2005 model
‘ALEX’ = Alexeev model
‘TS07’ = Tsyganenko and Sitnov 2007 model
Many of these models have limits placed on the valid range of input parameters, and outside these limits invalid (NaN) values will be returned.
MEAD : Mead & Fairfield [1975] (uses 0<=Kp<=9 - Valid for rGEO<=17. Re)
T87SHORT: Tsyganenko short [1987] (uses 0<=Kp<=9 - Valid for rGEO<=30. Re)
T87LONG : Tsyganenko long [1987] (uses 0<=Kp<=9 - Valid for rGEO<=70. Re)
T89 : Tsyganenko [1989] (uses 0<=Kp<=9 - Valid for rGEO<=70. Re)
OPQUIET : Olson & Pfitzer quiet [1977] (default - Valid for rGEO<=15. Re)
- OPDYNOlson & Pfitzer dynamic [1988] (uses 5.<=dens<=50., 300.<=velo<=500.,
-100.<=Dst<=20. - Valid for rGEO<=60. Re)
- T96Tsyganenko [1996] (uses -100.<=Dst (nT)<=20., 0.5<=Pdyn (nPa)<10.,
|ByIMF| (nT)<=10., |BzIMF| (nT)<=10. - Valid for rGEO<=40. Re)
- OSTAOstapenko & Maltsev [1997] (uses dst,Pdyn,BzIMF, Kp)
T01QUIET: Tsyganenko [2002a,b] (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)
- T01STORM: Tsyganenko, Singer & Kasper [2003] storm (uses Dst, Pdyn, ByIMF, BzIMF, G2, G3 -
there is no upper or lower limit for those inputs - Valid for xGSM>=-15. Re)
- T05Tsyganenko & Sitnov [2005] storm (uses Dst, Pdyn, ByIMF, BzIMF,
W1, W2, W3, W4, W5, W6 - no upper or lower limit for inputs - Valid for xGSM>=-15. Re)
TS07 : Tsyganenko and Sitnov [2007] model. Uses specially calculated coefficient files.
Functions
|
Find the equatorial pitch angle corresponding to a given second invariant K. |
|
thin layer to call coor_trans1 from irbem lib this will convert between systems GDZ, GEO, GSM, GSE, SM, GEI, MAG, SPH, RLL |
|
call find_mirror_point from irbem library and return a dictionary with values for Blocal, Bmirr and the GEO (cartesian) coordinates of the mirror point |
|
Find the last closed drift shell (LCDS) for a given equatorial pitch angle. |
|
Find the last closed drift shell (LCDS) for a given K (modified 2nd invariant). |
|
call find_foot_point1 from irbem library and return a dictionary with values for Bmin and the GEO (cartesian) coordinates of the magnetic equator |
|
call find_magequator from irbem library and return a dictionary with values for Bmin and the GEO (cartesian) coordinates of the magnetic equator |
|
will return the flux from the AE8-AP8 model |
|
Return magnetic field vector (in GEO) and magnitude |
|
Return the MacIlwain L value for a given location, time and model |
|
This will call make_lstar1 or make_lstar_shell_splitting_1 from the irbem library and will lookup omni values for given time if not provided (optional). |
|
will return the coordinate system type as string |
|
Converts |
|
Prepare inputs for direct IRBEM-LIB calls. |
|
Update coefficients for TS07 magnetic field model |
Classes
|
A class for performing dose calculations using Shieldose2 |
- spacepy.irbempy.AlphaOfK(ticks, loci, K, extMag='T01STORM', options=[0, 0, 3, 0, 0], omnivals=None)[source]¶
Find the equatorial pitch angle corresponding to a given second invariant K.
Uses the IRBEM library to determine K and searches via bisection to find Alpha(K).
- Parameters:
- - ticks (Ticktock class)containing time information for a single time
- - loci (Coords class)containing position information for a single point
- - K (numeric)value of second invariant (K) for search
- - extMag (string)optional; will choose the external magnetic field model
possible values [‘0’, ‘MEAD’, ‘T87SHORT’, ‘T87LONG’, ‘T89’, ‘OPQUIET’, ‘OPDYN’, ‘T96’, ‘OSTA’, ‘T01QUIET’, ‘T01STORM’, ‘T05’, ‘ALEX’, ‘TS07’]
- - options (optional list or array of integers length=5)explained in Lstar
- - omni values as dictionary (optional)if not provided, will use lookup table
- - (see Lstar documentation for further explanation)
- - tol (float)tolerance for search in radial distance
- - bracket (list): X-GSM coordinates to bracket bisection search
- Returns:
- AlphaEqequatorial pitch angle corresponding to K
See also
Examples
>>> t = spacepy.time.Ticktock(['2002-09-01T04:00:00'], 'ISO') >>> loci = spacepy.coordinates.Coords([-4,0,0], 'GSM', 'car', use_irbem=True) >>> spacepy.irbempy.AlphaOfK(t, loci, 0.11, extMag='T89') 48.984375
- spacepy.irbempy.coord_trans(loci, returntype, returncarsph)[source]¶
thin layer to call coor_trans1 from irbem lib this will convert between systems GDZ, GEO, GSM, GSE, SM, GEI, MAG, SPH, RLL
- Parameters:
- - loci (Coords instance)containing coordinate information, can contain n points
- - returntype (str)describing system as GDZ, GEO, GSM, GSE, SM, GEI, MAG, SPH, RLL
- - returncarsph (str)cartesian or spherical units ‘car’, ‘sph’
- Returns:
- - xout (ndarray)values after transformation in (n,3) dimensions
Examples
>>> c = Coords([[3,0,0],[2,0,0]], 'GEO', 'car', use_irbem=True) >>> c.ticks = Ticktock(['2002-02-02T12:00:00', '2002-02-02T12:10:00'], 'ISO') >>> coord_trans(c, 'GSM', 'car') array([[ 2.8639301 , -0.01848784, 0.89306361], [ 1.9124434 , 0.07209424, 0.58082929]])
- spacepy.irbempy.find_Bmirror(ticks, loci, alpha, extMag='T01STORM', options=[1, 0, 0, 0, 0], omnivals=None)[source]¶
call find_mirror_point from irbem library and return a dictionary with values for Blocal, Bmirr and the GEO (cartesian) coordinates of the mirror point
- Parameters:
- ticksTicktock class
containing time information
- lociCoords class
containing spatial information
- alphaarray-like
containing the pitch angles
- extMagstr
optional; will choose the external magnetic field model possible values [‘0’, ‘MEAD’, ‘T87SHORT’, ‘T87LONG’, ‘T89’, OPQUIET’, ‘OPDYN’, ‘T96’, ‘OSTA’, ‘T01QUIET’, ‘T01STORM’, ‘T05’, ‘ALEX’, ‘TS07’]
- optionsarray-like (optional)
length=5 : explained in Lstar
- omnivalsdict (optional)
if not provided, will use lookup table (see get_Lstar documentation for further explanation)
- Returns:
- resultsdictionary
containing keys: Blocal, Bmirr, GEOcar
See also
Notes
The resulting Bmirr is a one-dimensional array in time if alpha is a singleton, otherwise it is two-dimensional
Examples
>>> t = Ticktock(['2002-02-02T12:00:00', '2002-02-02T12:10:00'], 'ISO') >>> y = Coords([[3,0,0],[2,0,0]], 'GEO', 'car', use_irbem=True) >>> ib.find_Bmirror(t,y,[90,80,60,10]) {'Blocal': array([978.62167212, 3399.44983792]), 'Bmirr': array([[978.62167212, 1009.04166504, 1304.81635975, 32451.54718886], [3399.44983792, 3505.1329238, 4532.49001884, nan]]), 'loci': Coords([[0.6877979690266218, -0.11552153465118173, 0.9162488420030738], [nan, nan, nan]], 'GEO', 'car')}
- spacepy.irbempy.find_LCDS(ticks, alpha, extMag='T01STORM', options=[1, 0, 0, 0, 0], omnivals=None, tol=0.05, bracket=[3, 12], mlt=0, **kwargs)[source]¶
Find the last closed drift shell (LCDS) for a given equatorial pitch angle.
Uses the IRBEM library to determine L* and searches via bisection to find LCDS to a given tolerance in R (radial distance along GSM equator at local midnight). Note that this function is present to aid in reproducibility of older work, however it should be noted that drift shells are defined by constant K, not constant equatorial pitch angle. Therefore find_LCDS_K should be used for obtaining the LCDS.
- Parameters:
- - ticks (Ticktock class)containing time information for a single time
- - alpha (numeric)equatorial pitch angle for search
- - extMag (string)optional; will choose the external magnetic field model
possible values [‘0’, ‘MEAD’, ‘T87SHORT’, ‘T87LONG’, ‘T89’, ‘OPQUIET’, ‘OPDYN’, ‘T96’, ‘OSTA’, ‘T01QUIET’, ‘T01STORM’, ‘T05’, ‘ALEX’, ‘TS07’]
- - options (optional list or array of integers length=5)explained in Lstar
- - omni values as dictionary (optional)if not provided, will use lookup table
- - (see Lstar documentation for further explanation)
- - tol (float)tolerance for search in radial distance
- - bracket (list): X-GSM coordinates to bracket bisection search
- Returns:
- results (SpaceData, dictionary-like): contains keys LCDS, K, AlphaEq and UTC
See also
Examples
>>> t = spacepy.time.Ticktock(['2002-02-02T12:00:00'], 'ISO') >>> spacepy.irbempy.find_LCDS(t, 90, extMag='T89')
- spacepy.irbempy.find_LCDS_K(ticks, K, extMag='T01STORM', options=[1, 1, 3, 0, 0], omnivals=None, tol=0.05, bracket=[3, 12], mlt=0, **kwargs)[source]¶
Find the last closed drift shell (LCDS) for a given K (modified 2nd invariant).
Uses the IRBEM library to determine L* and searches via bisection to find LCDS to a given tolerance in R (radial distance along GSM equator at local midnight).
- Parameters:
- - ticks (Ticktock class)containing time information for a single time
- - K (numeric)K (modified 2nd adiabatic invariant) for search
- - extMag (string)optional; will choose the external magnetic field model
possible values [‘0’, ‘MEAD’, ‘T87SHORT’, ‘T87LONG’, ‘T89’, ‘OPQUIET’, ‘OPDYN’, ‘T96’, ‘OSTA’, ‘T01QUIET’, ‘T01STORM’, ‘T05’, ‘ALEX’, ‘TS07’]
- - options (optional list or array of integers length=5)explained in Lstar
- - omni values as dictionary (optional)if not provided, will use lookup table
- - (see Lstar documentation for further explanation)
- - tol (float)tolerance for search in radial distance
- - bracket (list): X-GSM coordinates to bracket bisection search
- Returns:
- results (SpaceData, dictionary-like): contains keys LCDS, K, AlphaEq and UTC
See also
Examples
>>> t = spacepy.time.Ticktock(['2002-02-02T12:00:00'], 'ISO') >>> spacepy.irbempy.find_LCDS_K(t, 0.2, extMag='T89')
- spacepy.irbempy.find_footpoint(ticks, loci, extMag='T01STORM', options=[1, 0, 3, 0, 0], hemi='same', alt=100, omnivals=None)[source]¶
call find_foot_point1 from irbem library and return a dictionary with values for Bmin and the GEO (cartesian) coordinates of the magnetic equator
- Parameters:
- - ticks (Ticktock class)containing time information
- - loci (Coords class)containing spatial information
- - extMag (string)optional; will choose the external magnetic field model
possible values [‘0’, ‘MEAD’, ‘T87SHORT’, ‘T87LONG’, ‘T89’, ‘OPQUIET’, ‘OPDYN’, ‘T96’, ‘OSTA’, ‘T01QUIET’, ‘T01STORM’, ‘T05’, ‘ALEX’, ‘TS07’]
- - options (optional list or array of integers length=5)explained in Lstar
- - omni values as dictionary (optional)if not provided, will use lookup table
- - (see Lstar documentation for further explanation)
- - hemi (string)optional (valid cases are ‘same’, ‘other’, ‘north’ or ‘south’)
will set the target hemisphere for tracing the footpoint
- - alt (numeric)optional keyword to set stop height [km] of fieldline trace (default 100km)
- Returns:
- - results (spacepy.datamodel.SpaceData)containing keys
Bfoot - Magnitude of B-field at footpoint [nT] loci - Coords instance with GDZ coordinates of the magnetic footpoint [alt, lat, lon] Bfootvec - Components of B-field at footpoint in cartesian GEO coordinates [nT]
See also
Examples
>>> t = Ticktock(['2002-02-02T12:00:00', '2002-02-02T12:10:00'], 'ISO') >>> y = Coords([[3,0,0],[3,0,0]], 'GEO', 'car', use_irbem=True) >>> spacepy.irbempy.find_footpoint(t, y) {'Bfoot': array([ 47559.04643444, 47542.84688657]), 'Bfootvec': array([[-38428.07217246, 4497.31549786, -27657.19291928], [-38419.08514332, 4501.45390964, -27641.14866517]]), 'loci': Coords( [[ 99.31443778 55.71415787 -10.21888955] [ 99.99397026 55.70716296 -10.22797462]] ), dtype=GDZ,sph, units=['km', 'deg', 'deg']}
- spacepy.irbempy.find_magequator(ticks, loci, extMag='T01STORM', options=[1, 0, 0, 0, 0], omnivals=None)[source]¶
call find_magequator from irbem library and return a dictionary with values for Bmin and the GEO (cartesian) coordinates of the magnetic equator
- Parameters:
- - ticks (Ticktock class)containing time information
- - loci (Coords class)containing spatial information
- - extMag (string)optional; will choose the external magnetic field model
possible values [‘0’, ‘MEAD’, ‘T87SHORT’, ‘T87LONG’, ‘T89’, ‘OPQUIET’, ‘OPDYN’, ‘T96’, ‘OSTA’, ‘T01QUIET’, ‘T01STORM’, ‘T05’, ‘ALEX’, ‘TS07’]
- - options (optional list or array of integers length=5)explained in Lstar
- - omni values as dictionary (optional)if not provided, will use lookup table
- - (see Lstar documentation for further explanation)
- Returns:
- - results (dictionary)containing keys: Bmin, Coords instance with GEO coordinates of
the magnetic equator
See also
Examples
>>> t = Ticktock(['2002-02-02T12:00:00', '2002-02-02T12:10:00'], 'ISO') >>> y = Coords([[3,0,0],[2,0,0]], 'GEO', 'car', use_irbem=True) >>> op.find_magequator(t,y) {'Bmin': array([ 945.63652413, 3373.64496167]), 'loci': Coords( [[ 2.99938371 0.00534151 -0.03213603] [ 2.00298822 -0.0073077 0.04584859]] ), dtype=GEO,car, units=['Re', 'Re', 'Re']}
- spacepy.irbempy.get_AEP8(energy, loci, model='min', fluxtype='diff', particles='e')[source]¶
will return the flux from the AE8-AP8 model
- Parameters:
- - energy (float)center energy in MeV; if fluxtype=RANGE, then needs to be a list [Emin, Emax]
- - loci (Coords)a Coords instance with the location inside the magnetosphere
optional instead of a Coords instance, one can also provide a list with [BBo, L] combination
- - model (str)MIN or MAX for solar cycle dependence
- - fluxtype (str)DIFF, RANGE, INT are possible types
- - particles (str): e or p or electrons or protons
- Returns:
- - floatflux from AE8/AP8 model
Examples
>>> spacepy.irbempy.get_aep8()
>>> import spacepy.time as spt >>> import spacepy.coordinates as spc >>> import spacepy.irbempy as ib >>> t = spt.Ticktock(['2017-02-02T12:00:00'], 'ISO') >>> y = spc.Coords([3,0,0], 'GEO', 'car', use_irbem=True) >>> y.ticks = t >>> energy = 1.0 #MeV >>> ib.get_AEP8(energy, y, model='max') 1932209.4427359989
- spacepy.irbempy.get_Bfield(ticks, loci, extMag='T01STORM', options=[1, 0, 0, 0, 0], omnivals=None)[source]¶
Return magnetic field vector (in GEO) and magnitude
Calls get_bfield from IRBEMlib and uses the underlying model implementations and coordinate transforms in IRBEMlib to obtain the result.
- Parameters:
- - ticks (Ticktock class)containing time information
- - loci (Coords class)containing spatial information
- - extMag (string)optional; will choose the external magnetic
field model possible values [‘0’, ‘MEAD’, ‘T87SHORT’, ‘T87LONG’, ‘T89’, ‘OPQUIET’, ‘OPDYN’, ‘T96’, ‘OSTA’, ‘T01QUIET’, ‘T01STORM’, ‘T05’, ‘ALEX’, ‘TS07’]
- - options (optional list or array of integers length=5)explained in Lstar
- - omni values as dictionary (optional)if not provided, will
use OMNI module to look up
- - (see Lstar documentation for further explanation)
- Returns:
- - results (dictionary)containing keys: Bvec, and Blocal
Bvec is specified in GEO coordinates
See also
Notes
Most parameterized external field models are subject to limits on the valid range of input parameters and will return NaN if evaluated outside the bounds.
Examples
>>> import spacepy.time as spt >>> import spacepy.coordinates as spc >>> import spacepy.irbempy as ib >>> t = spt.Ticktock(['2002-02-02T12:00:00', '2002-02-02T12:10:00'], 'ISO') >>> y = spc.Coords([[3,0,0],[2,0,0]], 'GEO', 'car', use_irbem=True) >>> ib.get_Bfield(t,y) {'Blocal': array([ 976.42565251, 3396.25991675]), 'Bvec': array([[ -5.01738885e-01, -1.65104338e+02, 9.62365503e+02], [ 3.33497974e+02, -5.42111173e+02, 3.33608693e+03]])}
- spacepy.irbempy.get_Lm(ticks, loci, alpha, extMag='T01STORM', intMag='IGRF', IGRFset=0, omnivals=None)[source]¶
Return the MacIlwain L value for a given location, time and model
- Parameters:
- - ticks (Ticktock class)containing time information
- - loci (Coords class)containing spatial information
- - alpha (list or ndarray)pitch angles in degrees
- - extMag (string)optional; will choose the external magnetic field model
possible values [‘0’, ‘MEAD’, ‘T87SHORT’, ‘T87LONG’, ‘T89’, ‘OPQUIET’, ‘OPDYN’, ‘T96’, ‘OSTA’, ‘T01QUIET’, ‘T01STORM’, ‘T05’, ‘ALEX’, ‘TS07’]
- - intMag (string)optional: select the internal field model
possible values [‘IGRF’,’EDIP’,’JC’,’GSFC’,’DUN’,’CDIP’] For full details see get_Lstar
- - omni values as dictionary (optional)if not provided, will use lookup table
- Returns:
- - results (dictionary)containing keys: Lm, Bmin, Blocal (or Bmirr), Xj, MLT
if pitch angles provided in “alpha” then drift shells are calculated and “Bmirr” is returned if not provided, then “Blocal” at spacecraft is returned. A negative value for Lm indicates the field line is closed but particles are lost to the atmosphere; the absolute value indicates the L value.
- spacepy.irbempy.get_Lstar(ticks, loci, alpha=90, extMag='T01STORM', options=[1, 0, 0, 0, 0], omnivals=None, landi2lstar=False)[source]¶
This will call make_lstar1 or make_lstar_shell_splitting_1 from the irbem library and will lookup omni values for given time if not provided (optional). If pitch angles are provided, drift shell splitting will be calculated and “Bmirr” will be returned. If they are not provided, then no drift shell splitting is calculated and “Blocal” is returned.
- Parameters:
- - ticks (Ticktock class)containing time information
- - loci (Coords class)containing spatial information
- - alpha (list or ndarray)optional pitch angles in degrees (default is 90);
if provided will calculate shell splitting; max 25 values
- - extMag (string)optional; will choose the external magnetic field model
possible values [‘0’, ‘MEAD’, ‘T87SHORT’, ‘T87LONG’, ‘T89’, ‘OPQUIET’, ‘OPDYN’, ‘T96’, ‘OSTA’, ‘T01QUIET’, ‘T01STORM’, ‘T05’, ‘ALEX’, ‘TS07’]
- - options (optional list or array of integers length=5)explained below
- - omni values as dictionary (optional)if not provided, will use lookup table
- - landi2lstarif True, will use the faster landi2lstar routine if possible. This
routine can only be used with OPQUIET+IGRF magnetic field models.
- Returns:
- - results (dictionary)containing keys: Lm, Lstar, Bmin, Blocal (or Bmirr), Xj, MLT
if pitch angles provided in “alpha” then drift shells are calculated and “Bmirr” is returned if not provided, then “Blocal” at spacecraft is returned. A negative value for Lm indicates the field line is closed but particles are lost to the atmosphere; the absolute value indicates the L value. A negative value for Lstar indicates the field line is closed but particles are lost to the atmosphere before completing a drift orbit; the absolute value indicates the drift shell.
Notes
- External Field
0 : no external field
MEAD : Mead & Fairfield [1975] (uses 0<=Kp<=9 - Valid for rGEO<=17. Re)
T87SHORT: Tsyganenko short [1987] (uses 0<=Kp<=9 - Valid for rGEO<=30. Re)
T87LONG : Tsyganenko long [1987] (uses 0<=Kp<=9 - Valid for rGEO<=70. Re)
T89 : Tsyganenko [1989] (uses 0<=Kp<=9 - Valid for rGEO<=70. Re)
OPQUIET : Olson & Pfitzer quiet [1977] (default - Valid for rGEO<=15. Re)
- OPDYNOlson & Pfitzer dynamic [1988] (uses 5.<=dens<=50., 300.<=velo<=500.,
-100.<=Dst<=20. - Valid for rGEO<=60. Re)
- T96Tsyganenko [1996] (uses -100.<=Dst (nT)<=20., 0.5<=Pdyn (nPa)<10.,
abs(ByIMF) (nT)<1=0., abs(BzIMF) (nT)<=10. - Valid for rGEO<=40. Re)
- OSTAOstapenko & Maltsev [1997] (uses dst,Pdyn,BzIMF, Kp)
T01QUIET: Tsyganenko [2002a,b] (uses -50.<Dst (nT)<20., 0.5<Pdyn (nPa)<=5., abs(ByIMF) (nT)<=5., abs(BzIMF) (nT)<=5., 0.<=G1<=10., 0.<=G2<=10. - Valid for xGSM>=-15. Re)
- T01STORM: Tsyganenko, Singer & Kasper [2003] storm (uses Dst, Pdyn, ByIMF, BzIMF, G2, G3 -
there is no upper or lower limit for those inputs - Valid for xGSM>=-15. Re)
- T05Tsyganenko & Sitnov [2005] storm (uses Dst, Pdyn, ByIMF, BzIMF,
W1, W2, W3, W4, W5, W6 - no upper or lower limit for inputs - Valid for xGSM>=-15. Re)
- OMNI contents
Kp: value of Kp as in OMNI2 files but has to be double instead of integer type
Dst: Dst index (nT)
dens: Solar Wind density (cm-3)
velo: Solar Wind velocity (km/s)
Pdyn: Solar Wind dynamic pressure (nPa)
ByIMF: GSM y component of IMF mag. field (nT)
BzIMF: GSM z component of IMF mag. field (nT)
- 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.
- G2: G2=< a*Vsw*Bs > where the <> mean an average over the previous 1 hour,
Vsw is solar wind speed, Bs=|IMF Bz| when IMF Bz < 0 and Bs=0 when IMF Bz > 0, a=0.005.
- 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.
W1 - W6: see definition in (JGR-A, v.110(A3), 2005.) (PDF 1.2MB)
AL: the auroral index
- Options
1st element: 0 - don’t compute L* or phi ; 1 - compute L*; 2- compute phi
- 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 …)
- 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 (theta) along the field line such as dtheta=(2pi)/(720*[options(3rd element)+1])
- 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 (phi) along the drift shell such as dphi=(2pi)/(25*[options(4th element)+1])
- 5th element: allows to select an internal magnetic field model (default is set to IGRF)
0 = IGRF
1 = Eccentric tilted dipole
2 = Jensen&Cain 1960
3 = GSFC 12/66 updated to 1970
4 = User-defined model (Default: Centred dipole + uniform [Dungey open model] )
5 = Centred dipole
Examples
>>> t = Ticktock(['2002-02-02T12:00:00', '2002-02-02T12:10:00'], 'ISO') >>> y = Coords([[3,0,0],[2,0,0]], 'GEO', 'car', use_irbem=True) >>> spacepy.irbempy.Lstar(t,y) {'Blocal': array([ 1020.40493286, 3446.08845227]), 'Bmin': array([ 1019.98404311, 3437.63865243]), 'Lm': array([ 3.08948304, 2.06022102]), 'Lstar': array([ 2.97684043, 1.97868577]), 'MLT': array([ 23.5728333 , 23.57287944]), 'Xj': array([ 0.00112884, 0.00286955])}
- spacepy.irbempy.get_dtype(sysaxes)[source]¶
will return the coordinate system type as string
- Parameters:
- - sysaxes (int)number according to the irbem, possible values: 0-8
- Returns:
- - dtype (str)coordinate system GDZ, GEO, GSM, GSE, SM, GEI, MAG, SPH, RLL
- - carsph (str)cartesian or spherical ‘car’, ‘sph’
Examples
>>> get_dtype(3) ('GSE', 'car')
- spacepy.irbempy.prep_ctypes(d)[source]¶
Converts
prep_irbem
output to correct types for ctypes call.This is primarily for internal use by
get_AEP8
andget_Lstar
.- Parameters:
- ddict
Output IRBEM parameters from
prep_irbem
- Returns:
- dict
d
with most values converted to ctypes objects or pointers.
- spacepy.irbempy.prep_irbem(ticks=None, loci=None, alpha=[], extMag='T01STORM', options=[1, 0, 0, 0, 0], omnivals=None)[source]¶
Prepare inputs for direct IRBEM-LIB calls. Not expected to be called by the user.
- spacepy.irbempy.updateTS07Coeffs(path=None, force=False, verbose=False, **kwargs)[source]¶
Update coefficients for TS07 magnetic field model
- Parameters:
- pathstr
Base path for TS07 dynamic coefficients. If None (default) then the TS07_DATA_PATH environment variable is used.
- forceboolean
If True (default is False) then missing paths will be created.
- verboseboolean
Print verbose output. Default is False.
- startdatetime.datetime or string
Start time for archive retrieval. If start is a string then it should ISO8601 format (YYYY-MM-DDTHH:mm:SS). Defaults to 2007-01-01.
- enddatetime.datetime or string
End time for archive retrieval. Required format same as start. Defaults to time of query (i.e. latest available).