Output statistics

Piff can also produce some statistics about the PSF model and residuals for the stars used to build the model. These are specified in a stats section of the config file.

class piff.Stats[source]

The base class for getting the statistics of a set of stars.

This is essentially an abstract base class intended to define the methods that should be implemented by any derived class.

The usual code flow for using a Stats instance is:

>>> stats = SomeKindofStats(...)
>>> stats.compute(psf, stars, logger)
>>> stats.write(file_name=file_name)

There is also a plot method if you want to make the matplot lib fig, ax and do something else with it besides just write it to a file.

compute(psf, stars, logger=None)[source]

Compute the given statistic for a PSF solution on a set of stars.

This needs to be done before the statistic is plotted or written to a file.

Parameters:
  • psf – A PSF Object

  • stars – A list of Star instances.

  • logger – A logger object for logging debug info. [default: None]

measureShapes(psf, stars, model_properties=None, fourth_order=False, raw_moments=False, logger=None)[source]

Compare PSF and true star shapes with HSM algorithm

Parameters:
  • psf – A PSF Object

  • stars – A list of Star instances.

  • model_properties – Optionally a dict of properties to use for the model rendering. The default behavior is to use the properties of the star itself for any properties that are not overridden by model_properties. [default: None]

  • fourth_order – Whether to include the fourth-order quantities as well in columns 7 through 11 of the output array. [default: False]

  • raw_moments – Whether to include the complete set of raw moments as calculated by calculate_moments. If requested, these come after the fourth-order quantities if those are being returned or after the regular quantities otherwise. Either way they are the last 18 columns in the output array. [default: False]

  • logger – A logger object for logging debug info. [default: None]

Returns:

positions of stars, shapes of stars, and shapes of models of stars (T, g1, g2)

classmethod parseKwargs(config_stats, logger=None)[source]

Parse the stats field of a configuration dict and return the kwargs to use for initializing an instance of the class.

The base class implementation just returns the kwargs as they are, but derived classes might want to override this if they need to do something more sophisticated with them.

Parameters:
  • config_stats – The stats field of the configuration dict, config[‘stats’]

  • logger – A logger object for logging debug info. [default: None]

Returns:

a kwargs dict to pass to the initializer

plot(logger=None, **kwargs)[source]

Make the plots for this statistic.

Parameters:
  • logger – A logger object for logging debug info. [default: None]

  • **kwargs – Optionally, provide extra kwargs for the matplotlib plot command.

Returns:

(fig, ax) The matplotlib figure and axis with the plot(s).

classmethod process(config_stats, logger=None)[source]

Parse the stats field of the config dict.

Parameters:
  • config_stats – The configuration dict for the stats field.

  • logger – A logger object for logging debug info. [default: None]

Returns:

a Stats instance

write(file_name=None, logger=None, **kwargs)[source]

Write plots to a file.

Parameters:
  • file_name – The name of the file to write to. [default: Use self.file_name, which is typically read from the config field.]

  • logger – A logger object for logging debug info. [default: None]

  • **kwargs – Optionally, provide extra kwargs for the matplotlib plot command.

class piff.ShapeHistStats(file_name=None, nbins=None, cut_frac=0.01, model_properties=None, logger=None)[source]

Stats class for calculating histograms of shape residuals

This will compute the size and shapes of the observed stars and the PSF models and make histograms of both the values and the residuals.

The plot will have 6 axes. The top row will have histograms of T, g1, g2, with the model and data color coded. The bottom row will have histograms of the differences.

After a call to compute(), the following attributes are accessible:

u:

The u positions in field coordinates.

v:

The v positions in field coordinates.

T:

The size (T = Iuu + Ivv) of the observed stars.

g1:

The g1 component of the shapes of the observed stars.

g2:

The g2 component of the shapes of the observed stars.

T_model:

The size of the PSF model at the same locations as the stars.

g1_model:

The g1 component of the PSF model at these locations.

g2_model:

The g2 component of the PSF model at these locations.

dT:

The size residual, T - T_model

dg1:

The g1 residual, g1 - g1_model

dg2:

The g2 residual, g2 - g2_model

Parameters:
  • file_name – Name of the file to output to. [default: None]

  • nbins – Number of bins to use. [default: sqrt(n_stars)]

  • cut_frac – Fraction to cut off from histograms at the high and low ends. [default: 0.01]

  • model_properties – Optionally a dict of properties to use for the model rendering. [default: None]

compute(psf, stars, logger=None)[source]
Parameters:
  • psf – A PSF Object

  • stars – A list of Star instances.

  • logger – A logger object for logging debug info. [default: None]

plot(logger=None, **kwargs)[source]

Make the plots.

Parameters:

logger – A logger object for logging debug info. [default: None]

Params **kwargs:

Any additional kwargs go into the matplotlib hist() function.

Returns:

fig, ax

class piff.RhoStats(min_sep=0.5, max_sep=300, bin_size=0.1, file_name=None, model_properties=None, logger=None, **kwargs)[source]

Stats class for calculating rho statistics.

This will plot the 5 rho statistics described in Jarvis et al, 2015, section 3.4.

e = e_psf; de = e_psf - e_model T is size; dT = T_psf - T_model

rho1 = < de* de > rho2 = < e* de > (in the rowe paper this is < e* de + de* e > rho3 = < (e* dT / T) (e dT / T) > rho4 = < de* (e dT / T) > rho5 = < e* (e dT / T) >

The plots for rho1, rho3, and rho4 will all be on the same axis (left), and the plots for rho2 and rho5 will be on the other axis (right).

Furthermore, these are technically complex quantities, but only the real parts are plotted, since the imaginary parts are uninteresting.

After a call to compute(), the following attributes are accessible:

rho1:

A TreeCorr GGCorrelation instance with the rho1 statistic.

rho2:

A TreeCorr GGCorrelation instance with the rho2 statistic.

rho3:

A TreeCorr GGCorrelation instance with the rho3 statistic.

rho4:

A TreeCorr GGCorrelation instance with the rho4 statistic.

rho5:

A TreeCorr GGCorrelation instance with the rho5 statistic.

The value of the canonical rho statistic is in the xip attribute of each of the above TreeCorr GGCorrelation instances. But there are other quantities that may be of interest in some cases, so we provide access to the full object.

Parameters:
  • min_sep – Minimum separation (in arcmin) for pairs. [default: 0.5]

  • max_sep – Maximum separation (in arcmin) for pairs. [default: 300]

  • bin_size – Size of bins in log(sep). [default 0.1]

  • file_name – Name of the file to output to. [default: None]

  • model_properties – Optionally a dict of properties to use for the model rendering. [default: None]

  • logger – A logger object for logging debug info. [default: None]

  • **kwargs – Any additional kwargs are passed on to TreeCorr.

compute(psf, stars, logger=None)[source]
Parameters:
  • psf – A PSF Object

  • stars – A list of Star instances.

  • logger – A logger object for logging debug info. [default: None]

plot(logger=None, **kwargs)[source]

Make the plots.

Parameters:

logger – A logger object for logging debug info. [default: None]

Params **kwargs:

Any additional kwargs go into the matplotlib plot() function. [ignored in this function]

Returns:

fig, ax

class piff.TwoDHistStats(nbins_u=20, nbins_v=20, reducing_function='np.median', file_name=None, model_properties=None, logger=None)[source]

Statistics class that can make pretty colormaps where each bin has some arbitrary function applied to it.

By default this will make a color map based on u and v coordinates of the input stars. The color scale is based on (by default) the median value of the objects in particular u-v voxel.

After a call to compute(), the following attributes are accessible:

twodhists:

A dictionary of two dimensional histograms, with keys ‘u’, ‘v’, ‘T’, ‘g1’, ‘g2’, ‘T_model’, ‘g1_model’, ‘g2_model’, ‘dT’, ‘dg1’, ‘dg2’

These histograms are two dimensional masked arrays where the value of the pixel corresponds to reducing_function([objects in u-v voxel])

Parameters:
  • nbins_u – Number of bins in u direction [default: 20]

  • nbins_v – Number of bins in v direction [default: 20]

  • reducing_function – Type of function to apply to grouped objects. numpy functions are prefixed by np. [default: ‘np.median’]

  • file_name – Name of the file to output to. [default: None]

  • model_properties – Optionally a dict of properties to use for the model rendering. [default: None]

  • logger – A logger object for logging debug info. [default: None]

classmethod _shiftedColorMap(cmap, start=0, midpoint=0.5, stop=1.0, name='shiftedcmap')[source]

Function to offset the “center” of a colormap. Useful for data with a negative min and positive max and you want the middle of the colormap’s dynamic range to be at zero

Taken from

https://github.com/olgabot/prettyplotlib/blob/master/prettyplotlib/colors.py

which makes beautiful plots by the way

Parameters:
  • cmap – The matplotlib colormap to be altered

  • start – Offset from lowest point in the colormap’s range. Defaults to 0.0 (no lower ofset). Should be between 0.0 and midpoint.

  • midpoint – The new center of the colormap. Defaults to 0.5 (no shift). Should be between 0.0 and 1.0. In general, this should be 1 - vmax/(vmax + abs(vmin)) For example if your data range from -15.0 to +5.0 and you want the center of the colormap at 0.0, midpoint should be set to 1 - 5/(5 + 15)) or 0.75

  • stop – Offset from highets point in the colormap’s range. Defaults to 1.0 (no upper ofset). Should be between midpoint and 1.0.

compute(psf, stars, logger=None)[source]
Parameters:
  • psf – A PSF Object

  • stars – A list of Star instances.

  • logger – A logger object for logging debug info. [default: None]

plot(logger=None, **kwargs)[source]

Make the plots.

Parameters:

logger – A logger object for logging debug info. [default: None]

Params **kwargs:

Any additional kwargs go into the matplotlib plot() function. [ignored in this function]

Returns:

fig, ax

class piff.WhiskerStats(file_name=None, nbins_u=20, nbins_v=20, reducing_function='np.median', scale=1, resid_scale=2, model_properties=None, logger=None)[source]

Statistics class that can make whiskerplots.

By default this will make a whisker plot based on u and v coordinates of the input stars. The whisker scale is based on (by default) the median value of the objects in a particular u-v voxel.

After a call to compute(), the following attributes are accessible:

twodhists:

A dictionary of two dimensional histograms, with keys ‘u’, ‘v’, ‘w1’, ‘w2’, ‘w1_model’, ‘w2_model’, ‘dw1’, ‘dw2’

These histograms are two dimensional masked arrays where the value of the pixel corresponds to reducing_function([objects in u-v voxel])

Note

There are a couple different ways to define your whiskers. Here we have taken the approach that the whisker represents the ellipticity as:

theta = arctan(e2, e1) / 2 r = sqrt(e1 ** 2 + e2 ** 2) w1 = r cos(theta) w2 = r sin(theta)

Because e1, e2 do not have units, w does not either.

Parameters:
  • file_name – Name of the file to output to. [default: None]

  • nbins_u – Number of bins in u direction [default: 20]

  • nbins_v – Number of bins in v direction [default: 20]

  • reducing_function – Type of function to apply to grouped objects. numpy functions are prefixed by np. [default: ‘np.median’]

  • scale – An overal scale factor by which to scale the size of the all whiskers. [default: 1]

  • resid_scale – An additional factor for the scale size of the residual whiskers only. [default: 2]

  • model_properties – Optionally a dict of properties to use for the model rendering. [default: None]

  • logger – A logger object for logging debug info. [default: None]

compute(psf, stars, logger=None)[source]
Parameters:
  • psf – A PSF Object

  • stars – A list of Star instances.

  • logger – A logger object for logging debug info. [default: None]

plot(logger=None, **kwargs)[source]

Make the plots.

Parameters:

logger – A logger object for logging debug info. [default: None]

Params **kwargs:

Any additional kwargs go into the matplotlib plot() function. [ignored in this function]

Returns:

fig, ax

class piff.StarStats(nplot=10, adjust_stars=False, include_reserve=True, only_reserve=False, include_flagged=False, include_ave=True, file_name=None, logger=None)[source]

This Stats class can take stars and make a set of plots of them as well as their models and residuals.

By default this will draw 5 random stars, make psf stars, and plot the residual of the two.

After a call to compute(), the following attributes are accessible:

stars:

List of stars used for plotting

models:

List of models of stars used for plotting

indices:

Indices of input stars that the plotting stars correspond to

Parameters:
  • nplot – Number of stars we wish to plot. If 0 or nplot > nstars in PSF, then we plot all stars. Otherwise, we draw nplot stars at random (without replacement). [default: 10]

  • adjust_stars – Boolean. If true, when computing, will also fit for best starfit center and flux to match observed star. [default: False]

  • include_reserve – Whether to inlude reserve stars. [default: True]

  • only_reserve – Whether to skip plotting non-reserve stars. [default: False]

  • include_flagged – Whether to include plotting flagged stars. [default: False]

  • include_ave – Whether to inlude the average image. [default: True]

  • file_name – Name of the file to output to. [default: None]

  • logger – A logger object for logging debug info. [default: None]

compute(psf, stars, logger=None)[source]
Parameters:
  • psf – A PSF Object

  • stars – A list of Star instances.

  • logger – A logger object for logging debug info. [default: None]

plot(logger=None, **kwargs)[source]

Make the plots.

Parameters:

logger – A logger object for logging debug info. [default: None]

Params **kwargs:

Any additional kwargs go into the matplotlib pcolor() function.

Returns:

fig, ax

class piff.HSMCatalogStats(file_name=None, model_properties=None, fourth_order=False, raw_moments=False, logger=None)[source]

Stats class for writing the shape information to an output file.

This will compute the size and shapes of the observed stars and the PSF models and write these data to a file.

The HSM adaptive momemnt measurements sometimes fail in various ways. When it does, we still output the information that we have for a star, but mark the failure with a flag: flag_data for errors in the data measurement, flag_model for errors in the model measurement. The meaning of these flags are (treated as a bit mask):

Flags:

  • 0 = Success

  • 1 = HSM returned a non-zero moments_status.

  • 2 = HSM returned a negative flux.

  • 4 = HSM’s centroid moved by more than 1 pixel from the input position.

The output file will include the following columns:

ra:

The RA of the star in degrees. (or 0 if the wcs is not a CelestialWCS)

dec:

The Declination of the star in degrees. (or 0 if the wcs is not a CelestialWCS)

u:

The u position in field coordinates.

v:

The v position in field coordinates.

x:

The x position in chip coordinates.

y:

The y position in chip coordinates.

T_data:

The size (T = Iuu + Ivv) of the observed star.

g1_data:

The g1 component of the shapes of the observed star.

g2_data:

The g2 component of the shapes of the observed star.

T_model:

The size of the PSF model at the same locations as the star.

g1_model:

The g1 component of the PSF model.

g2_model:

The g2 component of the PSF model.

reserve:

Whether the star was a reserve star.

flag_data:

0 where HSM succeeded on the observed star, >0 where it failed (see above).

flag_model:

0 where HSM succeeded on the PSF model, >0 where it failed (see above).

If fourth_order=True, then there are additional quantities calculated as well. We define the following notation:

\[\begin{split}T &= \int I(u,v) (u^2 + v^2) du dv \\ g &= \frac{\int I(u,v) (u + iv)^2 du dv}{T} \\ T^{(4)} &= \frac{\int I(u,v) (u^2 + v^2)^2 du dv}{T} \\ g^{(4)} &= \frac{\int I(u,v) (u^2 + v^2) (u + iv)^2 du dv}{T^2} - 3g \\ h^{(4)} &= \frac{\int I(u,v) (u + iv)^4 du dv}{T^2}\end{split}\]

I.e. \(T^{(4)}\) is a fourth-order spin-0 quantity, analogous to \(T\) at second order, \(g^{(4)}\) is a fourth-order spin-2 quantity, analogous to \(g\), and \(h^{(4)}\) is a spin-4 quantity. The denominators ensure that the units of \(T^{(4)}\) is \(\mathrm{arcsec}^2\), just like \(T\) and that \(g^{(4)}\) and \(h^{(4)}\) are dimensionless. And the \(-3g\) term for \(g^{(4)}\) subtracts off the dominant contribution to the fourth order quantity from the second order shape. For a pure elliptical Gaussian, this makes \(g^{(4)}\) come out very close to zero.

The output file contains the following additional columns:

T4_data:

The fourth-order “size”, \(T^{(4)}\), of the observed star.

g41_data:

The real component of \(g^{(4)}\) of the obserbed star.

g42_data:

The imaginary component of \(g^{(4)}\) of the obserbed star.

h41_data:

The real component of \(h^{(4)}\) of the obserbed star.

h42_data:

The imaginary component of \(h^{(4)}\) of the obserbed star.

T4_model:

The fourth-order “size”, \(T^{(4)}\), of the PSF model.

g41_model:

The real component of \(g^{(4)}\) of the PSF model.

g42_model:

The imaginary component of \(g^{(4)}\) of the PSF model.

h41_model:

The real component of \(h^{(4)}\) of the PSF model.

h42_model:

The imaginary component of \(h^{(4)}\) of the PSF model.

Parameters:
  • file_name – Name of the file to output to. [default: None]

  • model_properties – Optionally a dict of properties to use for the model rendering. The default behavior is to use the properties of the star itself for any properties that are not overridden by model_properties. [default: None]

  • fourth_order – Whether to include the additional fourth-order quantities described above. [default: False]

  • raw_moments – Whether to include the complete set of raw moments as calculated by piff.util.calculate_moments. [default: False]

compute(psf, stars, logger=None)[source]
Parameters:
  • psf – A PSF Object

  • stars – A list of Star instances.

  • logger – A logger object for logging debug info. [default: None]

write(file_name=None, logger=None)[source]

Write catalog to file.

Parameters:
  • file_name – The name of the file to write to. [default: Use self.file_name, which is typically read from the config field.]

  • logger – A logger object for logging debug info. [default: None]