lsst.pipe.tasks  21.0.0-119-g1a9a5132+45e0dac908
insertFakes.py
Go to the documentation of this file.
1 # This file is part of pipe tasks
2 #
3 # Developed for the LSST Data Management System.
4 # This product includes software developed by the LSST Project
5 # (http://www.lsst.org).
6 # See the COPYRIGHT file at the top-level directory of this distribution
7 # for details of code ownership.
8 #
9 # This program is free software: you can redistribute it and/or modify
10 # it under the terms of the GNU General Public License as published by
11 # the Free Software Foundation, either version 3 of the License, or
12 # (at your option) any later version.
13 #
14 # This program is distributed in the hope that it will be useful,
15 # but WITHOUT ANY WARRANTY; without even the implied warranty of
16 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 # GNU General Public License for more details.
18 #
19 # You should have received a copy of the GNU General Public License
20 # along with this program. If not, see <http://www.gnu.org/licenses/>.
21 
22 """
23 Insert fakes into deepCoadds
24 """
25 import galsim
26 from astropy.table import Table
27 import numpy as np
28 
29 import lsst.geom as geom
30 import lsst.afw.image as afwImage
31 import lsst.afw.math as afwMath
32 import lsst.pex.config as pexConfig
33 import lsst.pipe.base as pipeBase
34 
35 from lsst.pipe.base import CmdLineTask, PipelineTask, PipelineTaskConfig, PipelineTaskConnections
36 import lsst.pipe.base.connectionTypes as cT
37 from lsst.pex.exceptions import LogicError, InvalidParameterError
38 from lsst.coadd.utils.coaddDataIdContainer import ExistingCoaddDataIdContainer
39 from lsst.geom import SpherePoint, radians, Box2D, Point2D
40 
41 __all__ = ["InsertFakesConfig", "InsertFakesTask"]
42 
43 
44 def _add_fake_sources(exposure, objects, calibFluxRadius=12.0, logger=None):
45  """Add fake sources to the given exposure
46 
47  Parameters
48  ----------
49  exposure : `lsst.afw.image.exposure.exposure.ExposureF`
50  The exposure into which the fake sources should be added
51  objects : `typing.Iterator` [`tuple` ['lsst.geom.SpherePoint`, `galsim.GSObject`]]
52  An iterator of tuples that contains (or generates) locations and object
53  surface brightness profiles to inject.
54  calibFluxRadius : `float`, optional
55  Aperture radius (in pixels) used to define the calibration for this
56  exposure+catalog. This is used to produce the correct instrumental fluxes
57  within the radius. The value should match that of the field defined in
58  slot_CalibFlux_instFlux.
59  logger : `lsst.log.log.log.Log` or `logging.Logger`, optional
60  Logger.
61  """
62  exposure.mask.addMaskPlane("FAKE")
63  bitmask = exposure.mask.getPlaneBitMask("FAKE")
64  if logger:
65  logger.info(f"Adding mask plane with bitmask {bitmask}")
66 
67  wcs = exposure.getWcs()
68  psf = exposure.getPsf()
69 
70  bbox = exposure.getBBox()
71  fullBounds = galsim.BoundsI(bbox.minX, bbox.maxX, bbox.minY, bbox.maxY)
72  gsImg = galsim.Image(exposure.image.array, bounds=fullBounds)
73 
74  for spt, gsObj in objects:
75  pt = wcs.skyToPixel(spt)
76  posd = galsim.PositionD(pt.x, pt.y)
77  posi = galsim.PositionI(pt.x//1, pt.y//1)
78  if logger:
79  logger.debug(f"Adding fake source at {pt}")
80 
81  mat = wcs.linearizePixelToSky(spt, geom.arcseconds).getMatrix()
82  gsWCS = galsim.JacobianWCS(mat[0, 0], mat[0, 1], mat[1, 0], mat[1, 1])
83 
84  try:
85  psfArr = psf.computeKernelImage(pt).array
86  except InvalidParameterError:
87  # Try mapping to nearest point contained in bbox.
88  contained_pt = Point2D(
89  np.clip(pt.x, bbox.minX, bbox.maxX),
90  np.clip(pt.y, bbox.minY, bbox.maxY)
91  )
92  if pt == contained_pt: # no difference, so skip immediately
93  if logger:
94  logger.infof(
95  "Cannot compute Psf for object at {}; skipping",
96  pt
97  )
98  continue
99  # otherwise, try again with new point
100  try:
101  psfArr = psf.computeKernelImage(contained_pt).array
102  except InvalidParameterError:
103  if logger:
104  logger.infof(
105  "Cannot compute Psf for object at {}; skipping",
106  pt
107  )
108  continue
109  apCorr = psf.computeApertureFlux(calibFluxRadius)
110  psfArr /= apCorr
111  gsPSF = galsim.InterpolatedImage(galsim.Image(psfArr), wcs=gsWCS)
112 
113  conv = galsim.Convolve(gsObj, gsPSF)
114  stampSize = conv.getGoodImageSize(gsWCS.minLinearScale())
115  subBounds = galsim.BoundsI(posi).withBorder(stampSize//2)
116  subBounds &= fullBounds
117 
118  if subBounds.area() > 0:
119  subImg = gsImg[subBounds]
120  offset = posd - subBounds.true_center
121  # Note, for calexp injection, pixel is already part of the PSF and
122  # for coadd injection, it's incorrect to include the output pixel.
123  # So for both cases, we draw using method='no_pixel'.
124  conv.drawImage(
125  subImg,
126  add_to_image=True,
127  offset=offset,
128  wcs=gsWCS,
129  method='no_pixel'
130  )
131 
132  subBox = geom.Box2I(
133  geom.Point2I(subBounds.xmin, subBounds.ymin),
134  geom.Point2I(subBounds.xmax, subBounds.ymax)
135  )
136  exposure[subBox].mask.array |= bitmask
137 
138 
139 def _isWCSGalsimDefault(wcs, hdr):
140  """Decide if wcs = galsim.PixelScale(1.0) is explicitly present in header,
141  or if it's just the galsim default.
142 
143  Parameters
144  ----------
145  wcs : galsim.BaseWCS
146  Potentially default WCS.
147  hdr : galsim.fits.FitsHeader
148  Header as read in by galsim.
149 
150  Returns
151  -------
152  isDefault : bool
153  True if default, False if explicitly set in header.
154  """
155  if wcs != galsim.PixelScale(1.0):
156  return False
157  if hdr.get('GS_WCS') is not None:
158  return False
159  if hdr.get('CTYPE1', 'LINEAR') == 'LINEAR':
160  return not any(k in hdr for k in ['CD1_1', 'CDELT1'])
161  for wcs_type in galsim.fitswcs.fits_wcs_types:
162  # If one of these succeeds, then assume result is explicit
163  try:
164  wcs_type._readHeader(hdr)
165  return False
166  except Exception:
167  pass
168  else:
169  return not any(k in hdr for k in ['CD1_1', 'CDELT1'])
170 
171 
172 class InsertFakesConnections(PipelineTaskConnections,
173  defaultTemplates={"coaddName": "deep",
174  "fakesType": "fakes_"},
175  dimensions=("tract", "patch", "band", "skymap")):
176 
177  image = cT.Input(
178  doc="Image into which fakes are to be added.",
179  name="{coaddName}Coadd",
180  storageClass="ExposureF",
181  dimensions=("tract", "patch", "band", "skymap")
182  )
183 
184  fakeCat = cT.Input(
185  doc="Catalog of fake sources to draw inputs from.",
186  name="{fakesType}fakeSourceCat",
187  storageClass="DataFrame",
188  dimensions=("tract", "skymap")
189  )
190 
191  imageWithFakes = cT.Output(
192  doc="Image with fake sources added.",
193  name="{fakesType}{coaddName}Coadd",
194  storageClass="ExposureF",
195  dimensions=("tract", "patch", "band", "skymap")
196  )
197 
198 
199 class InsertFakesConfig(PipelineTaskConfig,
200  pipelineConnections=InsertFakesConnections):
201  """Config for inserting fake sources
202  """
203 
204  # Unchanged
205 
206  doCleanCat = pexConfig.Field(
207  doc="If true removes bad sources from the catalog.",
208  dtype=bool,
209  default=True,
210  )
211 
212  fakeType = pexConfig.Field(
213  doc="What type of fake catalog to use, snapshot (includes variability in the magnitudes calculated "
214  "from the MJD of the image), static (no variability) or filename for a user defined fits"
215  "catalog.",
216  dtype=str,
217  default="static",
218  )
219 
220  calibFluxRadius = pexConfig.Field(
221  doc="Aperture radius (in pixels) that was used to define the calibration for this image+catalog. "
222  "This will be used to produce the correct instrumental fluxes within the radius. "
223  "This value should match that of the field defined in slot_CalibFlux_instFlux.",
224  dtype=float,
225  default=12.0,
226  )
227 
228  coaddName = pexConfig.Field(
229  doc="The name of the type of coadd used",
230  dtype=str,
231  default="deep",
232  )
233 
234  doSubSelectSources = pexConfig.Field(
235  doc="Set to True if you wish to sub select sources to be input based on the value in the column"
236  "set in the sourceSelectionColName config option.",
237  dtype=bool,
238  default=False
239  )
240 
241  insertImages = pexConfig.Field(
242  doc="Insert images directly? True or False.",
243  dtype=bool,
244  default=False,
245  )
246 
247  doProcessAllDataIds = pexConfig.Field(
248  doc="If True, all input data IDs will be processed, even those containing no fake sources.",
249  dtype=bool,
250  default=False,
251  )
252 
253  trimBuffer = pexConfig.Field(
254  doc="Size of the pixel buffer surrounding the image. Only those fake sources with a centroid"
255  "falling within the image+buffer region will be considered for fake source injection.",
256  dtype=int,
257  default=100,
258  )
259 
260  sourceType = pexConfig.Field(
261  doc="The column name for the source type used in the fake source catalog.",
262  dtype=str,
263  default="sourceType",
264  )
265 
266  # New source catalog config variables
267 
268  ra_col = pexConfig.Field(
269  doc="Source catalog column name for RA (in radians).",
270  dtype=str,
271  default="ra",
272  )
273 
274  dec_col = pexConfig.Field(
275  doc="Source catalog column name for dec (in radians).",
276  dtype=str,
277  default="dec",
278  )
279 
280  bulge_semimajor_col = pexConfig.Field(
281  doc="Source catalog column name for the semimajor axis (in arcseconds) "
282  "of the bulge half-light ellipse.",
283  dtype=str,
284  default="bulge_semimajor",
285  )
286 
287  bulge_axis_ratio_col = pexConfig.Field(
288  doc="Source catalog column name for the axis ratio of the bulge "
289  "half-light ellipse.",
290  dtype=str,
291  default="bulge_axis_ratio",
292  )
293 
294  bulge_pa_col = pexConfig.Field(
295  doc="Source catalog column name for the position angle (measured from "
296  "North through East in degrees) of the semimajor axis of the bulge "
297  "half-light ellipse.",
298  dtype=str,
299  default="bulge_pa",
300  )
301 
302  bulge_n_col = pexConfig.Field(
303  doc="Source catalog column name for the Sersic index of the bulge.",
304  dtype=str,
305  default="bulge_n",
306  )
307 
308  disk_semimajor_col = pexConfig.Field(
309  doc="Source catalog column name for the semimajor axis (in arcseconds) "
310  "of the disk half-light ellipse.",
311  dtype=str,
312  default="disk_semimajor",
313  )
314 
315  disk_axis_ratio_col = pexConfig.Field(
316  doc="Source catalog column name for the axis ratio of the disk "
317  "half-light ellipse.",
318  dtype=str,
319  default="disk_axis_ratio",
320  )
321 
322  disk_pa_col = pexConfig.Field(
323  doc="Source catalog column name for the position angle (measured from "
324  "North through East in degrees) of the semimajor axis of the disk "
325  "half-light ellipse.",
326  dtype=str,
327  default="disk_pa",
328  )
329 
330  disk_n_col = pexConfig.Field(
331  doc="Source catalog column name for the Sersic index of the disk.",
332  dtype=str,
333  default="disk_n",
334  )
335 
336  bulge_disk_flux_ratio_col = pexConfig.Field(
337  doc="Source catalog column name for the bulge/disk flux ratio.",
338  dtype=str,
339  default="bulge_disk_flux_ratio",
340  )
341 
342  mag_col = pexConfig.Field(
343  doc="Source catalog column name template for magnitudes, in the format "
344  "``filter name``_mag_col. E.g., if this config variable is set to "
345  "``%s_mag``, then the i-band magnitude will be searched for in the "
346  "``i_mag`` column of the source catalog.",
347  dtype=str,
348  default="%s_mag"
349  )
350 
351  select_col = pexConfig.Field(
352  doc="Source catalog column name to be used to select which sources to "
353  "add.",
354  dtype=str,
355  default="select",
356  )
357 
358  # Deprecated config variables
359 
360  raColName = pexConfig.Field(
361  doc="RA column name used in the fake source catalog.",
362  dtype=str,
363  default="raJ2000",
364  deprecated="Use `ra_col` instead."
365  )
366 
367  decColName = pexConfig.Field(
368  doc="Dec. column name used in the fake source catalog.",
369  dtype=str,
370  default="decJ2000",
371  deprecated="Use `dec_col` instead."
372  )
373 
374  diskHLR = pexConfig.Field(
375  doc="Column name for the disk half light radius used in the fake source catalog.",
376  dtype=str,
377  default="DiskHalfLightRadius",
378  deprecated=(
379  "Use `disk_semimajor_col`, `disk_axis_ratio_col`, and `disk_pa_col`"
380  " to specify disk half-light ellipse."
381  )
382  )
383 
384  aDisk = pexConfig.Field(
385  doc="The column name for the semi major axis length of the disk component used in the fake source"
386  "catalog.",
387  dtype=str,
388  default="a_d",
389  deprecated=(
390  "Use `disk_semimajor_col`, `disk_axis_ratio_col`, and `disk_pa_col`"
391  " to specify disk half-light ellipse."
392  )
393  )
394 
395  bDisk = pexConfig.Field(
396  doc="The column name for the semi minor axis length of the disk component.",
397  dtype=str,
398  default="b_d",
399  deprecated=(
400  "Use `disk_semimajor_col`, `disk_axis_ratio_col`, and `disk_pa_col`"
401  " to specify disk half-light ellipse."
402  )
403  )
404 
405  paDisk = pexConfig.Field(
406  doc="The column name for the PA of the disk component used in the fake source catalog.",
407  dtype=str,
408  default="pa_disk",
409  deprecated=(
410  "Use `disk_semimajor_col`, `disk_axis_ratio_col`, and `disk_pa_col`"
411  " to specify disk half-light ellipse."
412  )
413  )
414 
415  nDisk = pexConfig.Field(
416  doc="The column name for the sersic index of the disk component used in the fake source catalog.",
417  dtype=str,
418  default="disk_n",
419  deprecated="Use `disk_n` instead."
420  )
421 
422  bulgeHLR = pexConfig.Field(
423  doc="Column name for the bulge half light radius used in the fake source catalog.",
424  dtype=str,
425  default="BulgeHalfLightRadius",
426  deprecated=(
427  "Use `bulge_semimajor_col`, `bulge_axis_ratio_col`, and "
428  "`bulge_pa_col` to specify disk half-light ellipse."
429  )
430  )
431 
432  aBulge = pexConfig.Field(
433  doc="The column name for the semi major axis length of the bulge component.",
434  dtype=str,
435  default="a_b",
436  deprecated=(
437  "Use `bulge_semimajor_col`, `bulge_axis_ratio_col`, and "
438  "`bulge_pa_col` to specify disk half-light ellipse."
439  )
440  )
441 
442  bBulge = pexConfig.Field(
443  doc="The column name for the semi minor axis length of the bulge component used in the fake source "
444  "catalog.",
445  dtype=str,
446  default="b_b",
447  deprecated=(
448  "Use `bulge_semimajor_col`, `bulge_axis_ratio_col`, and "
449  "`bulge_pa_col` to specify disk half-light ellipse."
450  )
451  )
452 
453  paBulge = pexConfig.Field(
454  doc="The column name for the PA of the bulge component used in the fake source catalog.",
455  dtype=str,
456  default="pa_bulge",
457  deprecated=(
458  "Use `bulge_semimajor_col`, `bulge_axis_ratio_col`, and "
459  "`bulge_pa_col` to specify disk half-light ellipse."
460  )
461  )
462 
463  nBulge = pexConfig.Field(
464  doc="The column name for the sersic index of the bulge component used in the fake source catalog.",
465  dtype=str,
466  default="bulge_n",
467  deprecated="Use `bulge_n` instead."
468  )
469 
470  magVar = pexConfig.Field(
471  doc="The column name for the magnitude calculated taking variability into account. In the format "
472  "``filter name``magVar, e.g. imagVar for the magnitude in the i band.",
473  dtype=str,
474  default="%smagVar",
475  deprecated="Use `mag_col` instead."
476  )
477 
478  sourceSelectionColName = pexConfig.Field(
479  doc="The name of the column in the input fakes catalogue to be used to determine which sources to"
480  "add, default is none and when this is used all sources are added.",
481  dtype=str,
482  default="templateSource",
483  deprecated="Use `select_col` instead."
484  )
485 
486 
487 class InsertFakesTask(PipelineTask, CmdLineTask):
488  """Insert fake objects into images.
489 
490  Add fake stars and galaxies to the given image, read in through the dataRef. Galaxy parameters are read in
491  from the specified file and then modelled using galsim.
492 
493  `InsertFakesTask` has five functions that make images of the fake sources and then add them to the
494  image.
495 
496  `addPixCoords`
497  Use the WCS information to add the pixel coordinates of each source.
498  `mkFakeGalsimGalaxies`
499  Use Galsim to make fake double sersic galaxies for each set of galaxy parameters in the input file.
500  `mkFakeStars`
501  Use the PSF information from the image to make a fake star using the magnitude information from the
502  input file.
503  `cleanCat`
504  Remove rows of the input fake catalog which have half light radius, of either the bulge or the disk,
505  that are 0. Also removes rows that have Sersic index outside of galsim's allowed paramters. If
506  the config option sourceSelectionColName is set then this function limits the catalog of input fakes
507  to only those which are True in this column.
508  `addFakeSources`
509  Add the fake sources to the image.
510 
511  """
512 
513  _DefaultName = "insertFakes"
514  ConfigClass = InsertFakesConfig
515 
516  def runDataRef(self, dataRef):
517  """Read in/write out the required data products and add fake sources to the deepCoadd.
518 
519  Parameters
520  ----------
521  dataRef : `lsst.daf.persistence.butlerSubset.ButlerDataRef`
522  Data reference defining the image to have fakes added to it
523  Used to access the following data products:
524  deepCoadd
525  """
526 
527  infoStr = "Adding fakes to: tract: %d, patch: %s, filter: %s" % (dataRef.dataId["tract"],
528  dataRef.dataId["patch"],
529  dataRef.dataId["filter"])
530  self.log.info(infoStr)
531 
532  # To do: should it warn when asked to insert variable sources into the coadd
533 
534  if self.config.fakeType == "static":
535  fakeCat = dataRef.get("deepCoadd_fakeSourceCat").toDataFrame()
536  # To do: DM-16254, the read and write of the fake catalogs will be changed once the new pipeline
537  # task structure for ref cats is in place.
538  self.fakeSourceCatType = "deepCoadd_fakeSourceCat"
539  else:
540  fakeCat = Table.read(self.config.fakeType).to_pandas()
541 
542  coadd = dataRef.get("deepCoadd")
543  wcs = coadd.getWcs()
544  photoCalib = coadd.getPhotoCalib()
545 
546  imageWithFakes = self.run(fakeCat, coadd, wcs, photoCalib)
547 
548  dataRef.put(imageWithFakes.imageWithFakes, "fakes_deepCoadd")
549 
550  def runQuantum(self, butlerQC, inputRefs, outputRefs):
551  inputs = butlerQC.get(inputRefs)
552  inputs["wcs"] = inputs["image"].getWcs()
553  inputs["photoCalib"] = inputs["image"].getPhotoCalib()
554 
555  outputs = self.run(**inputs)
556  butlerQC.put(outputs, outputRefs)
557 
558  @classmethod
559  def _makeArgumentParser(cls):
560  parser = pipeBase.ArgumentParser(name=cls._DefaultName)
561  parser.add_id_argument(name="--id", datasetType="deepCoadd",
562  help="data IDs for the deepCoadd, e.g. --id tract=12345 patch=1,2 filter=r",
563  ContainerClass=ExistingCoaddDataIdContainer)
564  return parser
565 
566  def run(self, fakeCat, image, wcs, photoCalib):
567  """Add fake sources to an image.
568 
569  Parameters
570  ----------
571  fakeCat : `pandas.core.frame.DataFrame`
572  The catalog of fake sources to be input
573  image : `lsst.afw.image.exposure.exposure.ExposureF`
574  The image into which the fake sources should be added
575  wcs : `lsst.afw.geom.SkyWcs`
576  WCS to use to add fake sources
577  photoCalib : `lsst.afw.image.photoCalib.PhotoCalib`
578  Photometric calibration to be used to calibrate the fake sources
579 
580  Returns
581  -------
582  resultStruct : `lsst.pipe.base.struct.Struct`
583  contains : image : `lsst.afw.image.exposure.exposure.ExposureF`
584 
585  Notes
586  -----
587  Adds pixel coordinates for each source to the fakeCat and removes objects with bulge or disk half
588  light radius = 0 (if ``config.doCleanCat = True``).
589 
590  Adds the ``Fake`` mask plane to the image which is then set by `addFakeSources` to mark where fake
591  sources have been added. Uses the information in the ``fakeCat`` to make fake galaxies (using galsim)
592  and fake stars, using the PSF models from the PSF information for the image. These are then added to
593  the image and the image with fakes included returned.
594 
595  The galsim galaxies are made using a double sersic profile, one for the bulge and one for the disk,
596  this is then convolved with the PSF at that point.
597  """
598  # Attach overriding wcs and photoCalib to image, but retain originals
599  # so we can reset at the end.
600  origWcs = image.getWcs()
601  origPhotoCalib = image.getPhotoCalib()
602  image.setWcs(wcs)
603  image.setPhotoCalib(photoCalib)
604 
605  band = image.getFilterLabel().bandLabel
606  fakeCat = self._standardizeColumns(fakeCat, band)
607 
608  fakeCat = self.addPixCoords(fakeCat, image)
609  fakeCat = self.trimFakeCat(fakeCat, image)
610 
611  if len(fakeCat) > 0:
612  if isinstance(fakeCat[self.config.sourceType].iloc[0], str):
613  galCheckVal = "galaxy"
614  starCheckVal = "star"
615  elif isinstance(fakeCat[self.config.sourceType].iloc[0], bytes):
616  galCheckVal = b"galaxy"
617  starCheckVal = b"star"
618  elif isinstance(fakeCat[self.config.sourceType].iloc[0], (int, float)):
619  galCheckVal = 1
620  starCheckVal = 0
621  else:
622  raise TypeError("sourceType column does not have required type, should be str, bytes or int")
623 
624  if not self.config.insertImages:
625  if self.config.doCleanCat:
626  fakeCat = self.cleanCat(fakeCat, starCheckVal)
627 
628  generator = self._generateGSObjectsFromCatalog(image, fakeCat, galCheckVal, starCheckVal)
629  else:
630  generator = self._generateGSObjectsFromImages(image, fakeCat)
631  _add_fake_sources(image, generator, calibFluxRadius=self.config.calibFluxRadius, logger=self.log)
632  elif len(fakeCat) == 0 and self.config.doProcessAllDataIds:
633  self.log.warn("No fakes found for this dataRef; processing anyway.")
634  image.mask.addMaskPlane("FAKE")
635  else:
636  raise RuntimeError("No fakes found for this dataRef.")
637 
638  # restore original exposure WCS and photoCalib
639  image.setWcs(origWcs)
640  image.setPhotoCalib(origPhotoCalib)
641 
642  resultStruct = pipeBase.Struct(imageWithFakes=image)
643 
644  return resultStruct
645 
646  def _standardizeColumns(self, fakeCat, band):
647  """Use config variables to 'standardize' the expected columns and column
648  names in the input catalog.
649 
650  Parameters
651  ----------
652  fakeCat : `pandas.core.frame.DataFrame`
653  The catalog of fake sources to be input
654  band : `str`
655  Label for the current band being processed.
656 
657  Returns
658  -------
659  outCat : `pandas.core.frame.DataFrame`
660  The standardized catalog of fake sources
661  """
662  cfg = self.config
663  replace_dict = {}
664 
665  # Prefer new config variables over deprecated config variables.
666  # The following are fairly simple to handle as they're just column name
667  # changes.
668  for new_name, depr_name, std_name in [
669  (cfg.ra_col, cfg.raColName, 'ra'),
670  (cfg.dec_col, cfg.decColName, 'dec'),
671  (cfg.bulge_n_col, cfg.nBulge, 'bulge_n'),
672  (cfg.bulge_pa_col, cfg.paBulge, 'bulge_pa'),
673  (cfg.disk_n_col, cfg.nDisk, 'disk_n'),
674  (cfg.disk_pa_col, cfg.paDisk, 'disk_pa'),
675  (cfg.mag_col%band, cfg.magVar%band, 'mag'),
676  (cfg.select_col, cfg.sourceSelectionColName, 'select')
677  ]:
678  # Only standardize "select" column if doSubSelectSources is True
679  if not cfg.doSubSelectSources and std_name == 'select':
680  continue
681  if new_name in fakeCat.columns:
682  replace_dict[new_name] = std_name
683  elif depr_name in fakeCat.columns:
684  replace_dict[depr_name] = std_name
685  else:
686  raise ValueError(f"Could not determine column for {std_name}.")
687  fakeCat = fakeCat.rename(columns=replace_dict, copy=False)
688 
689  # Handling the half-light radius and axis-ratio are trickier, since we
690  # moved from expecting (HLR, a, b) to expecting (semimajor, axis_ratio).
691  # Just handle these manually.
692  if (
693  cfg.bulge_semimajor_col in fakeCat.columns
694  and cfg.bulge_axis_ratio_col in fakeCat.columns
695  ):
696  fakeCat = fakeCat.rename(
697  columns={
698  cfg.bulge_semimajor_col: 'bulge_semimajor',
699  cfg.bulge_axis_ratio_col: 'bulge_axis_ratio',
700  cfg.disk_semimajor_col: 'disk_semimajor',
701  cfg.disk_axis_ratio_col: 'disk_axis_ratio',
702  },
703  copy=False
704  )
705  elif (
706  cfg.bulgeHLR in fakeCat.columns
707  and cfg.aBulge in fakeCat.columns
708  and cfg.bBulge in fakeCat.columns
709  ):
710  fakeCat['bulge_axis_ratio'] = (
711  fakeCat[cfg.bBulge]/fakeCat[cfg.aBulge]
712  )
713  fakeCat['bulge_semimajor'] = (
714  fakeCat[cfg.bulgeHLR]/np.sqrt(fakeCat['bulge_axis_ratio'])
715  )
716  fakeCat['disk_axis_ratio'] = (
717  fakeCat[cfg.bDisk]/fakeCat[cfg.aDisk]
718  )
719  fakeCat['disk_semimajor'] = (
720  fakeCat[cfg.diskHLR]/np.sqrt(fakeCat['disk_axis_ratio'])
721  )
722  else:
723  raise ValueError(
724  "Could not determine columns for half-light radius and axis "
725  "ratio."
726  )
727 
728  # Process the bulge/disk flux ratio if possible.
729  if cfg.bulge_disk_flux_ratio_col in fakeCat.columns:
730  fakeCat = fakeCat.rename(
731  columns={
732  cfg.bulge_disk_flux_ratio_col: 'bulge_disk_flux_ratio'
733  },
734  copy=False
735  )
736  else:
737  fakeCat['bulge_disk_flux_ratio'] = 1.0
738 
739  return fakeCat
740 
741  def _generateGSObjectsFromCatalog(self, exposure, fakeCat, galCheckVal, starCheckVal):
742  """Process catalog to generate `galsim.GSObject` s.
743 
744  Parameters
745  ----------
746  exposure : `lsst.afw.image.exposure.exposure.ExposureF`
747  The exposure into which the fake sources should be added
748  fakeCat : `pandas.core.frame.DataFrame`
749  The catalog of fake sources to be input
750  galCheckVal : `str`, `bytes` or `int`
751  The value that is set in the sourceType column to specifiy an object is a galaxy.
752  starCheckVal : `str`, `bytes` or `int`
753  The value that is set in the sourceType column to specifiy an object is a star.
754 
755  Yields
756  ------
757  gsObjects : `generator`
758  A generator of tuples of `lsst.geom.SpherePoint` and `galsim.GSObject`.
759  """
760  wcs = exposure.getWcs()
761  photoCalib = exposure.getPhotoCalib()
762 
763  self.log.info(f"Making {len(fakeCat)} objects for insertion")
764 
765  for (index, row) in fakeCat.iterrows():
766  ra = row['ra']
767  dec = row['dec']
768  skyCoord = SpherePoint(ra, dec, radians)
769  xy = wcs.skyToPixel(skyCoord)
770 
771  try:
772  flux = photoCalib.magnitudeToInstFlux(row['mag'], xy)
773  except LogicError:
774  continue
775 
776  sourceType = row[self.config.sourceType]
777  if sourceType == galCheckVal:
778  # GalSim convention: HLR = sqrt(a * b) = a * sqrt(b / a)
779  bulge_gs_HLR = row['bulge_semimajor']*np.sqrt(row['bulge_axis_ratio'])
780  bulge = galsim.Sersic(n=row['bulge_n'], half_light_radius=bulge_gs_HLR)
781  bulge = bulge.shear(q=row['bulge_axis_ratio'], beta=((90 - row['bulge_pa'])*galsim.degrees))
782 
783  disk_gs_HLR = row['disk_semimajor']*np.sqrt(row['disk_axis_ratio'])
784  disk = galsim.Sersic(n=row['disk_n'], half_light_radius=disk_gs_HLR)
785  disk = disk.shear(q=row['disk_axis_ratio'], beta=((90 - row['disk_pa'])*galsim.degrees))
786 
787  gal = bulge*row['bulge_disk_flux_ratio'] + disk
788  gal = gal.withFlux(flux)
789 
790  yield skyCoord, gal
791  elif sourceType == starCheckVal:
792  star = galsim.DeltaFunction()
793  star = star.withFlux(flux)
794  yield skyCoord, star
795  else:
796  raise TypeError(f"Unknown sourceType {sourceType}")
797 
798  def _generateGSObjectsFromImages(self, exposure, fakeCat):
799  """Process catalog to generate `galsim.GSObject` s.
800 
801  Parameters
802  ----------
803  exposure : `lsst.afw.image.exposure.exposure.ExposureF`
804  The exposure into which the fake sources should be added
805  fakeCat : `pandas.core.frame.DataFrame`
806  The catalog of fake sources to be input
807 
808  Yields
809  ------
810  gsObjects : `generator`
811  A generator of tuples of `lsst.geom.SpherePoint` and `galsim.GSObject`.
812  """
813  band = exposure.getFilterLabel().bandLabel
814  wcs = exposure.getWcs()
815  photoCalib = exposure.getPhotoCalib()
816 
817  self.log.info(f"Processing {len(fakeCat)} fake images")
818 
819  for (index, row) in fakeCat.iterrows():
820  ra = row['ra']
821  dec = row['dec']
822  skyCoord = SpherePoint(ra, dec, radians)
823  xy = wcs.skyToPixel(skyCoord)
824 
825  try:
826  flux = photoCalib.magnitudeToInstFlux(row['mag'], xy)
827  except LogicError:
828  continue
829 
830  imFile = row[band+"imFilename"]
831  try:
832  imFile = imFile.decode("utf-8")
833  except AttributeError:
834  pass
835  imFile = imFile.strip()
836  im = galsim.fits.read(imFile, read_header=True)
837 
838  # GalSim will always attach a WCS to the image read in as above. If
839  # it can't find a WCS in the header, then it defaults to scale = 1.0
840  # arcsec / pix. So if that's the scale, then we need to check if it
841  # was explicitly set or if it's just the default. If it's just the
842  # default then we should override with the pixel scale of the target
843  # image.
844  if _isWCSGalsimDefault(im.wcs, im.header):
845  im.wcs = galsim.PixelScale(
846  wcs.getPixelScale().asArcseconds()
847  )
848 
849  obj = galsim.InterpolatedImage(im)
850  obj = obj.withFlux(flux)
851  yield skyCoord, obj
852 
853  def processImagesForInsertion(self, fakeCat, wcs, psf, photoCalib, band, pixelScale):
854  """Process images from files into the format needed for insertion.
855 
856  Parameters
857  ----------
858  fakeCat : `pandas.core.frame.DataFrame`
859  The catalog of fake sources to be input
860  wcs : `lsst.afw.geom.skyWcs.skyWcs.SkyWc`
861  WCS to use to add fake sources
862  psf : `lsst.meas.algorithms.coaddPsf.coaddPsf.CoaddPsf` or
863  `lsst.meas.extensions.psfex.psfexPsf.PsfexPsf`
864  The PSF information to use to make the PSF images
865  photoCalib : `lsst.afw.image.photoCalib.PhotoCalib`
866  Photometric calibration to be used to calibrate the fake sources
867  band : `str`
868  The filter band that the observation was taken in.
869  pixelScale : `float`
870  The pixel scale of the image the sources are to be added to.
871 
872  Returns
873  -------
874  galImages : `list`
875  A list of tuples of `lsst.afw.image.exposure.exposure.ExposureF` and
876  `lsst.geom.Point2D` of their locations.
877  For sources labelled as galaxy.
878  starImages : `list`
879  A list of tuples of `lsst.afw.image.exposure.exposure.ExposureF` and
880  `lsst.geom.Point2D` of their locations.
881  For sources labelled as star.
882 
883  Notes
884  -----
885  The input fakes catalog needs to contain the absolute path to the image in the
886  band that is being used to add images to. It also needs to have the R.A. and
887  declination of the fake source in radians and the sourceType of the object.
888  """
889  galImages = []
890  starImages = []
891 
892  self.log.info("Processing %d fake images" % len(fakeCat))
893 
894  for (imFile, sourceType, mag, x, y) in zip(fakeCat[band + "imFilename"].array,
895  fakeCat["sourceType"].array,
896  fakeCat['mag'].array,
897  fakeCat["x"].array, fakeCat["y"].array):
898 
899  im = afwImage.ImageF.readFits(imFile)
900 
901  xy = geom.Point2D(x, y)
902 
903  # We put these two PSF calculations within this same try block so that we catch cases
904  # where the object's position is outside of the image.
905  try:
906  correctedFlux = psf.computeApertureFlux(self.config.calibFluxRadius, xy)
907  psfKernel = psf.computeKernelImage(xy).getArray()
908  psfKernel /= correctedFlux
909 
910  except InvalidParameterError:
911  self.log.info("%s at %0.4f, %0.4f outside of image" % (sourceType, x, y))
912  continue
913 
914  psfIm = galsim.InterpolatedImage(galsim.Image(psfKernel), scale=pixelScale)
915  galsimIm = galsim.InterpolatedImage(galsim.Image(im.array), scale=pixelScale)
916  convIm = galsim.Convolve([galsimIm, psfIm])
917 
918  try:
919  outIm = convIm.drawImage(scale=pixelScale, method="real_space").array
920  except (galsim.errors.GalSimFFTSizeError, MemoryError):
921  continue
922 
923  imSum = np.sum(outIm)
924  divIm = outIm/imSum
925 
926  try:
927  flux = photoCalib.magnitudeToInstFlux(mag, xy)
928  except LogicError:
929  flux = 0
930 
931  imWithFlux = flux*divIm
932 
933  if sourceType == b"galaxy":
934  galImages.append((afwImage.ImageF(imWithFlux), xy))
935  if sourceType == b"star":
936  starImages.append((afwImage.ImageF(imWithFlux), xy))
937 
938  return galImages, starImages
939 
940  def addPixCoords(self, fakeCat, image):
941 
942  """Add pixel coordinates to the catalog of fakes.
943 
944  Parameters
945  ----------
946  fakeCat : `pandas.core.frame.DataFrame`
947  The catalog of fake sources to be input
948  image : `lsst.afw.image.exposure.exposure.ExposureF`
949  The image into which the fake sources should be added
950 
951  Returns
952  -------
953  fakeCat : `pandas.core.frame.DataFrame`
954  """
955  wcs = image.getWcs()
956  ras = fakeCat['ra'].values
957  decs = fakeCat['dec'].values
958  xs, ys = wcs.skyToPixelArray(ras, decs)
959  fakeCat["x"] = xs
960  fakeCat["y"] = ys
961 
962  return fakeCat
963 
964  def trimFakeCat(self, fakeCat, image):
965  """Trim the fake cat to about the size of the input image.
966 
967  `fakeCat` must be processed with addPixCoords before using this method.
968 
969  Parameters
970  ----------
971  fakeCat : `pandas.core.frame.DataFrame`
972  The catalog of fake sources to be input
973  image : `lsst.afw.image.exposure.exposure.ExposureF`
974  The image into which the fake sources should be added
975 
976  Returns
977  -------
978  fakeCat : `pandas.core.frame.DataFrame`
979  The original fakeCat trimmed to the area of the image
980  """
981 
982  bbox = Box2D(image.getBBox()).dilatedBy(self.config.trimBuffer)
983  xs = fakeCat["x"].values
984  ys = fakeCat["y"].values
985 
986  isContained = xs >= bbox.minX
987  isContained &= xs <= bbox.maxX
988  isContained &= ys >= bbox.minY
989  isContained &= ys <= bbox.maxY
990 
991  return fakeCat[isContained]
992 
993  def mkFakeGalsimGalaxies(self, fakeCat, band, photoCalib, pixelScale, psf, image):
994  """Make images of fake galaxies using GalSim.
995 
996  Parameters
997  ----------
998  band : `str`
999  pixelScale : `float`
1000  psf : `lsst.meas.extensions.psfex.psfexPsf.PsfexPsf`
1001  The PSF information to use to make the PSF images
1002  fakeCat : `pandas.core.frame.DataFrame`
1003  The catalog of fake sources to be input
1004  photoCalib : `lsst.afw.image.photoCalib.PhotoCalib`
1005  Photometric calibration to be used to calibrate the fake sources
1006 
1007  Yields
1008  -------
1009  galImages : `generator`
1010  A generator of tuples of `lsst.afw.image.exposure.exposure.ExposureF` and
1011  `lsst.geom.Point2D` of their locations.
1012 
1013  Notes
1014  -----
1015 
1016  Fake galaxies are made by combining two sersic profiles, one for the bulge and one for the disk. Each
1017  component has an individual sersic index (n), a, b and position angle (PA). The combined profile is
1018  then convolved with the PSF at the specified x, y position on the image.
1019 
1020  The names of the columns in the ``fakeCat`` are configurable and are the column names from the
1021  University of Washington simulations database as default. For more information see the doc strings
1022  attached to the config options.
1023 
1024  See mkFakeStars doc string for an explanation of calibration to instrumental flux.
1025  """
1026 
1027  self.log.info("Making %d fake galaxy images" % len(fakeCat))
1028 
1029  for (index, row) in fakeCat.iterrows():
1030  xy = geom.Point2D(row["x"], row["y"])
1031 
1032  # We put these two PSF calculations within this same try block so that we catch cases
1033  # where the object's position is outside of the image.
1034  try:
1035  correctedFlux = psf.computeApertureFlux(self.config.calibFluxRadius, xy)
1036  psfKernel = psf.computeKernelImage(xy).getArray()
1037  psfKernel /= correctedFlux
1038 
1039  except InvalidParameterError:
1040  self.log.info("Galaxy at %0.4f, %0.4f outside of image" % (row["x"], row["y"]))
1041  continue
1042 
1043  try:
1044  flux = photoCalib.magnitudeToInstFlux(row['mag'], xy)
1045  except LogicError:
1046  flux = 0
1047 
1048  # GalSim convention: HLR = sqrt(a * b) = a * sqrt(b / a)
1049  bulge_gs_HLR = row['bulge_semimajor']*np.sqrt(row['bulge_axis_ratio'])
1050  bulge = galsim.Sersic(n=row['bulge_n'], half_light_radius=bulge_gs_HLR)
1051  bulge = bulge.shear(q=row['bulge_axis_ratio'], beta=((90 - row['bulge_pa'])*galsim.degrees))
1052 
1053  disk_gs_HLR = row['disk_semimajor']*np.sqrt(row['disk_axis_ratio'])
1054  disk = galsim.Sersic(n=row['disk_n'], half_light_radius=disk_gs_HLR)
1055  disk = disk.shear(q=row['disk_axis_ratio'], beta=((90 - row['disk_pa'])*galsim.degrees))
1056 
1057  gal = bulge*row['bulge_disk_flux_ratio'] + disk
1058  gal = gal.withFlux(flux)
1059 
1060  psfIm = galsim.InterpolatedImage(galsim.Image(psfKernel), scale=pixelScale)
1061  gal = galsim.Convolve([gal, psfIm])
1062  try:
1063  galIm = gal.drawImage(scale=pixelScale, method="real_space").array
1064  except (galsim.errors.GalSimFFTSizeError, MemoryError):
1065  continue
1066 
1067  yield (afwImage.ImageF(galIm), xy)
1068 
1069  def mkFakeStars(self, fakeCat, band, photoCalib, psf, image):
1070 
1071  """Make fake stars based off the properties in the fakeCat.
1072 
1073  Parameters
1074  ----------
1075  band : `str`
1076  psf : `lsst.meas.extensions.psfex.psfexPsf.PsfexPsf`
1077  The PSF information to use to make the PSF images
1078  fakeCat : `pandas.core.frame.DataFrame`
1079  The catalog of fake sources to be input
1080  image : `lsst.afw.image.exposure.exposure.ExposureF`
1081  The image into which the fake sources should be added
1082  photoCalib : `lsst.afw.image.photoCalib.PhotoCalib`
1083  Photometric calibration to be used to calibrate the fake sources
1084 
1085  Yields
1086  -------
1087  starImages : `generator`
1088  A generator of tuples of `lsst.afw.image.ImageF` of fake stars and
1089  `lsst.geom.Point2D` of their locations.
1090 
1091  Notes
1092  -----
1093  To take a given magnitude and translate to the number of counts in the image
1094  we use photoCalib.magnitudeToInstFlux, which returns the instrumental flux for the
1095  given calibration radius used in the photometric calibration step.
1096  Thus `calibFluxRadius` should be set to this same radius so that we can normalize
1097  the PSF model to the correct instrumental flux within calibFluxRadius.
1098  """
1099 
1100  self.log.info("Making %d fake star images" % len(fakeCat))
1101 
1102  for (index, row) in fakeCat.iterrows():
1103  xy = geom.Point2D(row["x"], row["y"])
1104 
1105  # We put these two PSF calculations within this same try block so that we catch cases
1106  # where the object's position is outside of the image.
1107  try:
1108  correctedFlux = psf.computeApertureFlux(self.config.calibFluxRadius, xy)
1109  starIm = psf.computeImage(xy)
1110  starIm /= correctedFlux
1111 
1112  except InvalidParameterError:
1113  self.log.info("Star at %0.4f, %0.4f outside of image" % (row["x"], row["y"]))
1114  continue
1115 
1116  try:
1117  flux = photoCalib.magnitudeToInstFlux(row['mag'], xy)
1118  except LogicError:
1119  flux = 0
1120 
1121  starIm *= flux
1122  yield ((starIm.convertF(), xy))
1123 
1124  def cleanCat(self, fakeCat, starCheckVal):
1125  """Remove rows from the fakes catalog which have HLR = 0 for either the buldge or disk component,
1126  also remove galaxies that have Sersic index outside the galsim min and max
1127  allowed (0.3 <= n <= 6.2).
1128 
1129  Parameters
1130  ----------
1131  fakeCat : `pandas.core.frame.DataFrame`
1132  The catalog of fake sources to be input
1133  starCheckVal : `str`, `bytes` or `int`
1134  The value that is set in the sourceType column to specifiy an object is a star.
1135 
1136  Returns
1137  -------
1138  fakeCat : `pandas.core.frame.DataFrame`
1139  The input catalog of fake sources but with the bad objects removed
1140  """
1141 
1142  rowsToKeep = (((fakeCat['bulge_semimajor'] != 0.0) & (fakeCat['disk_semimajor'] != 0.0))
1143  | (fakeCat[self.config.sourceType] == starCheckVal))
1144  numRowsNotUsed = len(fakeCat) - len(np.where(rowsToKeep)[0])
1145  self.log.info("Removing %d rows with HLR = 0 for either the bulge or disk" % numRowsNotUsed)
1146  fakeCat = fakeCat[rowsToKeep]
1147 
1148  minN = galsim.Sersic._minimum_n
1149  maxN = galsim.Sersic._maximum_n
1150  rowsWithGoodSersic = (((fakeCat['bulge_n'] >= minN) & (fakeCat['bulge_n'] <= maxN)
1151  & (fakeCat['disk_n'] >= minN) & (fakeCat['disk_n'] <= maxN))
1152  | (fakeCat[self.config.sourceType] == starCheckVal))
1153  numRowsNotUsed = len(fakeCat) - len(np.where(rowsWithGoodSersic)[0])
1154  self.log.info("Removing %d rows of galaxies with nBulge or nDisk outside of %0.2f <= n <= %0.2f" %
1155  (numRowsNotUsed, minN, maxN))
1156  fakeCat = fakeCat[rowsWithGoodSersic]
1157 
1158  if self.config.doSubSelectSources:
1159  numRowsNotUsed = len(fakeCat) - len(fakeCat['select'])
1160  self.log.info("Removing %d rows which were not designated as template sources" % numRowsNotUsed)
1161  fakeCat = fakeCat[fakeCat['select']]
1162 
1163  return fakeCat
1164 
1165  def addFakeSources(self, image, fakeImages, sourceType):
1166  """Add the fake sources to the given image
1167 
1168  Parameters
1169  ----------
1170  image : `lsst.afw.image.exposure.exposure.ExposureF`
1171  The image into which the fake sources should be added
1172  fakeImages : `typing.Iterator` [`tuple` ['lsst.afw.image.ImageF`, `lsst.geom.Point2d`]]
1173  An iterator of tuples that contains (or generates) images of fake sources,
1174  and the locations they are to be inserted at.
1175  sourceType : `str`
1176  The type (star/galaxy) of fake sources input
1177 
1178  Returns
1179  -------
1180  image : `lsst.afw.image.exposure.exposure.ExposureF`
1181 
1182  Notes
1183  -----
1184  Uses the x, y information in the ``fakeCat`` to position an image of the fake interpolated onto the
1185  pixel grid of the image. Sets the ``FAKE`` mask plane for the pixels added with the fake source.
1186  """
1187 
1188  imageBBox = image.getBBox()
1189  imageMI = image.maskedImage
1190 
1191  for (fakeImage, xy) in fakeImages:
1192  X0 = xy.getX() - fakeImage.getWidth()/2 + 0.5
1193  Y0 = xy.getY() - fakeImage.getHeight()/2 + 0.5
1194  self.log.debug("Adding fake source at %d, %d" % (xy.getX(), xy.getY()))
1195  if sourceType == "galaxy":
1196  interpFakeImage = afwMath.offsetImage(fakeImage, X0, Y0, "lanczos3")
1197  else:
1198  interpFakeImage = fakeImage
1199 
1200  interpFakeImBBox = interpFakeImage.getBBox()
1201  interpFakeImBBox.clip(imageBBox)
1202 
1203  if interpFakeImBBox.getArea() > 0:
1204  imageMIView = imageMI[interpFakeImBBox]
1205  clippedFakeImage = interpFakeImage[interpFakeImBBox]
1206  clippedFakeImageMI = afwImage.MaskedImageF(clippedFakeImage)
1207  clippedFakeImageMI.mask.set(self.bitmask)
1208  imageMIView += clippedFakeImageMI
1209 
1210  return image
1211 
1212  def _getMetadataName(self):
1213  """Disable metadata writing"""
1214  return None
def run(self, skyInfo, tempExpRefList, imageScalerList, weightList, altMaskList=None, mask=None, supplementaryData=None)
def addPixCoords(self, fakeCat, image)
Definition: insertFakes.py:940
def mkFakeStars(self, fakeCat, band, photoCalib, psf, image)
def mkFakeGalsimGalaxies(self, fakeCat, band, photoCalib, pixelScale, psf, image)
Definition: insertFakes.py:993
def cleanCat(self, fakeCat, starCheckVal)
def processImagesForInsertion(self, fakeCat, wcs, psf, photoCalib, band, pixelScale)
Definition: insertFakes.py:853
def trimFakeCat(self, fakeCat, image)
Definition: insertFakes.py:964
def addFakeSources(self, image, fakeImages, sourceType)