pygedm.ne2001_wrapper
Python wrapper for NE2001 model code.
References
[1] Cordes, J. M., & Lazio, T. J. W. (2002), NE2001.I. A New Model for the Galactic Distribution of Free Electrons and its Fluctuations, arXiv e-prints, astro-ph/0207156.
[2] Cordes, J. M., & Lazio, T. J. W. (2003), NE2001. II. Using Radio Propagation Data to Construct a Model for the Galactic Distribution of Free Electrons, arXiv e-prints, astro-ph/0301598.
- pygedm.ne2001_wrapper.TAUISS(d, sm, nu)
Convert a scattering measure (SM) to scattering timescale at given frequency.
Direct port from FORTRAN code scattering98.f
- Parameters:
d (float) – Distance in kpc
sm (float) – Scattering measure (kpc m^{-20/3})
nu (float) – Radio frequency in GHz
- Returns:
pulse broadening time (ms)
- Return type:
tauiss (float)
Fortran equiv:
REAL FUNCTION TAUISS(d, sm, nu) c c calculates the pulse broadening time in ms c from distance, scattering measure, and radio frequency c c input: d = pulsar distance (kpc) c sm = scattering measure (kpc m^{-20/3}) c nu = radio frequency (GHz) c output: tausis = pulse broadening time (ms) c implicit none real d, sm, nu tauiss = 1000. * (sm / 292.)**1.2 * d * nu**(-4.4) end
- pygedm.ne2001_wrapper.TRANSITION_FREQUENCY(d_eff: float, sm: float, xi: float = None) float
The transition frequency between weak and strong scattering
Python implementation of Equation 17 in C&L 2002 [1].
Computes nu = (318 GHz) * xi^(10/17) * sm^(6/17) * d_eff^(5/17)
- Parameters:
d_eff (float) – effective distance to the scattering medium (kpc)
sm (float) – Scattering measure (kpc m^{-20/3})
xi (float) – Fresnel scale scaling factor (~1). If not supplied, sqrt(2 pi) is used as preferred definition.
- Returns:
Transition frequency in GHz
- Return type:
nu (float)
- pygedm.ne2001_wrapper.calculate_electron_density_xyz(x, y, z)
Compute electron density at Galactocentric X, Y, Z coordinates
x,y,z are Galactocentric Cartesian coordinates, measured in kpc (NOT pc!) with the axes parallel to (l, b) = (90, 0), (180, 0), and (0, 90) degrees
- Parameters:
x (float) – Galactocentric coordinates in kpc
y (float) – Galactocentric coordinates in kpc
z (float) – Galactocentric coordinates in kpc
- Returns:
Electron density in cm-3
- Return type:
ne_out (astropy.Quantity)
- pygedm.ne2001_wrapper.dist_to_dm(l, b, dist, nu=1.0, full_output=False)
Convert distance to DM and compute scattering timescale
- Parameters:
l (float) – galactic longitude in degrees
b (float) – galactic latitude in degrees
dist (float) – Distance in kpc
nu (float in GHz or astropy.Quantity) – observing frequency (GHz)
full_output (bool) – Return full raw output (dict) from NE2001 if set to True
- Returns:
Dispersion measure (pc / cm3), scattering timescale at 1 GHz (s)
- Return type:
dm (astropy.Quantity), tau_sc (astropy.Quantity)
- pygedm.ne2001_wrapper.dm_to_dist(l, b, dm, nu=1.0, full_output=False)
Convert DM to distance and compute scattering timescale
- Parameters:
l (float) – galactic longitude in degrees
b (float) – galactic latitude in degrees
dm (floa) – Dispersion measure
nu (float in GHz or astropy.Quantity) – observing frequency (GHz)
full_output (bool) – Return full raw output (dict) from NE2001 if set to True
- Returns:
Distance (pc), scattering timescale at 1 GHz (s)
- Return type:
dist (astropy.Quantity), tau_sc (astropy.Quantity)
- pygedm.ne2001_wrapper.run_from_pkgdir(f)
Decorator function to chdir() into package directory when running
NE2001 code doesn’t know the relative path to its data files. This wraps the function call, changing into the right directory first, calling it, then changing back to original directory.