1D RM-synthesis

[1]:
from __future__ import annotations

import astropy.units as u
import matplotlib.pyplot as plt
import numpy as np
from astropy.visualization import quantity_support
from numpy.typing import NDArray
from rm_lite.tools_1d import rmsynth

plt.rcParams["figure.dpi"] = 150

_ = quantity_support()
rng = np.random.default_rng(42)

First generate some synthetic data

[2]:
from rm_lite.utils.synthesis import faraday_simple_spectrum, freq_to_lambda2


def faraday_slab_spectrum(
    freq_arr_hz: NDArray[np.float64],
    frac_pol: float,
    psi0_deg: float,
    rm_radm2: float,
    delta_rm_radm2: float,
) -> NDArray[np.complex128]:
    lambda_sq_arr_m2 = freq_to_lambda2(freq_arr_hz)

    return (
        frac_pol
        * np.exp(2j * (np.deg2rad(psi0_deg) + rm_radm2 * lambda_sq_arr_m2))
        * (
            np.sin(delta_rm_radm2 * lambda_sq_arr_m2)
            / (delta_rm_radm2 * lambda_sq_arr_m2)
        )
    )


def faraday_gaussian_spectrum(
    freq_arr_hz: NDArray[np.float64],
    frac_pol: float,
    psi0_deg: float,
    rm_radm2: float,
    sigma_rm_radm2: float,
):
    lambda_sq_arr_m2 = freq_to_lambda2(freq_arr_hz)
    rm_term = np.exp(2j * (np.deg2rad(psi0_deg) + rm_radm2 * lambda_sq_arr_m2))
    depol_term = np.exp(-2.0 * sigma_rm_radm2**2 * lambda_sq_arr_m2**2)
    return frac_pol * rm_term * depol_term

Here we’ll simulate RACS-all frequency coverage

[3]:
bw_low = 288
bw_mid = 144
bw_high = 288
low = np.linspace(943.5 - bw_low / 2, 943.5 + bw_low / 2, 36) * u.MHz
mid = np.linspace(1367.5 - bw_mid / 2, 1367.5 + bw_mid / 2, 9) * u.MHz
high = np.linspace(1655.5 - bw_high / 2, 1655.5 + bw_high / 2, 9) * u.MHz
freqs = np.concatenate([low, mid, high])
freq_hz = freqs.to(u.Hz).value

Now we make a Faraday simple spectrum with a single RM component. We will use the following parameters:

[4]:
delta_rm_radm2 = 30
rm_radm2 = 100
frac_pol = 0.5
psi0_deg = 10
complex_data_noiseless = faraday_simple_spectrum(
    freq_arr_hz=freq_hz,
    frac_pol=frac_pol,
    psi0_deg=psi0_deg,
    rm_radm2=rm_radm2,
)
[5]:
fig, ax = plt.subplots()
ax.plot(
    freq_hz, np.real(complex_data_noiseless), ".", label="Stokes Q", color="tab:red"
)
ax.plot(
    freq_hz, np.imag(complex_data_noiseless), ".", label="Stokes U", color="tab:blue"
)
ax.legend()
ax.set(
    xlabel=rf"$\nu$ / {u.Hz:latex_inline}",
    ylabel="Flux density",
    title="Stokes Q and U",
)
[5]:
[Text(0.5, 0, '$\\nu$ / $\\mathrm{Hz}$'),
 Text(0, 0.5, 'Flux density'),
 Text(0.5, 1.0, 'Stokes Q and U')]
../_images/examples_rmsynth_1d_8_1.png

Now we can run RM-synthesis by calling rmsynth.run_rmsynth

[6]:
help(rmsynth.run_rmsynth)
Help on function run_rmsynth in module rm_lite.tools_1d.rmsynth:

run_rmsynth(freq_arr_hz: 'NDArray[np.float64]', complex_pol_arr: 'NDArray[np.complex128]', complex_pol_error: 'NDArray[np.complex128]', stokes_i_arr: 'NDArray[np.float64] | None' = None, stokes_i_error_arr: 'NDArray[np.float64] | None' = None, stokes_i_model_arr: 'NDArray[np.float64] | None' = None, stokes_i_model_error: 'NDArray[np.float64] | None' = None, phi_max_radm2: 'float | None' = None, d_phi_radm2: 'float | None' = None, n_samples: 'float | None' = 10.0, weight_type: "Literal['variance', 'uniform']" = 'variance', do_fit_rmsf: 'bool' = False, do_fit_rmsf_real: 'bool' = False, fit_function: "Literal['log', 'linear']" = 'log', fit_order: 'int' = 2, ignore_stokes_i: 'bool' = False) -> 'RMSynth1DResults'
    Run RM-synthesis on 1D data

    Args:
        freq_arr_hz (NDArray[np.float64]): Frequencies in Hz
        complex_pol_arr (NDArray[np.complex128]): Complex polarisation values (Q + iU)
        complex_pol_error (NDArray[np.float64]): Complex polarisation errors (dQ + idU)
        stokes_i_arr (NDArray[np.float64] | None, optional): Total itensity values. Defaults to None.
        stokes_i_error_arr (NDArray[np.float64] | None, optional): Total intensity errors. Defaults to None.
        stokes_i_model_arr (NDArray[np.float64] | None, optional): Total intensity model array. Defaults to None.
        stokes_i_model_error (NDArray[np.float64] | None, optional): Total intensity model error. Defaults to None.
        phi_max_radm2 (float | None, optional): Maximum Faraday depth. Defaults to None.
        d_phi_radm2 (float | None, optional): Spacing in Faraday depth. Defaults to None.
        n_samples (float | None, optional): Number of samples across the RMSF. Defaults to 10.0.
        weight_type ("variance", "uniform", optional): Type of weighting. Defaults to "variance".
        do_fit_rmsf (bool, optional): Fit the RMSF main lobe. Defaults to False.
        do_fit_rmsf_real (bool, optional): The the real part of the RMSF. Defaults to False.
        fit_function ("log" | "linear", optional): _description_. Defaults to "log".
        fit_order (int, optional): Polynomial fit order. Defaults to 2. Negative values will iterate until the fit is good.

    Returns:
        RMSynth1DResults:
            fdf_parameters (pl.DataFrame): FDF parameters
            fdf_arrs (pl.DataFrame): RMSynth arrays
            rmsf_arrs (pl.DataFrame): RMSF arrays

[7]:
fdf_parameters, fdf_arrs, rmsf_arrs, stokes_i_arrs = rmsynth.run_rmsynth(
    freq_arr_hz=freq_hz,
    complex_pol_arr=complex_data_noiseless,
    complex_pol_error=np.zeros_like(complex_data_noiseless),
    do_fit_rmsf=True,
    n_samples=100,
)
WARNING rmsynth.run_rmsynth: Stokes I array/errors or model not provided. No fractional polarization will be calculated.
INFO synthesis.rmsynth_nufft: Running RM-synthesis using the NUFFTs over 2001 Faraday depth channels.
INFO synthesis.rmsynth_nufft: NUFFT complete in 0.00693 seconds.
INFO synthesis.get_rmsf_nufft: Fitting main lobe in each RMSF spectrum.
INFO rmsynth._run_rmsynth: RM-synthesis completed in 16.19ms.

The output values are Polars dataframes that can be inspected easily

[8]:
fdf_parameters
[8]:
shape: (1, 30)
fdf_error_madpeak_pi_fitpeak_pi_errorpeak_pi_fit_debiaspeak_pi_fit_snrpeak_pi_fit_indexpeak_rm_fitpeak_rm_fit_errorpeak_q_fitpeak_u_fitpeak_pa_fit_degpeak_pa_fit_deg_errorpeak_pa0_fit_degpeak_pa0_fit_deg_errorfit_functionlam_sq_0_m2ref_freq_hzfwhm_rmsf_radm2fdf_error_noisefdf_q_noisefdf_u_noisemin_freq_hzmax_freq_hzn_channelsmedian_d_freq_hzfrac_polfrac_pol_errorsigma_addsigma_add_minussigma_add_plus
f64f64f64f64f64i64f64f64f64f64f64f64f64f64strf64f64f64f64f64f64f64f64i64f64f64f64f64f64f64
NaN0.5002790.00.500279inf129699.9999040.0-0.202569-0.457124123.0500180.010.0004530.0"log"0.0825631.0433e934.5221610.00.00.07.995e81.7995e9548.2286e6NaNNaNNaNNaNNaN
[9]:
fdf_arrs
[9]:
shape: (2_001, 2)
phi_arr_radm2fdf_dirty_complex_arr
f64object
-336.725921(-0.012941228878331946+0.005575863548035223j)
-336.389195(-0.013583746733568334+0.005189841778271951j)
-336.052469(-0.01421246414961399+0.004794433774261817j)
-335.715743(-0.01482649113403289+0.004389815340038628j)
-335.379017(-0.015424947400465542+0.003976168650354153j)
335.379017(-0.0009008454185249535+0.00805335670839573j)
335.715743(-0.0009576577337929378+0.007885230167625068j)
336.052469(-0.0009822777769244494+0.007717713719139009j)
336.389195(-0.0009751524589886868+0.007551126341000855j)
336.725921(-0.0009367729099190582+0.007385776163138716j)
[10]:
rmsf_arrs
[10]:
shape: (4_003, 2)
phi2_arr_radm2rmsf_complex_arr
f64object
-673.788568(-0.03295766077582867-0.0360836275302219j)
-673.451842(-0.032467191441512396-0.03631484421234065j)
-673.115116(-0.031940535967838565-0.03653906199687933j)
-672.77839(-0.031377748766446875-0.03675625231592756j)
-672.441664(-0.03077893335526642-0.036966380167596694j)
672.441664(-0.03077893335568311+0.03696638016745783j)
672.77839(-0.03137774876683936+0.03675625231578416j)
673.115116(-0.031940535968206354+0.03653906199673062j)
673.451842(-0.03246719144185594+0.036314844212187315j)
673.788568(-0.03295765556152026+0.03608362656193195j)

Since we provided no Stokes \(I\) data, the stokes I model will just be unity with 0 error. The flag_arr array tells us which channels were not used in RM-synthesis or model fitting

[11]:
stokes_i_arrs
[11]:
shape: (54, 7)
freq_arr_hzlambda_sq_arr_m2stokes_i_model_arrstokes_i_model_errorflag_arrcomplex_pol_arrcomplex_pol_error
f64f64f64f64boolobjectobject
7.995e80.140606nullnulltrue(-0.49042944532403415-0.09735994638022451j)0j
8.0773e80.137756nullnulltrue(-0.4654241884816873+0.18270283187778674j)0j
8.1596e80.134992nullnulltrue(-0.3001379571186853+0.3998964949791661j)0j
8.2419e80.13231nullnulltrue(-0.053617325136419355+0.497116870006657j)0j
8.3241e80.129707nullnulltrue(0.20074123829142537+0.45793335240974226j)0j
1.6555e90.032793nullnulltrue(0.4056251419940925+0.29234952399871j)0j
1.6915e90.031412nullnulltrue(0.46997598588767076+0.17065336998990638j)0j
1.7275e90.030117nullnulltrue(0.49801253407054746+0.04453668048509923j)0j
1.7635e90.0289nullnulltrue(0.49406601726750204-0.07680345422849415j)0j
1.7995e90.027755nullnulltrue(0.463743096038396-0.1869287053309979j)0j

We can also easily visualise the data

[12]:
phi_arr_radm2 = fdf_arrs["phi_arr_radm2"].to_numpy()
fdf_dirty_arr = fdf_arrs["fdf_dirty_complex_arr"].to_numpy().astype(complex)

fig, ax = plt.subplots()
x1, x2, y1, y2 = 95, 105, 0.45, 0.55  # subregion of the original image
axins = ax.inset_axes(
    (0.9, 0.6, 0.4, 0.4), xlim=(x1, x2), ylim=(y1, y2), xticklabels=[], yticklabels=[]
)
for _ax in [ax, axins]:
    _ax.plot(
        phi_arr_radm2,
        fdf_dirty_arr.real,
        color="tab:red",
        label="Stokes Q",
    )
    _ax.plot(
        phi_arr_radm2,
        fdf_dirty_arr.imag,
        color="tab:blue",
        label="Stokes U",
    )
    _ax.plot(
        phi_arr_radm2,
        np.abs(fdf_dirty_arr),
        color="k",
        label="Polarized intensity",
    )

    _ax.errorbar(
        fdf_parameters["peak_rm_fit"],
        fdf_parameters["peak_pi_fit"],
        xerr=fdf_parameters["peak_rm_fit_error"],
        yerr=fdf_parameters["peak_pi_error"],
        fmt="o",
        lw=1,
        color="red",
        mfc="none",
        label="Fitted peak",
    )

ax.set(
    xlabel=rf"$\phi$ / {u.rad / u.m**2:latex_inline}",
    ylabel="Flux density",
    title="Dirty FDF",
    # xlim=[50, 150],
)
ax.indicate_inset_zoom(axins, edgecolor="black")
ax.legend()
[12]:
<matplotlib.legend.Legend at 0x7cea3666adb0>
../_images/examples_rmsynth_1d_19_1.png
[13]:
phi2_arr_radm2 = rmsf_arrs["phi2_arr_radm2"].to_numpy()
rmsf_arr = rmsf_arrs["rmsf_complex_arr"].to_numpy().astype(complex)

fig, ax = plt.subplots()
ax.plot(
    phi2_arr_radm2,
    rmsf_arr.real,
    color="tab:red",
    label="Stokes Q",
)
ax.plot(
    phi2_arr_radm2,
    rmsf_arr.imag,
    color="tab:blue",
    label="Stokes U",
)
ax.plot(
    phi2_arr_radm2,
    np.abs(rmsf_arr),
    color="k",
    label="Polarized intensity",
)
ax.legend()
ax.set(
    xlabel=rf"$\phi$ / {u.rad / u.m**2:latex_inline}",
    ylabel="RMSF",
    title="RMSF",
)
[13]:
[Text(0.5, 0, '$\\phi$ / $\\mathrm{rad\\,m^{-2}}$'),
 Text(0, 0.5, 'RMSF'),
 Text(0.5, 1.0, 'RMSF')]
../_images/examples_rmsynth_1d_20_1.png

Now lets do a more complex example. We’ll add noise and a Stokes \(I\) spectrum

[14]:
from rm_lite.utils.fitting import power_law
[15]:
delta_rm_radm2 = 30
rm_radm2 = 100
frac_pol = 0.5
psi0_deg = 10
complex_data_noiseless = faraday_slab_spectrum(
    freq_arr_hz=freq_hz,
    frac_pol=frac_pol,
    psi0_deg=psi0_deg,
    rm_radm2=rm_radm2,
    delta_rm_radm2=delta_rm_radm2,
)


stokes_i_flux = 1.0
spectral_index = -0.7
rms_noise = 0.1


stokes_i_model = power_law(order=1)
stokes_i_noiseless = stokes_i_model(
    freq_hz / (np.mean(freq_hz)), stokes_i_flux, spectral_index
)
stokes_i_noise = rng.normal(0, rms_noise, size=freq_hz.size)
stokes_i_noisy = stokes_i_noiseless + stokes_i_noise


stokes_q_noise = rng.normal(0, rms_noise, size=freq_hz.size)
stokes_u_noise = rng.normal(0, rms_noise, size=freq_hz.size)
complex_noise = stokes_q_noise + 1j * stokes_u_noise

complex_flux = complex_data_noiseless * stokes_i_noiseless
complex_data_noisy = complex_data_noiseless + complex_noise

Now we enable Stokes \(I\) model fitting through providing the data, and enabling fit_order. If fit_order<0 an iterative fit will be performed.

[16]:
fdf_parameters, fdf_arrs, rmsf_arrs, stokes_i_arrs = rmsynth.run_rmsynth(
    freq_arr_hz=freq_hz,
    complex_pol_arr=complex_data_noisy,
    complex_pol_error=np.ones_like(complex_data_noiseless)
    * (rms_noise + rms_noise * 1j),
    stokes_i_arr=stokes_i_noisy,
    stokes_i_error_arr=np.ones_like(stokes_i_noisy) * rms_noise,
    do_fit_rmsf=True,
    n_samples=100,
    fit_order=-3,
)
INFO synthesis.create_fractional_spectra: Fitting Stokes I model to calculate fractional spectra.
INFO fitting.dynamic_fit: Iteratively fitting Stokes I model of type log with max order 3.
INFO fitting.static_fit: Fitting Stokes I model of type log with order 0.
INFO fitting.static_fit: Fit results: ['1.04 ± 0.01']
INFO fitting.static_fit: Fitting Stokes I model of type log with order 1.
INFO fitting.static_fit: Fit results: ['1.07 ± 0.01', '-0.67 ± 0.06']
INFO fitting.static_fit: Fitting Stokes I model of type log with order 2.
INFO fitting.static_fit: Fit results: ['1.08 ± 0.02', '-0.62 ± 0.08', '-0.6 ± 0.7']
INFO fitting.static_fit: Fitting Stokes I model of type log with order 3.
INFO fitting.static_fit: Fit results: ['1.08 ± 0.02', '-0.6 ± 0.1', '-1 ± 1', '0 ± 7']
INFO fitting.dynamic_fit: Fit results for orders [0 1 2 3]:
INFO fitting.dynamic_fit: Best fit found with 2 parameters.
INFO synthesis.rmsynth_nufft: Running RM-synthesis using the NUFFTs over 2001 Faraday depth channels.
INFO synthesis.rmsynth_nufft: NUFFT complete in 0.00422 seconds.
INFO synthesis.get_rmsf_nufft: Fitting main lobe in each RMSF spectrum.
INFO rmsynth._run_rmsynth: RM-synthesis completed in 12.44ms.
[17]:
fdf_parameters
[17]:
shape: (1, 30)
fdf_error_madpeak_pi_fitpeak_pi_errorpeak_pi_fit_debiaspeak_pi_fit_snrpeak_pi_fit_indexpeak_rm_fitpeak_rm_fit_errorpeak_q_fitpeak_u_fitpeak_pa_fit_degpeak_pa_fit_deg_errorpeak_pa0_fit_degpeak_pa0_fit_deg_errorfit_functionlam_sq_0_m2ref_freq_hzfwhm_rmsf_radm2fdf_error_noisefdf_q_noisefdf_u_noisemin_freq_hzmax_freq_hzn_channelsmedian_d_freq_hzfrac_polfrac_pol_errorsigma_addsigma_add_minussigma_add_plus
f64f64f64f64f64i64f64f64f64f64f64f64f64f64strf64f64f64f64f64f64f64f64i64f64f64f64f64f64f64
NaN0.2048340.0145930.20363514.036319126589.3604821.229744-0.198321-0.05037497.1259112.04098334.40613311.814679"log"0.0825631.0433e934.5221610.0145930.0146140.0145727.995e81.7995e9548.2286e60.1909990.0136882.0502291.5832332.529693
[18]:
fdf_arrs
[18]:
shape: (2_001, 2)
phi_arr_radm2fdf_dirty_complex_arr
f64object
-336.725921(-0.014721594124084334+0.009366013126583498j)
-336.389195(-0.015252969370126173+0.00936499193961618j)
-336.052469(-0.01577723602372435+0.009356331489117206j)
-335.715743(-0.016293702916400685+0.009339655132559656j)
-335.379017(-0.01680167843043626+0.00931459759667603j)
335.379017(-0.020908748178999275+0.03203839573862219j)
335.715743(-0.02146098917615737+0.03131005570616383j)
336.052469(-0.02198910191668848+0.03058289777475265j)
336.389195(-0.02249294026843618+0.02985815339689203j)
336.725921(-0.02297238921613941+0.029137035926796093j)
[19]:
stokes_i_arrs
[19]:
shape: (54, 7)
freq_arr_hzlambda_sq_arr_m2stokes_i_model_arrstokes_i_model_errorflag_arrcomplex_pol_arrcomplex_pol_error
f64f64f64f64boolobjectobject
7.995e80.1406061.2742020.050221true(0.04341282114362869+0.08189446372194956j)(0.07849913233443775+0.07854682961227628j)
8.0773e80.1377561.2655350.048672true(0.023968486328812204-0.16571533762435478j)(0.07902336149652574+0.07927459640863825j)
8.1596e80.1349921.2569850.047169true(0.024593803941309264-0.027371875227656495j)(0.07956077634127369+0.0795620539228958j)
8.2419e80.132311.2485930.045739true(0.12769660914526085-0.054822075858302365j)(0.08022664072003177+0.08011532118598054j)
8.3241e80.1297071.2403250.044437true(-0.09814605418537536-0.11255065689955869j)(0.08070064971465857+0.08072478194165067j)
1.6555e90.0327930.7822510.049594true(0.533128538400172+0.4419237334961377j)(0.1322291190507992+0.1308705089614891j)
1.6915e90.0314120.7711080.050676true(0.4021200402362769+0.12522731387385583j)(0.13234869583719722+0.12994437124715538j)
1.7275e90.0301170.7603220.051727true(0.5424361210329705-0.1049219011507546j)(0.13660255418742617+0.13171683652493324j)
1.7635e90.02890.7498220.052746true(0.45271234902271473-0.21878346333295814j)(0.13711448607862017+0.1342500872933644j)
1.7995e90.0277550.7397510.053701true(0.5110934680301948-0.3225252980614816j)(0.140179786151518+0.13719325224540285j)
[20]:
fig, ax = plt.subplots()
ax.plot(freq_hz, stokes_i_noiseless, label="Input model")
ax.plot(freq_hz, stokes_i_noisy, ".", label="Noisy data")
ax.plot(
    stokes_i_arrs["freq_arr_hz"],
    stokes_i_arrs["stokes_i_model_arr"],
    "k--",
    label="Fitted model",
)
ax.fill_between(
    stokes_i_arrs["freq_arr_hz"],
    stokes_i_arrs["stokes_i_model_arr"] - stokes_i_arrs["stokes_i_model_error"],
    stokes_i_arrs["stokes_i_model_arr"] + stokes_i_arrs["stokes_i_model_error"],
    alpha=0.3,
    color="k",
    label="Fitted model error",
)
ax.legend()
ax.set(
    xlabel=rf"$\nu$ / {u.Hz:latex_inline}",
    ylabel="Flux density",
    title="Stokes I",
)
[20]:
[Text(0.5, 0, '$\\nu$ / $\\mathrm{Hz}$'),
 Text(0, 0.5, 'Flux density'),
 Text(0.5, 1.0, 'Stokes I')]
../_images/examples_rmsynth_1d_29_1.png
[21]:
phi_arr_radm2 = fdf_arrs["phi_arr_radm2"].to_numpy()
fdf_dirty_arr = fdf_arrs["fdf_dirty_complex_arr"].to_numpy().astype(complex)

fig, ax = plt.subplots()

ax.plot(
    phi_arr_radm2,
    fdf_dirty_arr.real,
    color="tab:red",
    label="Stokes Q",
)
ax.plot(
    phi_arr_radm2,
    fdf_dirty_arr.imag,
    color="tab:blue",
    label="Stokes U",
)
ax.plot(
    phi_arr_radm2,
    np.abs(fdf_dirty_arr),
    color="k",
    label="Polarized intensity",
)

ax.errorbar(
    fdf_parameters["peak_rm_fit"],
    fdf_parameters["peak_pi_fit"],
    xerr=fdf_parameters["peak_rm_fit_error"],
    yerr=fdf_parameters["peak_pi_error"],
    fmt="o",
    lw=1,
    color="red",
    mfc="none",
    label="Fitted peak",
)

ax.set(
    xlabel=rf"$\phi$ / {u.rad / u.m**2:latex_inline}",
    ylabel="Flux density",
    title="Dirty FDF",
)
ax.legend()
[21]:
<matplotlib.legend.Legend at 0x7cea38b64b60>
../_images/examples_rmsynth_1d_30_1.png
[22]:
phi2_arr_radm2 = rmsf_arrs["phi2_arr_radm2"].to_numpy()
rmsf_arr = rmsf_arrs["rmsf_complex_arr"].to_numpy().astype(complex)

fig, ax = plt.subplots()
ax.plot(
    phi2_arr_radm2,
    rmsf_arr.real,
    color="tab:red",
    label="Stokes Q",
)
ax.plot(
    phi2_arr_radm2,
    rmsf_arr.imag,
    color="tab:blue",
    label="Stokes U",
)
ax.plot(
    phi2_arr_radm2,
    np.abs(rmsf_arr),
    color="k",
    label="Polarized intensity",
)
ax.legend()
ax.set(
    xlabel=rf"$\phi$ / {u.rad / u.m**2:latex_inline}",
    ylabel="RMSF",
    title="RMSF",
)
[22]:
[Text(0.5, 0, '$\\phi$ / $\\mathrm{rad\\,m^{-2}}$'),
 Text(0, 0.5, 'RMSF'),
 Text(0.5, 1.0, 'RMSF')]
../_images/examples_rmsynth_1d_31_1.png