Reference on kuibit.gw_mismatch¶
The gw_mismatch
module functions to compute the mismatch between
two waves using a simple grid search for phase and time shifts. Since no
polarization shifts are performed, this is only relevant to the gravitational
wave 2,2 mode.
The two main interfaces are network_mismatch_from_psi4()
(when
computing the network mismatch starting from psi4 and the sky localization) and
mismatch_from_strains()
(when computing the mismatch from the
strains).
..warning:
Make sure to understand what is going on if you are using this module. You should read the code and comments in the code.
- kuibit.gw_mismatch.mismatch_from_strains(h1, h2, fmin=0, fmax=inf, noises=None, antenna_patterns=None, num_polarization_shifts=100, num_time_shifts=100, time_shift_start=-5, time_shift_end=5, force_numba=False)[source]¶
Compute the network-mismatch between
h1
andh2
by maximizing the overlap over time and polarization shifts.Network here means that the inner product is computed for N detectors, as provided by the lists antenna_patterns and noises. Noises and antenna patterns have to be properly ordered:
noises[i]
has to correspond toantenna_pattern[i]
.See Overlap and mismatch for formulas and details.
The mismatch is computed by maximizing over time and polarization shifts. Polarization shifts and are around the 2pi, time shifts are specified by time_shift_start and time_shift_end. If num_time_shifts is 1, then no time shift is performed. For times, we make sure that we always have to zero timeshift. All the transformations are done in h2.
This computation is a maximisation, which is very expensive. So, we have a very fast core function called _mismatch_core_numerical to do all the hard work. This function is compiled to native code by numba, resulting to enormous speed-up.There is an overhead in calling numba. Hence, by default we do not always use numba. We use it only when the number num_polarization_shifts * num_time_shifts is greater than 500*500. You can force using numba passing the keyword argument force_numba=True.
We do not perform phase shifts here, so this function makes sense only for the (2,2) mode.
h1 and h2 have to be already pre-processed for Fourier transform, so you should window them and zero pad as needed.
- Parameters:
h1 (
TimeSeries
) – First strain.h2 (
TimeSeries
) – Second strain (the one that will be modified).fmin (float) – Lower limit of the integration.
fmax (float) – Higher limit of the integration.
noises (list of
FrequencySeries
, or None) – Power spectral density of the noise for all the detectors. If None, a uniform noise is applied.antenna_patterns (list of tuples, or None) – Fc, Fp for all the detectors. It has to be ordered in the same way as noises. If None, a uniform antenna pattern is applied.
num_polarization_shifts (int) – How many points to divide the range (0, 2 pi) in the polarization shift.
num_time_shifts (int) – How many points to divide the range (time_shift_start, time_shift_end) in the time shift.
time_shift_start (float) – Minimum time shift applied. Search will be done linearly up to time_shift_end.
time_shift_end (float) – Largest value of time shift applied.
force_numba (bool) – Use numba irrespectively of the size of the input.
- kuibit.gw_mismatch.network_mismatch(h1, h2, right_ascension, declination, time_utc, fmin=0, fmax=inf, noises=None, num_polarization_shifts=200, num_time_shifts=1000, time_shift_start=-20, time_shift_end=20, force_numba=False)[source]¶
Compute network mismatch between strains h1 and h2.
This is a wrapper around
mismatch_from_strains()
(read the docstring for more information) that prepares the correct antenna patterns from the sky localization of the event. Moreover, it makes sure that noises (that has to be a gw_utils.Detectors object, or None) is correctly ordered.- Parameters:
h1 (
TimeSeries
) – First strain.h2 (
TimeSeries
) – Second strain (the one that will be modified).right_ascension (float) – Right ascension of the source in the sky.
declination (float) – Declination of the source in the sky.
time_utc (float) – Time UTC of the event.
fmin (float) – Lower limit of the integration.
fmax (float) – Higher limit of the integration.
noises (
Detector
, or None) – Power spectral density of the noise for all the detectors.num_polarization_shifts (int) – How many points to divide the range (0, 2 pi) in the polarization shift.
num_time_shifts (int) – How many points to divide the range (time_shift_start, time_shift_end) in the time shift.
time_shift_start (float) – Minimum time shift applied. Search will be done linearly up to time_shift_end.
time_shift_end (float) – Largest value of time shift applied.
force_numba (bool) – Use numba irrespectively of the size of the input.
- kuibit.gw_mismatch.network_mismatch_from_psi4(psi1, psi2, right_ascension, declination, time_utc, pcut1, pcut2, *args, window_function=None, align_at_peak=True, trim_ends=False, mass_scale1_msun=None, mass_scale2_msun=None, distance1=None, distance2=None, fmin=0, fmax=inf, noises=None, num_zero_pad=16384, num_time_shifts=100, num_polarization_shifts=100, time_shift_start=-50, time_shift_end=50, force_numba=False, time_to_keep_after_max=None, time_removed_beginning=None, **kwargs)[source]¶
Compute the mismatch between strains from psi1 and psi2.
This function is complex, read Overlap and mismatch for formulas and details.
- Parameters:
psi1 (
GravitationalWavesOneDet
) – \(\Psi_4\) for the first wave.psi2 (
GravitationalWavesOneDet
) – \(\Psi_4\) for the second wave (the one that will be modified).right_ascension (float) – Right ascension of the source in the sky.
declination (float) – Declination of the source in the sky.
mass_scale1_msun (float or None) – If not None, the signal h1 is converted from computational units to physical units assuming that
M = mass_scale1_msun
.mass_scale2_msun (float or None) – If not None, the signal h2 is converted from computational units to physical units assuming that M = mass_scale2_msun.
time_utc (float) – Time UTC of the event.
pcut1 (float) – Period associated with the threshold frequency
omega_0 = 2 * pi / pcut
for the fixed frequency integration of \(\Psi_4\)pcut2 (float) – Period associated with the threshold frequency
omega_0 = 2 * pi / pcut
for the fixed frequency integration of \(\Psi_4\)window_function (callable, str, or None) – If not None, apply window_function to the series before computing the strain.
trim_ends (bool) – If True, a portion of the resulting strain is removed at both the initial and final times. The amount removed is equal to pcut.
align_at_peak (bool) – Time-shifts the strain so that they have both the maximum amplitude at
t=0
.fmin (float) – Lower limit of the integration.
fmax (float) – Higher limit of the integration.
noises (
Detectors
, or None) – Power spectral density of the noise for all the detectors. If None, a uniform noise is applied.num_zero_pad (int) – How many points do the timeseries have to have before Fourier transforms are taken? This is not the number of zeros added (at the end) of the timeseries, this is the total number. If the series already have that length, no zeros will be added.
num_polarization_shifts (int) – How many points to divide the range (0, 2 pi) in the polarization shift.
num_time_shifts (int) – How many points to divide the range (time_shift_start, time_shift_end) in the time shift.
time_shift_start (float) – Minimum time shift applied. Search will be done linearly up to time_shift_end.
time_shift_end (float) – Largest value of time shift applied.
time_removed_beginning (float or None) – Remove this amount from the beginning of the strain signals before computing the mismatch. If None, nothing is removed.
time_to_keep_after_max (float or None) – If not None, remove all the signal that comes after t_max + time_to_keep_after_max, where t_max is the time at which the signal peaks.
force_numba (bool) – Use numba irrespectively of the size of the input.
args (anything) – All the other arguments are passed to the window function.
- kuibit.gw_mismatch.one_detector_mismatch_from_psi4(psi1, psi2, pcut1, pcut2, *args, window_function=None, align_at_peak=True, trim_ends=False, mass_scale1_msun=None, mass_scale2_msun=None, distance1=None, distance2=None, fmin=0, fmax=inf, noise=None, num_zero_pad=16384, num_time_shifts=100, num_polarization_shifts=100, time_shift_start=-50, time_shift_end=50, force_numba=False, time_to_keep_after_max=None, time_removed_beginning=None, **kwargs)[source]¶
Compute the mismatch between strains from psi1 and psi2.
This function is complex, read Overlap and mismatch for formulas and details.
- Parameters:
psi1 (
GravitationalWavesOneDet
) – \(\Psi_4\) for the first wave.psi2 (
GravitationalWavesOneDet
) – \(\Psi_4\) for the second wave (the one that will be modified).right_ascension (float) – Right ascension of the source in the sky.
declination (float) – Declination of the source in the sky.
mass_scale1_msun (float or None) – If not None, the signal h1 is converted from computational units to physical units assuming that M = mass_scale1_msun.
mass_scale2_msun (float or None) – If not None, the signal h2 is converted from computational units to physical units assuming that M = mass_scale2_msun.
time_utc (float) – Time UTC of the event.
pcut1 (float) – Period associated with the threshold frequency
omega_0 = 2 * pi / pcut
for the fixed frequency integration of \(\Psi_4\)pcut2 (float) – Period associated with the threshold frequency
omega_0 = 2 * pi / pcut
for the fixed frequency integration of \(\Psi_4\)window_function (callable, str, or None) – If not None, apply window_function to the series before computing the strain.
trim_ends (bool) – If True, a portion of the resulting strain is removed at both the initial and final times. The amount removed is equal to pcut.
align_at_peak (bool) – Time-shifts the strain so that they have both the maximum amplitude at
t=0
.fmin (float) – Lower limit of the integration.
fmax (float) – Higher limit of the integration.
noises (
Detectors
, or None) – Power spectral density of the noise for all the detectors. If None, a uniform noise is applied.num_zero_pad (int) – How many points do the timeseries have to have before Fourier transforms are taken? This is not the number of zeros added (at the end) of the timeseries, this is the total number. If the series already have that length, no zeros will be added.
num_polarization_shifts (int) – How many points to divide the range (0, 2 pi) in the polarization shift.
num_time_shifts (int) – How many points to divide the range (time_shift_start, time_shift_end) in the time shift.
time_shift_start (float) – Minimum time shift applied. Search will be done linearly up to time_shift_end.
time_shift_end (float) – Largest value of time shift applied.
time_removed_beginning (float or None) – Remove this amount from the beginning of the strain signals before computing the mismatch. If None, nothing is removed.
time_to_keep_after_max (float or None) – If not None, remove all the signal that comes after t_max + time_to_keep_after_max, where t_max is the time at which the signal peaks.
force_numba (bool) – Use numba irrespectively of the size of the input.
args (anything) – All the other arguments are passed to the window function.