rm_lite.utils.synthesis¶
RM-synthesis utils
Attributes¶
Classes¶
Options for RM-synthesis |
|
RM spread function parameters |
|
Results of the RMSF calculation |
|
Parameters for RM-synthesis calculation |
|
Results of the RM-synthesis calculation |
|
Sigma_add complexity metrics |
|
Stokes parameters and errors |
|
Stokes Sigma_add complexity metrics |
|
Theoretical noise of the FDF |
Functions¶
|
Calculate the 2nd moment of the polarised intensity FDF. Can be applied to |
|
Calculate the most likely value of additional scatter, assuming the |
|
|
|
Return the value at a given percentile of a cumulative distribution function |
|
|
|
Calculate the parameters for RM-synthesis. |
|
|
|
|
Create a simple Faraday spectrum with a single component. |
|
|
Convert frequency to lambda^2. |
|
Measure standard parameters from a complex Faraday Dispersion Function. |
|
Calculate the FWHM of the RMSF. |
|
|
|
Compute the RMSF for a given set of lambda^2 values. |
Inverse RM-synthesis - FDF to Stokes Q and U in wavelength^2 space. |
|
|
Convert lambda^2 to frequency. |
|
|
|
Construct a Faraday depth array. |
|
|
|
|
|
Run RM-synthesis on a cube of Stokes Q and U data using the NUFFT method. |
Module Contents¶
- class rm_lite.utils.synthesis.FDFOptions¶
Bases:
NamedTupleOptions for RM-synthesis
- weight_type: Literal['variance', 'uniform'] = 'variance'¶
Weight type
- class rm_lite.utils.synthesis.FWHM¶
Bases:
NamedTuple
- class rm_lite.utils.synthesis.FractionalSpectra¶
Bases:
NamedTuple- fit_result: rm_lite.utils.fitting.FitResult | None¶
- no_nan_idx: numpy.typing.NDArray[numpy.bool_]¶
- stokes_data: StokesData¶
- class rm_lite.utils.synthesis.RMSFParams¶
Bases:
NamedTupleRM spread function parameters
- class rm_lite.utils.synthesis.RMSFResults¶
Bases:
NamedTupleResults of the RMSF calculation
- fit_status_arr: numpy.typing.NDArray[numpy.float64]¶
The status of the RMSF fit
- fwhm_rmsf_arr: numpy.typing.NDArray[numpy.float64]¶
The FWHM of the RMSF main lobe
- phi_double_arr_radm2: numpy.typing.NDArray[numpy.float64]¶
The (double length) Faraday depth array
- rmsf_cube: numpy.typing.NDArray[numpy.float64]¶
The RMSF cube
- class rm_lite.utils.synthesis.RMSynthParams¶
Bases:
NamedTupleParameters for RM-synthesis calculation
- lambda_sq_arr_m2: numpy.typing.NDArray[numpy.float64]¶
Wavelength^2 values in m^2
- phi_arr_radm2: numpy.typing.NDArray[numpy.float64]¶
Faraday depth values in rad/m^2
- weight_arr: numpy.typing.NDArray[numpy.float64]¶
Weight array
- class rm_lite.utils.synthesis.RMsynthResults¶
Bases:
NamedTupleResults of the RM-synthesis calculation
- fdf_dirty_cube: numpy.typing.NDArray[numpy.float64]¶
The Faraday dispersion function cube
- class rm_lite.utils.synthesis.SigmaAdd¶
Bases:
NamedTupleSigma_add complexity metrics
- sigma_add_arrays: SigmaAddArrays¶
Sigma_add arrays
- class rm_lite.utils.synthesis.SigmaAddArrays¶
Bases:
NamedTuple- cdf: numpy.typing.NDArray[numpy.float64]¶
CDF array of the additional noise term
- pdf: numpy.typing.NDArray[numpy.float64]¶
PDF array of the additional noise term
- sigma_add_arr: numpy.typing.NDArray[numpy.float64]¶
Array of additional noise values
- class rm_lite.utils.synthesis.StokesData¶
Bases:
NamedTupleStokes parameters and errors
- with_options(**kwargs)¶
Create a new StokesData instance with keywords updated
- Returns:
New StokesData instance with updated attributes
- Return type:
- complex_pol_arr: numpy.typing.NDArray[numpy.complex128]¶
Stokes Q and U array
- complex_pol_error: numpy.typing.NDArray[numpy.complex128]¶
Stokes Q and U error array
- freq_arr_hz: numpy.typing.NDArray[numpy.float64]¶
Frequency array in Hz
- class rm_lite.utils.synthesis.StokesSigmaAdd¶
Bases:
NamedTupleStokes Sigma_add complexity metrics
- class rm_lite.utils.synthesis.TheoreticalNoise¶
Bases:
NamedTupleTheoretical noise of the FDF
- with_options(**kwargs)¶
Create a new TheoreticalNoise instance with keywords updated
- Returns:
New TheoreticalNoise instance with updated attributes
- Return type:
- rm_lite.utils.synthesis.calc_mom2_fdf(complex_fdf_arr: numpy.typing.NDArray[numpy.complex128], phi_arr_radm2: numpy.typing.NDArray[numpy.float64]) float¶
Calculate the 2nd moment of the polarised intensity FDF. Can be applied to a clean component spectrum or a standard FDF
- rm_lite.utils.synthesis.calculate_sigma_add(y_arr: numpy.typing.NDArray[numpy.float64], dy_arr: numpy.typing.NDArray[numpy.float64], median: float | None = None, noise: float | None = None, n_samples: int = 1000) SigmaAdd¶
Calculate the most likely value of additional scatter, assuming the input data is drawn from a normal distribution. The total uncertainty on each data point Y_i is modelled as dYtot_i**2 = dY_i**2 + dYadd**2.
- rm_lite.utils.synthesis.calculate_sigma_add_arr(y_arr: numpy.typing.NDArray[numpy.float64], dy_arr: numpy.typing.NDArray[numpy.float64], median: float | None = None, noise: float | None = None, n_samples: int = 1000) SigmaAddArrays¶
- rm_lite.utils.synthesis.cdf_percentile(values: numpy.typing.NDArray[numpy.float64], cdf: numpy.typing.NDArray[numpy.float64], q=50.0) float¶
Return the value at a given percentile of a cumulative distribution function
- rm_lite.utils.synthesis.compute_rmsf_params(freq_arr_hz: numpy.typing.NDArray[numpy.float64], weight_arr: numpy.typing.NDArray[numpy.float64]) RMSFParams¶
- rm_lite.utils.synthesis.compute_rmsynth_params(freq_arr_hz: numpy.typing.NDArray[numpy.float64], complex_pol_arr: numpy.typing.NDArray[numpy.complex128], complex_pol_error: numpy.typing.NDArray[numpy.complex128], fdf_options: FDFOptions) RMSynthParams¶
Calculate the parameters for RM-synthesis.
- Parameters:
freq_arr_hz (NDArray[np.float64]) – Frequency array in Hz
pol_arr (NDArray[np.complex128]) – Complex polarisation array
real_qu_error (NDArray[np.float64 | np.float32]) – Error in Stokes Q and U (real)
fdf_options (FDFOptions) – Options for RM-synthesis
- Raises:
ValueError – If d_phi_radm2 is not provided and n_samples is None.
- Returns:
Wavelength^2 values, reference wavelength^2, Faraday depth values, weight array
- Return type:
- rm_lite.utils.synthesis.compute_theoretical_noise(complex_pol_error: numpy.typing.NDArray[numpy.complex128], weight_arr: numpy.typing.NDArray[numpy.float64]) TheoreticalNoise¶
- rm_lite.utils.synthesis.create_fractional_spectra(stokes_data: StokesData, ref_freq_hz: float, fit_order: int = 2, fit_function: Literal['log', 'linear'] = 'log', n_error_samples: int = 10000) FractionalSpectra | None¶
- rm_lite.utils.synthesis.faraday_simple_spectrum(freq_arr_hz: numpy.typing.NDArray[numpy.float64], frac_pol: float, psi0_deg: float, rm_radm2: float) numpy.typing.NDArray[numpy.complex128]¶
Create a simple Faraday spectrum with a single component.
- rm_lite.utils.synthesis.freq_to_lambda2(freq_hz: T) T¶
Convert frequency to lambda^2.
- rm_lite.utils.synthesis.get_fdf_parameters(fdf_arr: numpy.typing.NDArray[numpy.complex128], phi_arr_radm2: numpy.typing.NDArray[numpy.float64], fwhm_rmsf_radm2: float, freq_arr_hz: numpy.typing.NDArray[numpy.float64], complex_pol_arr: numpy.typing.NDArray[numpy.complex128], complex_pol_error: numpy.typing.NDArray[numpy.complex128], lambda_sq_arr_m2: numpy.typing.NDArray[numpy.float64], lam_sq_0_m2: float, stokes_i_reference_flux: float, theoretical_noise: TheoreticalNoise, fit_function: Literal['log', 'linear'], bias_correction_snr: float = 5.0) polars.DataFrame¶
Measure standard parameters from a complex Faraday Dispersion Function. Currently this function assumes that the noise levels in the Stokes Q and U spectra are the same. Returns a dictionary containing measured parameters.
- rm_lite.utils.synthesis.get_fwhm_rmsf(lambda_sq_arr_m2: numpy.typing.NDArray[numpy.float64]) FWHM¶
Calculate the FWHM of the RMSF.
- Parameters:
lambda_sq_arr_m2 (NDArray[np.float64]) – Wavelength^2 values in m^2
super_resolution (bool, optional) – Use Cotton+Rudnick superresolution. Defaults to False.
- Returns:
FWHM of the RMSF main lobe, maximum difference in lambda^2 values, range of lambda^2 values
- Return type:
fwhm_rmsf_arr
- rm_lite.utils.synthesis.get_mask_index(stokes_data: StokesData) numpy.typing.NDArray[numpy.bool_]¶
- rm_lite.utils.synthesis.get_rmsf_nufft(lambda_sq_arr_m2: numpy.typing.NDArray[numpy.float64], phi_arr_radm2: numpy.typing.NDArray[numpy.float64], weight_arr: numpy.typing.NDArray[numpy.float64], lam_sq_0_m2: float, mask_arr: numpy.typing.NDArray[numpy.bool_] | None = None, do_fit_rmsf: bool = False, do_fit_rmsf_real=False, eps: float = 1e-06) RMSFResults¶
Compute the RMSF for a given set of lambda^2 values.
- Parameters:
lambda_sq_arr_m2 (NDArray[np.float64]) – Wavelength^2 values in m^2
phi_arr_radm2 (NDArray[np.float64]) – Faraday depth values in rad/m^2
weight_arr (NDArray[np.float64]) – Weight array
lam_sq_0_m2 (float) – Reference wavelength^2 value
super_resolution (bool, optional) – Use superresolution. Defaults to False.
mask_arr (Optional[NDArray[np.float64]], optional) – Mask array. Defaults to None.
do_fit_rmsf (bool, optional) – Fit the RMSF with a Gaussian. Defaults to False.
do_fit_rmsf_real (bool, optional) – Fit the real part of the. Defaults to False.
eps (float, optional) – NUFFT tolerance. Defaults to 1e-6.
- Raises:
ValueError – If the wavelength^2 and weight arrays are not the same shape.
ValueError – If the mask dimensions are > 3.
ValueError – If the mask depth does not match the lambda^2 vector.
- Returns:
rmsf_cube, phi_double_arr_radm2, fwhm_rmsf_arr, fit_status_arr
- Return type:
- rm_lite.utils.synthesis.inverse_rmsynth_nufft(complex_fdf_arr: numpy.typing.NDArray[numpy.complex128], lambda_sq_arr_m2: numpy.typing.NDArray[numpy.float64], phi_arr_radm2: numpy.typing.NDArray[numpy.float64], lam_sq_0_m2: float, eps: float = 1e-06) numpy.typing.NDArray[numpy.complex128]¶
Inverse RM-synthesis - FDF to Stokes Q and U in wavelength^2 space.
- Parameters:
complex_fdf_arr (NDArray[np.complex128]) – Complex polarisation array in Faraday depth space
lambda_sq_arr_m2 (NDArray[np.float64]) – Wavelength^2 values in m^2
phi_arr_radm2 (NDArray[np.float64]) – Faraday depth values in rad/m^2
lam_sq_0_m2 (float) – Reference wavelength^2 value
eps (float, optional) – NUFFT tolerance. Defaults to 1e-6.
- Raises:
ValueError – If the Stokes Q and U data arrays are not the same shape.
ValueError – If the data dimensions are > 3.
ValueError – If the data depth does not match the lambda^2 vector.
- Returns:
Complex polarisation array in wavelength^2 space
- Return type:
NDArray[np.float64]
- rm_lite.utils.synthesis.lambda2_to_freq(lambda_sq_m2: T) T¶
Convert lambda^2 to frequency.
- Parameters:
lambda_sq_m2 (NDArray[np.float64]) – Wavelength^2 in m^2
- Returns:
Frequency in Hz
- Return type:
NDArray[np.float64]
- rm_lite.utils.synthesis.make_double_phi_arr(phi_arr_radm2: numpy.typing.NDArray[numpy.float64]) numpy.typing.NDArray[numpy.float64]¶
- rm_lite.utils.synthesis.make_phi_arr(phi_max_radm2: float, d_phi_radm2: float) numpy.typing.NDArray[numpy.float64]¶
Construct a Faraday depth array.
- rm_lite.utils.synthesis.measure_fdf_complexity(phi_arr_radm2: numpy.typing.NDArray[numpy.float64], complex_fdf_arr: numpy.typing.NDArray[numpy.complex128]) float¶
- rm_lite.utils.synthesis.measure_qu_complexity(freq_arr_hz: numpy.typing.NDArray[numpy.float64], complex_pol_arr: numpy.typing.NDArray[numpy.complex128], complex_pol_error: numpy.typing.NDArray[numpy.complex128], frac_pol: float, psi0_deg: float, rm_radm2: float) StokesSigmaAdd¶
- rm_lite.utils.synthesis.rmsynth_nufft(complex_pol_arr: numpy.typing.NDArray[numpy.complex128], lambda_sq_arr_m2: numpy.typing.NDArray[numpy.float64], phi_arr_radm2: numpy.typing.NDArray[numpy.float64], weight_arr: numpy.typing.NDArray[numpy.float64], lam_sq_0_m2: float, eps: float = 1e-06) numpy.typing.NDArray[numpy.complex128]¶
Run RM-synthesis on a cube of Stokes Q and U data using the NUFFT method.
- Parameters:
complex_pol_arr (NDArray[np.complex128]) – Complex polarisation values (Q + iU)
lambda_sq_arr_m2 (NDArray[np.float64]) – Wavelength^2 values in m^2
phi_arr_radm2 (NDArray[np.float64]) – Faraday depth values in rad/m^2
weight_arr (NDArray[np.float64]) – Weight array
lam_sq_0_m2 (Optional[float], optional) – Reference wavelength^2 in m^2. Defaults to None.
eps (float, optional) – NUFFT tolerance. Defaults to 1e-6.
- Raises:
ValueError – If the weight and lambda^2 arrays are not the same shape.
ValueError – If the Stokes Q and U data arrays are not the same shape.
ValueError – If the data dimensions are > 3.
ValueError – If the data depth does not match the lambda^2 vector.
- Returns:
Dirty Faraday dispersion function cube
- Return type:
NDArray[np.float64]
- rm_lite.utils.synthesis.T¶
- rm_lite.utils.synthesis.fdf_params_schema¶
- rm_lite.utils.synthesis.fdf_params_schema_df¶