irf.psf

ximpol.irf.psf.gauss_king(r, W, sigma, N, r_c, eta)[source]

Functional representation of the Gaussian plus King PSF profile described in Fabiani et al., 2014, equation (2):

\text{PSF}(r) = W \exp^{-(\frac{r^2}{2\sigma^2})} +
N\left( 1 + \left( \frac{r}{r_c} \right)^2 \right)^{-\eta}

Parameters:
  • r (float or array) – The radial distance from the true source position is arcsec.
  • W (float) – Normalization of the Gaussian component.
  • sigma (float) – Width of the Gaussian component.
  • N (float) – Normalization of the King component.
  • r_c (float) – Characteristic radius of the King component.
  • eta (float) – Exponent of the King component.
ximpol.irf.psf.gauss_king_eef_at_infinity(W, sigma, N, r_c, eta)[source]

Return the value of the Encircled Energy Fraction (EEF) at infinity, given the parameters of the functional representation, see equation (4) of Fabiani et al., 2014.

\text{EEF}(\infty) = 2\pi W\sigma^2 +
\pi\frac{r_c^2 N}{\eta - 1}

Parameters:
  • r (float or array) – The radial distance from the true source position is arcsec.
  • W (float) – Normalization of the Gaussian component.
  • sigma (float) – Width of the Gaussian component.
  • N (float) – Normalization of the King component.
  • r_c (float) – Characteristic radius of the King component.
  • eta (float) – Exponent of the King component.
class ximpol.irf.psf.xBinTableHDUPSF(data=None, keywords=[], comments=[])[source]

Binary table for the PSF extension of a psf file.

class ximpol.irf.psf.xPointSpreadFunction(psf_file_path)[source]

Class describing a (simplified, energy independent) PSF.

The effective area is essentially a linear spline, with built-in facilities for evaluation and plotting.

Parameters:
  • psf_file_path (str) – The path to the .psf FITS file containing the effective area table.
  • rmax (float) – The maximum radial distance (in arcsec) for the PSF.

Example

>>> from ximpol import XIMPOL_IRF
>>> from ximpol.utils.matplotlib_ import pyplot as plt
>>>
>>> file_path = os.path.join(XIMPOL_IRF,'fits','xipe_baseline.psf')
>>> psf = xPointSpreadFunction(file_path)
>>> print(psf.rvs(10))
>>> ra, dec = 5.0, 12.3
>>> print(psf.smear_single(ra, dec, 10))

Note

The parametrization is taken from Fabiani et al., 2014, table 2. The PSF is technically energy dependent, but the dependence is not wild and for the moment we stick with the values at 4.51 keV in the table.

build_eef()[source]

Build the Encircled Energy Fraction (EEF) as a function of r.

And, while we’re at it, we also calculate and cache the HEW.

delta(size=1)[source]

Return an array of random offset (in ra, dec or L, B) due to the PSF.

Note the output is converted in degrees.

draw_psf_circle(image, x, y, text='PSF', color='white', lw=2, number=False)[source]

Add the PSF circle to the image with labels. This function must be called after the (possible) image recenter.

Note the x and y are coordinates relative to the figure axes (0.0 is left or bottom and 1.0 is right or top).

rvs(size)[source]

Extract values of the radial distance according to the PSF shape.

smear(ra, dec)[source]

Smear a pair of arrays of coordinates.

smear_single(ra, dec, num_times=1)[source]

Smear a pair of coordinates for an arbitrary number of times.

view(show=True)[source]

Overloaded plot method (with default log scale on the y-axis).