Point-spread functions

There are a number of profiles that are designed to be used for point-spread functions (PSFs). There is nothing in GalSim restricting these classes to be used only for PSFs or stars, but that was the use case they were designed for.

Airy Profile

class galsim.Airy(lam_over_diam=None, lam=None, diam=None, obscuration=0.0, flux=1.0, scale_unit=None, gsparams=None)[source]

Bases: GSObject

A class describing the surface brightness profile for an Airy disk (perfect diffraction-limited PSF for a circular aperture), with an optional central obscuration.

For more information, refer to

http://en.wikipedia.org/wiki/Airy_disc

The Airy profile is defined in terms of the diffraction angle, which is a function of the ratio lambda / D, where lambda is the wavelength of the light (say in the middle of the bandpass you are using) and D is the diameter of the telescope.

The natural units for this value is radians, which is not normally a convenient unit to use for other GSObject dimensions. Assuming that the other sky coordinates you are using are all in arcsec (e.g. the pixel scale when you draw the image, the size of the galaxy, etc.), then you should convert this to arcsec as well:

>>> lam = 700  # nm
>>> diam = 4.0    # meters
>>> lam_over_diam = (lam * 1.e-9) / diam  # radians
>>> lam_over_diam *= 206265  # Convert to arcsec
>>> airy = galsim.Airy(lam_over_diam)

To make this process a bit simpler, we recommend instead providing the wavelength and diameter separately using the parameters lam (in nm) and diam (in m). GalSim will then convert this to any of the normal kinds of angular units using the scale_unit parameter:

>>> airy = galsim.Airy(lam=lam, diam=diam, scale_unit=galsim.arcsec)

When drawing images, the scale_unit should match the unit used for the pixel scale or the WCS. e.g. in this case, a pixel scale of 0.2 arcsec/pixel would be specified as pixel_scale=0.2.

Parameters:
  • lam_over_diam – The parameter that governs the scale size of the profile. See above for details about calculating it.

  • lam – Lambda (wavelength) in units of nanometers. Must be supplied with diam, and in this case, image scales (scale) should be specified in units of scale_unit.

  • diam – Telescope diameter in units of meters. Must be supplied with lam, and in this case, image scales (scale) should be specified in units of scale_unit.

  • obscuration – The linear dimension of a central obscuration as a fraction of the pupil dimension. [default: 0]

  • flux – The flux (in photons/cm^2/s) of the profile. [default: 1]

  • scale_unit – Units to use for the sky coordinates when calculating lam/diam if these are supplied separately. Note that the results of using properties like fwhm will be returned in units of scale_unit as well. Should be either a galsim.AngleUnit or a string that can be used to construct one (e.g., ‘arcsec’, ‘radians’, etc.). [default: galsim.arcsec]

  • gsparams – An optional GSParams argument. [default: None]

property fwhm

The FWHM of this Airy profile (only supported for obscuration = 0.).

property half_light_radius

The half light radius of this Airy profile (only supported for obscuration = 0.).

property lam_over_diam

The input lambda/diam value.

property obscuration

The input obscuration.

withFlux(flux)[source]

Create a version of the current object with a different flux.

This function is equivalent to obj.withScaledFlux(flux / obj.flux).

It creates a new object that has the same profile as the original, but with the surface brightness at every location rescaled such that the total flux will be the given value. Note that if flux is an SED, the return value will be a ChromaticObject with specified SED.

Parameters:

flux – The new flux for the object.

Returns:

the object with the new flux

Moffat Profile

class galsim.Moffat(beta, scale_radius=None, half_light_radius=None, fwhm=None, trunc=0.0, flux=1.0, gsparams=None)[source]

Bases: GSObject

A class describing a Moffat surface brightness profile.

The Moffat surface brightness profile is

\[I(R) \sim \left(1 + (r/r_0)^2\right)^{-\beta}\]

where \(r_0\) is scale_radius.

The GalSim representation of a Moffat profile also includes an optional truncation beyond a given radius.

For more information, refer to

A Moffat can be initialized using one (and only one) of three possible size parameters: scale_radius, fwhm, or half_light_radius. Exactly one of these three is required.

Parameters:
  • beta – The beta parameter of the profile.

  • scale_radius – The scale radius of the profile. Typically given in arcsec. [One of scale_radius, fwhm, or half_light_radius is required.]

  • half_light_radius – The half-light radius of the profile. Typically given in arcsec. [One of scale_radius, fwhm, or half_light_radius is required.]

  • fwhm – The full-width-half-max of the profile. Typically given in arcsec. [One of scale_radius, fwhm, or half_light_radius is required.]

  • trunc – An optional truncation radius at which the profile is made to drop to zero, in the same units as the size parameter. [default: 0, indicating no truncation]

  • flux – The flux (in photons/cm^2/s) of the profile. [default: 1]

  • gsparams – An optional GSParams argument. [default: None]

property beta

The beta parameter of this Moffat profile.

property half_light_radius

The half-light radius of this Moffat profile.

property scale_radius

The scale radius of this Moffat profile.

property trunc

The truncation radius (if any) of this Moffat profile.

withFlux(flux)[source]

Create a version of the current object with a different flux.

This function is equivalent to obj.withScaledFlux(flux / obj.flux).

It creates a new object that has the same profile as the original, but with the surface brightness at every location rescaled such that the total flux will be the given value. Note that if flux is an SED, the return value will be a ChromaticObject with specified SED.

Parameters:

flux – The new flux for the object.

Returns:

the object with the new flux

Kolmogorov Profile

class galsim.Kolmogorov(lam_over_r0=None, fwhm=None, half_light_radius=None, lam=None, r0=None, r0_500=None, flux=1.0, scale_unit=None, gsparams=None)[source]

Bases: GSObject

A class describing a Kolmogorov surface brightness profile, which represents a long exposure atmospheric PSF.

The Kolmogorov profile is defined in Fourier space. Its transfer function is:

\[T(k) \sim \exp(-D(k)/2)\]

where

\[D(k) = 6.8839 \left(\frac{\lambda k}{2\pi r_0} \right)^{5/3},\]

\(\lambda\) is the wavelength of the light (say in the middle of the bandpass you are using), and \(r_0\) is the Fried parameter. Typical values for the Fried parameter are on the order of 0.1 m for most observatories and up to 0.2 m for excellent sites. The values are usually quoted at \(\lambda\) = 500nm and \(r_0\) depends on wavelength as \(r_0 \sim \lambda^{6/5}\).

For more information, refer to

The Kolmogorov profile is normally defined in terms of the ratio \(\lambda / r_0\). The natural units for this ratio is radians, which is not normally a convenient unit to use for other GSObject dimensions. Assuming that the other sky coordinates you are using are all in arcsec (e.g. the pixel scale when you draw the image, the size of the galaxy, etc.), then you should convert this to arcsec as well:

>>> lam = 700  # nm
>>> r0 = 0.15 * (lam/500)**1.2  # meters
>>> lam_over_r0 = (lam * 1.e-9) / r0  # radians
>>> lam_over_r0 *= 206265  # Convert to arcsec
>>> psf = galsim.Kolmogorov(lam_over_r0)

To make this process a bit simpler, we recommend instead providing the wavelength and Fried parameter separately using the parameters lam (in nm) and either r0 or r0_500 (in m). GalSim will then convert this to any of the normal kinds of angular units using the scale_unit parameter:

>>> psf = galsim.Kolmogorov(lam=lam, r0=r0, scale_unit=galsim.arcsec)

or:

>>> psf = galsim.Kolmogorov(lam=lam, r0_500=0.15, scale_unit=galsim.arcsec)

When drawing images, the scale_unit should match the unit used for the pixel scale or the WCS. e.g. in this case, a pixel scale of 0.2 arcsec/pixel would be specified as pixel_scale=0.2.

A Kolmogorov object may also be initialized using fwhm or half_light_radius. Exactly one of these four ways to define the size is required.

The FWHM of the Kolmogorov PSF is \(\sim 0.976 \lambda/r_0\) arcsec. (e.g., Racine 1996, PASP 699, 108).

Parameters:
  • lam_over_r0 – The parameter that governs the scale size of the profile. See above for details about calculating it. [One of lam_over_r0, fwhm, half_light_radius, or lam (along with either r0 or r0_500) is required.]

  • fwhm – The full-width-half-max of the profile. Typically given in arcsec. [One of lam_over_r0, fwhm, half_light_radius, or lam (along with either r0 or r0_500) is required.]

  • half_light_radius – The half-light radius of the profile. Typically given in arcsec. [One of lam_over_r0, fwhm, half_light_radius, or lam (along with either r0 or r0_500) is required.]

  • lam – Lambda (wavelength) in units of nanometers. Must be supplied with either r0 or r0_500, and in this case, image scales (scale) should be specified in units of scale_unit.

  • r0 – The Fried parameter in units of meters. Must be supplied with lam, and in this case, image scales (scale) should be specified in units of scale_unit.

  • r0_500 – The Fried parameter in units of meters at 500 nm. The Fried parameter at the given wavelength, lam, will be computed using the standard relation r0 = r0_500 * (lam/500)**1.2.

  • flux – The flux (in photons/cm^2/s) of the profile. [default: 1]

  • scale_unit – Units to use for the sky coordinates when calculating lam/r0 if these are supplied separately. Note that the results of using properties like fwhm will be returned in units of scale_unit as well. Should be either a galsim.AngleUnit or a string that can be used to construct one (e.g., ‘arcsec’, ‘radians’, etc.). [default: galsim.arcsec]

  • gsparams – An optional GSParams argument. [default: None]

In addition to the usual GSObject methods, Kolmogorov has the following access properties:

Attributes:
  • lam_over_r0 – The input ratio lambda / r0

  • fwhm – The full-width half-max size

  • half_light_radius – The half-light radius

property fwhm

The FWHM of this Kolmogorov profile.

property half_light_radius

The half light radius of this Kolmogorov profile.

property lam_over_r0

The input lambda / r0 for this Kolmogorov profile.

withFlux(flux)[source]

Create a version of the current object with a different flux.

This function is equivalent to obj.withScaledFlux(flux / obj.flux).

It creates a new object that has the same profile as the original, but with the surface brightness at every location rescaled such that the total flux will be the given value. Note that if flux is an SED, the return value will be a ChromaticObject with specified SED.

Parameters:

flux – The new flux for the object.

Returns:

the object with the new flux

Von Karman Profile

class galsim.VonKarman(lam, r0=None, r0_500=None, L0=25.0, flux=1, scale_unit=coord.arcsec, force_stepk=0.0, do_delta=False, suppress_warning=False, gsparams=None)[source]

Bases: GSObject

Class describing a von Karman surface brightness profile, which represents a long exposure atmospheric PSF. The difference between the von Karman profile and the related Kolmogorov profile is that the von Karman profile includes a parameter for the outer scale of atmospheric turbulence, which is a physical scale beyond which fluctuations in the refractive index stop growing, typically between 10 and 100 meters. Quantitatively, the von Karman phase fluctuation power spectrum at spatial frequency f is proportional to

\[P(f) = r_0^{-5/3} \left(f^2 + L_0^{-2}\right)^{-11/6}\]

where \(r_0\) is the Fried parameter which sets the overall turbulence amplitude and \(L_0\) is the outer scale in meters. The Kolmogorov power spectrum proportional to \(r_0^{-5/3} f^{-11/3}\) is recovered as \(L_0 \rightarrow \infty\).

For more information, we recommend the following references:

Martinez et al. 2010 A&A vol. 516 Conan 2008 JOSA vol. 25

Note

If one blindly follows the math for converting the von Karman power spectrum into a PSF, one finds that the PSF contains a delta-function at the origin with fractional flux of:

\[F_\mathrm{delta} = e^{-0.086 (r_0/L_0)^{-5/3}}\]

In almost all cases of interest this evaluates to something tiny, often on the order of \(10^{-100}\) or smaller. By default, GalSim will ignore this delta function entirely since it usually doesn’t make any difference, but can complicate some calculations like drawing using method=’real_space’ or by formally requiring huge Fourier transforms for drawing using method=’fft’. If for some reason you want to keep the delta function though, then you can pass the do_delta=True argument to the VonKarman initializer.

Parameters:
  • lam – Wavelength in nanometers.

  • r0 – Fried parameter at specified wavelength lam in meters. Exactly one of r0 and r0_500 should be specified.

  • r0_500 – Fried parameter at 500 nm in meters. Exactly one of r0 and r0_500 should be specified.

  • L0 – Outer scale in meters. [default: 25.0]

  • flux – The flux (in photons/cm^2/s) of the profile. [default: 1]

  • scale_unit – Units assumed when drawing this profile or evaluating xValue, kValue, etc. Should be a galsim.AngleUnit or a string that can be used to construct one (e.g., ‘arcsec’, ‘radians’, etc.). [default: galsim.arcsec]

  • force_stepk – By default, VonKarman will derive a value of stepk from a computed real-space surface brightness profile and gsparams settings. Although this profile is cached for future instantiations of identical VonKarman objects, it is relatively slow to compute for the first instance and can dominate the compute time when drawing many VonKaman’s with different parameters using method ‘fft’, ‘sb’, or ‘no_pixel’, a situation that may occur, e.g., in a fitting context. This keyword enables a user to bypass the real-space profile computation by directly specifying a stepk value. Note that if the .half_light_radius property is queried, or the object is drawn using method ‘phot’ or ‘real_space’, then the real-space profile calculation is performed (if not cached) at that point. [default: 0.0, which means do not force a value of stepk]

  • do_delta – Include delta-function at origin? (not recommended; see above). [default: False]

  • suppress_warning – For some combinations of r0 and L0, it may become impossible to satisfy the gsparams.maxk_threshold criterion (see above). In that case, the code will emit a warning alerting the user that they may have entered a non-physical regime. However, this warning can be suppressed with this keyword. [default: False]

  • gsparams – An optional GSParams argument. [default: None]

property L0

The input L0 value.

property delta_amplitude

The amplitude of the delta function at the center.

property do_delta

Whether to include the delta function at the center.

property force_stepk

Forced value of stepk or 0.0.

property half_light_radius

The half-light radius.

property lam

The input lam value.

property r0

The input r0 value.

property r0_500

The input r0_500 value.

property scale_unit

The input scale_units.

withFlux(flux)[source]

Create a version of the current object with a different flux.

This function is equivalent to obj.withScaledFlux(flux / obj.flux).

It creates a new object that has the same profile as the original, but with the surface brightness at every location rescaled such that the total flux will be the given value. Note that if flux is an SED, the return value will be a ChromaticObject with specified SED.

Parameters:

flux – The new flux for the object.

Returns:

the object with the new flux

Optical PSF

class galsim.OpticalPSF(lam_over_diam=None, lam=None, diam=None, tip=0.0, tilt=0.0, defocus=0.0, astig1=0.0, astig2=0.0, coma1=0.0, coma2=0.0, trefoil1=0.0, trefoil2=0.0, spher=0.0, aberrations=None, annular_zernike=False, aper=None, circular_pupil=True, obscuration=0.0, interpolant=None, oversampling=1.5, pad_factor=1.5, ii_pad_factor=None, flux=1.0, nstruts=0, strut_thick=0.05, strut_angle=coord.Angle(0.0, coord.radians), pupil_plane_im=None, pupil_plane_scale=None, pupil_plane_size=None, pupil_angle=coord.Angle(0.0, coord.radians), scale_unit=coord.arcsec, fft_sign='+', gsparams=None, _force_stepk=0.0, _force_maxk=0.0, suppress_warning=False, geometric_shooting=False)[source]

Bases: GSObject

A class describing aberrated PSFs due to telescope optics. Its underlying implementation uses an InterpolatedImage to characterize the profile.

The diffraction effects are characterized by the diffraction angle, which is a function of the ratio lambda / D, where lambda is the wavelength of the light and D is the diameter of the telescope. The natural unit for this value is radians, which is not normally a convenient unit to use for other GSObject dimensions. Assuming that the other sky coordinates you are using are all in arcsec (e.g. the pixel scale when you draw the image, the size of the galaxy, etc.), then you should convert this to arcsec as well:

>>> lam = 700  # nm
>>> diam = 4.0    # meters
>>> lam_over_diam = (lam * 1.e-9) / diam  # radians
>>> lam_over_diam *= 206265  # Convert to arcsec
>>> psf = galsim.OpticalPSF(lam_over_diam, ...)

To make this process a bit simpler, we recommend instead providing the wavelength and diameter separately using the parameters lam (in nm) and diam (in m). GalSim will then convert this to any of the normal kinds of angular units using the scale_unit parameter:

>>> psf = galsim.OpticalPSF(lam=lam, diam=diam, scale_unit=galsim.arcsec, ...)

When drawing images, the scale_unit should match the unit used for the pixel scale or the WCS. e.g. in this case, a pixel scale of 0.2 arcsec/pixel would be specified as pixel_scale=0.2.

Input aberration coefficients are assumed to be supplied in units of wavelength, and correspond to the Zernike polynomials in the Noll convention defined in Noll, J. Opt. Soc. Am. 66, 207-211(1976). For a brief summary of the polynomials, refer to http://en.wikipedia.org/wiki/Zernike_polynomials#Zernike_polynomials. By default, the aberration coefficients indicate the amplitudes of _circular_ Zernike polynomials, which are orthogonal over a circle. If you would like the aberration coefficients to instead be interpretted as the amplitudes of _annular_ Zernike polynomials, which are orthogonal over an annulus (see Mahajan, J. Opt. Soc. Am. 71, 1 (1981)), set the annular_zernike keyword argument to True.

There are two ways to specify the geometry of the pupil plane, i.e., the obscuration disk size and the areas that will be illuminated outside of it. The first way is to use keywords that specify the size of the obscuration, and the nature of the support struts holding up the secondary mirror (or prime focus cage, etc.). These are taken to be rectangular obscurations extending from the outer edge of the pupil to the outer edge of the obscuration disk (or the pupil center if obscuration = 0.). You can specify how many struts there are (evenly spaced in angle), how thick they are as a fraction of the pupil diameter, and what angle they start at relative to the positive y direction.

The second way to specify the pupil plane configuration is by passing in an image of it. This can be useful for example if the struts are not evenly spaced or are not radially directed, as is assumed by the simple model for struts described above. In this case, keywords related to struts are ignored; moreover, the obscuration keyword is used to ensure that the images are properly sampled (so it is still needed), but the keyword is then ignored when using the supplied image of the pupil plane. Note that for complicated pupil configurations, it may be desireable to increase pad_factor for more fidelity at the expense of slower running time. The pupil_plane_im that is passed in can be rotated during internal calculations by specifying a pupil_angle keyword.

If you choose to pass in a pupil plane image, it must be a square array in which the image of the pupil is centered. The areas that are illuminated should have some value >0, and the other areas should have a value of precisely zero. Based on what the OpticalPSF class thinks is the required sampling to make the PSF image, the image that is passed in of the pupil plane might be zero-padded during internal calculations. The pixel scale of the pupil plane can be specified in one of three ways. In descending order of priority, these are:

  1. The pupil_plane_scale keyword argument (units are meters).

  2. The pupil_plane_im.scale attribute (units are meters).

  3. If (1) and (2) are both None, then the scale will be inferred by assuming that the illuminated pixel farthest from the image center is at a physical distance of self.diam/2.

Note that if the scale is specified by either (1) or (2) above (which always includes specifying the pupil_plane_im as a filename, since the default scale then will be 1.0), then the lam_over_diam keyword must not be used, but rather the lam and diam keywords are required separately. Finally, to ensure accuracy of calculations using a pupil plane image, we recommend sampling it as finely as possible.

As described above, either specify the lam/diam ratio directly in arbitrary units:

>>> optical_psf = galsim.OpticalPSF(lam_over_diam=lam_over_diam, defocus=0., ...)

or, use separate keywords for the telescope diameter and wavelength in meters and nanometers, respectively:

>>> optical_psf = galsim.OpticalPSF(lam=lam, diam=diam, defocus=0., ...)

Either of these options initializes optical_psf as an OpticalPSF instance.

Parameters:
  • lam_over_diam – Lambda / telescope diameter in the physical units adopted for scale (user responsible for consistency). Either lam_over_diam, or lam and diam, must be supplied.

  • lam – Lambda (wavelength) in units of nanometers. Must be supplied with diam, and in this case, image scales (scale) should be specified in units of scale_unit.

  • diam – Telescope diameter in units of meters. Must be supplied with lam, and in this case, image scales (scale) should be specified in units of scale_unit.

  • tip – Tip in units of incident light wavelength. [default: 0]

  • tilt – Tilt in units of incident light wavelength. [default: 0]

  • defocus – Defocus in units of incident light wavelength. [default: 0]

  • astig1 – Astigmatism (like e2) in units of incident light wavelength. [default: 0]

  • astig2 – Astigmatism (like e1) in units of incident light wavelength. [default: 0]

  • coma1 – Coma along y in units of incident light wavelength. [default: 0]

  • coma2 – Coma along x in units of incident light wavelength. [default: 0]

  • trefoil1 – Trefoil (one of the arrows along y) in units of incident light wavelength. [default: 0]

  • trefoil2 – Trefoil (one of the arrows along x) in units of incident light wavelength. [default: 0]

  • spher – Spherical aberration in units of incident light wavelength. [default: 0]

  • aberrations – Optional keyword, to pass in a list, tuple, or NumPy array of aberrations in units of reference wavelength (ordered according to the Noll convention), rather than passing in individual values for each individual aberration. Note that aberrations[1] is piston (and not aberrations[0], which is unused.) This list can be arbitrarily long to handle Zernike polynomial aberrations of arbitrary order.

  • annular_zernike – Boolean indicating that aberrations specify the amplitudes of annular Zernike polynomials instead of circular Zernike polynomials. [default: False]

  • aperAperture object to use when creating PSF. [default: None]

  • circular_pupil – Adopt a circular pupil? [default: True]

  • obscuration – Linear dimension of central obscuration as fraction of pupil linear dimension, [0., 1.). This should be specified even if you are providing a pupil_plane_im, since we need an initial value of obscuration to use to figure out the necessary image sampling. [default: 0]

  • interpolant – Either an Interpolant instance or a string indicating which interpolant should be used. Options are ‘nearest’, ‘sinc’, ‘linear’, ‘cubic’, ‘quintic’, or ‘lanczosN’ where N should be the integer order to use. [default: galsim.Quintic()]

  • oversampling – Optional oversampling factor for the InterpolatedImage. Setting oversampling < 1 will produce aliasing in the PSF (not good). Usually oversampling should be somewhat larger than 1. 1.5 is usually a safe choice. [default: 1.5]

  • pad_factor – Additional multiple by which to zero-pad the PSF image to avoid folding compared to what would be employed for a simple Airy. Note that pad_factor may need to be increased for stronger aberrations, i.e. those larger than order unity. [default: 1.5]

  • ii_pad_factor – Zero-padding factor by which to extend the image of the PSF when creating the InterpolatedImage. See the InterpolatedImage docstring for more details. [default: 1.5]

  • suppress_warning – If pad_factor is too small, the code will emit a warning telling you its best guess about how high you might want to raise it. However, you can suppress this warning by using suppress_warning=True. [default: False]

  • geometric_shooting – If True, then when drawing using photon shooting, use geometric optics approximation where the photon angles are derived from the phase screen gradient. If False, then first draw using Fourier optics and then shoot from the derived InterpolatedImage. [default: False]

  • flux – Total flux of the profile. [default: 1.]

  • 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_angleAngle made between the vertical and the strut starting closest to it, defined to be positive in the counter-clockwise direction; must be an Angle instance. [default: 0. * galsim.degrees]

  • pupil_plane_im – The GalSim.Image, NumPy array, or name of file containing the pupil plane image, to be used instead of generating one based on the obscuration and strut parameters. [default: None]

  • pupil_angle – If pupil_plane_im is not None, rotation angle for the pupil plane (positive in the counter-clockwise direction). Must be an Angle instance. [default: 0. * galsim.degrees]

  • pupil_plane_scale – Sampling interval in meters to use for the pupil plane array. In most cases, it’s a good idea to leave this as None, in which case GalSim will attempt to find a good value automatically. The exception is when specifying the pupil arrangement via an image, in which case this keyword can be used to indicate the sampling of that image. See also pad_factor for adjusting the pupil sampling scale. [default: None]

  • pupil_plane_size – Size in meters to use for the pupil plane array. In most cases, it’s a good idea to leave this as None, in which case GalSim will attempt to find a good value automatically. See also oversampling for adjusting the pupil size. [default: None]

  • scale_unit – Units to use for the sky coordinates when calculating lam/diam if these are supplied separately. Should be either a galsim.AngleUnit or a string that can be used to construct one (e.g., ‘arcsec’, ‘radians’, etc.). [default: galsim.arcsec]

  • fft_sign – The sign (+/-) to use in the exponent of the Fourier kernel when evaluating the Fourier optics PSF. As of version 2.3, GalSim uses a plus sign by default, which we believe to be consistent with, for example, how Zemax computes a Fourier optics PSF on DECam. Before version 2.3, the default was a negative sign. Input should be either the string ‘+’ or the string ‘-’. [default: ‘+’]

  • gsparams – An optional GSParams argument. [default: None]

withFlux(flux)[source]

Create a version of the current object with a different flux.

This function is equivalent to obj.withScaledFlux(flux / obj.flux).

It creates a new object that has the same profile as the original, but with the surface brightness at every location rescaled such that the total flux will be the given value. Note that if flux is an SED, the return value will be a ChromaticObject with specified SED.

Parameters:

flux – The new flux for the object.

Returns:

the object with the new flux