Galaxies
There are a number of profiles that are designed to be used for galaxy profiles. Again though, there is nothing restricting these classes to be used only for that purpose if you have another use case for which one would be relevant.
Exponenatial Profile
- class galsim.Exponential(half_light_radius=None, scale_radius=None, flux=1.0, gsparams=None)[source]
Bases:
GSObject
A class describing an exponential profile.
Surface brightness profile with
\[I(r) \sim e^{-r/r_0}\]where \(r_0\) is the
scale_radius
. This is a special case of theSersic
profile, but is given a separate class since the Fourier transform has closed form and can be generated without lookup tables.An Exponential can be initialized using one (and only one) of two possible size parameters:
scale_radius
orhalf_light_radius
. Exactly one of these two is required.- Parameters:
half_light_radius – The half-light radius of the profile. Typically given in arcsec. [One of
scale_radius
orhalf_light_radius
is required.]scale_radius – The scale radius of the profile. Typically given in arcsec. [One of
scale_radius
orhalf_light_radius
is required.]flux – The flux (in photons/cm^2/s) of the profile. [default: 1]
gsparams – An optional
GSParams
argument. [default: None]
- property half_light_radius
The half-light radius of the profile.
- property scale_radius
The scale radius of the 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 anSED
, the return value will be aChromaticObject
with specifiedSED
.- Parameters:
flux – The new flux for the object.
- Returns:
the object with the new flux
De Vaucouleurs Profile
- class galsim.DeVaucouleurs(half_light_radius=None, scale_radius=None, flux=1.0, trunc=0.0, flux_untruncated=False, gsparams=None)[source]
Bases:
Sersic
A class describing DeVaucouleurs profile objects.
Surface brightness profile with
\[I(r) \sim e^{-(r/r_0)^{1/4}}\]where \(r_0\) is the
scale_radius
. This is completely equivalent to a Sersic with n=4.For more information, refer to
http://en.wikipedia.org/wiki/De_Vaucouleurs’_law
A DeVaucouleurs can be initialized using one (and only one) of two possible size parameters:
scale_radius
orhalf_light_radius
. Exactly one of these two is required.- Parameters:
scale_radius – The value of scale radius of the profile. Typically given in arcsec. [One of
scale_radius
orhalf_light_radius
is required.]half_light_radius – The half-light radius of the profile. Typically given in arcsec. [One of
scale_radius
orhalf_light_radius
is required.]flux – The flux (in photons/cm^2/s) of the profile. [default: 1]
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_untruncated – Should the provided
flux
andhalf_light_radius
refer to the untruncated profile? See the docstring for Sersic for more details. [default: False]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 anSED
, the return value will be aChromaticObject
with specifiedSED
.- Parameters:
flux – The new flux for the object.
- Returns:
the object with the new flux
Sersic Profile
- class galsim.Sersic(n, half_light_radius=None, scale_radius=None, flux=1.0, trunc=0.0, flux_untruncated=False, gsparams=None)[source]
Bases:
GSObject
A class describing a Sersic profile.
The Sersic surface brightness profile is characterized by three properties: its Sersic index
n
, itsflux
, and either thehalf_light_radius
orscale_radius
. Given these properties, the surface brightness profile scales as\[I(r) \sim e^{-(r/r_0)^{1/n}}\]where \(r_0\) is the
scale_radius
, or\[I(r) \sim e^{-b (r/r_e)^{1/n}}\]where \(r_e\) is the
half_light_radius
and \(b\) is calculated to give the right half-light radius.For more information, refer to
http://en.wikipedia.org/wiki/Sersic_profile
The allowed range of values for the
n
parameter is 0.3 <= n <= 6.2. An exception will be thrown if you provide a value outside that range. Below n=0.3, there are severe numerical problems. Above n=6.2, we found that the code begins to be inaccurate when sheared or magnified (at the level of upcoming shear surveys), so we do not recommend extending beyond this. See Issues #325 and #450 for more details.Sersic profile calculations take advantage of Hankel transform tables that are precomputed for a given value of n when the Sersic profile is initialized. Making additional objects with the same n can therefore be many times faster than making objects with different values of n that have not been used before. Moreover, these Hankel transforms are only cached for a maximum of 100 different n values at a time. For this reason, for large sets of simulations, it is worth considering the use of only discrete n values rather than allowing it to vary continuously. For more details, see https://github.com/GalSim-developers/GalSim/issues/566.
Note that if you are building many Sersic profiles using truncation, the code will be more efficient if the truncation is always the same multiple of
scale_radius
, since it caches many calculations that depend on the ratiotrunc/scale_radius
.A Sersic can be initialized using one (and only one) of two possible size parameters:
scale_radius
orhalf_light_radius
. Exactly one of these two is required.Flux of a truncated profile:
If you are truncating the profile, the optional parameter,
flux_untruncated
, specifies whether theflux
andhalf_light_radius
specifications correspond to the untruncated profile (True
) or to the truncated profile (False
, default). The impact of this parameter is a little subtle, so we’ll go through a few examples to show how it works.First, let’s examine the case where we specify the size according to the half-light radius. If
flux_untruncated
is True (andtrunc > 0
), then the profile will be identical to the version without truncation up to the truncation radius, beyond which it drops to 0. In this case, the actual half-light radius will be different from the specified half-light radius. The half_light_radius property will return the true half-light radius. Similarly, the actual flux will not be the same as the specified value; the true flux is also returned by the flux property.Example:
>>> sersic_obj1 = galsim.Sersic(n=3.5, half_light_radius=2.5, flux=40.) >>> sersic_obj2 = galsim.Sersic(n=3.5, half_light_radius=2.5, flux=40., trunc=10.) >>> sersic_obj3 = galsim.Sersic(n=3.5, half_light_radius=2.5, flux=40., trunc=10., \\ flux_untruncated=True) >>> sersic_obj1.xValue(galsim.PositionD(0.,0.)) 237.3094228615618 >>> sersic_obj2.xValue(galsim.PositionD(0.,0.)) 142.54505376530574 # Normalization and scale radius adjusted (same half-light radius) >>> sersic_obj3.xValue(galsim.PositionD(0.,0.)) 237.30942286156187 >>> sersic_obj1.xValue(galsim.PositionD(10.0001,0.)) 0.011776164687304694 >>> sersic_obj2.xValue(galsim.PositionD(10.0001,0.)) 0.0 >>> sersic_obj3.xValue(galsim.PositionD(10.0001,0.)) 0.0 >>> sersic_obj1.half_light_radius 2.5 >>> sersic_obj2.half_light_radius 2.5 >>> sersic_obj3.half_light_radius 1.9795101383056892 # The true half-light radius is smaller than the specified value >>> sersic_obj1.flux 40.0 >>> sersic_obj2.flux 40.0 >>> sersic_obj3.flux 34.56595186009519 # Flux is missing due to truncation >>> sersic_obj1.scale_radius 0.003262738739834598 >>> sersic_obj2.scale_radius 0.004754602453641744 # the scale radius needed adjustment to accommodate HLR >>> sersic_obj3.scale_radius 0.003262738739834598 # the scale radius is still identical to the untruncated case
When the truncated Sersic scale is specified with
scale_radius
, the behavior between the three cases (untruncated,flux_untruncated=True
andflux_untruncated=False
) will be somewhat different from above. Since it is the scale radius that is being specified, and since truncation does not change the scale radius the way it can change the half-light radius, the scale radius will remain unchanged in all cases. This also results in the half-light radius being the same between the two truncated cases (although different from the untruncated case). The flux normalization is the only difference betweenflux_untruncated=True
andflux_untruncated=False
in this case.Example:
>>> sersic_obj1 = galsim.Sersic(n=3.5, scale_radius=0.05, flux=40.) >>> sersic_obj2 = galsim.Sersic(n=3.5, scale_radius=0.05, flux=40., trunc=10.) >>> sersic_obj3 = galsim.Sersic(n=3.5, scale_radius=0.05, flux=40., trunc=10., \\ flux_untruncated=True) >>> sersic_obj1.xValue(galsim.PositionD(0.,0.)) 1.010507575186637 >>> sersic_obj2.xValue(galsim.PositionD(0.,0.)) 5.786692612210923 # Normalization adjusted to accomodate the flux within trunc radius >>> sersic_obj3.xValue(galsim.PositionD(0.,0.)) 1.010507575186637 >>> sersic_obj1.half_light_radius 38.311372735390016 >>> sersic_obj2.half_light_radius 5.160062547614234 >>> sersic_obj3.half_light_radius 5.160062547614234 # For the truncated cases, the half-light radii are the same >>> sersic_obj1.flux 40.0 >>> sersic_obj2.flux 40.0 >>> sersic_obj3.flux 6.985044085834393 # Flux is missing due to truncation >>> sersic_obj1.scale_radius 0.05 >>> sersic_obj2.scale_radius 0.05 >>> sersic_obj3.scale_radius 0.05 half_light_radius: The half-light radius
- Parameters:
n – The Sersic index, n.
half_light_radius – The half-light radius of the profile. Typically given in arcsec. [One of
scale_radius
orhalf_light_radius
is required.]scale_radius – The scale radius of the profile. Typically given in arcsec. [One of
scale_radius
orhalf_light_radius
is required.]flux – The flux (in photons/cm^2/s) of the profile. [default: 1]
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_untruncated – Should the provided
flux
andhalf_light_radius
refer to the untruncated profile? See below for more details. [default: False]gsparams – An optional
GSParams
argument. [default: None]
- calculateIntegratedFlux(r)[source]
Return the fraction of the total flux enclosed within a given radius, r
- property half_light_radius
The half-light radius.
- property n
The Sersic parameter n.
- property scale_radius
The scale radius.
- property trunc
The truncation radius (if any).
- 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 anSED
, the return value will be aChromaticObject
with specifiedSED
.- Parameters:
flux – The new flux for the object.
- Returns:
the object with the new flux
Inclined Exponential Profile
- class galsim.InclinedExponential(inclination, half_light_radius=None, scale_radius=None, scale_height=None, scale_h_over_r=None, flux=1.0, gsparams=None)[source]
Bases:
GSObject
A class describing an inclined exponential profile.
The Inclined Exponential surface brightness profile is characterized by three properties: its inclination angle (where 0 degrees = face-on and 90 degrees = edge-on), its scale radius, and its scale height. The 3D light distribution function is:
\[I(R,z) \sim \mathrm{sech}^2 (z/h_s) \, \exp\left(-R/R_s\right)\]where \(z\) is the distance along the minor axis, \(R\) is the radial distance from the minor axis, \(R_s\) is the scale radius of the disk, and \(h_s\) is the scale height of the disk. The 2D light distribution function is then determined from the scale height and radius here, along with the inclination angle.
In this implementation, the profile is inclined along the y-axis. This means that it will likely need to be rotated in most circumstances.
At present, this profile is not enabled for photon-shooting.
A profile can be initialized using one (and only one) of two possible size parameters:
scale_radius
orhalf_light_radius
. Exactly one of these two is required. Similarly, at most one ofscale_height
andscale_h_over_r
is required; if neither is given, the default of scale_h_over_r = 0.1 will be used. Note that if half_light_radius and scale_h_over_r are supplied (or the default value of scale_h_over_r is used), scale_h_over_r will be assumed to refer to the scale radius, not the half-light radius.- Parameters:
inclination – The inclination angle, which must be a
galsim.Angle
instancescale_radius – The scale radius of the exponential disk. Typically given in arcsec. This can be compared to the ‘scale_radius’ parameter of the
galsim.Exponential
class, and in the face-on case, the same scale radius will result in the same 2D light distribution as with that class.half_light_radius – The half-light radius of the exponential disk, as an alternative to the scale radius.
scale_height – The scale height of the exponential disk. Typically given in arcsec. [default: None]
scale_h_over_r – In lieu of the scale height, you may also specify the ratio of the scale height to the scale radius. [default: 0.1]
flux – The flux (in photons) of the profile. [default: 1]
gsparams – An optional
GSParams
argument. [default: None]
- property disk_half_light_radius
The half-light radius of the exponential disk.
- property inclination
The inclination angle.
- property scale_h_over_r
The ratio scale_height / scale_radius
- property scale_height
The scale height of the disk.
- property scale_radius
The scale radius of the exponential disk.
- 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 anSED
, the return value will be aChromaticObject
with specifiedSED
.- Parameters:
flux – The new flux for the object.
- Returns:
the object with the new flux
Inclined Sersic Profile
- class galsim.InclinedSersic(n, inclination, half_light_radius=None, scale_radius=None, scale_height=None, scale_h_over_r=None, flux=1.0, trunc=0.0, flux_untruncated=False, gsparams=None)[source]
Bases:
GSObject
A class describing an inclined sersic profile. This class is general, and so for certain special cases, more specialized classes will be more efficient. For the case where n==1 with no truncation, the
InclinedExponential
class will be much more efficient. For the case where the inclination angle is zero (face-on), theSersic
class will be slightly more efficient.The InclinedSersic surface brightness profile is characterized by four properties: its Sersic index
n
, its inclination angle (where 0 degrees = face-on and 90 degrees = edge-on), its scale radius, and its scale height. The 3D light distribution function is:\[I(R,z) \sim \mathrm{sech}^2 (z/h_s) \, \exp\left(-(R/R_s)^{1/n}\right)\]where \(z\) is the distance along the minor axis, \(R\) is the radial distance from the minor axis, \(R_s\) is the scale radius of the disk, and \(h_s\) is the scale height of the disk. The 2D light distribution function is then determined from the scale height and radius here, along with the inclination angle.
In this implementation, the profile is inclined along the y-axis. This means that it will likely need to be rotated in most circumstances.
At present, this profile is not enabled for photon-shooting.
The allowed range of values for the
n
parameter is 0.3 <= n <= 6.2. An exception will be thrown if you provide a value outside that range, matching the range of theSersic
profile.This class shares the caching of Hankel transformations with the
Sersic
class; see that class for documentation on efficiency considerations with regards to caching.A profile can be initialized using one (and only one) of two possible size parameters:
scale_radius
orhalf_light_radius
. Exactly one of these two is required. Similarly, at most one ofscale_height
andscale_h_over_r
is required; if neither is given, the default of scale_h_over_r = 0.1 will be used. Note that if half_light_radius and scale_h_over_r are supplied (or the default value of scale_h_over_r is used), scale_h_over_r will be assumed to refer to the scale radius, not the half-light radius.- Parameters:
n – The Sersic index, n.
inclination – The inclination angle, which must be a
galsim.Angle
instancescale_radius – The scale radius of the disk. Typically given in arcsec. This can be compared to the ‘scale_radius’ parameter of the
galsim.Sersic
class, and in the face-on case, the same scale radius will result in the same 2D light distribution as with that class. Exactly one of this and half_light_radius must be provided.half_light_radius – The half-light radius of disk when seen face-on. Exactly one of this and scale_radius must be provided.
scale_height – The scale height of the exponential disk. Typically given in arcsec. [default: None]
scale_h_over_r – In lieu of the scale height, you may specify the ratio of the scale height to the scale radius. [default: 0.1]
flux – The flux (in photons) of the profile. [default: 1]
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_untruncated – Should the provided
flux
andhalf_light_radius
refer to the untruncated profile? See the documentation of theSersic
class for more details. [default: False]gsparams – An optional
GSParams
argument. [default: None]
- property disk_half_light_radius
The half-light radius of the exponential disk.
- property inclination
The inclination angle.
- property n
The Sersic parameter n.
- property scale_h_over_r
The ratio scale_height / scale_radius.
- property scale_height
The scale height of the disk.
- property scale_radius
The scale radius of the exponential disk.
- property trunc
The truncation radius (if any).
- 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 anSED
, the return value will be aChromaticObject
with specifiedSED
.- Parameters:
flux – The new flux for the object.
- Returns:
the object with the new flux
Spergel Profile
- class galsim.Spergel(nu, half_light_radius=None, scale_radius=None, flux=1.0, gsparams=None)[source]
Bases:
GSObject
A class describing a Spergel profile.
The Spergel surface brightness profile is characterized by three properties: its Spergel index
nu
, itsflux
, and either thehalf_light_radius
orscale_radius
. Given these properties, the surface brightness profile scales as\[I(r) \sim \left(\frac{r}{r_0}\right)^\nu K_\nu\left(\frac{r}{r_0}\right)\]where \(r_0\) is the
scale_radius
and \(K_\nu\) is the modified Bessel function of the second kind.The Spergel profile is intended as a generic galaxy profile, somewhat like a
Sersic
profile, but with the advantage of being analytic in both real space and Fourier space. The Spergel index \(\nu\) plays a similar role to the Sersic index \(n\), in that it adjusts the relative peakiness of the profile core and the relative prominence of the profile wings. At \(\nu = 0.5\), the Spergel profile is equivalent to anExponential
profile (or alternatively an :math`n = 1`Sersic
profile). At \(\nu = -0.6\) (and in the radial range near the half-light radius), the Spergel profile is similar to aDeVaucouleurs
profile or \(n = 4\)Sersic
profile.Note that for \(\nu <= 0\), the Spergel profile surface brightness diverges at the origin. This may lead to rendering problems if the profile is not convolved by either a PSF or a pixel and the profile center is precisely on a pixel center.
Due to its analytic Fourier transform and depending on the indices \(n\) and \(\nu\), the Spergel profile can be considerably faster to draw than the roughly equivalent
Sersic
profile. For example, the \(\nu = -0.6\) Spergel profile is roughly 3x faster to draw than an \(n = 4\)Sersic
profile once theSersic
profile cache has been set up. However, if not taking advantage of the cache, for example, if drawingSersic
profiles with \(n\) continuously varying near 4.0 and Spergel profiles with \(\nu\) continuously varying near -0.6, then the Spergel profiles are about 50x faster to draw. At the other end of the galaxy profile spectrum, the \(\nu = 0.5\) Spergel profile, \(n = 1\)Sersic
profile, and theExponential
profile all take about the same amount of time to draw if cached, and the Spergel profile is about 2x faster than theSersic
profile if uncached.For more information, refer to
D. N. Spergel, “ANALYTICAL GALAXY PROFILES FOR PHOTOMETRIC AND LENSING ANALYSIS,” ASTROPHYS J SUPPL S 191(1), 58-65 (2010) [doi:10.1088/0067-0049/191/1/58].
The allowed range of values for the
nu
parameter is -0.85 <=nu
<= 4. An exception will be thrown if you provide a value outside that range. The lower limit is set above the theoretical lower limit of -1 due to numerical difficulties integrating the very peakynu
< -0.85 profiles. The upper limit is set to avoid numerical difficulties evaluating the modified Bessel function of the second kind.A Spergel profile can be initialized using one (and only one) of two possible size parameters:
scale_radius
orhalf_light_radius
. Exactly one of these two is required.- Parameters:
nu – The Spergel index, nu.
half_light_radius – The half-light radius of the profile. Typically given in arcsec. [One of
scale_radius
orhalf_light_radius
is required.]scale_radius – The scale radius of the profile. Typically given in arcsec. [One of
scale_radius
orhalf_light_radius
is required.]flux – The flux (in photons/cm^2/s) of the profile. [default: 1]
gsparams – An optional
GSParams
argument. [default: None]
- property half_light_radius
The half-light radius
- property nu
The Spergel index, nu
- property scale_radius
The scale radius
- 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 anSED
, the return value will be aChromaticObject
with specifiedSED
.- Parameters:
flux – The new flux for the object.
- Returns:
the object with the new flux
Knots of Star Formation
- class galsim.RandomKnots(npoints, half_light_radius=None, flux=None, profile=None, rng=None, gsparams=None)[source]
Bases:
GSObject
A class for generating a set of point sources, following either a
Gaussian
profile or a specified input profile.Uses of this profile include representing an “irregular” galaxy, or adding this profile to an Exponential to represent knots of star formation.
RandomKnots profiles have “shape noise” that depends on the number of point sources used. For example, using the default Gaussian distribution, with 100 points, the shape noise is g~0.05, and this will decrease as more points are added. The profile can be sheared to give additional ellipticity, for example to follow that of an associated disk.
The requested half light radius (hlr) should be thought of as a rough value. With a finite number point sources the actual realized hlr will be noisy.
Note
If providing an input
profile
object, it must be “shoot-able”. Objects that cannot be drawn withmethod='phot'
cannot be used as theprofile
parameter here.- Parameters:
npoints – Number of point sources to generate.
half_light_radius – Optional half light radius of the distribution of points. This value is used for a Gaussian distribution if an explicit profile is not sent. This is the mean half light radius produced by an infinite number of points. A single instance will be noisy. [default None]
flux – Optional total flux in all point sources. This value is used for a Gaussian distribution if an explicit profile is not sent. Defaults to None if profile is sent, otherwise 1. [default: None]
profile – Optional profile to use for drawing points. If a profile is sent, the half_light_radius and flux keywords are invalid. [default: None]
rng – Optional random number generator. Can be any
galsim.BaseDeviate
. If None, the rng is created internally. [default: None]gsparams – Optional
GSParams
for the objects representing each point source. [default: None]
- Attributes:
npoints – The number of points to use as knots
input_half_light_radius – The input half_light_radius
flux – The flux
points – The array of x,y offsets used to create the point sources
Note
The algorithm was originally a modified version of that presented in https://arxiv.org/abs/1312.5514v3. However, we now use the GalSim photon shooting mechanism, which allows the knots to trace any profile, not just a
Gaussian
.- dilate(scale)[source]
Dilate the linear size of the profile by the given
scale
factor, while preserving flux.e.g.
half_light_radius
<–half_light_radius * scale
See expand() and magnify() for versions that preserve surface brightness, and thus changes the flux.
- Parameters:
scale – The linear rescaling factor to apply.
- Returns:
the dilated object.
- expand(scale)[source]
Expand the linear size of the profile by the given
scale
factor, while preserving surface brightness.e.g.
half_light_radius
<–half_light_radius * scale
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.
- Returns:
the expanded object.
- property input_half_light_radius
Get the input half light radius (HLR).
Note the input HLR is not necessarily the realized HLR, due to the finite number of points used in the profile.
If a profile is sent, and that profile is a Transformation object (e.g. it has been sheared, its flux set, etc), then this value will be None.
You can get the calculated half light radius using the calculateHLR method. That value will be valid in all cases.
- property npoints
The number of point sources.
- rotate(theta)[source]
Rotate this object by an
Angle
theta
.- Parameters:
theta – Rotation angle (
Angle
object, positive means anticlockwise).- Returns:
the rotated object.
- shear(*args, **kwargs)[source]
Create a version of the current object with an area-preserving shear applied to it.
The arguments may be either a
Shear
instance or arguments to 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.
- Parameters:
shear – The
Shear
to be applied. Or, as described above, you may instead supply parameters do construct a shear directly. eg.obj.shear(g1=g1,g2=g2)
.- Returns:
the sheared object.
- shift(*args, **kwargs)[source]
Create a version of the current object shifted by some amount in real space.
After this call, the caller’s type will be a
GSObject
. This means that if the caller was a derived type that had extra methods or properties beyond those defined inGSObject
(e.g.Gaussian.sigma
), then these methods are no longer available.Note: in addition to the dx,dy parameter names, you may also supply dx,dy as a tuple, or as a
Position
object.The shift coordinates here are sky coordinates.
GSObject
profiles are always defined in sky coordinates and only later (when they are drawn) is the connection to pixel coordinates established (via a pixel_scale or WCS). So a shift of dx moves the object horizontally in the sky (e.g. west in the local tangent plane of the observation), and dy moves the object vertically (north in the local tangent plane).The units are typically arcsec, but we don’t enforce that anywhere. The units here just need to be consistent with the units used for any size values used by the
GSObject
. The connection of these units to the eventual image pixels is defined by either thepixel_scale
or thewcs
parameter ofGSObject.drawImage
.Note: if you want to shift the object by a set number (or fraction) of pixels in the drawn image, you probably want to use the
offset
parameter ofGSObject.drawImage
rather than this method.- Parameters:
dx – Horizontal shift to apply.
dy – Vertical shift to apply.
Alternatively, you may supply a single parameter as a
Position
instance, rather than the two components separately if that is more convenient.- Parameter:
offset: The shift to apply, given as PositionD(dx,dy) or PositionI(dx,dy)
- Returns:
the shifted object.
- transform(dudx, dudy, dvdx, dvdy)[source]
Create a version of the current object with an arbitrary Jacobian matrix transformation applied to it.
This applies a Jacobian matrix to the coordinate system in which this object is defined. It changes a profile defined in terms of (x,y) to one defined in terms of (u,v) where:
u = dudx x + dudy y v = dvdx x + dvdy y
That is, an arbitrary affine transform, but without the translation (which is easily effected via the
shift
method).Note that this function is similar to expand in that it preserves surface brightness, not flux. If you want to preserve flux, you should also do:
>>> prof *= 1./abs(dudx*dvdy - dudy*dvdx)
- 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(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 anSED
, the return value will be aChromaticObject
with specifiedSED
.- Parameters:
flux – The new flux for the object.
- Returns:
the object with the new flux
- withScaledFlux(flux_ratio)[source]
Create a version of the current object with the flux scaled by the given
flux_ratio
.This function is equivalent to
obj.withFlux(flux_ratio * obj.flux)
. Indeed, withFlux() is implemented in terms of this one.It creates a new object that has the same profile as the original, but with the surface brightness at every location scaled by the given amount. If
flux_ratio
is anSED
, then the returned object is aChromaticObject
with theSED
multiplied by its currentflux
.Note that in this case the
flux
attribute of theGSObject
being scaled gets interpreted as being dimensionless, instead of having its normal units of [photons/s/cm^2]. The photons/s/cm^2 units are (optionally) carried by theSED
instead, or even left out entirely if theSED
is dimensionless itself (see discussion in theChromaticObject
docstring). TheGSObject
flux
attribute does still contribute to theChromaticObject
normalization, though. For example, the following are equivalent:>>> chrom_obj = gsobj.withScaledFlux(sed * 3.0) >>> chrom_obj2 = (gsobj * 3.0).withScaledFlux(sed)
An equivalent, and usually simpler, way to effect this scaling is:
>>> obj = obj * flux_ratio
- Parameters:
flux_ratio – The ratio by which to rescale the flux of the object when creating a new one.
- Returns:
the object with the new flux.