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.

_finish_write(writer)[source]

Finish the writing process with any class-specific steps.

The base class implementation doesn’t do anything, which is often appropriate, but this hook exists in case any Model classes need to write extra information to the fits file.

Parameters:

writer – A writer object that encapsulates the serialization format.

_finish_read(reader)[source]

Finish the reading process with any class-specific steps.

The base class implementation doesn’t do anything, which is often appropriate, but this hook exists in case any Model classes need to read extra information from the fits file.

Parameters:

reader – A reader object that encapsulates the serialization format.

classmethod _fix_kwargs(kwargs)[source]

Fix the kwargs read in from an input file.

This is intended to make it easier to preserve backwards compatibility if a class has changed something about the kwargs, this provides a way for old parameter names or defaults to be updated for a newer version of Piff than the one that wrong them.

Usually, this is a no op.

Parameters:

kwargs – The old kwargs read in from a previous version Piff output file.

draw(star, copy_image=True)[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]

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, default_init=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]

  • default_init – The default initilization method if the user doesn’t specify one. [default: None]

Returns:

Star instance with the appropriate initial fit values

initialize_iteration()[source]

Do any required initialization at the start of an iteration of fitting.

Usually a no op.

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(reader, name)[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:
  • reader – A reader object that encapsulates the serialization format.

  • name – Name associated with this model in the serialized output.

Returns:

a model built with a information in the FITS file

set_num(num)[source]

If there are multiple components involved in the fit, set the number to use for this model.

write(writer, name)[source]

Write a Model via a Writer object.

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:
  • writer – A writer object that encapsulates the serialization format.

  • name – A name to associate with this model in the serialized output.

Models based on GalSim objects

class piff.GSObjectModel(gsobj, fastfit=False, centered=True, include_pixel=True, init=None, fit_flux=False, 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.

The following initialization methods are available for the init parameter.

  • hsm Start with flux and size values that match the hsm moments of the star.

  • zero Start with flux = 1.e-6 x the hsm flux.

  • delta Start with size = 1.e-6 x the hsm size.

  • (flux,size) (Note: flux and size here are numeric values.) Start with flux and size

    scaled by the given values relative to hsm.

All initialization methods start with zero shear and zero centroid offset.

Use type name “GSObject” in a config field to use this model.

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]

  • init – Initialization method. [default: None, which uses hsm unless a PSF class specifies a different default.]

  • fit_flux – If True, the PSF model will include the flux value. This is useful when this model is an element of a Sum composite PSF. [default: False, unless init==’zero’, in which case it is automatically True.]

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

  • logger

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

    _resid(params, star, convert_func, draw_method)[source]

    Residual function to use with least_squares.

    Essentially “chi” from “chisq”, but not summed over pixels yet.

    Parameters:
    • params – A numpy array of model 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.

    • draw_method – The method to use with the GalSim drawImage command to determine the residuals. [if draw_method is None, use self._method]

    Returns:

    ”chi” as a flattened numpy array.

    _get_params(star)[source]

    Generate an array of model parameters.

    Parameters:

    star – A Star from which to initialize parameter values.

    Returns:

    a numpy array

fit(star, fastfit=None, logger=None, convert_func=None, draw_method=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; if provided, fastfit is ignored and taken to be False.]

  • draw_method – The method to use with the GalSim drawImage command to determine the residuals. [default: None, which uses ‘auto’ if include_pixel = True, else ‘no_pixel’]

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, default_init=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]

  • default_init – The default initilization method if the user doesn’t specify one. [default: None]

Returns:

a new initialized Star.

least_squares_fit(star, logger=None, convert_func=None, draw_method=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]

  • draw_method – The method to use with galsim.DrawImage to determine the residuals. [default: ‘auto’ if include_pixel = True, else ‘no_pixel’]

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(**kwargs)[source]

Model PSFs as elliptical Gaussians.

Use type name “Gaussian” in a config field to use this model.

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]

  • init – Initialization method. [default: None, which means to start with size=1, flux=1, and not shift or shear]

  • fit_flux – If True, the PSF model will include the flux value. This is useful when this model is an element of a Sum composite PSF. [default: False]

  • 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(**kwargs)[source]

Model PSFs as elliptical Kolmogorovs.

Use type name “Kolmogorov” in a config field to use this model.

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]

  • init – Initialization method. [default: None, which means to start with size=1, flux=1, and not shift or shear]

  • fit_flux – If True, the PSF model will include the flux value. This is useful when this model is an element of a Sum composite PSF. [default: False]

  • 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, **kwargs)[source]

Model PSFs as elliptical Moffats.

Use type name “Moffat” in a config field to use this model.

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]

  • init – Initialization method. [default: None, which means to start with size=1, flux=1, and not shift or shear]

  • fit_flux – If True, the PSF model will include the flux value. This is useful when this model is an element of a Sum composite PSF. [default: False]

  • 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, init=None, fit_flux=False, 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)

The following initialization methods are available for the init parameter.

  • hsm Start with flux and size values that match the hsm moments of the star.

  • zero Start with flux = 1.e-6 x the hsm flux.

  • delta Start with size = 1.e-6 x the hsm size.

All initialization methods start with zero shear and zero centroid offset.

Use type name “PixelGrid” in a config field to use this model.

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]

  • init – Initialization method. [default: None, which uses hsm unless a PSF class specifies a different default.]

  • fit_flux – If True, the PSF model will include the flux value. This is useful when this model is an element of a Sum composite PSF. [default: False, unless init==’zero’, in which case it is automatically True.]

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

chisq(star, logger=None, convert_func=None, draw_method=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]

  • draw_method – The method to use with the GalSim drawImage command to determine the residuals. [PixelGrid always uses ‘no_pixel’; this parameter is only present for API compatibility. It must be either None or ‘no_pixel’.]

Returns:

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

fit(star, logger=None, convert_func=None, draw_method=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]

  • draw_method – The method to use with the GalSim drawImage command to determine the residuals. [PixelGrid always uses ‘no_pixel’; this parameter is only present for API compatibility. It must be either None or ‘no_pixel’.]

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, default_init=None)[source]

Initialize a star to work with the current model.

Parameters:
  • star – The Star to initialize.

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

  • default_init – The default initilization method if the user doesn’t specify one. [default: None]

Returns:

a star instance with the appropriate initial fit values

initialize_iteration()[source]

Do any required initialization at the start of an iteration of fitting.

Usually a no op, but in PixelGrid, it resets the cached stepk/maxk, so they can potentially adapt to a new size for this iteration.

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, gsparams=None, atmo_type='VonKarman', logger=None, **kwargs)[source]

Initialize the Optical+Atmosphere Model

Use type name “Optical” in a config field to use this model.

There are four components to this model that are convolved together.

First, there is an optical component, which uses a galsim.OpticalPSF to model the profile. 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]

  • [default (oversampling Numerical factor by which to oversample the FFT) – 1]

  • [default – 1]

  • 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]

  • base_aberrations – A list of aberrations to which any fitted aberrations are added. [default: None]

The next two parameters are specific to our treatment of the Optical PSF. If given, they are used to create an additional galsim.PhaseScreen describing the mirror’s own contribution to the wavefront.

Parameters:
  • mirror_figure_im – The name of a file containing an additional phase contribution, for example from the figure of the primary mirror. [default: None]

  • mirror_figure_scale – The pixel scale to use in mirror_figure_im. If omittited, Piff will use the pixel scale in the image itself. [default: None]

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

:param atmo_type The name of the Atmospheric kernel. [default ‘VonKarman’]

Finally, there is both an additional Gaussian component to describe CCD diffusion and separately an applied shear.

:param sigma Gaussian CCD diffusion in arcsec [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_templates to use for setting values of these aperture parameters. [default: None]

  • gsparams – A key word in the dict piff.optical_model.gsparams_templates to use for setting values of these Galsim 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.Optical(template='des', lam=1000)

Note that the Zernike coeffients (zernike_coeff), Atmospheric parameters (ie. r0,L0) and Shear (g1,g2) are fitted parameters, passed via the arguments to getProfile.

fit(star, logger=None, convert_func=None, draw_method=None)[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

getAtmosphere(r0, L0=None, g1=0.0, g2=0.0)[source]

Get the Atmospheric PSF component.

Parameters:
  • r0 – Fried parameter

  • L0 – Outer scale

  • g1 – Shear g1 component

  • g2 – Shear g2 component

Returns:

a galsim.GSObject instance

getOptics(zernike_coeff)[source]

Get the strictly Optics PSF component. Use two phase screens, one from mirror_figure_im, and the other from zernike aberrations

Parameters:

zernike_coeff – An array or list with Zernike Coefficients, in units of waves. Array indexed as [0, 0, z2, z3, z4,…z11…]

Returns:

a galsim.GSObject instance

getProfile(params=None, zernike_coeff=None, r0=None, L0=None, g1=0.0, g2=0.0)[source]

Get a version of the model as a GalSim GSObject

:param params An ndarray with the parameters ordered via idx_PARAM data members :param zernike_coeff A list of the Zernike coefficients in units of [waves] :param r0 Atmospheric Fried parameter [meters] :param L0 Atmospheric Outer scale [meters] :param g1 Atmospheric g1 Shear :param g2 Atmospheric g2 Shear

Returns:

a galsim.GSObject instance of the combined optics, atmosphere and diffusion PSF

initialize(star, logger=None, default_init=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]

  • default_init – Ignored. [default: None]

Returns:

a new initialized Star.

kwargs_to_params(zernike_coeff=[0, 0, 0, 0, 0], r0=0.15, L0=10.0, g1=0.0, g2=0.0)[source]

Fill params ndarray from kwargs

Parameters:
  • zernike_coeff – A ndarray or list with [0, 0, z2, z3, z4, z5, z6…z11..z37] [default=[0,0,0,0,0]]

  • r0 – Value of r0 for Atmospheric Kernel. [default=0.15]

  • L0 – Value of L0 for Atmospheric Kernel. [default=10.0]

  • g1 – Value of g1 to Shear PSF [default=0.]

  • g2 – Value of g2 to Shear PSF [default=0.]

Returns:

an ndarray of parameters

params_to_kwargs(params)[source]

Fill kwarg from params ndarray or list

:param params An ndarray or list with parameters in order

Returns:

a dictionary of named parameters