Diagnostics

lisatools.diagnostic.inner_product(sig1: ndarray | list | DataResidualArray, sig2: ndarray | list | DataResidualArray, dt: float | None = None, df: float | None = None, f_arr: float | None = None, psd: str | None | ndarray | SensitivityMatrix = 'LISASens', psd_args: tuple | None = (), psd_kwargs: dict | None = {}, normalize: bool | str | None = False, complex: bool | None = False) float | complex

Compute the inner product between two signals weighted by a psd.

The inner product between time series \(a(t)\) and \(b(t)\) is

\[\langle a | b \rangle = 2\int_{f_\text{min}}^{f_\text{max}} \frac{\tilde{a}(f)^*\tilde{b}(f) + \tilde{a}(f)\tilde{b}(f)^*}{S_n(f)} df\ \ ,\]

where \(\tilde{a}(f)\) is the Fourier transform of \(a(t)\) and \(S_n(f)\) is the one-sided Power Spectral Density of the noise.

The inner product can be left complex using the complex kwarg.

GPU Capability: Pass CuPy arrays rather than NumPy arrays.

Parameters:
  • sig1 – First signal to use for the inner product. Can be time-domain or frequency-domain. Must be 1D np.ndarray, list of 1D np.ndarray``s, or 2D ``np.ndarray across channels with shape (nchannels, data length).

  • sig2 – Second signal to use for the inner product. Can be time-domain or frequency-domain. Must be 1D np.ndarray, list of 1D np.ndarray``s, or 2D ``np.ndarray across channels with shape (nchannels, data length).

  • dt – Time step in seconds. If provided, assumes time-domain signals.

  • df – Constant frequency spacing. This will assume a frequency domain signal with constant frequency spacing.

  • f_arr – Array of specific frequencies at which the signal is given.

  • psd – Indicator of what psd to use. If a str, this will be passed as the sens_fn kwarg to get_sensitivity(). If None, it will be an array of ones. Or, you can pass a 1D np.ndarray of psd values that must be the same length as the frequency domain signals.

  • psd_args – Arguments to pass to the psd function if type(psd) == str.

  • psd_kwargs – Keyword arguments to pass to the psd function if type(psd) == str.

  • normalize – Normalize the inner product. If True, it will normalize the square root of the product of individual signal inner products. You can also pass "sig1" or "sig2" to normalize with respect to one signal.

  • complex – If True, return the complex value of the inner product rather than just its real-valued part.

Returns:

Inner product value.

lisatools.diagnostic.residual_source_likelihood_term(data_res_arr: DataResidualArray, **kwargs: dict) float | complex

Calculate the source term in the Likelihood for a data residual (d - h).

The source term in the likelihood is given by,

\[\log{\mathcal{L}}_\text{src} = -\frac{1}{2}\langle \vec{d} - \vec{h} | \vec{d} - \vec{h}\rangle.\]
Parameters:
  • data_res_arr – Data residual.

  • **kwargs – Keyword arguments to pass to inner_product().

Returns:

Source term Likelihood value.

lisatools.diagnostic.noise_likelihood_term(psd: SensitivityMatrix) float

Calculate the noise term in the Likelihood.

The noise term in the likelihood is given by,

\[\log{\mathcal{L}}_n = -\sum \log{\vec{S}_n}.\]
Parameters:

psd – Sensitivity information.

Returns:

Noise term Likelihood value.

lisatools.diagnostic.residual_full_source_and_noise_likelihood(data_res_arr: DataResidualArray, psd: str | None | ndarray | SensitivityMatrix, **kwargs: dict) float | complex

Calculate the full Likelihood including noise and source terms.

The noise term is calculated with noise_likelihood_term().

The source term is calcualted with residual_source_likelihood_term().

Parameters:
  • data_res_arr – Data residual.

  • psd – Sensitivity information.

  • **kwargs – Keyword arguments to pass to inner_product().

Returns:

Full Likelihood value.

lisatools.diagnostic.data_signal_source_likelihood_term(data_arr: DataResidualArray, sig_arr: DataResidualArray, **kwargs: dict) float | complex

Calculate the source term in the Likelihood for separate signal and data.

The source term in the likelihood is given by,

\[\log{\mathcal{L}}_\text{src} = -\frac{1}{2}\left(\langle \vec{d} | \vec{d}\rangle + \langle \vec{h} | \vec{h}\rangle - 2\langle \vec{d} | \vec{h}\rangle \right)\ \ .\]
Parameters:
  • data_arr – Data.

  • sig_arr – Signal.

  • **kwargs – Keyword arguments to pass to inner_product().

Returns:

Source term Likelihood value.

lisatools.diagnostic.data_signal_full_source_and_noise_likelihood(data_arr: DataResidualArray, sig_arr: DataResidualArray, psd: str | None | ndarray | SensitivityMatrix, **kwargs: dict) float | complex

Calculate the full Likelihood including noise and source terms.

Here, the signal is treated separate from the data.

The noise term is calculated with noise_likelihood_term().

The source term is calcualted with data_signal_source_likelihood_term().

Parameters:
  • data_arr – Data.

  • sig_arr – Signal.

  • psd – Sensitivity information.

  • **kwargs – Keyword arguments to pass to inner_product().

Returns:

Full Likelihood value.

lisatools.diagnostic.snr(sig1: ndarray | list | DataResidualArray, *args: Any, data: ndarray | list | DataResidualArray | None = None, **kwargs: Any) float

Compute the snr between two signals weighted by a psd.

The signal-to-noise ratio of a signal is \(\sqrt{\langle a|a\rangle}\).

This will be the optimal SNR if data==None. If a data array is given, it will be the observed SNR: \(\langle d|a\rangle/\sqrt{\langle a|a\rangle}\).

GPU Capability: Pass CuPy arrays rather than NumPy arrays.

Parameters:
  • sig1 – Signal to use as the templatefor the SNR. Can be time-domain or frequency-domain. Must be 1D np.ndarray, list of 1D np.ndarray``s, or 2D ``np.ndarray across channels with shape (nchannels, data length).

  • *args – Arguments to pass to inner_product().

  • data – Data becomes the sig2 argument to inner_product().

  • **kwargs – Keyword arguments to pass to inner_product().

Returns:

Optimal or detected SNR value (depending on data kwarg).

lisatools.diagnostic.h_var_p_eps(step: float, waveform_model: callable, params: ndarray | list, index: int, parameter_transforms: TransformContainer | None = None, waveform_args: tuple | None = (), waveform_kwargs: dict | None = {}) ndarray

Calculate the waveform with a perturbation step of the variable V[i]

Parameters:
  • waveform_model – Callable function to the waveform generator with signature (*params, **waveform_kwargs).

  • params – Source parameters that are over derivatives (not in fill dict of parameter transforms)

  • step – Absolute step size for variable of interest.

  • index – Index to parameter of interest.

  • parameter_transformsTransformContainer object to transform from the derivative parameter basis to the waveform parameter basis. This class can also fill in fixed parameters where the derivatives are not being taken.

  • waveform_args – args (beyond parameters) for the waveform generator.

  • waveform_kwargs – kwargs for the waveform generation.

Returns:

Perturbation to the waveform in the given parameter. Will always be 2D array with shape (num channels, data length)

lisatools.diagnostic.dh_dlambda(eps: float, *args: tuple, more_accurate: bool | None = True, **kwargs: dict) ndarray

Derivative of the waveform

Calculate the derivative of the waveform with precision of order (step^4) with respect to the variable V in the i direction.

Parameters:
  • eps – Absolute derivative step size for variable of interest.

  • *args – Arguments passed to h_var_p_eps().

  • more_accurate – If True, run a more accurate derivate requiring 2x more waveform generations.

  • **kwargs – Keyword arguments passed to h_var_p_eps().

Returns:

Numerical derivative of the waveform with respect to a varibale of interest. Will be 2D array with shape (num channels, data length).

lisatools.diagnostic.info_matrix(eps: float | ndarray, waveform_model: callable, params: ndarray | list, deriv_inds: ndarray | list | None = None, inner_product_kwargs: dict | None = {}, return_derivs: bool | None = False, **kwargs: dict) ndarray | Tuple[ndarray, list]

Calculate Information Matrix.

This calculates the information matrix for a given waveform model at a given set of parameters. The inverse of the information matrix gives the covariance matrix.

This is also referred to as the Fisher information matrix, but @MichaelKatz has chosen to leave out the name because of this.

The info matrix is given by:

\[M_{ij} = \langle h_i | h_j \rangle \text{ with } h_i = \frac{\partial h}{\partial \lambda_i}.\]
Parameters:
  • eps – Absolute derivative step size for variable of interest. Can be provided as a float value that applies to all variables or an array, one for each parameter being evaluated in the information matrix.

  • waveform_model – Callable function to the waveform generator with signature (*params, **waveform_kwargs).

  • params – Source parameters.

  • deriv_inds – Subset of parameters of interest for which to calculate the information matrix, by index. If None, it will be np.arange(len(params)).

  • inner_product_kwargs – Keyword arguments for the inner product function.

  • return_derivs – If True, also returns computed numerical derivatives.

  • **kwargs – Keyword arguments passed to dh_dlambda().

Returns:

If not return_derivs, this will be the information matrix as a numpy array. If return_derivs is True, it will be a tuple with the first entry as the information matrix and the second entry as the partial derivatives.

lisatools.diagnostic.covariance(*info_mat_args: tuple, info_mat: ndarray | None = None, diagonalize: bool | None = False, return_info_mat: bool | None = False, precision: bool | None = False, **info_mat_kwargs: dict) ndarray | list

Calculate covariance matrix for a set of EMRI parameters, computing the information matrix if not supplied.

Parameters:
  • *info_mat_args – Set of arguments to pass to info_matrix(). Not required if info_mat is not None.

  • info_mat – Pre-computed information matrix. If supplied, this matrix will be inverted.

  • diagonalize – If True, diagonalizes the covariance matrix.

  • return_info_mat – If True, also returns the computed information matrix.

  • precision – If True, uses 500-dps precision to compute the information matrix inverse (requires mpmath). This is typically a good idea as the information matrix can be highly ill-conditioned.

  • **info_mat_kwargs – Keyword arguments to pass to info_matrix().

Returns:

Covariance matrix. If return_info_mat is True. A list will be returned with the covariance as the first entry and the information matrix as the second entry. If return_derivs is True (keyword argument to info_matrix()), then another entry will be added to the list for the derivatives.

lisatools.diagnostic.plot_covariance_corner(params: ndarray, cov: ndarray, nsamp: int | None = 25000, fig: Figure | None = None, **kwargs: dict) Figure

Construct a corner plot for a given covariance matrix.

The corner module is required for this.

Parameters:
  • params – The set of parameters used for the event (the mean vector of the covariance matrix).

  • cov – Covariance matrix from which to construct the corner plot.

  • nsamp – Number of samples to draw from the multivariate distribution.

  • fig – Matplotlib plt.Figure object. Use this if passing an existing corner plot figure.

  • **kwargs – Keyword arguments for the corner plot - see the module documentation for more info.

Returns:

The corner plot figure.

lisatools.diagnostic.plot_covariance_contour(params: ndarray, cov: ndarray, horizontal_index: int, vertical_index: int, nsamp: int | None = 25000, ax: Axes | None = None, **kwargs: dict) Axes | Tuple[Figure, Axes]

Construct a contour plot for a given covariance matrix on a single axis object.

The corner module is required for this.

Parameters:
  • params – The set of parameters used for the event (the mean vector of the covariance matrix).

  • cov – Covariance matrix from which to construct the corner plot.

  • horizontal_index – Parameter index to plot along the horizontal axis of the contour plot.

  • vertical_index – Parameter index to plot along the vertical axis of the contour plot.

  • nsamp – Number of samples to draw from the multivariate distribution.

  • fig – Matplotlib plt.Figure object. Use this if passing an existing corner plot figure.

  • **kwargs – Keyword arguments for the corner plot - see the module documentation for more info.

Returns:

If ax is provided, the return will be that ax object. If it is not provided, then a Matplotlib Figure and Axes obejct is created and returned as a tuple: (plt.Figure, plt.Axes).

lisatools.diagnostic.get_eigeninfo(arr: ndarray, high_precision: bool | None = False) Tuple[ndarray, ndarray]

Performs eigenvalue decomposition and returns the eigenvalues and right-eigenvectors for the supplied fisher/covariance matrix.

Parameters:
  • arr – Input matrix for which to perform eigenvalue decomposition.

  • high_precision

    If True, use 500-dps precision to ensure accurate eigenvalue decomposition (requires mpmath to be installed). Defaults to False.

Returns:

,k] corresponds to evals[k].

Return type:

Tuple containing Eigenvalues and right-Eigenvectors for the supplied array, constructed such that evects[

lisatools.diagnostic.cutler_vallisneri_bias(waveform_model_true: callable, waveform_model_approx: callable, params: ndarray, eps: float | ndarray, input_diagnostics: dict | None = None, info_mat: ndarray | None = None, deriv_inds: ndarray | list | None = None, return_derivs: bool | None = False, return_cov: bool | None = False, parameter_transforms: TransformContainer | None = None, waveform_true_args: tuple | None = (), waveform_true_kwargs: dict | None = {}, waveform_approx_args: tuple | None = (), waveform_approx_kwargs: dict | None = {}, inner_product_kwargs: dict | None = {}) list

Calculate the Cutler-Vallisneri bias.

# TODO: add basic math

Parameters:
  • waveform_model_true – Callable function to the true waveform generator with signature (*params, **waveform_kwargs).

  • waveform_model_approx – Callable function to the approximate waveform generator with signature (*params, **waveform_kwargs).

  • params – Source parameters.

  • eps – Absolute derivative step size. See info_matrix().

  • input_diagnostics – Dictionary including the diagnostic information if it is precomputed. Dictionary must include keys "cov" (covariance matrix, output of covariance()), "h_true" (the true waveform), and "dh" (derivatives of the waveforms, list of outputs from dh_dlambda()).

  • info_mat – Pre-computed information matrix. If supplied, this matrix will be inverted to find the covariance.

  • deriv_inds – Subset of parameters of interest. See info_matrix().

  • return_derivs – If True, also returns computed numerical derivatives.

  • return_cov – If True, also returns computed covariance matrix.

  • parameter_transforms

    TransformContainer object. See info_matrix().

  • waveform_true_args – Arguments for the true waveform generator.

  • waveform_true_kwargs – Keyword arguments for the true waveform generator.

  • waveform_approx_args – Arguments for the approximate waveform generator.

  • waveform_approx_kwargs – Keyword arguments for the approximate waveform generator.

  • inner_product_kwargs – Keyword arguments for the inner product function.

Returns:

List of return information. By default, it is [systematic error, bias]. If return_derivs or return_cov are True, they will be added to the list with derivs added before covs.

lisatools.diagnostic.scale_to_snr(target_snr: float, sig: ndarray | list, *snr_args: tuple, return_orig_snr: bool | None = False, **snr_kwargs: dict) ndarray | list | Tuple[ndarray | list, float]

Calculate the SNR and scale a signal.

Parameters:
  • target_snr – Desired SNR value for the injected signal.

  • sig – Signal to adjust. A copy will be made. Can be time-domain or frequency-domain. Must be 1D np.ndarray, list of 1D np.ndarray``s, or 2D ``np.ndarray across channels with shape (nchannels, data length).

  • *snr_args – Arguments to pass to snr().

  • return_orig_snr – If True, return the original SNR in addition to the adjusted data.

  • **snr_kwargs – Keyword arguments to pass to snr().

Returns:

Returns the copied input signal adjusted to the target SNR. If return_orig_snr is True, the original SNR is added as the second entry of a tuple with the adjusted signal (as the first entry in the tuple).