aboutsummaryrefslogtreecommitdiffstats
path: root/fg21sim/utils/convert.py
blob: 9eb906121190b4f03f09b76c588010949d89f586 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
# Copyright (c) 2016-2017 Weitian LI <weitian@aaronly.me>
# MIT license

"""
Utilities for conversion among common astronomical quantities.
"""

import astropy.units as au

from .units import (UnitConversions as AUC, Constants as AC)


def Fnu_to_Tb(Fnu, omega, freq):
    """
    Convert flux density to brightness temperature, using the
    Rayleigh-Jeans law, in the Rayleigh-Jeans limit.

    Parameters
    ----------
    Fnu : `~astropy.units.Quantity`
        Input flux density, e.g., `1.0*au.Jy`
    omega : `~astropy.units.Quantity`
        Source angular size/area, e.g., `100*au.arcsec**2`
    freq : `~astropy.units.Quantity`
        Frequency where the flux density measured, e.g., `120*au.MHz`

    Returns
    -------
    Tb : `~astropy.units.Quantity`
        Brightness temperature, with default unit `au.K`

    References
    ----------
    - Brightness and Flux
      http://www.cv.nrao.edu/course/astr534/Brightness.html
    - Wikipedia: Brightness Temperature
      https://en.wikipedia.org/wiki/Brightness_temperature
    - NJIT: Physics 728: Introduction to Radio Astronomy: Lecture #1
      https://web.njit.edu/~gary/728/Lecture1.html
    - Astropy: Equivalencies: Brightness Temperature / Flux Density
      http://docs.astropy.org/en/stable/units/equivalencies.html
    """
    equiv = au.brightness_temperature(omega, freq)
    Tb = Fnu.to(au.K, equivalencies=equiv)
    return Tb


def Sb_to_Tb(Sb, freq):
    """
    Convert surface brightness to brightness temperature, using the
    Rayleigh-Jeans law, in the Rayleigh-Jeans limit.

    Parameters
    ----------
    Sb : `~astropy.units.Quantity`
        Input surface brightness, e.g., `1.0*(au.Jy/au.arcsec**2)`
    freq : `~astropy.units.Quantity`
        Frequency where the flux density measured, e.g., `120*au.MHz`

    Returns
    -------
    Tb : `~astropy.units.Quantity`
        Brightness temperature, with default unit `au.K`
    """
    omega = 1.0 * au.arcsec**2
    Fnu = (Sb * omega).to(au.Jy)  # [Jy]
    return Fnu_to_Tb(Fnu, omega, freq)


def Sb_to_Tb_fast(Sb, freq):
    """
    Convert surface brightness to brightness temperature, using the
    Rayleigh-Jeans law, in the Rayleigh-Jeans limit.

    Avoid using `astropy.units` to optimize the speed.

        Tb = Sb * c^2 / (2 * k_B * nu^2)

    where `Sb` is the surface brightness density measured at a certain
    frequency, in units of [Jy/arcsec^2].

    1 [Jy] = 1e-23 [erg/s/cm^2/Hz] = 1e-26 [W/m^2/Hz]

    Parameters
    ----------
    Sb : float
        Input surface brightness
        Unit: [Jy/arcsec^2]
    freq : float
        Frequency where the flux density measured
        Unit: [MHz]

    Returns
    -------
    Tb : float
        Calculated brightness temperature
        Unit: [K]
    """
    # NOTE: [rad] & [sr] are dimensionless
    arcsec2 = AUC.arcsec2rad ** 2  # [sr]
    Sb /= arcsec2  # [Jy/arcsec^2] -> [Jy/sr]
    coef = 1e-35  # unit conversion coefficient
    Tb = coef * (Sb * AC.c**2) / (2*AC.k_B * freq**2)  # [K]
    return Tb


def Fnu_to_Tb_fast(Fnu, omega, freq):
    """
    Convert flux density to brightness temperature, using the
    Rayleigh-Jeans law, in the Rayleigh-Jeans limit.

    Avoid using `astropy.units` to optimize the speed.

    Parameters
    ----------
    Fnu : float
        Input flux density
        Unit: [Jy] = 1e-23 [erg/s/cm^2/Hz] = 1e-26 [W/m^2/Hz]
    omega : float
        Source angular size/area
        Unit: [arcsec^2]
    freq : float
        Frequency where the flux density measured
        Unit: [MHz]

    Returns
    -------
    Tb : float
        Calculated brightness temperature
        Unit: [K]
    """
    Sb = Fnu / omega  # [Jy/arcsec^2]
    return Sb_to_Tb_fast(Sb, freq)


def JyPerPix_to_K(freq, pixelsize):
    """
    The factor that converts [Jy/pixel] to [K] (brightness temperature).

    Parameters
    ----------
    freq : float
        The frequency where the flux density measured.
        Unit: [Jy]
    pixelsize : float
        The pixel size.
        Unit: [arcsec]
    """
    factor = Fnu_to_Tb_fast(Fnu=1.0, omega=pixelsize**2, freq=freq)
    return factor