Chromatic Profiles

The ChromaticObject class and its various subclasses Define wavelength-dependent surface brightness profiles.

Implementation is done by constructing GSObject instances as functions of wavelength. The ChromaticObject.drawImage method then integrates over wavelength while also multiplying by a throughput function (a galsim.Bandpass instance).

These classes can be used to simulate galaxies with color gradients, observe a given galaxy through different filters, or implement wavelength-dependent point spread functions.

So far photon-shooting a ChromaticObject is not implemented, but there is ongoing work to include this option in GalSim, as the current FFT drawing method is very slow. So these are not yet particularly useful for large image simulations, especially ones including many faint sources.

class galsim.ChromaticObject(obj)[source]

Base class for defining wavelength-dependent objects.

This class primarily serves as the base class for chromatic subclasses. See the docstrings for subclasses for more details.

A ChromaticObject can be instantiated directly from an existing GSObject. In this case, the newly created ChromaticObject will act in nearly the same way as the original GSObject works, except that it has access to the ChromaticObject transformation methods described below (e.g., expand(), dilate(), shift(), withFlux(), …) These can all take functions as arguments to describe wavelength-dependent transformations. E.g.,:

>>> gsobj = galsim.Gaussian(fwhm=1)
>>> chrom_obj = galsim.ChromaticObject(gsobj).dilate(lambda wave: (wave/500.)**(-0.2))

In this and similar cases, the argument to the transformation method should be a python callable that accepts wavelength in nanometers and returns whatever type the transformation method normally accepts (so an int or float above).

One caveat to creating a ChromaticObject directly from a GSObject like this is that even though the source GSObject instance has flux units in photons/s/cm^2, the newly formed ChromaticObject will be interpreted as dimensionless, i.e., it will have a dimensionless SED (and have its .dimensionless attribute set to True). See below for more discussion about the dimensions of ChromaticObjects.

Another way to instantiate a ChromaticObject from a GSObject is to multiply by an SED. This can be useful to consistently generate the same galaxy observed through different filters, or, with ChromaticSum, to construct multi-component galaxies, each component with a different SED. For example, a bulge+disk galaxy could be constructed:

>>> bulge_sed = user_function_to_get_bulge_spectrum()
>>> disk_sed = user_function_to_get_disk_spectrum()
>>> bulge_mono = galsim.DeVaucouleurs(half_light_radius=1.0)
>>> disk_mono = galsim.Exponential(half_light_radius=2.0)
>>> bulge = bulge_mono * bulge_sed
>>> disk = disk_mono * disk_sed
>>> gal = bulge + disk

The sed instances above describe the flux density in photons/nm/cm^2/s of an object, possibly normalized with either the SED.withFlux or SED.withMagnitude methods (see their docstrings for details about these and other normalization options). Note that for dimensional consistency, in this case, the flux attribute of the multiplied GSObject is interpreted as being dimensionless instead of in its normal units of [photons/s/cm^2]. The photons/s/cm^2 units are (optionally) carried by the SED instead, or even left out entirely if the SED is dimensionless itself (see discussion on ChromaticObject dimensions below). The GSObject flux attribute does still contribute to the ChromaticObject normalization, though.

For example, the following are equivalent:

>>> chrom_obj = (sed * 3.0) * gsobj
>>> chrom_obj2 = sed * (gsobj * 3.0)

Subclasses that instantiate a ChromaticObject directly, such as ChromaticAtmosphere, also exist. Even in this case, however, the underlying implementation always eventually wraps one or more GSObject instances.

Dimensions:

ChromaticObjects can generally be sorted into two distinct types: those that represent galaxies or stars and have dimensions of [photons/wavelength-interval/area/time/solid-angle], and those that represent other types of wavelength dependence besides flux, like chromatic PSFs (these have dimensions of [1/solid-angle]). The former category of ChromaticObjects will have their .spectral attribute set to True, while the latter category of ChromaticObjects will have their .dimensionless attribute set to True. These two classes of ChromaticObjects have different restrictions associated with them. For example, only spectral ChromaticObjects can be drawn using drawImage, only ChromaticObjects of the same type can be added together, and at most one spectral ChromaticObject can be part of a ChromaticConvolution.

Multiplying a dimensionless ChromaticObject a spectral SED produces a spectral ChromaticObject (though note that the new object’s SED may not be equal to the SED being multiplied by since the original ChromaticObject may not have had unit normalization.)

Methods:

evaluateAtWavelength returns the monochromatic surface brightness profile (as a GSObject) at a given wavelength (in nanometers).

interpolate can be used for non-separable ChromaticObjects to expedite the image rendering process. See the docstring of that method for more details and discussion of when this is a useful tool (and the interplay between interpolation, object transformations, and convolutions).

Also, ChromaticObject has most of the same methods as GSObject with the following exceptions:

The GSObject access methods (e.g. GSObject.xValue, GSObject.maxk, etc.) are not available. Instead, you would need to evaluate the profile at a particular wavelength and access what you want from that.

The withFlux, withFluxDensity, and withMagnitude methods will return a new chromatic object with the appropriate spatially integrated flux, flux density, or magnitude.

The drawImage method draws the object as observed through a particular bandpass, so the arguments are somewhat different. See the docstring of drawImage for more details.

__mul__(flux_ratio)[source]

Scale the flux of the object by the given flux ratio, which may be an SED, a float, or a univariate callable function (of wavelength in nanometers) that returns a float.

The normalization of a ChromaticObject is tracked through its .sed attribute, which may have dimensions of either [photons/wavelength-interval/area/time/solid-angle] or [1/solid-angle].

If flux_ratio is a spectral SED (i.e., flux_ratio.spectral==True), then self.sed must be dimensionless for dimensional consistency. The returned object will have a spectral sed attribute. On the other hand, if flux_ratio is a dimensionless SED, float, or univariate callable function, then the returned object will have .spectral and .dimensionless matching self.spectral and self.dimensionless.

Parameters:

flux_ratio – The factor by which to scale the normalization of the object. flux_ratio may be a float, univariate callable function, in which case the argument should be wavelength in nanometers and return value the flux ratio for that wavelength, or an SED.

Returns:

a new object with scaled flux.

applyTo(photon_array, local_wcs=None, rng=None)[source]

Apply the chromatic profile as a convolution to an existing photon array.

This method allows instances of this class to duck type as a PhotonOp, so one can apply it in a photon_ops list.

Parameters:
  • photon_array – A PhotonArray to apply the operator to.

  • local_wcs – A LocalWCS instance defining the local WCS for the current photon bundle in case the operator needs this information. [default: None]

  • rng – A random number generator to use to effect the convolution. [default: None]

atRedshift(redshift)[source]

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

This will both adjust the SED to have the given redshift and set a redshift attribute with the given value.

Returns:

the object with the new redshift

calculateCentroid(bandpass)[source]

Determine the centroid of the wavelength-integrated surface brightness profile.

Parameters:

bandpass – The bandpass through which the observation is made.

Returns:

the centroid of the integrated surface brightness profile, as a PositionD.

calculateFlux(bandpass)[source]

Return the flux (photons/cm^2/s) of the ChromaticObject through a Bandpass bandpass.

Parameters:

bandpass – A Bandpass object representing a filter, or None to compute the bolometric flux. For the bolometric flux the integration limits will be set to (0, infinity) unless overridden by non-None SED attributes blue_limit or red_limit. Note that an SED defined using a LookupTable automatically has blue_limit and red_limit set.

Returns:

the flux through the bandpass.

calculateMagnitude(bandpass)[source]

Return the ChromaticObject magnitude through a Bandpass bandpass.

Note that this requires bandpass to have been assigned a zeropoint using Bandpass.withZeropoint.

Parameters:

bandpass – A Bandpass object representing a filter, or None to compute the bolometric magnitude. For the bolometric magnitude the integration limits will be set to (0, infinity) unless overridden by non-None SED attributes blue_limit or red_limit. Note that an SED defined using a LookupTable automatically has blue_limit and red_limit set.

Returns:

the bandpass magnitude.

dilate(scale)[source]

Dilate the linear size of the profile by the given (possibly wavelength-dependent) scale, while preserving flux.

e.g. half_light_radius <– half_light_radius * scale

See expand() and magnify() for versions that preserve surface brightness, and thus change the flux.

Parameters:

scale – The linear rescaling factor to apply. In addition, scale may be a callable function, in which case the argument should be wavelength in nanometers and the return value the scale for that wavelength.

Returns:

the dilated object.

property dimensionless

Boolean indicating if the ChromaticObject is dimensionless.

drawImage(bandpass, image=None, integrator='quadratic', **kwargs)[source]

Base implementation for drawing an image of a ChromaticObject.

Some subclasses may choose to override this for specific efficiency gains. For instance, most GalSim use cases will probably finish with a convolution, in which case ChromaticConvolution.drawImage will be used.

The task of drawImage() in a chromatic context is to integrate a chromatic surface brightness profile multiplied by the throughput of bandpass, over the wavelength interval indicated by bandpass.

Several integrators are available in galsim.integ to do this integration when using the first method (non-interpolated integration). By default, galsim.integ.SampleIntegrator will be used if either bandpass.wave_list or self.wave_list have len() > 0.

If lengths of both are zero, which may happen if both the bandpass throughput and the SED associated with self are analytic python functions for example, then galsim.integ.ContinuousIntegrator will be used instead. This latter case by default will evaluate the integrand at 250 equally-spaced wavelengths between bandpass.blue_limit and bandpass.red_limit.

By default, the above two integrators will use the rule galsim.integ.quadRule for integration. The midpoint rule for integration can be specified instead by passing an integrator that has been initialized with the rule set to galsim.integ.midptRule. When creating a ContinuousIntegrator, the number of samples N is also an argument. For example:

>>> integrator = galsim.integ.ContinuousIntegrator(rule=galsim.integ.midptRule, N=100)
>>> image = chromatic_obj.drawImage(bandpass, integrator=integrator)

Finally, this method uses a cache to avoid recomputing the integral over the product of the bandpass and object SED when possible (i.e., for separable profiles). Because the cache size is finite, users may find that it is more efficient when drawing many images to group images using the same SEDs and bandpasses together in order to hit the cache more often. The default cache size is 10, but may be resized using the ChromaticObject.resize_multiplier_cache method.

Parameters:
  • bandpass – A Bandpass object representing the filter against which to integrate.

  • image – Optionally, the Image to draw onto. (See GSObject.drawImage for details.) [default: None]

  • integrator – When doing the exact evaluation of the profile, this argument should be one of the image integrators from galsim.integ, or a string ‘trapezoidal’, ‘midpoint’, or ‘quadratic’, in which case the routine will use a SampleIntegrator or ContinuousIntegrator depending on whether or not the object has a wave_list. [default: ‘quadratic’, which will try to select an appropriate integrator using the quadratic integration rule automatically.]

  • **kwargs – For all other kwarg options, see GSObject.drawImage

Returns:

the drawn Image.

drawKImage(bandpass, image=None, integrator='quadratic', **kwargs)[source]

Base implementation for drawing the Fourier transform of a ChromaticObject.

The task of drawKImage() in a chromatic context is exactly analogous to the task of drawImage in a chromatic context: to integrate the sed * bandpass weighted Fourier profiles over wavelength.

See drawImage for details on integration options.

Parameters:
  • bandpass – A Bandpass object representing the filter against which to integrate.

  • image – If provided, this will be the complex Image onto which to draw the k-space image. If image is None, then an automatically-sized image will be created. If image is given, but its bounds are undefined, then it will be resized appropriately based on the profile’s size. [default: None]

  • integrator – When doing the exact evaluation of the profile, this argument should be one of the image integrators from galsim.integ, or a string ‘trapezoidal’, ‘midpoint’, or ‘quadratic’, in which case the routine will use a SampleIntegrator or ContinuousIntegrator depending on whether or not the object has a wave_list. [default: ‘quadratic’, which will try to select an appropriate integrator using the quadratic integration rule automatically.]

  • **kwargs – For all other kwarg options, see GSObject.drawKImage.

Returns:

a complex Image instance (created if necessary)

evaluateAtWavelength(wave)[source]

Evaluate this chromatic object at a particular wavelength.

Parameters:

wave – Wavelength in nanometers.

Returns:

the monochromatic object at the given wavelength.

expand(scale)[source]

Expand the linear size of the profile by the given (possibly wavelength-dependent) scale factor scale, while preserving surface brightness.

This doesn’t correspond to either of the normal operations one would typically want to do to a galaxy. The functions dilate() and magnify() are the more typical usage. But this function is conceptually simple. It rescales the linear dimension of the profile, while preserving surface brightness. As a result, the flux will necessarily change as well.

See dilate() for a version that applies a linear scale factor while preserving flux.

See magnify() for a version that applies a scale factor to the area while preserving surface brightness.

Parameters:

scale – The factor by which to scale the linear dimension of the object. In addition, scale may be a callable function, in which case the argument should be wavelength in nanometers and the return value the scale for that wavelength.

Returns:

the expanded object

property gsparams

The GSParams for this object.

interpolate(waves, **kwargs)[source]

Build interpolation images to (possibly) speed up subsequent drawImage calls.

This method is used as a pre-processing step that can expedite image rendering using objects that have to be built up as sums of GSObject instances with different parameters at each wavelength, by interpolating between images at each wavelength instead of making a more costly instantiation of the relevant GSObject at each value of wavelength at which the bandpass is defined.

This routine does a costly initialization process to build up a grid of Image instances to be used for the interpolation later on. However, the object can get reused with different bandpasses, so there should not be any need to make many versions of this object, and there is a significant savings each time it is drawn into an image.

As a general rule of thumb, chromatic objects that are separable do not benefit from this particular optimization, whereas those that involve making GSObject instances with wavelength-dependent keywords or transformations do benefit from it.

Note that the interpolation scheme is simple linear interpolation in wavelength, and no extrapolation beyond the originally-provided range of wavelengths is permitted. However, the overall flux at each wavelength will use the exact SED at that wavelength to give more accurate final flux values. You can disable this feature by setting use_exact_sed = False.

The speedup involved in using interpolation depends in part on the bandpass used for rendering (since that determines how many full profile evaluations are involved in rendering the image). However, for ChromaticAtmosphere with simple profiles like Kolmogorov, the speedup in some simple example cases is roughly a factor of three, whereas for more expensive to render profiles like the ChromaticOpticalPSF, the speedup is more typically a factor of 10-50.

Achromatic transformations can be applied either before or after setting up interpolation, with the best option depending on the application. For example, when rendering many times with the same achromatic transformation applied, it is typically advantageous to apply the transformation before setting up the interpolation. But there is no value in this when applying different achromatic transformation to each object. Chromatic transformations should be applied before setting up interpolation; attempts to render images of ChromaticObject instances with interpolation followed by a chromatic transformation will result in the interpolation being unset and the full calculation being done.

Because of the clever way that the ChromaticConvolution routine works, convolutions of separable chromatic objects with non-separable ones that use interpolation will still benefit from these optimizations. For example, a non-separable chromatic PSF that uses interpolation, when convolved with a sum of two separable galaxy components each with their own SED, will be able to take advantage of this optimization. In contrast, when convolving two non-separable profiles that already have interpolation set up, there is no way to take advantage of that interpolation optimization, so it will be ignored and the full calculation will be done. However, interpolation can be set up for the convolution of two non-separable profiles, after the convolution step. This could be beneficial for example when convolving a chromatic optical PSF and chromatic atmosphere, before convolving with multiple galaxy profiles.

For use cases requiring a high level of precision, we recommend a comparison between the interpolated and the more accurate calculation for at least one case, to ensure that the required precision has been reached.

The input parameter waves determines the input grid on which images are precomputed. It is difficult to give completely general guidance as to how many wavelengths to choose or how they should be spaced; some experimentation compared with the exact calculation is warranted for each particular application. The best choice of settings might depend on how strongly the parameters of the object depend on wavelength.

Parameters:
  • waves – The list, tuple, or NumPy array of wavelengths to be used when building up the grid of images for interpolation. The wavelengths should be given in nanometers, and they should span the full range of wavelengths covered by any bandpass to be used for drawing an Image (i.e., this class will not extrapolate beyond the given range of wavelengths). They can be spaced any way the user likes, not necessarily linearly, though interpolation will be linear in wavelength between the specified wavelengths.

  • oversample_fac – Factor by which to oversample the stored profiles compared to the default, which is to sample them at the Nyquist frequency for whichever wavelength has the highest Nyquist frequency. oversample_fac>1 results in higher accuracy but costlier pre-computations (more memory and time). [default: 1]

  • use_exact_sed – If true, then rescale the interpolated image for a given wavelength by the ratio of the exact SED at that wavelength to the linearly interpolated SED at that wavelength. Thus, the flux of the interpolated object should be correct, at the possible expense of other features. [default: True]

Returns:

the version of the Chromatic object that uses interpolation (This will be an InterpolatedChromaticObject instance.)

lens(g1, g2, mu)[source]

Apply a lensing shear and magnification to this object.

This ChromaticObject method applies a lensing (reduced) shear and magnification. The shear must be specified using the g1, g2 definition of shear (see Shear for details). This is the same definition as the outputs of the PowerSpectrum and NFWHalo classes, which compute shears according to some lensing power spectrum or lensing by an NFW dark matter halo. The magnification determines the rescaling factor for the object area and flux, preserving surface brightness.

While gravitational lensing is achromatic, we do allow the parameters g1, g2, and mu to be callable functions to be parallel to all the other transformations of chromatic objects. In this case, the functions should take the wavelength in nanometers as the argument, and the return values are the corresponding value at that wavelength.

Parameters:
  • g1 – First component of lensing (reduced) shear to apply to the object.

  • g2 – Second component of lensing (reduced) shear to apply to the object.

  • mu – Lensing magnification to apply to the object. This is the factor by which the solid angle subtended by the object is magnified, preserving surface brightness.

Returns:

the lensed object.

magnify(mu)[source]

Apply a lensing magnification, scaling the area and flux by mu at fixed surface brightness.

This process applies a lensing magnification mu, which scales the linear dimensions of the image by the factor sqrt(mu), i.e., half_light_radius <– half_light_radius * sqrt(mu) while increasing the flux by a factor of mu. Thus, magnify() preserves surface brightness.

See dilate() for a version that applies a linear scale factor while preserving flux.

Parameters:

mu – The lensing magnification to apply. In addition, mu may be a callable function, in which case the argument should be wavelength in nanometers and the return value the magnification for that wavelength.

Returns:

the magnified object.

property redshift

The redshift of the object.

static resize_multiplier_cache(maxsize)[source]

Resize the cache (default size=10) containing the integral over the product of an SED and a Bandpass, which is used by ChromaticObject.drawImage.

Parameters:

maxsize – The new number of products to cache.

rotate(theta)[source]

Rotate this object by an Angle theta.

Parameters:

theta – Rotation angle (Angle object, +ve anticlockwise). In addition, theta may be a callable function, in which case the argument should be wavelength in nanometers and the return value the rotation angle for that wavelength, returned as a galsim.Angle instance.

Returns:

the rotated object.

shear(*args, **kwargs)[source]

Apply an area-preserving shear to this object, where arguments are either a Shear, or arguments that will be used to initialize one.

For more details about the allowed keyword arguments, see the Shear docstring.

The shear() method precisely preserves the area. To include a lensing distortion with the appropriate change in area, either use shear() with magnify(), or use lens(), which combines both operations.

Note that, while gravitational shear is monochromatic, the shear method may be used for many other use cases including some which may be wavelength-dependent, such as intrinsic galaxy shape, telescope dilation, atmospheric PSF shape, etc. Thus, the shear argument is allowed to be a function of wavelength like other transformations.

Parameters:

shear – The shear to be applied. Or, as described above, you may instead supply parameters to construct a Shear directly. eg. obj.shear(g1=g1,g2=g2). In addition, the shear parameter may be a callable function, in which case the argument should be wavelength in nanometers and the return value the shear for that wavelength, returned as a galsim.Shear instance.

Returns:

the sheared object.

shift(*args, **kwargs)[source]

Apply a (possibly wavelength-dependent) (dx, dy) shift to this chromatic object.

For a wavelength-independent shift, you may supply dx,dy as either two arguments, as a tuple, or as a PositionD or PositionI object.

For a wavelength-dependent shift, you may supply two functions of wavelength in nanometers which will be interpreted as dx(wave) and dy(wave), or a single function of wavelength in nanometers that returns either a 2-tuple, PositionD, or PositionI.

Parameters:
  • dx – Horizontal shift to apply (float or function).

  • dy – Vertical shift to apply (float or function).

Returns:

the shifted object.

property spectral

Boolean indicating if the ChromaticObject has units compatible with a spectral density.

transform(dudx, dudy, dvdx, dvdy)[source]

Apply a transformation to this object defined by an arbitrary Jacobian matrix.

This works the same as GSObject.transform, so see that method’s docstring for more details.

As with the other more specific chromatic trasnformations, dudx, dudy, dvdx, and dvdy may be callable functions, in which case the argument should be wavelength in nanometers and the return value the appropriate value for that wavelength.

Parameters:
  • dudx – du/dx, where (x,y) are the current coords, and (u,v) are the new coords.

  • dudy – du/dy, where (x,y) are the current coords, and (u,v) are the new coords.

  • dvdx – dv/dx, where (x,y) are the current coords, and (u,v) are the new coords.

  • dvdy – dv/dy, where (x,y) are the current coords, and (u,v) are the new coords.

Returns:

the transformed object.

withFlux(target_flux, bandpass)[source]

Return a new ChromaticObject with flux through the Bandpass bandpass set to target_flux.

Parameters:
  • target_flux – The desired flux normalization of the ChromaticObject.

  • bandpass – A Bandpass object defining a filter bandpass.

Returns:

the new normalized ChromaticObject.

withFluxDensity(target_flux_density, wavelength)[source]

Return a new ChromaticObject with flux density set to target_flux_density at wavelength wavelength.

Parameters:
  • target_flux_density – The target normalization in photons/nm/cm^2/s.

  • wavelength – The wavelength, in nm, at which the flux density will be set.

Returns:

the new normalized SED.

withGSParams(gsparams=None, **kwargs)[source]

Create a version of the current object with the given gsparams

Note: if this object wraps other objects (e.g. Convolution, Sum, Transformation, etc.) those component objects will also have their gsparams updated to the new value.

withMagnitude(target_magnitude, bandpass)[source]

Return a new ChromaticObject with magnitude through bandpass set to target_magnitude. Note that this requires bandpass to have been assigned a zeropoint using Bandpass.withZeropoint.

Parameters:
  • target_magnitude – The desired magnitude of the ChromaticObject.

  • bandpass – A Bandpass object defining a filter bandpass.

Returns:

the new normalized ChromaticObject.

withScaledFlux(flux_ratio)[source]

Multiply the flux of the object by flux_ratio

Parameters:

flux_ratio – The factor by which to scale the normalization of the object. flux_ratio may be a float, univariate callable function, in which case the argument should be wavelength in nanometers and return value the flux ratio for that wavelength, or an SED.

Returns:

a new object with scaled flux.

class galsim.ChromaticAtmosphere(base_obj, base_wavelength, scale_unit=None, **kwargs)[source]

Bases: ChromaticObject

A ChromaticObject implementing two atmospheric chromatic effects: differential chromatic refraction (DCR) and wavelength-dependent seeing.

Due to DCR, blue photons land closer to the zenith than red photons. Kolmogorov turbulence also predicts that blue photons get spread out more by the atmosphere than red photons, specifically FWHM is proportional to wavelength^(-0.2). Both of these effects can be implemented by wavelength-dependent shifts and dilations.

Since DCR depends on the zenith angle and the parallactic angle (which is the position angle of the zenith measured from North through East) of the object being drawn, these must be specified via keywords. There are four ways to specify these values:

  1. explicitly provide zenith_angle as a keyword of type Angle, and parallactic_angle will be assumed to be 0 by default.

  2. explicitly provide both zenith_angle and parallactic_angle as keywords of type Angle.

  3. provide the coordinates of the object obj_coord and the coordinates of the zenith zenith_coord as keywords of type CelestialCoord.

  4. provide the coordinates of the object obj_coord as a CelestialCoord, the hour angle of the object HA as an Angle, and the latitude of the observer latitude as an Angle.

DCR also depends on temperature, pressure and water vapor pressure of the atmosphere. The default values for these are expected to be appropriate for LSST at Cerro Pachon, Chile, but they are broadly reasonable for most observatories.

Note that a ChromaticAtmosphere by itself is NOT the correct thing to use to draw an image of a star. Stars (and galaxies too, of course) have an SED that is not flat. To draw a real star, you should either multiply the ChromaticAtmosphere object by an SED, or convolve it with a point source multiplied by an SED:

>>> psf = galsim.ChromaticAtmosphere(...)
>>> star = galsim.DeltaFunction() * psf_sed
>>> final_star = galsim.Convolve( [psf, star] )
>>> final_star.drawImage(bandpass = bp, ...)
Parameters:
  • base_obj – Fiducial PSF, equal to the monochromatic PSF at base_wavelength

  • base_wavelength – Wavelength represented by the fiducial PSF, in nanometers.

  • scale_unit – Units used by base_obj for its linear dimensions. [default: galsim.arcsec]

  • alpha – Power law index for wavelength-dependent seeing. [default: -0.2, the prediction for Kolmogorov turbulence]

  • zenith_angleAngle from object to zenith [default: 0]

  • parallactic_angle – Parallactic angle, i.e. the position angle of the zenith, measured from North through East. [default: 0]

  • obj_coord – Celestial coordinates of the object being drawn as a CelestialCoord. [default: None]

  • zenith_coord – Celestial coordinates of the zenith as a CelestialCoord. [default: None]

  • HA – Hour angle of the object as an Angle. [default: None]

  • latitude – Latitude of the observer as an Angle. [default: None]

  • pressure – Air pressure in kiloPascals. [default: 69.328 kPa]

  • temperature – Temperature in Kelvins. [default: 293.15 K]

  • H2O_pressure – Water vapor pressure in kiloPascals. [default: 1.067 kPa]

build_obj()[source]

Build a ChromaticTransformation object for this ChromaticAtmosphere.

We don’t do this right away to help make ChromaticAtmosphere objects be picklable. Building this is quite fast, so we do it on the fly in evaluateAtWavelength and ChromaticObject.drawImage.

evaluateAtWavelength(wave)[source]

Evaluate this chromatic object at a particular wavelength.

Parameters:

wave – Wavelength in nanometers.

Returns:

the monochromatic object at the given wavelength.

property gsparams

The GSParams for this object.

withGSParams(gsparams=None, **kwargs)[source]

Create a version of the current object with the given gsparams

Note: if this object wraps other objects (e.g. Convolution, Sum, Transformation, etc.) those component objects will also have their gsparams updated to the new value.

class galsim.ChromaticOpticalPSF(lam, diam=None, lam_over_diam=None, aberrations=None, scale_unit=None, gsparams=None, **kwargs)[source]

Bases: ChromaticObject

A subclass of ChromaticObject meant to represent chromatic optical PSFs.

Chromaticity plays two roles in optical PSFs. First, it determines the diffraction limit, via the wavelength/diameter factor. Second, aberrations such as defocus, coma, etc. are typically defined in physical distances, but their impact on the PSF depends on their size in units of wavelength. Other aspects of the optical PSF do not require explicit specification of their chromaticity, e.g., once the obscuration and struts are specified in units of the aperture diameter, their chromatic dependence gets taken care of automatically. Note that the ChromaticOpticalPSF implicitly defines diffraction limits in units of scale_units, which by default are arcsec, but can in principle be set to any of our GalSim angle units.

When using interpolation to speed up image rendering (see the ChromaticObject.interpolate method for details), the ideal number of wavelengths to use across a given bandpass depends on the application and accuracy requirements. In general it will be necessary to do a test in comparison with a more exact calculation to ensure convergence. However, a typical calculation might use ~10-15 samples across a typical optical bandpass, with oversample_fac in the range 1.5-2; for moderate accuracy, ~5 samples across the bandpass and oversample_fac=1 may suffice. All of these statements assume that aberrations are not very large (typically <~0.25 waves, which is commonly satisfied by space telescopes); if they are larger than that, then more stringent settings are required.

Note that a ChromaticOpticalPSF by itself is NOT the correct thing to use to draw an image of a star. Stars (and galaxies too, of course) have an SED that is not flat. To draw a real star, you should either multiply the ChromaticOpticalPSF object by an SED, or convolve it with a point source multiplied by an SED:

>>> psf = galsim.ChromaticOpticalPSF(...)
>>> star = galsim.DeltaFunction() * psf_sed
>>> final_star = galsim.Convolve( [psf, star] )
>>> final_star.drawImage(bandpass = bp, ...)

Note

When geometric_shooting is False (the default), the photon shooting implementation is only approximately correct with respect to the wavelength dependence. It is also not particularly fast, since it generates three optical screens to span the wavelength range and shoots from these (with a subsequent adjustment to improve the accuracy of this approximation). We expect that most users who want to use photon shooting in conjunction with this class will prefer to make an InterpolatedChromaticObject (by calling psf.interpolate(...)), especially if it is a good approximation to use the same optical PSF for a whole exposure or CCD image, so the setup time for the interpolation is able to be amortized for many objects.

Parameters:
  • lam – Fiducial wavelength for which diffraction limit and aberrations are initially defined, in nanometers.

  • diam – Telescope diameter in meters. Either diam or lam_over_diam must be specified.

  • lam_over_diam – Ratio of (fiducial wavelength) / telescope diameter in units of scale_unit. Either diam or lam_over_diam must be specified.

  • aberrations – An array of aberrations, in units of fiducial wavelength lam. The size and format of this array is described in the OpticalPSF docstring.

  • scale_unit – Units used to define the diffraction limit and draw images. [default: galsim.arcsec]

  • gsparams – An optional GSParams argument. See the docstring for GSParams for details. [default: None]

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

  • **kwargs – Any other keyword arguments to be passed to OpticalPSF, for example, related to struts, obscuration, oversampling, etc. See OpticalPSF docstring for a complete list of options.

evaluateAtWavelength(wave)[source]

Method to directly instantiate a monochromatic instance of this object.

Parameters:

wave – Wavelength in nanometers.

property gsparams

The GSParams for this object.

withGSParams(gsparams=None, **kwargs)[source]

Create a version of the current object with the given gsparams

Note: if this object wraps other objects (e.g. Convolution, Sum, Transformation, etc.) those component objects will also have their gsparams updated to the new value.

class galsim.ChromaticAiry(lam, diam=None, lam_over_diam=None, scale_unit=None, gsparams=None, **kwargs)[source]

Bases: ChromaticObject

A subclass of ChromaticObject meant to represent chromatic Airy profiles.

For more information about the basics of Airy profiles, please see galsim.Airy.

This class is a chromatic representation of Airy profiles, including the wavelength-dependent diffraction limit. One can also get this functionality using the ChromaticOpticalPSF class, but that class includes additional complications beyond a simple Airy profile, and thus has a more complicated internal representation. For users who only want a (possibly obscured) Airy profile, the ChromaticAiry class is likely to be a less computationally expensive and more accurate option.

Parameters:
  • lam – Fiducial wavelength for which diffraction limit is initially defined, in nanometers.

  • diam – Telescope diameter in meters. Either diam or lam_over_diam must be specified.

  • lam_over_diam – Ratio of (fiducial wavelength) / telescope diameter in units of scale_unit. Either diam or lam_over_diam must be specified.

  • scale_unit – Units used to define the diffraction limit and draw images. [default: galsim.arcsec]

  • gsparams – An optional GSParams argument. See the docstring for GSParams for details. [default: None]

  • **kwargs – Any other keyword arguments to be passed to Airy: either flux, or gsparams. See galsim.Airy docstring for a complete description of these options.

evaluateAtWavelength(wave)[source]

Method to directly instantiate a monochromatic instance of this object.

Parameters:

wave – Wavelength in nanometers.

property gsparams

The GSParams for this object.

withGSParams(gsparams=None, **kwargs)[source]

Create a version of the current object with the given gsparams

Note: if this object wraps other objects (e.g. Convolution, Sum, Transformation, etc.) those component objects will also have their gsparams updated to the new value.

class galsim.ChromaticRealGalaxy(real_galaxy_catalogs, index=None, id=None, random=False, rng=None, gsparams=None, logger=None, **kwargs)[source]

Bases: ChromaticSum

A class describing real galaxies over multiple wavelengths, using some multi-band training dataset. The underlying implementation models multi-band images of individual galaxies as chromatic PSF convolutions (and integrations over wavelength) with a sum of profiles separable into spatial and spectral components. The spectral components are specified by the user, and the spatial components are determined one Fourier mode at a time by the class. This decomposition can be thought of as a constrained chromatic deconvolution of the multi-band images by the associated PSFs, similar in spirit to RealGalaxy.

Because ChromaticRealGalaxy involves an InterpolatedKImage, method = 'phot' is unavailable for the ChromaticObject.drawImage function.

Fundamentally, the required inputs for this class are:

  1. a series of high resolution input Image instances of a single galaxy in different bands,

  2. a list of Bandpass corresponding to those images,

  3. the PSFs of those images as either GSObject or ChromaticObject instances, and

  4. the noise properties of the input images as BaseCorrelatedNoise instances.

If you want to specify these inputs directly, that is possible via the makeFromImages factory method of this class:

>>> crg = galsim.ChromaticRealGalaxy.makeFromImages(imgs, bands, PSFs, xis, ...)

Alternatively, you may create a ChromaticRealGalaxy via a list of RealGalaxyCatalog that correspond to a set of galaxies observed in different bands:

>>> crg = galsim.ChromaticRealGalaxy(real_galaxy_catalogs, index=0, ...)

The above will use the 1st object in the catalogs, which should be the same galaxy, just observed in different bands. Note that there are multiple keywords for choosing a galaxy from a catalog; exactly one must be set. In the future we may add more such options, e.g., to choose at random but accounting for the non-constant weight factors (probabilities for objects to make it into the training sample).

The flux normalization of the returned object will by default match the original data, scaled to correspond to a 1 second HST exposure (though see the area_norm parameter). If you want a flux appropriate for a longer exposure or telescope with different collecting area, you can use the ChromaticObject.withScaledFlux method on the returned object, or use the exptime and area keywords to ChromaticObject.drawImage.

Note that while you can also use ChromaticObject.withFlux, ChromaticObject.withMagnitude, and ChromaticObject.withFluxDensity to set the absolute normalization, these methods technically adjust the flux of the entire postage stamp image (including noise!) and not necessarily the flux of the galaxy itself. (These two fluxes will be strongly correlated for high signal-to-noise ratio galaxies, but may be considerably different at low signal-to-noise ratio.)

Note that ChromaticRealGalaxy objects use arcsec for the units of their linear dimension. If you are using a different unit for other things (the PSF, WCS, etc.), then you should dilate the resulting object with gal.dilate(galsim.arcsec / scale_unit).

Noise from the original images is propagated by this class, though certain restrictions apply to when and how that noise is made available. The propagated noise depends on which Bandpass the ChromaticRealGalaxy is being imaged through, so the noise is only available after the ChromaticObject.drawImage method has been called. Also, since ChromaticRealGalaxy will only produce reasonable images when convolved with a (suitably wide) PSF, the noise attribute is attached to the ChromaticConvolution (or ChromaticTransformation of the ChromaticConvolution) which holds as one of its convolutants the ChromaticRealGalaxy.:

>>> crg = galsim.ChromaticRealGalaxy(...)
>>> psf = ...
>>> obj = galsim.Convolve(crg, psf)
>>> bandpass = galsim.Bandpass(...)
>>> assert not hasattr(obj, 'noise')
>>> image = obj.drawImage(bandpass)
>>> assert hasattr(obj, 'noise')
>>> noise1 = obj.noise

Note that the noise attribute is only associated with the most recently used bandpass. If you draw another image of the same object using a different bandpass, the noise object will be replaced.:

>>> bandpass2 = galsim.Bandpass(...)
>>> image2 = obj.drawImage(bandpass2)
>>> assert noise1 != obj.noise
Parameters:
  • real_galaxy_catalogs – A list of RealGalaxyCatalog objects from which to create ChromaticRealGalaxy objects. Each catalog should represent the same set of galaxies, and in the same order, just imaged through different filters.

  • index – Index of the desired galaxy in the catalog. [One of index, id, or random is required.]

  • id – Object ID for the desired galaxy in the catalog. [One of index, id, or random is required.]

  • random – If True, then just select a completely random galaxy from the catalog. [One of index, id, or random is required.]

  • rng – A random number generator to use for selecting a random galaxy (may be any kind of BaseDeviate or None) and to use in generating any noise field when padding.

  • SEDs – An optional list of SED instances to use when representing real galaxies as sums of separable profiles. By default, it will use len(real_galaxy_catalogs) SEDs that are polynomials in wavelength. Note that if given, len(SEDs) must equal len(real_galaxy_catalogs). [default: None]

  • k_interpolant – Either an Interpolant instance or a string indicating which k-space interpolant should be used. Options are ‘nearest’, ‘sinc’, ‘linear’, ‘cubic’, ‘quintic’, or ‘lanczosN’ where N should be the integer order to use. We strongly recommend leaving this parameter at its default value; see text above for details. [default: galsim.Quintic()]

  • maxk – Optional maxk argument. If you know you will be convolving the resulting ChromaticRealGalaxy with a “fat” PSF in a subsequent step, then it can be more efficient to limit the range of Fourier modes used when solving for the sum of separable profiles below. [default: None]

  • pad_factor – Factor by which to internally oversample the Fourier-space images that represent the ChromaticRealGalaxy (equivalent to zero-padding the real-space profiles). We strongly recommend leaving this parameter at its default value; see text in Realgalaxy docstring for details. [default: 4]

  • noise_pad_size – If provided, the image will be padded out to this size (in arcsec) with the noise specified in the real galaxy catalog. This is important if you are planning to whiten the resulting image. You should make sure that the padded image is larger than the postage stamp onto which you are drawing this object. [default: None]

  • area_norm – Area in cm^2 by which to normalize the flux of the returned object. When area_norm=1 (the default), using exptime=1 and area=1 arguments in ChromaticObject.drawImage (also the default) will simulate an image with the appropriate number of counts for a 1 second exposure with the original telescope/camera (e.g., with HST when using the COSMOS catalog). If you would rather explicitly specify the collecting area of the telescope when using ChromaticObject.drawImage with a ChromaticRealGalaxy, then you should set area_norm equal to the collecting area of the source catalog telescope when creating the ChromaticRealGalaxy (e.g., area_norm=45238.93416 for HST). [default: 1]

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

  • logger – A logger object for output of progress statements if the user wants them. [default: None]

classmethod makeFromImages(images, bands, PSFs, xis, **kwargs)[source]

Create a ChromaticRealGalaxy directly from images, bandpasses, PSFs, and noise descriptions. See the ChromaticRealGalaxy docstring for more information.

Parameters:
  • images – An iterable of high resolution Image instances of a galaxy through different bandpasses.

  • bands – An iterable of Bandpass objects corresponding to the input images.

  • PSFs – Either an iterable of GSObject or ChromaticObject indicating the PSFs of the different input images, or potentially a single GSObject or ChromaticObject that will be used as the PSF for all images.

  • xis – An iterable of BaseCorrelatedNoise objects characterizing the noise in the input images.

  • SEDs – An optional list of SED instances to use when representing real galaxies as sums of separable profiles. By default, it will use len(images) SEDs that are polynomials in wavelength. Note that if given, len(SEDs) must equal len(images). [default: None]

  • k_interpolant – Either an Interpolant instance or a string indicating which k-space interpolant should be used. Options are ‘nearest’, ‘sinc’, ‘linear’, ‘cubic’, ‘quintic’, or ‘lanczosN’ where N should be the integer order to use. We strongly recommend leaving this parameter at its default value; see text above for details. [default: galsim.Quintic()]

  • maxk – Optional maxk argument. If you know you will be convolving the resulting ChromaticRealGalaxy with a “fat” PSF in a subsequent step, then it can be more efficient to limit the range of Fourier modes used when solving for the sum of separable profiles below. [default: None]

  • pad_factor – Factor by which to internally oversample the Fourier-space images that represent the ChromaticRealGalaxy (equivalent to zero-padding the real-space profiles). We strongly recommend leaving this parameter at its default value; see text in Realgalaxy docstring for details. [default: 4]

  • noise_pad_size – If provided, the image will be padded out to this size (in arcsec) with the noise specified in the real galaxy catalog. This is important if you are planning to whiten the resulting image. You should make sure that the padded image is larger than the postage stamp onto which you are drawing this object. [default: None]

  • area_norm – Area in cm^2 by which to normalize the flux of the returned object. When area_norm=1 (the default), using exptime=1 and area=1 arguments in ChromaticObject.drawImage (also the default) will simulate an image with the appropriate number of counts for a 1 second exposure with the original telescope/camera (e.g., with HST when using the COSMOS catalog). If you would rather explicitly specify the collecting area of the telescope when using ChromaticObject.drawImage with a ChromaticRealGalaxy, then you should set area_norm equal to the collecting area of the source catalog telescope when creating the ChromaticRealGalaxy (e.g., area_norm=45238.93416 for HST). [default: 1]

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

  • logger – A logger object for output of progress statements if the user wants them. [default: None]

class galsim.InterpolatedChromaticObject(original, waves, oversample_fac=1.0, use_exact_sed=True, use_exact_SED=None)[source]

Bases: ChromaticObject

A ChromaticObject that uses interpolation of predrawn images to speed up subsequent rendering.

This class wraps another ChromaticObject, which is stored in the attribute deinterpolated. Any ChromaticObject can be used, although the interpolation procedure is most effective for non-separable objects, which can sometimes be very slow to render.

Normally, you would not create an InterpolatedChromaticObject directly. It is the return type from ChromaticObject.interpolate. See the description of that function for more details.

Parameters:
  • original – The ChromaticObject to be interpolated.

  • waves – The list, tuple, or NumPy array of wavelengths to be used when building up the grid of images for interpolation. The wavelengths should be given in nanometers, and they should span the full range of wavelengths covered by any bandpass to be used for drawing an Image (i.e., this class will not extrapolate beyond the given range of wavelengths). They can be spaced any way the user likes, not necessarily linearly, though interpolation will be linear in wavelength between the specified wavelengths.

  • oversample_fac – Factor by which to oversample the stored profiles compared to the default, which is to sample them at the Nyquist frequency for whichever wavelength has the highest Nyquist frequency. oversample_fac>1 results in higher accuracy but costlier pre-computations (more memory and time). [default: 1]

  • use_exact_sed – If true, then rescale the interpolated image for a given wavelength by the ratio of the exact SED at that wavelength to the linearly interpolated SED at that wavelength. Thus, the flux of the interpolated object should be correct, at the possible expense of other features. [default: True]

drawImage(bandpass, image=None, integrator='quadratic', **kwargs)[source]

Draw an image as seen through a particular bandpass using the stored interpolated images at the specified wavelengths.

This integration will take place using interpolation between stored images that were setup when the object was constructed. (See interpolate() for more details.)

Parameters:
  • bandpass – A Bandpass object representing the filter against which to integrate.

  • image – Optionally, the Image to draw onto. (See GSObject.drawImage for details.) [default: None]

  • integrator – The integration algorithm to use, given as a string. Either ‘midpoint’, ‘trapezoidal’, or ‘quadratic’ is allowed. [default: ‘quadratic’]

  • **kwargs – For all other kwarg options, see GSObject.drawImage.

Returns:

the drawn Image.

evaluateAtWavelength(wave)[source]

Evaluate this ChromaticObject at a particular wavelength using interpolation.

Parameters:

wave – Wavelength in nanometers.

Returns:

the monochromatic object at the given wavelength, as a GSObject.

property gsparams

The GSParams for this object.

withGSParams(gsparams=None, **kwargs)[source]

Create a version of the current object with the given gsparams

Note: if this object wraps other objects (e.g. Convolution, Sum, Transformation, etc.) those component objects will also have their gsparams updated to the new value.

class galsim.ChromaticSum(*args, **kwargs)[source]

Bases: ChromaticObject

A sum of several ChromaticObject and/or GSObject instances.

Any GSObject in the sum is assumed to have a flat SED with spectral density of 1 photon/s/cm**2/nm.

This is the type returned from galsim.Add if any of the objects are a ChromaticObject.

Typically, you do not need to construct a ChromaticSum object explicitly. Normally, you would just use the + operator, which returns a ChromaticSum when used with chromatic objects:

>>> bulge = galsim.Sersic(n=3, half_light_radius=0.8) * bulge_sed
>>> disk = galsim.Exponential(half_light_radius=1.4) * disk_sed
>>> gal = bulge + disk

You can also use the Add factory function, which returns a ChromaticSum object if any of the individual objects are chromatic:

>>> gal = galsim.Add([bulge,disk])
Parameters:
  • args – Unnamed args should be a list of objects to add.

  • gsparams – An optional GSParams argument. See the docstring for GSParams for details. [default: None]

  • propagate_gsparams – Whether to propagate gsparams to each of the components. This is normally a good idea, but there may be use cases where one would not want to do this. [default: True]

atRedshift(redshift)[source]

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

This will both adjust the SED to have the given redshift and set a redshift attribute with the given value.

Returns:

the object with the new redshift

drawImage(bandpass, image=None, integrator='quadratic', **kwargs)[source]

Slightly optimized draw method for ChromaticSum instances.

Draws each summand individually and add resulting images together. This might waste time if two or more summands are separable and have the same SED, and another summand with a different SED is also added, in which case the summands should be added together first and the resulting Sum object can then be chromaticized. In general, however, drawing individual sums independently can help with speed by identifying chromatic profiles that are separable into spectral and spatial factors.

Parameters:
  • bandpass – A Bandpass object representing the filter against which to integrate.

  • image – Optionally, the Image to draw onto. (See GSObject.drawImage for details.) [default: None]

  • integrator – When doing the exact evaluation of the profile, this argument should be one of the image integrators from galsim.integ, or a string ‘trapezoidal’, ‘midpoint’, ‘quadratic’, in which case the routine will use a SampleIntegrator or ContinuousIntegrator depending on whether or not the object has a wave_list. [default: ‘quadratic’, which will try to select an appropriate integrator using the quadratic integration rule automatically.]

  • **kwargs – For all other kwarg options, see GSObject.drawImage.

Returns:

the drawn Image.

evaluateAtWavelength(wave)[source]

Evaluate this chromatic object at a particular wavelength wave.

Parameters:

wave – Wavelength in nanometers.

Returns:

the monochromatic object at the given wavelength.

property gsparams

The GSParams for this object.

property obj_list

The list of objects being added.

withGSParams(gsparams=None, **kwargs)[source]

Create a version of the current object with the given gsparams

Note: if this object wraps other objects (e.g. Convolution, Sum, Transformation, etc.) those component objects will also have their gsparams updated to the new value.

withScaledFlux(flux_ratio)[source]

Multiply the flux of the object by flux_ratio

Parameters:

flux_ratio – The factor by which to scale the flux.

Returns:

the object with the new flux.

class galsim.ChromaticConvolution(*args, **kwargs)[source]

Bases: ChromaticObject

A convolution of several ChromaticObject and/or GSObject instances.

Any GSObject in the convolution is assumed to have a flat SED with spectral density of 1 photon/s/cm**2/nm.

This is the type returned from galsim.Convolve if any of the objects is a ChromaticObject.

The normal way to use this class is to use the Convolve factory function:

>>> gal = galsim.Sersic(n, half_light_radius) * galsim.SED(sed_file, 'nm', 'flambda')
>>> psf = galsim.ChromaticAtmosphere(...)
>>> final = galsim.Convolve([gal, psf])

The objects to be convolved may be provided either as multiple unnamed arguments (e.g. Convolve(psf, gal, pix)) or as a list (e.g. Convolve([psf, gal, pix])). Any number of objects may be provided using either syntax. (Well, the list has to include at least 1 item.)

Parameters:
  • args – Unnamed args should be a list of objects to convolve.

  • real_space – Whether to use real space convolution. [default: None, which means to automatically decide this according to whether the objects have hard edges.]

  • gsparams – An optional GSParams argument. See the docstring for GSParams for details. [default: None]

  • propagate_gsparams – Whether to propagate gsparams to each of the components. This is normally a good idea, but there may be use cases where one would not want to do this. [default: True]

drawImage(bandpass, image=None, integrator='quadratic', iimult=None, **kwargs)[source]

Optimized draw method for the ChromaticConvolution class.

Works by finding sums of profiles which include separable portions, which can then be integrated before doing any convolutions, which are pushed to the end.

This method uses a cache to avoid recomputing ‘effective’ profiles, which are the wavelength-integrated products of inseparable profiles, the spectral components of separable profiles, and the bandpass. Because the cache size is finite, users may find that it is more efficient when drawing many images to group images using the same SEDs, bandpasses, and inseparable profiles (generally PSFs) together in order to hit the cache more often. The default cache size is 10, but may be resized using the ChromaticConvolution.resize_effective_prof_cache method.

Parameters:
  • bandpass – A Bandpass object representing the filter against which to integrate.

  • image – Optionally, the Image to draw onto. (See GSObject.drawImage for details.) [default: None]

  • integrator – When doing the exact evaluation of the profile, this argument should be one of the image integrators from galsim.integ, or a string ‘trapezoidal’, ‘midpoint’, or ‘quadratic’, in which case the routine will use a SampleIntegrator or ContinuousIntegrator depending on whether or not the object has a wave_list. [default: ‘quadratic’, which will try to select an appropriate integrator using the quadratic integration rule automatically.]

  • iimult – Oversample any intermediate InterpolatedImage created to hold effective profiles by this amount. [default: None]

  • **kwargs – For all other kwarg options, see GSObject.drawImage.

Returns:

the drawn Image.

evaluateAtWavelength(wave)[source]

Evaluate this chromatic object at a particular wavelength wave.

Parameters:

wave – Wavelength in nanometers.

Returns:

the monochromatic object at the given wavelength.

property gsparams

The GSParams for this object.

property obj_list

The list of objects being convolved.

static resize_effective_prof_cache(maxsize)[source]

Resize the cache containing effective profiles.

These are wavelength-integrated products of separable profile SEDs, inseparable profiles, and Bandpasses) used by ChromaticConvolution.drawImage.

Parameters:

maxsize – The new number of effective profiles to cache.

withGSParams(gsparams=None, **kwargs)[source]

Create a version of the current object with the given gsparams

Note: if this object wraps other objects (e.g. Convolution, Sum, Transformation, etc.) those component objects will also have their gsparams updated to the new value.

class galsim.ChromaticDeconvolution(obj, gsparams=None, propagate_gsparams=True)[source]

Bases: ChromaticObject

A class for deconvolving a ChromaticObject.

The ChromaticDeconvolution class represents a wavelength-dependent deconvolution kernel.

You may also specify a gsparams argument. See the docstring for GSParams for more information about this option. Note: if gsparams is unspecified (or None), then the ChromaticDeconvolution instance inherits the same GSParams as the object being deconvolved.

This is the type returned from galsim.Deconvolve if the argument is a ChromaticObject. This is the normal way to construct this class.

Parameters:
  • obj – The object to deconvolve.

  • gsparams – An optional GSParams argument. See the docstring for GSParams for details. [default: None]

  • propagate_gsparams – Whether to propagate gsparams to each of the components. This is normally a good idea, but there may be use cases where one would not want to do this. [default: True]

evaluateAtWavelength(wave)[source]

Evaluate this chromatic object at a particular wavelength wave.

Parameters:

wave – Wavelength in nanometers.

Returns:

the monochromatic object at the given wavelength.

property gsparams

The GSParams for this object.

withGSParams(gsparams=None, **kwargs)[source]

Create a version of the current object with the given gsparams

Note: if this object wraps other objects (e.g. Convolution, Sum, Transformation, etc.) those component objects will also have their gsparams updated to the new value.

class galsim.ChromaticAutoConvolution(obj, real_space=None, gsparams=None, propagate_gsparams=True)[source]

Bases: ChromaticObject

A special class for convolving a ChromaticObject with itself.

It is equivalent in functionality to galsim.Convolve([obj,obj]), but takes advantage of the fact that the two profiles are the same for some efficiency gains.

This is the type returned from galsim.AutoConvolve if the argument is a ChromaticObject. This is the normal way to construct this class.

Parameters:
  • obj – The object to be convolved with itself.

  • real_space – Whether to use real space convolution. [default: None, which means to automatically decide this according to whether the objects have hard edges.]

  • gsparams – An optional GSParams argument. See the docstring for GSParams for details. [default: None]

  • propagate_gsparams – Whether to propagate gsparams to each of the components. This is normally a good idea, but there may be use cases where one would not want to do this. [default: True]

evaluateAtWavelength(wave)[source]

Evaluate this chromatic object at a particular wavelength wave.

Parameters:

wave – Wavelength in nanometers.

Returns:

the monochromatic object at the given wavelength.

property gsparams

The GSParams for this object.

withGSParams(gsparams=None, **kwargs)[source]

Create a version of the current object with the given gsparams

Note: if this object wraps other objects (e.g. Convolution, Sum, Transformation, etc.) those component objects will also have their gsparams updated to the new value.

class galsim.ChromaticAutoCorrelation(obj, real_space=None, gsparams=None, propagate_gsparams=True)[source]

Bases: ChromaticObject

A special class for correlating a ChromaticObject with itself.

It is equivalent in functionality to:

galsim.Convolve([obj,obj.rotate(180.*galsim.degrees)])

but takes advantage of the fact that the two profiles are the same for some efficiency gains.

This is the type returned from galsim.AutoCorrelate if the argument is a ChromaticObject. This is the normal way to construct this class.

Parameters:
  • obj – The object to be convolved with itself.

  • real_space – Whether to use real space convolution. [default: None, which means to automatically decide this according to whether the objects have hard edges.]

  • gsparams – An optional GSParams argument. See the docstring for GSParams for details. [default: None]

  • propagate_gsparams – Whether to propagate gsparams to each of the components. This is normally a good idea, but there may be use cases where one would not want to do this. [default: True]

evaluateAtWavelength(wave)[source]

Evaluate this chromatic object at a particular wavelength wave.

Parameters:

wave – Wavelength in nanometers.

Returns:

the monochromatic object at the given wavelength.

property gsparams

The GSParams for this object.

withGSParams(gsparams=None, **kwargs)[source]

Create a version of the current object with the given gsparams

Note: if this object wraps other objects (e.g. Convolution, Sum, Transformation, etc.) those component objects will also have their gsparams updated to the new value.

class galsim.ChromaticTransformation(obj, jac=None, offset=(0.0, 0.0), flux_ratio=1.0, redshift=None, gsparams=None, propagate_gsparams=True)[source]

Bases: ChromaticObject

A class for modeling a wavelength-dependent affine transformation of a ChromaticObject instance.

Typically, you do not need to construct a ChromaticTransformation object explicitly. This is the type returned by the various transformation methods of ChromaticObject such as ChromaticObject.shear, ChromaticObject.rotate, ChromaticObject.shift, etc.

All the various transformations can be described as a combination of a jacobian matrix (i.e. ChromaticObject.transform) and a translation (ChromaticObject.shift), which are described by (dudx,dudy,dvdx,dvdy) and (dx,dy) respectively.

Parameters:
  • obj – The object to be transformed.

  • jac – A list or tuple (dudx, dudy, dvdx, dvdy), or a numpy.array object [[dudx, dudy], [dvdx, dvdy]] describing the Jacobian to apply. May also be a function of wavelength returning a numpy array. Use None to indicate that the Jacobian is the 2x2 unit matrix. [default: None]

  • offset – A galsim.PositionD or list or tuple or numpy array giving the offset (dx,dy) by which to shift the profile. May also be a function of wavelength returning a numpy array. [default: None]

  • flux_ratio – A factor by which to multiply the flux of the object. [default: 1]

  • redshift – A redshift to apply to the wavelength when evaluating. [default: None]

  • gsparams – An optional GSParams argument. See the docstring for GSParams for details. [default: None]

  • propagate_gsparams – Whether to propagate gsparams to each of the components. This is normally a good idea, but there may be use cases where one would not want to do this. [default: True]

drawImage(bandpass, image=None, integrator='quadratic', **kwargs)[source]

See ChromaticObject.drawImage for a full description.

This version usually just calls that one, but if the transformed object (self.original) is an InterpolatedChromaticObject, and the transformation is achromatic, then it will still be able to use the interpolation.

Parameters:
  • bandpass – A Bandpass object representing the filter against which to integrate.

  • image – Optionally, the Image to draw onto. (See GSObject.drawImage for details.) [default: None]

  • integrator – When doing the exact evaluation of the profile, this argument should be one of the image integrators from galsim.integ, or a string ‘trapezoidal’, ‘midpoint’, ‘quadratic’, in which case the routine will use a SampleIntegrator or ContinuousIntegrator depending on whether or not the object has a wave_list. [default: ‘quadratic’, which will try to select an appropriate integrator using the quadratic integration rule automatically.] If the object being transformed is an InterpolatedChromaticObject, then integrator can only be a string, either ‘midpoint’, ‘trapezoidal’, or ‘quadratic’.

  • **kwargs – For all other kwarg options, see GSObject.drawImage.

Returns:

the drawn Image.

evaluateAtWavelength(wave)[source]

Evaluate this chromatic object at a particular wavelength.

Parameters:

wave – Wavelength in nanometers.

Returns:

the monochromatic object at the given wavelength.

property gsparams

The GSParams for this object.

property original

The original object that was transformed.

withGSParams(gsparams=None, **kwargs)[source]

Create a version of the current object with the given gsparams

Note: if this object wraps other objects (e.g. Convolution, Sum, Transformation, etc.) those component objects will also have their gsparams updated to the new value.

class galsim.ChromaticFourierSqrtProfile(obj, gsparams=None, propagate_gsparams=True)[source]

Bases: ChromaticObject

A class for computing the Fourier-space square root of a ChromaticObject.

The ChromaticFourierSqrtProfile class represents a wavelength-dependent Fourier-space square root of a profile.

You may also specify a gsparams argument. See the docstring for GSParams for more information about this option. Note: if gsparams is unspecified (or None), then the ChromaticFourierSqrtProfile inherits the same GSParams as the object being operated on.

The normal way to use this class is to use the FourierSqrt factory function:

>>> fourier_sqrt = galsim.FourierSqrt(chromatic_obj)

If chromatic_obj is indeed a ChromaticObject, then that function will create a ChromaticFourierSqrtProfile object.

Parameters:
  • obj – The object to compute the Fourier-space square root of.

  • gsparams – An optional GSParams argument. See the docstring for GSParams for details. [default: None]

  • propagate_gsparams – Whether to propagate gsparams to each of the components. This is normally a good idea, but there may be use cases where one would not want to do this. [default: True]

evaluateAtWavelength(wave)[source]

Evaluate this chromatic object at a particular wavelength wave.

Parameters:

wave – Wavelength in nanometers.

Returns:

the monochromatic object at the given wavelength.

property gsparams

The GSParams for this object.

withGSParams(gsparams=None, **kwargs)[source]

Create a version of the current object with the given gsparams

Note: if this object wraps other objects (e.g. Convolution, Sum, Transformation, etc.) those component objects will also have their gsparams updated to the new value.