LSST Applications  21.0.0+04719a4bac,21.0.0-1-ga51b5d4+f5e6047307,21.0.0-11-g2b59f77+a9c1acf22d,21.0.0-11-ga42c5b2+86977b0b17,21.0.0-12-gf4ce030+76814010d2,21.0.0-13-g1721dae+760e7a6536,21.0.0-13-g3a573fe+768d78a30a,21.0.0-15-g5a7caf0+f21cbc5713,21.0.0-16-g0fb55c1+b60e2d390c,21.0.0-19-g4cded4ca+71a93a33c0,21.0.0-2-g103fe59+bb20972958,21.0.0-2-g45278ab+04719a4bac,21.0.0-2-g5242d73+3ad5d60fb1,21.0.0-2-g7f82c8f+8babb168e8,21.0.0-2-g8f08a60+06509c8b61,21.0.0-2-g8faa9b5+616205b9df,21.0.0-2-ga326454+8babb168e8,21.0.0-2-gde069b7+5e4aea9c2f,21.0.0-2-gecfae73+1d3a86e577,21.0.0-2-gfc62afb+3ad5d60fb1,21.0.0-25-g1d57be3cd+e73869a214,21.0.0-3-g357aad2+ed88757d29,21.0.0-3-g4a4ce7f+3ad5d60fb1,21.0.0-3-g4be5c26+3ad5d60fb1,21.0.0-3-g65f322c+e0b24896a3,21.0.0-3-g7d9da8d+616205b9df,21.0.0-3-ge02ed75+a9c1acf22d,21.0.0-4-g591bb35+a9c1acf22d,21.0.0-4-g65b4814+b60e2d390c,21.0.0-4-gccdca77+0de219a2bc,21.0.0-4-ge8a399c+6c55c39e83,21.0.0-5-gd00fb1e+05fce91b99,21.0.0-6-gc675373+3ad5d60fb1,21.0.0-64-g1122c245+4fb2b8f86e,21.0.0-7-g04766d7+cd19d05db2,21.0.0-7-gdf92d54+04719a4bac,21.0.0-8-g5674e7b+d1bd76f71f,master-gac4afde19b+a9c1acf22d,w.2021.13
LSST Data Management Base Package
Classes | Functions | Variables
lsst.pipe.tasks.insertFakes Namespace Reference

Classes

class  InsertFakesConnections
 

Functions

def addPixCoords (self, fakeCat, wcs)
 
def trimFakeCat (self, fakeCat, image, wcs)
 
def mkFakeGalsimGalaxies (self, fakeCat, band, photoCalib, pixelScale, psf, image)
 
def mkFakeStars (self, fakeCat, band, photoCalib, psf, image)
 
def cleanCat (self, fakeCat, starCheckVal)
 
def addFakeSources (self, image, fakeImages, sourceType)
 

Variables

 correctedFlux = psf.computeApertureFlux(self.config.calibFluxRadius, xy)
 
 psfKernel = psf.computeKernelImage(xy).getArray()
 
 psfIm = galsim.InterpolatedImage(galsim.Image(psfKernel), scale=pixelScale)
 
 galsimIm = galsim.InterpolatedImage(galsim.Image(im.array), scale=pixelScale)
 
 convIm = galsim.Convolve([galsimIm, psfIm])
 
 outIm = convIm.drawImage(scale=pixelScale, method="real_space").array
 
 imSum = np.sum(outIm)
 
 divIm = outIm/imSum
 
 flux = photoCalib.magnitudeToInstFlux(mag, xy)
 
 imWithFlux = flux*divIm
 

Function Documentation

◆ addFakeSources()

def lsst.pipe.tasks.insertFakes.addFakeSources (   self,
  image,
  fakeImages,
  sourceType 
)
Add the fake sources to the given image

Parameters
----------
image : `lsst.afw.image.exposure.exposure.ExposureF`
            The image into which the fake sources should be added
fakeImages : `typing.Iterator` [`tuple` ['lsst.afw.image.ImageF`, `lsst.geom.Point2d`]]
            An iterator of tuples that contains (or generates) images of fake sources,
            and the locations they are to be inserted at.
sourceType : `str`
            The type (star/galaxy) of fake sources input

Returns
-------
image : `lsst.afw.image.exposure.exposure.ExposureF`

Notes
-----
Uses the x, y information in the ``fakeCat`` to position an image of the fake interpolated onto the
pixel grid of the image. Sets the ``FAKE`` mask plane for the pixels added with the fake source.

Definition at line 716 of file insertFakes.py.

716  def addFakeSources(self, image, fakeImages, sourceType):
717  """Add the fake sources to the given image
718 
719  Parameters
720  ----------
721  image : `lsst.afw.image.exposure.exposure.ExposureF`
722  The image into which the fake sources should be added
723  fakeImages : `typing.Iterator` [`tuple` ['lsst.afw.image.ImageF`, `lsst.geom.Point2d`]]
724  An iterator of tuples that contains (or generates) images of fake sources,
725  and the locations they are to be inserted at.
726  sourceType : `str`
727  The type (star/galaxy) of fake sources input
728 
729  Returns
730  -------
731  image : `lsst.afw.image.exposure.exposure.ExposureF`
732 
733  Notes
734  -----
735  Uses the x, y information in the ``fakeCat`` to position an image of the fake interpolated onto the
736  pixel grid of the image. Sets the ``FAKE`` mask plane for the pixels added with the fake source.
737  """
738 
739  imageBBox = image.getBBox()
740  imageMI = image.maskedImage
741 
742  for (fakeImage, xy) in fakeImages:
743  X0 = xy.getX() - fakeImage.getWidth()/2 + 0.5
744  Y0 = xy.getY() - fakeImage.getHeight()/2 + 0.5
745  self.log.debug("Adding fake source at %d, %d" % (xy.getX(), xy.getY()))
746  if sourceType == "galaxy":
747  interpFakeImage = afwMath.offsetImage(fakeImage, X0, Y0, "lanczos3")
748  else:
749  interpFakeImage = fakeImage
750 
751  interpFakeImBBox = interpFakeImage.getBBox()
752  interpFakeImBBox.clip(imageBBox)
753 
754  if interpFakeImBBox.getArea() > 0:
755  imageMIView = imageMI[interpFakeImBBox]
756  clippedFakeImage = interpFakeImage[interpFakeImBBox]
757  clippedFakeImageMI = afwImage.MaskedImageF(clippedFakeImage)
758  clippedFakeImageMI.mask.set(self.bitmask)
759  imageMIView += clippedFakeImageMI
760 
761  return image
762 
std::shared_ptr< ImageT > offsetImage(ImageT const &image, float dx, float dy, std::string const &algorithmName="lanczos5", unsigned int buffer=0)
Return an image offset by (dx, dy) using the specified algorithm.
Definition: offsetImage.cc:41
def addFakeSources(self, image, fakeImages, sourceType)
Definition: insertFakes.py:716

◆ addPixCoords()

def lsst.pipe.tasks.insertFakes.addPixCoords (   self,
  fakeCat,
  wcs 
)
Add pixel coordinates to the catalog of fakes.

Parameters
----------
fakeCat : `pandas.core.frame.DataFrame`
            The catalog of fake sources to be input
wcs : `lsst.afw.geom.SkyWcs`
            WCS to use to add fake sources

Returns
-------
fakeCat : `pandas.core.frame.DataFrame`

Notes
-----
The default option is to use the WCS information from the image. If the ``useUpdatedCalibs`` config
option is set then it will use the updated WCS from jointCal.

Definition at line 476 of file insertFakes.py.

476  def addPixCoords(self, fakeCat, wcs):
477 
478  """Add pixel coordinates to the catalog of fakes.
479 
480  Parameters
481  ----------
482  fakeCat : `pandas.core.frame.DataFrame`
483  The catalog of fake sources to be input
484  wcs : `lsst.afw.geom.SkyWcs`
485  WCS to use to add fake sources
486 
487  Returns
488  -------
489  fakeCat : `pandas.core.frame.DataFrame`
490 
491  Notes
492  -----
493  The default option is to use the WCS information from the image. If the ``useUpdatedCalibs`` config
494  option is set then it will use the updated WCS from jointCal.
495  """
496 
497  ras = fakeCat[self.config.raColName].values
498  decs = fakeCat[self.config.decColName].values
499  skyCoords = [SpherePoint(ra, dec, radians) for (ra, dec) in zip(ras, decs)]
500  pixCoords = wcs.skyToPixel(skyCoords)
501  xs = [coord.getX() for coord in pixCoords]
502  ys = [coord.getY() for coord in pixCoords]
503  fakeCat["x"] = xs
504  fakeCat["y"] = ys
505 
506  return fakeCat
507 
lsst::geom::SpherePoint SpherePoint
Definition: misc.h:35
def addPixCoords(self, fakeCat, wcs)
Definition: insertFakes.py:476

◆ cleanCat()

def lsst.pipe.tasks.insertFakes.cleanCat (   self,
  fakeCat,
  starCheckVal 
)
Remove rows from the fakes catalog which have HLR = 0 for either the buldge or disk component,
   also remove galaxies that have Sersic index outside the galsim min and max
   allowed (0.3 <= n <= 6.2).

Parameters
----------
fakeCat : `pandas.core.frame.DataFrame`
            The catalog of fake sources to be input
starCheckVal : `str`, `bytes` or `int`
            The value that is set in the sourceType column to specifiy an object is a star.

Returns
-------
fakeCat : `pandas.core.frame.DataFrame`
            The input catalog of fake sources but with the bad objects removed

Notes
-----
If the config option sourceSelectionColName is set then only objects with this column set to True
will be added.

Definition at line 665 of file insertFakes.py.

665  def cleanCat(self, fakeCat, starCheckVal):
666  """Remove rows from the fakes catalog which have HLR = 0 for either the buldge or disk component,
667  also remove galaxies that have Sersic index outside the galsim min and max
668  allowed (0.3 <= n <= 6.2).
669 
670  Parameters
671  ----------
672  fakeCat : `pandas.core.frame.DataFrame`
673  The catalog of fake sources to be input
674  starCheckVal : `str`, `bytes` or `int`
675  The value that is set in the sourceType column to specifiy an object is a star.
676 
677  Returns
678  -------
679  fakeCat : `pandas.core.frame.DataFrame`
680  The input catalog of fake sources but with the bad objects removed
681 
682  Notes
683  -----
684  If the config option sourceSelectionColName is set then only objects with this column set to True
685  will be added.
686  """
687 
688  rowsToKeep = (((fakeCat[self.config.bulgeHLR] != 0.0) & (fakeCat[self.config.diskHLR] != 0.0))
689  | (fakeCat[self.config.sourceType] == starCheckVal))
690  numRowsNotUsed = len(fakeCat) - len(np.where(rowsToKeep)[0])
691  self.log.info("Removing %d rows with HLR = 0 for either the bulge or disk" % numRowsNotUsed)
692  fakeCat = fakeCat[rowsToKeep]
693 
694  minN = galsim.Sersic._minimum_n
695  maxN = galsim.Sersic._maximum_n
696  rowsWithGoodSersic = (((fakeCat[self.config.nBulge] >= minN) & (fakeCat[self.config.nBulge] <= maxN)
697  & (fakeCat[self.config.nDisk] >= minN) & (fakeCat[self.config.nDisk] <= maxN))
698  | (fakeCat[self.config.sourceType] == starCheckVal))
699  numRowsNotUsed = len(fakeCat) - len(np.where(rowsWithGoodSersic)[0])
700  self.log.info("Removing %d rows of galaxies with nBulge or nDisk outside of %0.2f <= n <= %0.2f" %
701  (numRowsNotUsed, minN, maxN))
702  fakeCat = fakeCat[rowsWithGoodSersic]
703 
704  if self.config.doSubSelectSources:
705  try:
706  rowsSelected = (fakeCat[self.config.sourceSelectionColName])
707  except KeyError:
708  raise KeyError("Given column, %s, for source selection not found." %
709  self.config.sourceSelectionColName)
710  numRowsNotUsed = len(fakeCat) - len(rowsSelected)
711  self.log.info("Removing %d rows which were not designated as template sources" % numRowsNotUsed)
712  fakeCat = fakeCat[rowsSelected]
713 
714  return fakeCat
715 
def cleanCat(self, fakeCat, starCheckVal)
Definition: insertFakes.py:665

◆ mkFakeGalsimGalaxies()

def lsst.pipe.tasks.insertFakes.mkFakeGalsimGalaxies (   self,
  fakeCat,
  band,
  photoCalib,
  pixelScale,
  psf,
  image 
)
Make images of fake galaxies using GalSim.

Parameters
----------
band : `str`
pixelScale : `float`
psf : `lsst.meas.extensions.psfex.psfexPsf.PsfexPsf`
            The PSF information to use to make the PSF images
fakeCat : `pandas.core.frame.DataFrame`
            The catalog of fake sources to be input
photoCalib : `lsst.afw.image.photoCalib.PhotoCalib`
            Photometric calibration to be used to calibrate the fake sources

Yields
-------
galImages : `generator`
            A generator of tuples of `lsst.afw.image.exposure.exposure.ExposureF` and
            `lsst.geom.Point2D` of their locations.

Notes
-----

Fake galaxies are made by combining two sersic profiles, one for the bulge and one for the disk. Each
component has an individual sersic index (n), a, b and position angle (PA). The combined profile is
then convolved with the PSF at the specified x, y position on the image.

The names of the columns in the ``fakeCat`` are configurable and are the column names from the
University of Washington simulations database as default. For more information see the doc strings
attached to the config options.

See mkFakeStars doc string for an explanation of calibration to instrumental flux.

Definition at line 535 of file insertFakes.py.

535  def mkFakeGalsimGalaxies(self, fakeCat, band, photoCalib, pixelScale, psf, image):
536  """Make images of fake galaxies using GalSim.
537 
538  Parameters
539  ----------
540  band : `str`
541  pixelScale : `float`
542  psf : `lsst.meas.extensions.psfex.psfexPsf.PsfexPsf`
543  The PSF information to use to make the PSF images
544  fakeCat : `pandas.core.frame.DataFrame`
545  The catalog of fake sources to be input
546  photoCalib : `lsst.afw.image.photoCalib.PhotoCalib`
547  Photometric calibration to be used to calibrate the fake sources
548 
549  Yields
550  -------
551  galImages : `generator`
552  A generator of tuples of `lsst.afw.image.exposure.exposure.ExposureF` and
553  `lsst.geom.Point2D` of their locations.
554 
555  Notes
556  -----
557 
558  Fake galaxies are made by combining two sersic profiles, one for the bulge and one for the disk. Each
559  component has an individual sersic index (n), a, b and position angle (PA). The combined profile is
560  then convolved with the PSF at the specified x, y position on the image.
561 
562  The names of the columns in the ``fakeCat`` are configurable and are the column names from the
563  University of Washington simulations database as default. For more information see the doc strings
564  attached to the config options.
565 
566  See mkFakeStars doc string for an explanation of calibration to instrumental flux.
567  """
568 
569  self.log.info("Making %d fake galaxy images" % len(fakeCat))
570 
571  for (index, row) in fakeCat.iterrows():
572  xy = geom.Point2D(row["x"], row["y"])
573 
574  # We put these two PSF calculations within this same try block so that we catch cases
575  # where the object's position is outside of the image.
576  try:
577  correctedFlux = psf.computeApertureFlux(self.config.calibFluxRadius, xy)
578  psfKernel = psf.computeKernelImage(xy).getArray()
579  psfKernel /= correctedFlux
580 
581  except InvalidParameterError:
582  self.log.info("Galaxy at %0.4f, %0.4f outside of image" % (row["x"], row["y"]))
583  continue
584 
585  try:
586  flux = photoCalib.magnitudeToInstFlux(row[self.config.magVar % band], xy)
587  except LogicError:
588  flux = 0
589 
590  bulge = galsim.Sersic(row[self.config.nBulge], half_light_radius=row[self.config.bulgeHLR])
591  axisRatioBulge = row[self.config.bBulge]/row[self.config.aBulge]
592  bulge = bulge.shear(q=axisRatioBulge, beta=((90 - row[self.config.paBulge])*galsim.degrees))
593 
594  disk = galsim.Sersic(row[self.config.nDisk], half_light_radius=row[self.config.diskHLR])
595  axisRatioDisk = row[self.config.bDisk]/row[self.config.aDisk]
596  disk = disk.shear(q=axisRatioDisk, beta=((90 - row[self.config.paDisk])*galsim.degrees))
597 
598  gal = disk + bulge
599  gal = gal.withFlux(flux)
600 
601  psfIm = galsim.InterpolatedImage(galsim.Image(psfKernel), scale=pixelScale)
602  gal = galsim.Convolve([gal, psfIm])
603  try:
604  galIm = gal.drawImage(scale=pixelScale, method="real_space").array
605  except (galsim.errors.GalSimFFTSizeError, MemoryError):
606  continue
607 
608  yield (afwImage.ImageF(galIm), xy)
609 
def mkFakeGalsimGalaxies(self, fakeCat, band, photoCalib, pixelScale, psf, image)
Definition: insertFakes.py:535

◆ mkFakeStars()

def lsst.pipe.tasks.insertFakes.mkFakeStars (   self,
  fakeCat,
  band,
  photoCalib,
  psf,
  image 
)
Make fake stars based off the properties in the fakeCat.

Parameters
----------
band : `str`
psf : `lsst.meas.extensions.psfex.psfexPsf.PsfexPsf`
            The PSF information to use to make the PSF images
fakeCat : `pandas.core.frame.DataFrame`
            The catalog of fake sources to be input
image : `lsst.afw.image.exposure.exposure.ExposureF`
            The image into which the fake sources should be added
photoCalib : `lsst.afw.image.photoCalib.PhotoCalib`
            Photometric calibration to be used to calibrate the fake sources

Yields
-------
starImages : `generator`
            A generator of tuples of `lsst.afw.image.ImageF` of fake stars and
            `lsst.geom.Point2D` of their locations.

Notes
-----
To take a given magnitude and translate to the number of counts in the image
we use photoCalib.magnitudeToInstFlux, which returns the instrumental flux for the
given calibration radius used in the photometric calibration step.
Thus `calibFluxRadius` should be set to this same radius so that we can normalize
the PSF model to the correct instrumental flux within calibFluxRadius.

Definition at line 610 of file insertFakes.py.

610  def mkFakeStars(self, fakeCat, band, photoCalib, psf, image):
611 
612  """Make fake stars based off the properties in the fakeCat.
613 
614  Parameters
615  ----------
616  band : `str`
617  psf : `lsst.meas.extensions.psfex.psfexPsf.PsfexPsf`
618  The PSF information to use to make the PSF images
619  fakeCat : `pandas.core.frame.DataFrame`
620  The catalog of fake sources to be input
621  image : `lsst.afw.image.exposure.exposure.ExposureF`
622  The image into which the fake sources should be added
623  photoCalib : `lsst.afw.image.photoCalib.PhotoCalib`
624  Photometric calibration to be used to calibrate the fake sources
625 
626  Yields
627  -------
628  starImages : `generator`
629  A generator of tuples of `lsst.afw.image.ImageF` of fake stars and
630  `lsst.geom.Point2D` of their locations.
631 
632  Notes
633  -----
634  To take a given magnitude and translate to the number of counts in the image
635  we use photoCalib.magnitudeToInstFlux, which returns the instrumental flux for the
636  given calibration radius used in the photometric calibration step.
637  Thus `calibFluxRadius` should be set to this same radius so that we can normalize
638  the PSF model to the correct instrumental flux within calibFluxRadius.
639  """
640 
641  self.log.info("Making %d fake star images" % len(fakeCat))
642 
643  for (index, row) in fakeCat.iterrows():
644  xy = geom.Point2D(row["x"], row["y"])
645 
646  # We put these two PSF calculations within this same try block so that we catch cases
647  # where the object's position is outside of the image.
648  try:
649  correctedFlux = psf.computeApertureFlux(self.config.calibFluxRadius, xy)
650  starIm = psf.computeImage(xy)
651  starIm /= correctedFlux
652 
653  except InvalidParameterError:
654  self.log.info("Star at %0.4f, %0.4f outside of image" % (row["x"], row["y"]))
655  continue
656 
657  try:
658  flux = photoCalib.magnitudeToInstFlux(row[self.config.magVar % band], xy)
659  except LogicError:
660  flux = 0
661 
662  starIm *= flux
663  yield ((starIm.convertF(), xy))
664 
def mkFakeStars(self, fakeCat, band, photoCalib, psf, image)
Definition: insertFakes.py:610

◆ trimFakeCat()

def lsst.pipe.tasks.insertFakes.trimFakeCat (   self,
  fakeCat,
  image,
  wcs 
)
Trim the fake cat to about the size of the input image.

`fakeCat` must be processed with addPixCoords before using this method.

Parameters
----------
fakeCat : `pandas.core.frame.DataFrame`
            The catalog of fake sources to be input
image : `lsst.afw.image.exposure.exposure.ExposureF`
            The image into which the fake sources should be added
wcs : `lsst.afw.geom.SkyWcs`
            WCS to use to add fake sources

Returns
-------
fakeCat : `pandas.core.frame.DataFrame`
            The original fakeCat trimmed to the area of the image

Definition at line 508 of file insertFakes.py.

508  def trimFakeCat(self, fakeCat, image, wcs):
509  """Trim the fake cat to about the size of the input image.
510 
511  `fakeCat` must be processed with addPixCoords before using this method.
512 
513  Parameters
514  ----------
515  fakeCat : `pandas.core.frame.DataFrame`
516  The catalog of fake sources to be input
517  image : `lsst.afw.image.exposure.exposure.ExposureF`
518  The image into which the fake sources should be added
519  wcs : `lsst.afw.geom.SkyWcs`
520  WCS to use to add fake sources
521 
522  Returns
523  -------
524  fakeCat : `pandas.core.frame.DataFrame`
525  The original fakeCat trimmed to the area of the image
526  """
527 
528  bbox = Box2D(image.getBBox())
529 
530  def trim(row):
531  return bbox.dilatedBy(self.config.trimBuffer).contains(row["x"], row["y"])
532 
533  return fakeCat[fakeCat.apply(trim, axis=1)]
534 
def trimFakeCat(self, fakeCat, image, wcs)
Definition: insertFakes.py:508

Variable Documentation

◆ convIm

lsst.pipe.tasks.insertFakes.convIm = galsim.Convolve([galsimIm, psfIm])

Definition at line 452 of file insertFakes.py.

◆ correctedFlux

lsst.pipe.tasks.insertFakes.correctedFlux = psf.computeApertureFlux(self.config.calibFluxRadius, xy)
raColName = pexConfig.Field(
    doc="RA column name used in the fake source catalog.",
    dtype=str,
    default="raJ2000",
)

decColName = pexConfig.Field(
    doc="Dec. column name used in the fake source catalog.",
    dtype=str,
    default="decJ2000",
)

doCleanCat = pexConfig.Field(
    doc="If true removes bad sources from the catalog.",
    dtype=bool,
    default=True,
)

diskHLR = pexConfig.Field(
    doc="Column name for the disk half light radius used in the fake source catalog.",
    dtype=str,
    default="DiskHalfLightRadius",
)

bulgeHLR = pexConfig.Field(
    doc="Column name for the bulge half light radius used in the fake source catalog.",
    dtype=str,
    default="BulgeHalfLightRadius",
)

magVar = pexConfig.Field(
    doc="The column name for the magnitude calculated taking variability into account. In the format "
        "``filter name``magVar, e.g. imagVar for the magnitude in the i band.",
    dtype=str,
    default="%smagVar",
)

nDisk = pexConfig.Field(
    doc="The column name for the sersic index of the disk component used in the fake source catalog.",
    dtype=str,
    default="disk_n",
)

nBulge = pexConfig.Field(
    doc="The column name for the sersic index of the bulge component used in the fake source catalog.",
    dtype=str,
    default="bulge_n",
)

aDisk = pexConfig.Field(
    doc="The column name for the semi major axis length of the disk component used in the fake source"
        "catalog.",
    dtype=str,
    default="a_d",
)

aBulge = pexConfig.Field(
    doc="The column name for the semi major axis length of the bulge component.",
    dtype=str,
    default="a_b",
)

bDisk = pexConfig.Field(
    doc="The column name for the semi minor axis length of the disk component.",
    dtype=str,
    default="b_d",
)

bBulge = pexConfig.Field(
    doc="The column name for the semi minor axis length of the bulge component used in the fake source "
        "catalog.",
    dtype=str,
    default="b_b",
)

paDisk = pexConfig.Field(
    doc="The column name for the PA of the disk component used in the fake source catalog.",
    dtype=str,
    default="pa_disk",
)

paBulge = pexConfig.Field(
    doc="The column name for the PA of the bulge component used in the fake source catalog.",
    dtype=str,
    default="pa_bulge",
)

sourceType = pexConfig.Field(
    doc="The column name for the source type used in the fake source catalog.",
    dtype=str,
    default="sourceType",
)

fakeType = pexConfig.Field(
    doc="What type of fake catalog to use, snapshot (includes variability in the magnitudes calculated "
        "from the MJD of the image), static (no variability) or filename for a user defined fits"
        "catalog.",
    dtype=str,
    default="static",
)

calibFluxRadius = pexConfig.Field(
    doc="Aperture radius (in pixels) that was used to define the calibration for this image+catalog. "
    "This will be used to produce the correct instrumental fluxes within the radius. "
    "This value should match that of the field defined in slot_CalibFlux_instFlux.",
    dtype=float,
    default=12.0,
)

coaddName = pexConfig.Field(
    doc="The name of the type of coadd used",
    dtype=str,
    default="deep",
)

doSubSelectSources = pexConfig.Field(
    doc="Set to True if you wish to sub select sources to be input based on the value in the column"
        "set in the sourceSelectionColName config option.",
    dtype=bool,
    default=False
)

sourceSelectionColName = pexConfig.Field(
    doc="The name of the column in the input fakes catalogue to be used to determine which sources to"
        "add, default is none and when this is used all sources are added.",
    dtype=str,
    default="templateSource"
)

insertImages = pexConfig.Field(
    doc="Insert images directly? True or False.",
    dtype=bool,
    default=False,
)

doProcessAllDataIds = pexConfig.Field(
    doc="If True, all input data IDs will be processed, even those containing no fake sources.",
    dtype=bool,
    default=False,
)

trimBuffer = pexConfig.Field(
    doc="Size of the pixel buffer surrounding the image. Only those fake sources with a centroid"
    "falling within the image+buffer region will be considered for fake source injection.",
    dtype=int,
    default=100,
)


class InsertFakesTask(PipelineTask, CmdLineTask):

Definition at line 442 of file insertFakes.py.

◆ divIm

lsst.pipe.tasks.insertFakes.divIm = outIm/imSum

Definition at line 460 of file insertFakes.py.

◆ flux

int lsst.pipe.tasks.insertFakes.flux = photoCalib.magnitudeToInstFlux(mag, xy)

Definition at line 463 of file insertFakes.py.

◆ galsimIm

lsst.pipe.tasks.insertFakes.galsimIm = galsim.InterpolatedImage(galsim.Image(im.array), scale=pixelScale)

Definition at line 451 of file insertFakes.py.

◆ imSum

lsst.pipe.tasks.insertFakes.imSum = np.sum(outIm)

Definition at line 459 of file insertFakes.py.

◆ imWithFlux

lsst.pipe.tasks.insertFakes.imWithFlux = flux*divIm

Definition at line 467 of file insertFakes.py.

◆ outIm

lsst.pipe.tasks.insertFakes.outIm = convIm.drawImage(scale=pixelScale, method="real_space").array

Definition at line 455 of file insertFakes.py.

◆ psfIm

lsst.pipe.tasks.insertFakes.psfIm = galsim.InterpolatedImage(galsim.Image(psfKernel), scale=pixelScale)

Definition at line 450 of file insertFakes.py.

◆ psfKernel

lsst.pipe.tasks.insertFakes.psfKernel = psf.computeKernelImage(xy).getArray()

Definition at line 443 of file insertFakes.py.