PSF classes

The PSF classes define the interaction between the models and the interpolants. The simplest version SimplePSF has one model and one interpolant, but it is possible to have more complicated combinations of these.

class piff.PSF[source]

The base class for describing a PSF model across a field of view.

The usual way to create a PSF is through one of the two factory functions:

>>> psf = piff.PSF.process(config, logger)
>>> psf = piff.PSF.read(file_name, logger)

The first is used to build a PSF model from the data according to a config dict. The second is used to read in a PSF model from disk.

draw(x, y, chipnum=None, flux=1.0, center=None, offset=None, stamp_size=48, image=None, logger=None, **kwargs)[source]

Draws an image of the PSF at a given location.

The normal usage would be to specify (chipnum, x, y), in which case Piff will use the stored wcs information for that chip to interpolate to the given position and draw an image of the PSF:

>>> image = psf.draw(chipnum=4, x=103.3, y=592.0, stamp_size=48)

However, if the PSF interpolation used extra properties for the interpolation (cf. psf.interp_property_names), you need to provide them as additional kwargs.

>>> print(psf.interp_property_names)
('u','v','ri_color')
>>> image = psf.draw(chipnum=4, x=103.3, y=592.0, ri_color=0.23, stamp_size=48)

Normally, the image is constructed automatically based on stamp_size, in which case the WCS will be taken to be the local Jacobian at this location on the original image. However, if you provide your own image using the :image: argument, then whatever WCS is present in that image will be respected. E.g. if you want an image of the PSF in sky coordinates rather than image coordinates, you can provide an image with just a pixel scale for the WCS.

When drawing the PSF, there are a few options regarding how the profile will be centered on the image.

  1. The default behavior (center==None) is to draw the profile centered at the same (x,y) as you requested for the location of the PSF in the original image coordinates. The returned image will not (normally) be as large as the full image – it will just be a postage stamp centered roughly around (x,y). The image.bounds give the bounding box of this stamp, and within this, the PSF will be centered at position (x,y).

  2. If you want to center the profile at some other arbitrary position, you may provide a center parameter, which should be a tuple (xc,yc) giving the location at which you want the PSF to be centered. The bounding box will still be around the nominal (x,y) position, so this should only be used for small adjustments to the (x,y) value if you want it centered at a slightly different location.

  3. If you provide your own image with the image parameter, then you may set the center to any location in this box (or technically off it – it doesn’t check that the center is actually inside the bounding box). This may be useful if you want to draw on an image with origin at (0,0) or (1,1) and just put the PSF at the location you want.

  4. If you want the PSf centered exactly in the center of the image, then you can use center=True. This will work for either an automatically built image or one that you provide.

  5. With any of the above options you may additionally supply an offset parameter, which will apply a slight offset to the calculated center. This is probably only useful in conjunction with the default center=None or center=True.

Parameters
  • x – The x position of the desired PSF in the original image coordinates.

  • y – The y position of the desired PSF in the original image coordinates.

  • chipnum – Which chip to use for WCS information. [required if the psf model covers more than a single chip]

  • flux – Flux of PSF to be drawn [default: 1.0]

  • center – (xc,yc) tuple giving the location on the image where you want the nominal center of the profile to be drawn. Also allowed is the value center=True to place in the center of the image. [default: None, which means draw at the position (x,y) of the PSF.]

  • offset – Optional (dx,dy) tuple giving an additional offset relative to the center. [default: None]

  • stamp_size – The size of the image to construct if no image is provided. [default: 48]

  • image – An existing image on which to draw, if desired. [default: None]

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

  • **kwargs – Any additional properties required for the interpolation.

Returns

A GalSim Image of the PSF

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

Generate PSF image for a given star.

Note

If the star already has the fit parameters calculated, then this will trust those values and not redo the interpolation. If this might be a concern, you can force the interpolation to be redone by running

>>> star = psf.interpolateList(star)

before running drawStar.

Parameters
  • star – Star instance holding information needed for interpolation as well as an image/WCS into which PSF will be rendered.

  • 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

Star instance with its image filled with rendered PSF

drawStarList(stars, copy_image=True)[source]

Generate PSF images for given stars. Takes advantage of interpolateList for significant speedup with some interpolators.

Note

If the stars already have the fit parameters calculated, then this will trust those values and not redo the interpolation. If this might be a concern, you can force the interpolation to be redone by running

>>> stars = psf.interpolateList(stars)

before running drawStarList.

Parameters
  • stars – List of Star instances holding information needed for interpolation as well as an image/WCS into which PSF will be rendered.

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

Returns

List of Star instances with its image filled with rendered PSF

get_profile(x, y, chipnum=None, flux=1.0, logger=None, **kwargs)[source]

Get the PSF profile at the given position as a GalSim GSObject.

The normal usage would be to specify (chipnum, x, y), in which case Piff will use the stored wcs information for that chip to interpolate to the given position and draw an image of the PSF:

>>> prof, method = psf.get_profile(chipnum=4, x=103.3, y=592.0)

The first return value, prof, is the GSObject describing the PSF profile. The second one, method, is the method parameter that should be used when drawing the profile using prof.drawImage(..., method=method). This may be either ‘no_pixel’ or ‘auto’ depending on whether the PSF model already includes the pixel response or not. Some underlying models includ the pixel response, and some don’t, so this difference needs to be accounted for properly when drawing. This method is also appropriate if you first convolve the PSF by some other (e.g. galaxy) profile and then draw that.

If the PSF interpolation used extra properties for the interpolation (cf. psf.extra_interp_properties), you need to provide them as additional kwargs.

>>> print(psf.extra_interp_properties)
('ri_color',)
>>> prof, method = psf.get_profile(chipnum=4, x=103.3, y=592.0, ri_color=0.23)
Parameters
  • x – The x position of the desired PSF in the original image coordinates.

  • y – The y position of the desired PSF in the original image coordinates.

  • chipnum – Which chip to use for WCS information. [required if the psf model covers more than a single chip]

  • flux – Flux of PSF model [default: 1.0]

  • **kwargs – Any additional properties required for the interpolation.

Returns

(profile, method) profile = A GalSim GSObject of the PSF method = either ‘no_pixel’ or ‘auto’ indicating which method to use when drawing the profile on an image.

interpolateStar(star)[source]

Update the star to have the current interpolated fit parameters according to the current PSF model.

Parameters

star – Star instance to update.

Returns

Star instance with its fit parameters updated.

interpolateStarList(stars)[source]

Update the stars to have the current interpolated fit parameters according to the current PSF model.

Parameters

stars – List of Star instances to update.

Returns

List of Star instances with their fit parameters updated.

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

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

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

Returns

a kwargs dict to pass to the initializer

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

Process the config dict and return a PSF instance.

As the PSF class is an abstract base class, the returned type will in fact be some subclass of PSF according to the contents of the config dict.

The provided config dict is typically the ‘psf’ field in the base config dict in a YAML file, although for compound PSF types, it may be the field for just one of several components.

This function merely creates a “blank” PSF object. It does not actually do any part of the solution yet. Typically this will be followed by fit:

>>> psf = piff.PSF.process(config['psf'])
>>> stars, wcs, pointing = piff.Input.process(config['input'])
>>> psf.fit(stars, wcs, pointing)

at which point, the psf instance would have a solution to the PSF model.

Parameters
  • config_psf – A dict specifying what type of PSF to build along with the appropriate kwargs for building it.

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

Returns

a PSF instance of the appropriate type.

classmethod read(file_name, logger=None)[source]

Read a PSF object from a file.

Parameters
  • file_name – The name of the file to read.

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

Returns

a PSF instance

classmethod readWCS(fits, extname, logger)[source]

Read the WCS information from a FITS file.

Parameters
  • fits – An open fitsio.FITS object

  • extname – The name of the extension to read from

  • logger – A logger object for logging debug info.

Returns

wcs, pointing where wcs is a dict of galsim.BaseWCS instances and pointing is a galsim.CelestialCoord instance

write(file_name, logger=None)[source]

Write a PSF object to a file.

Parameters
  • file_name – The name of the file to write to.

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

writeWCS(fits, extname, logger)[source]

Write the WCS information to a FITS file.

Parameters
  • fits – An open fitsio.FITS object

  • extname – The name of the extension to write to

  • logger – A logger object for logging debug info.

The simple case of one model/interp

class piff.SimplePSF(model, interp, outliers=None, chisq_thresh=0.1, max_iter=30)[source]

A PSF class that uses a single model and interpolator.

A SimplePSF is built from a Model and an Interp object. The model defines the functional form of the surface brightness profile, and the interpolator defines how the parameters of the model vary across the field of view.

Parameters
  • model – A Model instance used for modeling the surface brightness profile.

  • interp – An Interp instance used to interpolate across the field of view.

  • outliers – Optionally, an Outliers instance used to remove outliers. [default: None]

  • chisq_thresh – Change in reduced chisq at which iteration will terminate. [default: 0.1]

  • max_iter – Maximum number of iterations to try. [default: 30]

fit(stars, wcs, pointing, logger=None, convert_func=None)[source]

Fit interpolated PSF model to star data using standard sequence of operations.

Parameters
  • stars – A list of Star instances.

  • wcs – A dict of WCS solutions indexed by chipnum.

  • pointing – A galsim.CelestialCoord object giving the telescope pointing. [Note: pointing should be None if the WCS is not a CelestialWCS]

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

interpolateStar(star)[source]

Update the star to have the current interpolated fit parameters according to the current PSF model.

Parameters

star – Star instance to update.

Returns

Star instance with its fit parameters updated.

interpolateStarList(stars)[source]

Update the stars to have the current interpolated fit parameters according to the current PSF model.

Parameters

stars – List of Star instances to update.

Returns

List of Star instances with their fit parameters updated.

classmethod parseKwargs(config_psf, logger)[source]

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

Parameters
  • config_psf – The psf field of the configuration dict, config[‘psf’]

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

Returns

a kwargs dict to pass to the initializer

Using a separate solution for each chip

class piff.SingleChipPSF(single_psf, nproc=1)[source]

A PSF class that uses a separate PSF solution for each chip

Parameters
  • single_psf – A PSF instance to use for the PSF solution on each chip. (This will be turned into nchips copies of the provided object.)

  • nproc – How many multiprocessing processes to use for running multiple chips at once. [default: 1]

fit(stars, wcs, pointing, logger=None, convert_func=None)[source]

Fit interpolated PSF model to star data using standard sequence of operations.

Parameters
  • stars – A list of Star instances.

  • wcs – A dict of WCS solutions indexed by chipnum.

  • pointing – A galsim.CelestialCoord object giving the telescope pointing. [Note: pointing should be None if the WCS is not a CelestialWCS]

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

interpolateStar(star)[source]

Update the star to have the current interpolated fit parameters according to the current PSF model.

Parameters

star – Star instance to update.

Returns

Star instance with its fit parameters updated.

classmethod parseKwargs(config_psf, logger)[source]

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

Parameters
  • config_psf – The psf field of the configuration dict, config[‘psf’]

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

Returns

a kwargs dict to pass to the initializer