Models

Models govern what functional form to use for the PSF at a single location.

The Model base class

class piff.Model[source]

The base class for modeling a single PSF (i.e. no interpolation yet)

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

draw(star, copy_image=True, center=None)[source]

Draw the model on the given image.

Parameters
  • star – A Star instance with the fitted parameters to use for drawing and a data field that acts as a template image for the drawn model.

  • copy_image – If False, will use the same image object. If True, will copy the image and then overwrite it. [default: True]

  • center – An optional tuple (x,y) location for where to center the drawn profile in the image. [default: None, which draws at the star’s location.]

Returns

a new Star instance with the data field having an image of the drawn model.

fit(star, convert_func=None)[source]

Fit the Model to the star’s data to yield iterative improvement on its PSF parameters, their uncertainties, and flux (and center, if free). The returned star.fit.alpha will be inverse covariance of solution if it is estimated, else is None.

Parameters
  • star – A Star instance

  • convert_func – An optional function to apply to the profile being fit before drawing it onto the image. This is used by composite PSFs to isolate the effect of just this model component. [default: None]

Returns

New Star instance with updated fit information

initialize(star, logger=None)[source]

Initialize a star to work with the current model.

Parameters
  • star – A Star instance with the raw data.

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

Returns

Star instance with the appropriate initial fit values

normalize(star)[source]

Make sure star.fit.params are normalized properly.

Note: This modifies the input star in place.

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

Parse the model 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_model – The model field of the configuration dict, config[‘model’]

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

Returns

a kwargs dict to pass to the initializer

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

Parse the model field of the config dict.

Parameters
  • config_model – The configuration dict for the model field.

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

Returns

a Model instance

classmethod read(fits, extname)[source]

Read a Model from a FITS file.

Note: the returned Model will not have its parameters set. This just initializes a fresh model that can be used to interpret interpolated vectors.

Parameters
  • fits – An open fitsio.FITS object

  • extname – The name of the extension with the model information.

Returns

a model built with a information in the FITS file.

reflux(star, fit_center=True, logger=None)[source]

Fit the Model to the star’s data, varying only the flux (and center, if it is free). Flux and center are updated in the Star’s attributes. This is a single-step solution if only solving for flux, otherwise an iterative operation. DOF in the result assume only flux (& center) are free parameters.

Parameters
  • star – A Star instance

  • fit_center – If False, disable any motion of center

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

Returns

New Star instance, with updated flux, center, chisq, dof

write(fits, extname)[source]

Write a Model to a FITS file.

Note: this only writes the initialization kwargs to the fits extension, not the parameters.

The base class implemenation works if the class has a self.kwargs attribute and these are all simple values (str, float, or int)

Parameters
  • fits – An open fitsio.FITS object

  • extname – The name of the extension to write the model information.

Models based on GalSim objects

class piff.GSObjectModel(gsobj, fastfit=False, centered=True, include_pixel=True, scipy_kwargs=None, logger=None)[source]

Model that takes a fiducial GalSim.GSObject and dilates, shifts, and shears it to get a good match to stars.

Parameters
  • gsobj – GSObject to use as fiducial profile.

  • fastfit – Use HSM moments for fitting. Approximate, but fast. [default: False]

  • centered – If True, PSF model centroid is forced to be (0,0), and the PSF fitting will marginalize over stellar position. If False, stellar position is fixed at input value and the fitted PSF may be off-center. [default: True]

  • include_pixel – Include integration over pixel when drawing? [default: True]

  • scipy_kwargs – Optional kwargs to pass to scipy.optimize.least_squares [default: None]

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

fit(star, fastfit=None, logger=None, convert_func=None)[source]

Fit the image either using HSM or least-squares minimization.

If fastfit is True, then the galsim.hsm module will be used to estimate the transformation parameters that take the fiducial moments into the data moments. If fastfit is False, then the Levenberg-Marquardt minimization algorithm will be used instead. The latter should generally be more accurate, but slower due to the need to iteratively propose model improvements.

Parameters
  • star – A Star to fit.

  • fastfit – Use fast HSM moments to fit? [default: None, which means use fitting mode specified in the constructor.]

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

  • convert_func – An optional function to apply to the profile being fit before drawing it onto the image. This is used by composite PSFs to isolate the effect of just this model component. [default: None]

Returns

a new Star with the fitted parameters in star.fit

getProfile(params)[source]

Get a version of the model as a GalSim GSObject

Parameters

params – A numpy array with either [ size, g1, g2 ] or [ cenu, cenv, size, g1, g2 ] depending on if the center of the model is being forced to (0.0, 0.0) or not.

Returns

a galsim.GSObject instance

initialize(star, logger=None)[source]

Initialize the given star’s fit parameters.

Parameters
  • star – The Star to initialize.

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

Returns

a new initialized Star.

least_squares_fit(star, logger=None, convert_func=None)[source]

Fit parameters of the given star using least-squares minimization.

Parameters
  • star – A Star to fit.

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

  • convert_func – An optional function to apply to the profile being fit before drawing it onto the image. This is used by composite PSFs to isolate the effect of just this model component. [default: None]

Returns

(flux, dx, dy, scale, g1, g2, flag)

moment_fit(star, logger=None)[source]

Estimate transformations needed to bring self.gsobj towards given star.

class piff.Gaussian(fastfit=False, centered=True, include_pixel=True, scipy_kwargs=None, logger=None)[source]

Model PSFs as elliptical Gaussians.

Parameters
  • fastfit – Use HSM moments for fitting. Approximate, but fast. [default: False]

  • centered – If True, PSF model centroid is forced to be (0,0), and the PSF fitting will marginalize over stellar position. If False, stellar position is fixed at input value and the fitted PSF may be off-center. [default: True]

  • include_pixel – Include integration over pixel when drawing? [default: True]

  • scipy_kwargs – Optional kwargs to pass to scipy.optimize.least_squares [default: None]

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

class piff.Kolmogorov(fastfit=False, centered=True, include_pixel=True, scipy_kwargs=None, logger=None)[source]

Model PSFs as elliptical Kolmogorovs.

Parameters
  • fastfit – Use HSM moments for fitting. Approximate, but fast. [default: False]

  • centered – If True, PSF model centroid is forced to be (0,0), and the PSF fitting will marginalize over stellar position. If False, stellar position is fixed at input value and the fitted PSF may be off-center. [default: True]

  • include_pixel – Include integration over pixel when drawing? [default: True]

  • scipy_kwargs – Optional kwargs to pass to scipy.optimize.least_squares [default: None]

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

class piff.Moffat(beta, trunc=0.0, fastfit=False, centered=True, include_pixel=True, scipy_kwargs=None, logger=None)[source]

Model PSFs as elliptical Moffats.

Parameters
  • beta – Moffat shape parameter.

  • trunc – Optional truncation radius at which profile drops to zero. Measured in half light radii. [default: 0, indicating no truncation]

  • fastfit – Use HSM moments for fitting. Approximate, but fast. [default: False]

  • centered – If True, PSF model centroid is forced to be (0,0), and the PSF fitting will marginalize over stellar position. If False, stellar position is fixed at input value and the fitted PSF may be off-center. [default: True]

  • include_pixel – Include integration over pixel when drawing? [default: True]

  • scipy_kwargs – Optional kwargs to pass to scipy.optimize.least_squares [default: None]

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

Pixel grid model

class piff.PixelGrid(scale, size, interp=None, centered=True, logger=None)[source]

A PSF modeled as interpolation between a grid of points.

The parameters of the model are the values at the grid points, although the sum of the values is constrained to be 1/scale**2, to give it unit total flux. The grid is in uv space, with the scale and size specified on construction. Interpolation will always assume values of zero outside of grid.

PixelGrid also needs to specify an interpolant to define how to values between grid points are determined from the pixelated model. Any galsim.Interpolant type is allowed. The default interpolant is galsim.Lanczos(7)

Parameters
  • scale – Pixel scale of the PSF model (in arcsec)

  • size – Number of pixels on each side of square grid.

  • interp – An Interpolant to be used [default: Lanczos(7)]

  • centered – If True, PSF model centroid is forced to be (0,0), and the PSF fitting will marginalize over stellar position. If False, stellar position is fixed at input value and the fitted PSF may be off-center. [default: True]

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

chisq(star, logger=None, convert_func=None)[source]

Calculate dependence of chi^2 = -2 log L(D|p) on PSF parameters for single star. as a quadratic form chi^2 = dp^T AT A dp - 2 bT A dp + chisq, where dp is the shift from current parameter values. Returned Star instance has the resultant (A, b, chisq) attributes, but params vector has not have been updated yet (could be degenerate).

Parameters
  • star – A Star instance

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

  • convert_func – An optional function to apply to the profile being fit before drawing it onto the image. This is used by composite PSFs to isolate the effect of just this model component. [default: None]

Returns

a new Star instance with updated fit parameters. (esp. A,b)

fit(star, logger=None, convert_func=None)[source]

Fit the Model to the star’s data to yield iterative improvement on its PSF parameters and uncertainties.

Parameters
  • star – A Star instance

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

  • convert_func – An optional function to apply to the profile being fit before drawing it onto the image. This is used by composite PSFs to isolate the effect of just this model component. [default: None]

Returns

a new Star instance with updated fit information

getProfile(params)[source]

Get a version of the model as a GalSim GSObject

Parameters

params – The fit parameters for a given star.

Returns

a galsim.GSObject instance

initialize(star, logger=None)[source]

Initialize a star to work with the current model.

Parameters
  • star – A Star instance with the raw data.

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

Returns

a star instance with the appropriate initial fit values

normalize(star)[source]

Make sure star.fit.params are normalized properly.

Note: This modifies the input star in place.

Optical model

class piff.Optical(template=None, logger=None, **kwargs)[source]

Initialize the Optical Model

There are potentially three components to this model that are convolved together.

First, there is an optical component, which uses a galsim.OpticalPSF to model the profile. The aberrations are considered fitted parameters, but the other attributes are fixed and are given at initialization. These parameters are passed to GalSim, so they have the same definitions as used there.

Parameters
  • diam – Diameter of telescope aperture in meters. [required (but cf. template option)]

  • lam – Wavelength of observations in nanometers. [required (but cf. template option)]

  • obscuration – Linear dimension of central obscuration as fraction of pupil linear dimension, [0., 1.). [default: 0]

  • nstruts – Number of radial support struts to add to the central obscuration. [default: 0]

  • strut_thick – Thickness of support struts as a fraction of pupil diameter. [default: 0.05]

  • strut_angle – Angle made between the vertical and the strut starting closest to it, defined to be positive in the counter-clockwise direction. [default: 0. * galsim.degrees]

  • pupil_plane_im – The name of a file containing the pupil plane image to use instead of creating one from obscuration, struts, etc. [default: None]

Second, there may be an atmospheric component, which uses either a galsim.Kolmogorov or galsim.VonKarman to model the profile.

Parameters
  • fwhm – The full-width half-max of the atmospheric part of the PSF. [default: None]

  • r0 – The Fried parameter in units of meters to use for the Kolmogorov profile. [default: None]

  • L0 – The outer scale in units of meters if desired, in which case the atmospheric part will be a VonKarman. [default: None]

Finally, there is allowed to be a final Gaussian component and an applied shear.

Parameters
  • sigma – Convolve with gaussian of size sigma. [default: 0]

  • g2 (g1,) – Shear to apply to final image. Simulates vibrational modes. [default: 0]

Since there are a lot of parameters here, we provide the option of setting many of them from a template value. e.g. template = ‘des’ will use the values stored in the dict piff.optical_model.optical_templates[‘des’].

Parameters

template – A key word in the dict piff.optical_model.optical_template to use for setting values of these parameters. [default: None]

If you use a template as well as other specific parameters, the specific parameters will override the values from the template. e.g. to simulate what DES would be like at lambda=1000 nm (the default is 700), you could do:

>>> model = piff.OpticalModel(template='des', lam=1000)
fit(star)[source]

Warning: This method just updates the fit with the chisq and dof!

Parameters

star – A Star instance

Returns

a new Star with the fitted parameters in star.fit

getProfile(params)[source]

Get a version of the model as a GalSim GSObject

Parameters

params – A np array with [z4, z5, z6…z11]

Returns

a galsim.GSObject instance