nenupy.astro.astro_tools
Astronomical tools
- class nenupy.astro.astro_tools.AstroObject[source]
Bases:
ABC
Abstract base class for all astronomy related classes.
- coordinates
- property custom_ho_coordinates
Allows to modify horizontal coordinates without messing up with the actual coordinates object.
- property frequency
test de doc
- property ground_projection
- property horizontal_coordinates
- observer
- polarization
- time
- value
- nenupy.astro.astro_tools.altaz_to_radec(altaz, fast_compute=False)[source]
Converts a celestial object horizontal coordinates to equatorial coordinates.
If
fast_compute=True
is selected, the computation is accelerated using Local Sidereal Time approximation (seelocal_sidereal_time()
). The right ascension \(\alpha\) and declination \(\delta\) are computed as follows:\[\begin{split}\cases{ \delta = \sin^(-1) \left( \sin l \sin \theta + \cos l \cos \theta \cos \varphi \right)\\ h = \cos^{-1} \left( \frac{\sin \theta - \sin l \sin \delta}{\cos l \cos \delta} \right)\\ \alpha = t_{\rm{LST}} - h }\end{split}\]with \(\theta\) the object’s elevation, \(\varphi\) the azimuth, \(l\) the
observer
’s latitude and \(h\) the Local Hour Angle (seehour_angle()
). If \(\sin(h^{\prime}) \geq 0\), then \(h = - (h^{\prime} - \pi)\). Otherwise,transform_to()
is used.- Parameters:
altaz – Celestial object horizontal coordinates.
fast_compute (
bool
) – If set toTrue
, it enables faster computation time for the conversion, mainly relying on an approximation of the local sidereal time. All other values would lead to accurate coordinates computation. Differences in coordinates values are of the order of \(10^{-2}\) degrees or less.
- Returns:
Celestial object’s equatorial coordinates.
- Return type:
- Example:
from nenupy.astro import altaz_to_radec from nenupy import nenufar_position from astropy.time import Time from astropy.coordinates import SkyCoord, AltAz radec = altaz_to_radec( altaz=SkyCoord( 300, 45, unit="deg", frame=AltAz( obstime=Time("2022-01-01T12:00:00"), location=nenufar_position ) ), fast_compute=True )
- nenupy.astro.astro_tools.dispersion_delay(frequency, dispersion_measure)[source]
Dispersion delay induced to a radio wave of
frequency
(\(\nu\)) propagating through an electron plasma of uniform density \(n_e\).The pulse travel time \(\Delta t_p\) emitted at a distance \(d\) is:
\[\Delta t_p = \frac{d}{c} + \frac{e^2}{2\pi m_e c} \frac{\int_0^d n_e\, dl}{\nu^2}\]where \(\mathcal{D}\mathcal{M} = \int_0^d n_e\, dl\) is the Dispersion Measure (
dispersion_measure
). Therefore, the time delay \(\Delta t_d\) due to the dispersion is:\[\Delta t_d = \frac{e^2}{2 \pi m_e c} \frac{\mathcal{D}\mathcal{M}}{\nu^2}\]and computed as:
\[\Delta t_d = 4140 \left( \frac{\mathcal{D}\mathcal{M}}{\rm{pc}\,\rm{cm}^{-3}} \right) \left( \frac{\nu}{1\, \rm{MHz}} \right)^{-2}\, \rm{sec}\]- Parameters:
- Returns:
Dispersion delay in seconds.
- Return type:
- Example:
from nenupy.astro import dispersion_delay import astropy.units as u dispersion_delay( frequency=50*u.MHz, dispersion_measure=12.4*u.pc/(u.cm**3) )
- nenupy.astro.astro_tools.etrs_to_enu(positions, location=<EarthLocation (4323936.68522791, 165534.49991696, 4670345.36540385) m>)[source]
Local east, north, up (ENU) coordinates centered on the position
location
.The conversion from cartesian coordinates \((x, y, z)\) to ENU \((e, n, u)\) is done as follows:
\[\begin{split}\pmatrix{ e \\ n \\ u } = \pmatrix{ -\sin(b) & \cos(l) & 0\\ -\sin(l) \cos(b) & -\sin(l) \sin(b) & \cos(l)\\ \cos(l)\cos(b) & \cos(l) \sin(b) & \sin(l) } \pmatrix{ \delta x\\ \delta y\\ \delta z }\end{split}\]where \(b\) is the longitude, \(l\) is the latitude and \((\delta x, \delta y, \delta z)\) are the cartesian coordinates with respect to the center
location
.- Parameters:
positions (
ndarray
) – ETRS positionslocation (
EarthLocation
) – Center of ENU frame. Default is NenuFAR’s location.
- Returns:
Wavelength in meters, same shape as
frequency
.- Return type:
- Example:
from nenupy import nenufar_position from nenupy.astro import etrs_to_enu etrs_positions = np.array([ [4323934.57369062, 165585.71569665, 4670345.01314493], [4323949.24009871, 165567.70236494, 4670332.18016874] ]) enu = etrs_to_enu( positions=etrs_positions, location=nenufar_position )
- nenupy.astro.astro_tools.faraday_angle(frequency, rotation_measure, inverse=False)[source]
Computes the Faraday rotation angle.
- Parameters:
- Returns:
Faraday rotation angle in radians.
- Return type:
- nenupy.astro.astro_tools.geo_to_etrs(location=<EarthLocation (4323936.68522791, 165534.49991696, 4670345.36540385) m>)[source]
Transforms geographic coordinates to ETRS (European Terrestrial Reference System).
- Parameters:
location (
EarthLocation
) – Location to convert.- Returns:
ETRS positions.
- Return type:
- Example:
from nenupy.astro import geo_to_etrs from astropy.coordinates import EarthLocation import astropy.units as u locs = EarthLocation( lat=[30, 40] * u.deg, lon=[0, 10] * u.deg, height=[100, 200] * u.m ) geo_to_etrs(location=locs)
- nenupy.astro.astro_tools.geo_to_l93(location=<EarthLocation (4323936.68522791, 165534.49991696, 4670345.36540385) m>)[source]
Convert geographic coordinates to RGF93 coordinates. This function takes in a location in geographic coordinates (EPSG:4326) and converts it to RGF93 coordinates (EPSG:2154).
- Parameters:
location (
EarthLocation
) – The location to be converted. Defaults tonenufar_position
.- Returns:
The location in RGF93 coordinates, as an array of 3 elements.
- Return type:
- nenupy.astro.astro_tools.hour_angle(radec, time, observer=<EarthLocation (4323936.68522791, 165534.49991696, 4670345.36540385) m>, fast_compute=True)[source]
Local Hour Angle of an object in the
observer
’s sky. It is defined as the angular distance on the celestial sphere measured westward along the celestial equator from the meridian to the hour circle passing through a point.The local hour angle \(h\) is computed with respect to the local sidereal time \(t_{\rm LST}\) and the astronomical source (defined as
radec
) right ascension \(\alpha\):\[h = t_{\rm LST} - \alpha\]with the rule that if \(h < 0\), a \(2\pi\) angle is added or if \(h > 2\pi\), a \(2\pi\) angle is subtracted.
- Parameters:
radec (
SkyCoord
) – Sky coordinates to convert to Local Hour Angles.time (
Time
) – UTC time a which the hour angle is computed.observer (
EarthLocation
) – Earth location where the observer is at. Default is NenuFAR’s position.fast_compute (
bool
) – If set toTrue
, an approximation is made while computing the local sidereal time. Default isTrue
.
- Returns:
Local hour angle.
- Return type:
- Example:
from nenupy.astro import hour_angle from astropy.coordinates import SkyCoord from astropy.time import Time ha = hour_angle( radec=SkyCoord(300, 45, unit="deg"), time=Time("2022-01-01T12:00:00"), fast_compute=True )
- nenupy.astro.astro_tools.l93_to_etrs(positions)[source]
Transforms Lambert93 coordinates to ETRS (European Terrestrial Reference System).
- Parameters:
positions (
ndarray
) – Lambert-93 positions (in meters).- Returns:
ETRS positions.
- Return type:
- Example:
from nenupy.astro import l93_to_etrs import numpy as np l93 = np.array([ [6.39113316e+05, 6.69766347e+06, 1.81735000e+02], [6.39094578e+05, 6.69764471e+06, 1.81750000e+02] ]) 93_to_etrs(positions=l93)
- nenupy.astro.astro_tools.l93_to_geo(positions)[source]
Convert RGF93 to geographic coordinates. This function takes in positions in RGF93 coordinates (EPSG:2154) and converts it to geographic coordinates (EPSG:4326).
- Parameters:
positions (
ndarray
) – The positions to be converted.- Returns:
Geographic coordinates.
- Return type:
- nenupy.astro.astro_tools.local_sidereal_time(time, observer=<EarthLocation (4323936.68522791, 165534.49991696, 4670345.36540385) m>, fast_compute=True)[source]
Computes the Local Sidereal Time. Viewed from
observer
, a fixed celestial object seen at one position in the sky will be seen at the same position on another night at the same sidereal time. LST angle indicates the Right Ascension on the sky that is currently crossing the Local Meridian.- Parameters:
time (
Time
) – UT Time to be converted.observer (
EarthLocation
) – Earth location of the observer.fast_compute (
bool
) – If set toTrue
, this computes an approximation of the local sidereal time, otherwisesidereal_time()
is used.
- Returns:
Local Sidereal Time angle
- Return type:
- Example:
from nenupy.astro import local_sidereal_time from astropy.time import Time lst = local_sidereal_time( time=Time("2022-01-01T12:00:00"), fast_compute=True )
- nenupy.astro.astro_tools.parallactic_angle(coordinates, time, observer=<EarthLocation (4323936.68522791, 165534.49991696, 4670345.36540385) m>)[source]
TO DO/Done?
\[q = - {\rm atan2 }( -\sin(h) , \cos(\delta) \tan(l) - \sin(\delta)\cos(h))\]where \(q\) is the parallactic angle, \(h\) is the hour angle, \(\delta\) is the declination and \(l\) is the observers’s latitude.
- Parameters:
coordinates (
SkyCoord
) – Coordinates at which the parallactic angle is evaluated.time (:class:~`astropy.time.Time`) – UT time(s) at which the parallactic angle is evaluated.
observer (
EarthLocation
) – Observer location on the Earth.
- Returns:
The parallactic angle.
- Return type:
- Example:
from nenupy.astro.astro_tools import parallactic_angle from nenupy.astro.target import FixedTarget from astropy.time import TimeDelta cyga = FixedTarget.from_name('Cyg A') transit = cyga.next_meridian_transit() dt = TimeDelta(5*3600, format='sec') steps = 20 times = transit - dt + np.arange(steps)*(dt*2/steps) pa = parallactic_angle( ra_j2000=cyga.coordinates.ra, dec_j2000=cyga.coordinates.dec, time=times )
- nenupy.astro.astro_tools.radec_to_altaz(radec, time, observer=<EarthLocation (4323936.68522791, 165534.49991696, 4670345.36540385) m>, fast_compute=True)[source]
Converts a celestial object equatorial coordinates to horizontal coordinates.
If
fast_compute=True
is selected, the computation is accelerated using Local Sidereal Time approximation (seelocal_sidereal_time()
). The altitude \(\theta\) and azimuth \(\varphi\) are computed as follows:\[\begin{split}\cases{ \sin(\theta) = \sin(\delta) \sin(l) + \cos(\delta) \cos(l) \cos(h)\\ \cos(\varphi) = \frac{\sin(\delta) - \sin(l) \sin(\theta)}{\cos(l)\cos(\varphi)} }\end{split}\]with \(\delta\) the object’s declination, \(l\) the
observer
’s latitude and \(h\) the Local Hour Angle (seehour_angle()
). If \(\sin(h) \geq 0\), then \(\varphi = 2 \pi - \varphi\). Otherwise,transform_to()
is used.- Parameters:
radec (
SkyCoord
) – Celestial object equatorial coordinates.time (
Time
) – Coordinated universal time.observer (
EarthLocation
) – Earth location where the observer is at. Default is NenuFAR’s position.fast_compute (
bool
) – If set toTrue
, it enables faster computation time for the conversion, mainly relying on an approximation of the local sidereal time. All other values would lead to accurate coordinates computation. Differences in coordinates values are of the order of \(10^{-2}\) degrees or less.
- Returns:
Celestial object’s horizontal coordinates. If either
radec
ortime
are 1D arrays, the resulting object will be of shape(time, positions)
.- Return type:
- Example:
from nenupy.astro import radec_to_altaz from astropy.time import Time from astropy.coordinates import SkyCoord altaz = radec_to_altaz( radec=SkyCoord([300, 200], [45, 45], unit="deg"), time=Time("2022-01-01T12:00:00"), fast_compute=True )
- nenupy.astro.astro_tools.sky_temperature(frequency=<Quantity 50. MHz>)[source]
Sky temperature at a given frequency
frequency
(strongly dominated by Galactic emission).\[T_{\rm sky} = T_0 \lambda^{2.55}\]with \(T_0 = 60 \pm 20\,\rm{K}\) for Galactic latitudes between 10 and 90 degrees and \(\lambda\) the wavelength.
- Parameters:
frequency (
Quantity
) – Frequency at which computing the sky temperature. Default is50 MHz
.- Returns:
Sky temperature in Kelvins
- Return type:
- Example:
from nenupy.astro import sky_temperature import astropy.units as u temp = sky_temperature(frequency=50*u.MHz)
See also
LOFAR website, Haslam et al. (1982) and Mozdzen et al. (2017, 2019)
- nenupy.astro.astro_tools.solar_system_source(name, time=<Time object: scale='utc' format='datetime' value=2024-05-03 08:29:39.106431>, observer=<EarthLocation (4323936.68522791, 165534.49991696, 4670345.36540385) m>)[source]
Returns a Solar System body in the ICRS reference system.
- Parameters:
name (
str
) – Name of the Solar System object (a call is made toastropy.coordinates.get_body()
).time (
Time
) – Time at which the Solar System object is observed.observer (
EarthLocation
) – Earth location from which the source is observed.
- Returns:
Coordinates object in ICRS reference frame of the Solar System object.
- Return type:
- Example:
from astropy.time import Time from nenupy.astro import solar_system_source sun = solar_system_source( name="Sun", time=Time("2022-01-01T12:00:00") )