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.process(config, logger) >>> psf = piff.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.
- _write(writer, name, logger)[source]
This is the function that actually does the work for the write function. Composite PSF classes that need to iterate can call this multiple times as needed.
- Parameters:
writer – A writer object that encapsulates the serialization format.
name – A name to associate with this PSF in the serialized output.
logger – A logger object for logging debug info.
- classmethod _read(reader, name, logger)[source]
This is the function that actually does the work for the read function. Composite PSF classes that need to iterate can call this multiple times as needed.
- Parameters:
reader – A reader object that encapsulates the serialization format.
name – Name associated with this PSF in the serialized output.
logger – A logger object for logging debug info.
- 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.
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).If you want to center the profile at some other arbitrary position, you may provide a
centerparameter, 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.If you provide your own image with the
imageparameter, then you may set thecenterto 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.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.With any of the above options you may additionally supply an
offsetparameter, which will apply a slight offset to the calculated center. This is probably only useful in conjunction with the defaultcenter=Noneorcenter=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)[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]
- 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.
- Returns:
List of Star instances with its image filled with rendered PSF
- fit(stars, wcs, pointing, logger=None, convert_funcs=None, draw_method=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_funcs – An optional list of function to apply to the profiles being fit before drawing it onto the image. This is used by composite PSFs to isolate the effect of just this model component. If provided, it should be the same length as stars. [default: None]
draw_method – The method to use with the GalSim drawImage command. If not given, use the default method for the PSF model being fit. [default: None]
- property fit_center
Whether to fit the center of the star in reflux.
This is generally set in the model specifications. If all component models includes a shift, then this is False. Otherwise it is True.
- 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.interp_property_names), you need to provide them as additional kwargs.
>>> print(psf.interp_property_names) ('u','v','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.
- property include_model_centroid
Whether a model that we want to center can have a non-zero centroid during iterations.
- initialize_flux_center(stars, logger=None)[source]
Initialize the flux and center of the stars.
The flux is just a simple sum of unmasked pixels. The center is a simple centroid relative to the nominal position. It is only updated if the PSF model is centered. (I.e. if self.fit_center is True.)
- Parameters:
stars – The initial list of Star instances that will be used to constrain the PSF.
logger – A logger object for logging progress. [default: None]
- Returns:
the initialized stars
- initialize_params(stars, logger=None, default_init=None)[source]
Initialize the psf solver to begin an iterative solution.
- Parameters:
stars – The initial list of Star instances that will be used to constrain the PSF.
logger – A logger object for logging progress. [default: None]
default_init – The default initilization method if the user doesn’t specify one. [default: None]
- Returns:
the initialized stars, nremoved
- interpolateStar(star, inplace=False)[source]
Update the star to have the current interpolated fit parameters according to the current PSF model.
- Parameters:
star – Star instance to update.
inplace – Whether to update the parameters in place, in which case the returned star is the same object as the input star. [default: False]
- Returns:
Star instance with its fit parameters updated.
- interpolateStarList(stars, inplace=False)[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.
inplace – Whether to update the parameters in place, in which case the returned stars are the same objects as the input stars. [default: False]
- 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
psfinstance 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
- reflux(star, logger=None)[source]
Fit the PSF 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
logger – A logger object for logging debug info. [default: None]
- Returns:
New Star instance, with updated flux, center, chisq, dof
- remove_outliers(stars, iteration, logger)[source]
Look for and flag outliers from the list of stars
- Parameters:
stars – The complete list of stars to consider
iteration – The number of the iteration that was just completed.
logger – A logger object for logging progress.
- Returns:
new_stars, nremoved
- set_num(num)[source]
If there are multiple components involved in the fit, set the number to use for this model.
- single_iteration(stars, logger, convert_funcs, draw_method)[source]
Perform a single iteration of the solver.
Note that some object might fail at some point in the fitting, so some objects can be flagged as bad during this step, prior to the outlier rejection step. This information is reported in the return tuple as nremoved.
- Parameters:
stars – The list of stars to use for constraining the PSF.
logger – A logger object for logging progress.
convert_funcs – An optional list of function to apply to the profiles being fit before drawing it onto the image. If not None, it should be the same length as stars.
draw_method – The method to use with the GalSim drawImage command. If None, use the default method for the PSF model being fit.
- Returns:
an updated list of all_stars, nremoved
The simple case of one model/interp
- class piff.SimplePSF(model, interp, outliers=None, chisq_thresh=0.1, min_iter=2, 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.
Use type name “Simple” in a config field to use this psf type, or leave off the type name, as this is the default PSF type.
- 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]
min_iter – Minimum number of iterations to try. [default: 2]
max_iter – Maximum number of iterations to try. [default: 30]
- _finish_write(writer, logger)[source]
Finish the writing process with any class-specific steps.
- Parameters:
writer – A writer object that encapsulates the serialization format.
logger – A logger object for logging debug info.
- _finish_read(reader, logger)[source]
Finish the reading process with any class-specific steps.
- Parameters:
reader – A reader object that encapsulates the serialization format.
name – Name associated with this PSF in the serialized output.
logger – A logger object for logging debug info.
- property fit_center
Whether to fit the center of the star in reflux.
This is generally set in the model specifications. If all component models includes a shift, then this is False. Otherwise it is True.
- property include_model_centroid
Whether a model that we want to center can have a non-zero centroid during iterations.
- initialize_params(stars, logger=None, default_init=None)[source]
Initialize the psf solver to begin an iterative solution.
- Parameters:
stars – The initial list of Star instances that will be used to constrain the PSF.
logger – A logger object for logging progress. [default: None]
default_init – The default initilization method if the user doesn’t specify one. [default: None]
- Returns:
the initialized stars, nremoved
- interpolateStar(star, inplace=False)[source]
Update the star to have the current interpolated fit parameters according to the current PSF model.
- Parameters:
star – Star instance to update.
inplace – Whether to update the parameters in place, in which case the returned star is the same object as the input star. [default: False]
- Returns:
Star instance with its fit parameters updated.
- interpolateStarList(stars, inplace=False)[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.
inplace – Whether to update the parameters in place, in which case the returned stars are the same objects as the input stars. [default: False]
- 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
- set_num(num)[source]
If there are multiple components involved in the fit, set the number to use for this model.
- single_iteration(stars, logger, convert_funcs, draw_method)[source]
Perform a single iteration of the solver.
Note that some object might fail at some point in the fitting, so some objects can be flagged as bad during this step, prior to the outlier rejection step. This information is reported in the return tuple as nremoved.
- Parameters:
stars – The list of stars to use for constraining the PSF.
logger – A logger object for logging progress.
convert_funcs – An optional list of function to apply to the profiles being fit before drawing it onto the image. If not None, it should be the same length as stars.
draw_method – The method to use with the GalSim drawImage command. If None, use the default method for the PSF model being fit.
- Returns:
an updated list of all_stars, nremoved
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
Use type name “SingleChip” in a config field to use this psf type.
- 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]
- _finish_write(writer, logger)[source]
Finish the writing process with any class-specific steps.
- Parameters:
writer – A writer object that encapsulates the serialization format.
logger – A logger object for logging debug info.
- _finish_read(reader, logger)[source]
Finish the reading process with any class-specific steps.
- Parameters:
reader – A reader object that encapsulates the serialization format.
logger – A logger object for logging debug info.
- fit(stars, wcs, pointing, logger=None, convert_funcs=None, draw_method=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_funcs – An optional list of function to apply to the profiles being fit before drawing it onto the image. This is used by composite PSFs to isolate the effect of just this model component. If provided, it should be the same length as stars. [default: None]
draw_method – The method to use with the GalSim drawImage command. If not given, use the default method for the PSF model being fit. [default: None]
- property fit_center
Whether to fit the center of the star in reflux.
This is generally set in the model specifications. If all component models includes a shift, then this is False. Otherwise it is True.
- property include_model_centroid
Whether a model that we want to center can have a non-zero centroid during iterations.
- interpolateStar(star, inplace=False)[source]
Update the star to have the current interpolated fit parameters according to the current PSF model.
- Parameters:
star – Star instance to update.
inplace – Whether to update the parameters in place, in which case the returned star is the same object as the input star. [default: False]
- 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
Composite PSF types
- class piff.SumPSF(components, outliers=None, chisq_thresh=0.1, min_iter=2, max_iter=30)[source]
A PSF class that is the Sum of two or more other PSF types.
A SumPSF is built from an ordered list of other PSFs.
When fitting the Sum, the default pattern is that all components after the first one are initialized to zero, and the first component is fit as usual, but just using a single iteration of the fit. Then the residuals of this model are fit using the second component, and so on. Once all components are fit, outliers may be rejected, and then the process is iterated.
This pattern can be tweaked somewhat using the initialization options available to PSF models. If a component should be initialized to something other than a zero model, then one should explicitly set it.
Use type name “Sum” in a config field to use this psf type.
- Parameters:
components – A list of PSF instances defining the components to be summed.
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]
min_iter – Minimum number of iterations to try. [default: 2]
max_iter – Maximum number of iterations to try. [default: 30]
- _finish_write(writer, logger)[source]
Finish the writing process with any class-specific steps.
- Parameters:
writer – A writer object that encapsulates the serialization format.
logger – A logger object for logging debug info.
- _finish_read(reader, logger)[source]
Finish the reading process with any class-specific steps.
- Parameters:
reader – A reader object that encapsulates the serialization format.
logger – A logger object for logging debug info.
- property fit_center
Whether to fit the center of the star in reflux.
This is generally set in the model specifications. If all component models includes a shift, then this is False. Otherwise it is True.
- property include_model_centroid
Whether a model that we want to center can have a non-zero centroid during iterations.
- initialize_params(stars, logger, default_init=None)[source]
Initialize the psf solver to begin an iterative solution.
- Parameters:
stars – The initial list of Star instances that will be used to constrain the PSF.
logger – A logger object for logging progress. [default: None]
default_init – The default initilization method if the user doesn’t specify one. [default: None]
- Returns:
the initialized stars, nremoved
- interpolateStar(star, inplace=False)[source]
Update the star to have the current interpolated fit parameters according to the current PSF model.
- Parameters:
star – Star instance to update.
inplace – Whether to update the parameters in place, in which case the returned star is the same object as the input star. [default: False]
- Returns:
Star instance with its fit parameters updated.
- interpolateStarList(stars, inplace=False)[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.
inplace – Whether to update the parameters in place, in which case the returned stars are the same objects as the input stars. [default: False]
- 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
- set_num(num)[source]
If there are multiple components involved in the fit, set the number to use for this model.
- single_iteration(stars, logger, convert_funcs, draw_method)[source]
Perform a single iteration of the solver.
Note that some object might fail at some point in the fitting, so some objects can be flagged as bad during this step, prior to the outlier rejection step. This information is reported in the return tuple as nremoved.
- Parameters:
stars – The list of stars to use for constraining the PSF.
logger – A logger object for logging progress.
convert_funcs – An optional list of function to apply to the profiles being fit before drawing it onto the image. If not None, it should be the same length as stars.
draw_method – The method to use with the GalSim drawImage command. If None, use the default method for the PSF model being fit.
- Returns:
an updated list of all_stars, nremoved
- class piff.ConvolvePSF(components, outliers=None, chisq_thresh=0.1, min_iter=2, max_iter=30)[source]
A PSF class that is the Convolution of two or more other PSF types.
A ConvolvePSF is built from an ordered list of other PSFs.
When fitting the convolution, the default pattern is that all components after the first one are initialized as (approximately) a delta function, and the first component is fit as usual, but just using a single iteration of the fit. Then the residuals of this model are fit using the second component, and so on. Once all components are fit, outliers may be rejected, and then the process is iterated.
This pattern can be tweaked somewhat using the initialization options available to PSF models. If a component should be initialized to something other than a delta-function. then one should explicitly set it.
Use type name “Convolve” in a config field to use this psf type.
- Parameters:
components – A list of PSF instances defining the components to be convolved.
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]
min_iter – Minimum number of iterations to try. [default: 2]
max_iter – Maximum number of iterations to try. [default: 30]
- _finish_write(writer, logger)[source]
Finish the writing process with any class-specific steps.
- Parameters:
writer – A writer object that encapsulates the serialization format.
logger – A logger object for logging debug info.
- _finish_read(reader, logger)[source]
Finish the reading process with any class-specific steps.
- Parameters:
reader – A reader object that encapsulates the serialization format.
logger – A logger object for logging debug info.
- property fit_center
Whether to fit the center of the star in reflux.
This is generally set in the model specifications. If all component models includes a shift, then this is False. Otherwise it is True.
- property include_model_centroid
Whether a model that we want to center can have a non-zero centroid during iterations.
- initialize_params(stars, logger, default_init=None)[source]
Initialize the psf solver to begin an iterative solution.
- Parameters:
stars – The initial list of Star instances that will be used to constrain the PSF.
logger – A logger object for logging progress. [default: None]
default_init – The default initilization method if the user doesn’t specify one. [default: None]
- Returns:
the initialized stars, nremoved
- interpolateStar(star, inplace=False)[source]
Update the star to have the current interpolated fit parameters according to the current PSF model.
- Parameters:
star – Star instance to update.
inplace – Whether to update the parameters in place, in which case the returned star is the same object as the input star. [default: False]
- Returns:
Star instance with its fit parameters updated.
- interpolateStarList(stars, inplace=False)[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.
inplace – Whether to update the parameters in place, in which case the returned stars are the same objects as the input stars. [default: False]
- 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
- set_num(num)[source]
If there are multiple components involved in the fit, set the number to use for this model.
- single_iteration(stars, logger, convert_funcs, draw_method)[source]
Perform a single iteration of the solver.
Note that some object might fail at some point in the fitting, so some objects can be flagged as bad during this step, prior to the outlier rejection step. This information is reported in the return tuple as nremoved.
- Parameters:
stars – The list of stars to use for constraining the PSF.
logger – A logger object for logging progress.
convert_funcs – An optional list of function to apply to the profiles being fit before drawing it onto the image. If not None, it should be the same length as stars.
draw_method – The method to use with the GalSim drawImage command. If None, use the default method for the PSF model being fit.
- Returns:
an updated list of all_stars, nremoved