lsst.fgcmcal  21.0.0-11-gd78879e+d13a45ff19
fgcmOutputProducts.py
Go to the documentation of this file.
1 # See COPYRIGHT file at the top of the source tree.
2 #
3 # This file is part of fgcmcal.
4 #
5 # Developed for the LSST Data Management System.
6 # This product includes software developed by the LSST Project
7 # (https://www.lsst.org).
8 # See the COPYRIGHT file at the top-level directory of this distribution
9 # for details of code ownership.
10 #
11 # This program is free software: you can redistribute it and/or modify
12 # it under the terms of the GNU General Public License as published by
13 # the Free Software Foundation, either version 3 of the License, or
14 # (at your option) any later version.
15 #
16 # This program is distributed in the hope that it will be useful,
17 # but WITHOUT ANY WARRANTY; without even the implied warranty of
18 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19 # GNU General Public License for more details.
20 #
21 # You should have received a copy of the GNU General Public License
22 # along with this program. If not, see <https://www.gnu.org/licenses/>.
23 """Make the final fgcmcal output products.
24 
25 This task takes the final output from fgcmFitCycle and produces the following
26 outputs for use in the DM stack: the FGCM standard stars in a reference
27 catalog format; the model atmospheres in "transmission_atmosphere_fgcm"
28 format; and the zeropoints in "fgcm_photoCalib" format. Optionally, the
29 task can transfer the 'absolute' calibration from a reference catalog
30 to put the fgcm standard stars in units of Jansky. This is accomplished
31 by matching stars in a sample of healpix pixels, and applying the median
32 offset per band.
33 """
34 import sys
35 import traceback
36 import copy
37 
38 import numpy as np
39 import healpy as hp
40 import esutil
41 from astropy import units
42 
43 import lsst.daf.base as dafBase
44 import lsst.pex.config as pexConfig
45 import lsst.pipe.base as pipeBase
46 from lsst.pipe.base import connectionTypes
47 from lsst.afw.image import TransmissionCurve
48 from lsst.meas.algorithms import LoadIndexedReferenceObjectsTask
49 from lsst.meas.algorithms import ReferenceObjectLoader
50 from lsst.pipe.tasks.photoCal import PhotoCalTask
51 import lsst.geom
52 import lsst.afw.image as afwImage
53 import lsst.afw.math as afwMath
54 import lsst.afw.table as afwTable
55 from lsst.meas.algorithms import IndexerRegistry
56 from lsst.meas.algorithms import DatasetConfig
57 from lsst.meas.algorithms.ingestIndexReferenceTask import addRefCatMetadata
58 
59 from .utilities import computeApproxPixelAreaFields
60 from .utilities import lookupStaticCalibrations
61 
62 import fgcm
63 
64 __all__ = ['FgcmOutputProductsConfig', 'FgcmOutputProductsTask', 'FgcmOutputProductsRunner']
65 
66 
67 class FgcmOutputProductsConnections(pipeBase.PipelineTaskConnections,
68  dimensions=("instrument",),
69  defaultTemplates={"cycleNumber": "0"}):
70  camera = connectionTypes.PrerequisiteInput(
71  doc="Camera instrument",
72  name="camera",
73  storageClass="Camera",
74  dimensions=("instrument",),
75  lookupFunction=lookupStaticCalibrations,
76  isCalibration=True,
77  )
78 
79  fgcmLookUpTable = connectionTypes.PrerequisiteInput(
80  doc=("Atmosphere + instrument look-up-table for FGCM throughput and "
81  "chromatic corrections."),
82  name="fgcmLookUpTable",
83  storageClass="Catalog",
84  dimensions=("instrument",),
85  deferLoad=True,
86  )
87 
88  fgcmVisitCatalog = connectionTypes.PrerequisiteInput(
89  doc="Catalog of visit information for fgcm",
90  name="fgcmVisitCatalog",
91  storageClass="Catalog",
92  dimensions=("instrument",),
93  deferLoad=True,
94  )
95 
96  fgcmStandardStars = connectionTypes.PrerequisiteInput(
97  doc="Catalog of standard star data from fgcm fit",
98  name="fgcmStandardStars{cycleNumber}",
99  storageClass="SimpleCatalog",
100  dimensions=("instrument",),
101  deferLoad=True,
102  )
103 
104  fgcmZeropoints = connectionTypes.PrerequisiteInput(
105  doc="Catalog of zeropoints from fgcm fit",
106  name="fgcmZeropoints{cycleNumber}",
107  storageClass="Catalog",
108  dimensions=("instrument",),
109  deferLoad=True,
110  )
111 
112  fgcmAtmosphereParameters = connectionTypes.PrerequisiteInput(
113  doc="Catalog of atmosphere parameters from fgcm fit",
114  name="fgcmAtmosphereParameters{cycleNumber}",
115  storageClass="Catalog",
116  dimensions=("instrument",),
117  deferLoad=True,
118  )
119 
120  refCat = connectionTypes.PrerequisiteInput(
121  doc="Reference catalog to use for photometric calibration",
122  name="cal_ref_cat",
123  storageClass="SimpleCatalog",
124  dimensions=("skypix",),
125  deferLoad=True,
126  multiple=True,
127  )
128 
129  fgcmBuildStarsTableConfig = connectionTypes.PrerequisiteInput(
130  doc="Config used to build FGCM input stars",
131  name="fgcmBuildStarsTable_config",
132  storageClass="Config",
133  )
134 
135  fgcmPhotoCalib = connectionTypes.Output(
136  doc=("Per-visit photometric calibrations derived from fgcm calibration. "
137  "These catalogs use detector id for the id and are sorted for "
138  "fast lookups of a detector."),
139  name="fgcmPhotoCalibCatalog",
140  storageClass="ExposureCatalog",
141  dimensions=("instrument", "visit",),
142  multiple=True,
143  )
144 
145  fgcmTransmissionAtmosphere = connectionTypes.Output(
146  doc="Per-visit atmosphere transmission files produced from fgcm calibration",
147  name="transmission_atmosphere_fgcm",
148  storageClass="TransmissionCurve",
149  dimensions=("instrument",
150  "visit",),
151  multiple=True,
152  )
153 
154  fgcmOffsets = connectionTypes.Output(
155  doc="Per-band offsets computed from doReferenceCalibration",
156  name="fgcmReferenceCalibrationOffsets",
157  storageClass="Catalog",
158  dimensions=("instrument",),
159  multiple=False,
160  )
161 
162  def __init__(self, *, config=None):
163  super().__init__(config=config)
164 
165  if str(int(config.connections.cycleNumber)) != config.connections.cycleNumber:
166  raise ValueError("cycleNumber must be of integer format")
167  if config.connections.refCat != config.refObjLoader.ref_dataset_name:
168  raise ValueError("connections.refCat must be the same as refObjLoader.ref_dataset_name")
169 
170  if config.doRefcatOutput:
171  raise ValueError("FgcmOutputProductsTask (Gen3) does not support doRefcatOutput")
172 
173  if not config.doReferenceCalibration:
174  self.prerequisiteInputs.remove("refCat")
175  if not config.doAtmosphereOutput:
176  self.prerequisiteInputs.remove("fgcmAtmosphereParameters")
177  if not config.doZeropointOutput:
178  self.prerequisiteInputs.remove("fgcmZeropoints")
179  if not config.doReferenceCalibration:
180  self.outputs.remove("fgcmOffsets")
181 
182 
183 class FgcmOutputProductsConfig(pipeBase.PipelineTaskConfig,
184  pipelineConnections=FgcmOutputProductsConnections):
185  """Config for FgcmOutputProductsTask"""
186 
187  cycleNumber = pexConfig.Field(
188  doc="Final fit cycle from FGCM fit",
189  dtype=int,
190  default=None,
191  )
192 
193  # The following fields refer to calibrating from a reference
194  # catalog, but in the future this might need to be expanded
195  doReferenceCalibration = pexConfig.Field(
196  doc=("Transfer 'absolute' calibration from reference catalog? "
197  "This afterburner step is unnecessary if reference stars "
198  "were used in the full fit in FgcmFitCycleTask."),
199  dtype=bool,
200  default=False,
201  )
202  doRefcatOutput = pexConfig.Field(
203  doc="Output standard stars in reference catalog format",
204  dtype=bool,
205  default=True,
206  )
207  doAtmosphereOutput = pexConfig.Field(
208  doc="Output atmospheres in transmission_atmosphere_fgcm format",
209  dtype=bool,
210  default=True,
211  )
212  doZeropointOutput = pexConfig.Field(
213  doc="Output zeropoints in fgcm_photoCalib format",
214  dtype=bool,
215  default=True,
216  )
217  doComposeWcsJacobian = pexConfig.Field(
218  doc="Compose Jacobian of WCS with fgcm calibration for output photoCalib?",
219  dtype=bool,
220  default=True,
221  )
222  doApplyMeanChromaticCorrection = pexConfig.Field(
223  doc="Apply the mean chromatic correction to the zeropoints?",
224  dtype=bool,
225  default=True,
226  )
227  refObjLoader = pexConfig.ConfigurableField(
228  target=LoadIndexedReferenceObjectsTask,
229  doc="reference object loader for 'absolute' photometric calibration",
230  )
231  photoCal = pexConfig.ConfigurableField(
232  target=PhotoCalTask,
233  doc="task to perform 'absolute' calibration",
234  )
235  referencePixelizationNside = pexConfig.Field(
236  doc="Healpix nside to pixelize catalog to compare to reference catalog",
237  dtype=int,
238  default=64,
239  )
240  referencePixelizationMinStars = pexConfig.Field(
241  doc=("Minimum number of stars per healpix pixel to select for comparison"
242  "to the specified reference catalog"),
243  dtype=int,
244  default=200,
245  )
246  referenceMinMatch = pexConfig.Field(
247  doc="Minimum number of stars matched to reference catalog to be used in statistics",
248  dtype=int,
249  default=50,
250  )
251  referencePixelizationNPixels = pexConfig.Field(
252  doc=("Number of healpix pixels to sample to do comparison. "
253  "Doing too many will take a long time and not yield any more "
254  "precise results because the final number is the median offset "
255  "(per band) from the set of pixels."),
256  dtype=int,
257  default=100,
258  )
259  datasetConfig = pexConfig.ConfigField(
260  dtype=DatasetConfig,
261  doc="Configuration for writing/reading ingested catalog",
262  )
263 
264  def setDefaults(self):
265  pexConfig.Config.setDefaults(self)
266 
267  # In order to transfer the "absolute" calibration from a reference
268  # catalog to the relatively calibrated FGCM standard stars (one number
269  # per band), we use the PhotoCalTask to match stars in a sample of healpix
270  # pixels. These basic settings ensure that only well-measured, good stars
271  # from the source and reference catalogs are used for the matching.
272 
273  # applyColorTerms needs to be False if doReferenceCalibration is False,
274  # as is the new default after DM-16702
275  self.photoCal.applyColorTerms = False
276  self.photoCal.fluxField = 'instFlux'
277  self.photoCal.magErrFloor = 0.003
278  self.photoCal.match.referenceSelection.doSignalToNoise = True
279  self.photoCal.match.referenceSelection.signalToNoise.minimum = 10.0
280  self.photoCal.match.sourceSelection.doSignalToNoise = True
281  self.photoCal.match.sourceSelection.signalToNoise.minimum = 10.0
282  self.photoCal.match.sourceSelection.signalToNoise.fluxField = 'instFlux'
283  self.photoCal.match.sourceSelection.signalToNoise.errField = 'instFluxErr'
284  self.photoCal.match.sourceSelection.doFlags = True
285  self.photoCal.match.sourceSelection.flags.good = []
286  self.photoCal.match.sourceSelection.flags.bad = ['flag_badStar']
287  self.photoCal.match.sourceSelection.doUnresolved = False
288  self.datasetConfig.ref_dataset_name = 'fgcm_stars'
289  self.datasetConfig.format_version = 1
290 
291  def validate(self):
292  super().validate()
293 
294  # Force the connections to conform with cycleNumber
295  self.connections.cycleNumber = str(self.cycleNumber)
296 
297 
298 class FgcmOutputProductsRunner(pipeBase.ButlerInitializedTaskRunner):
299  """Subclass of TaskRunner for fgcmOutputProductsTask
300 
301  fgcmOutputProductsTask.run() takes one argument, the butler, and
302  does not run on any data in the repository.
303  This runner does not use any parallelization.
304  """
305 
306  @staticmethod
307  def getTargetList(parsedCmd):
308  """
309  Return a list with one element, the butler.
310  """
311  return [parsedCmd.butler]
312 
313  def __call__(self, butler):
314  """
315  Parameters
316  ----------
317  butler: `lsst.daf.persistence.Butler`
318 
319  Returns
320  -------
321  exitStatus: `list` with `pipeBase.Struct`
322  exitStatus (0: success; 1: failure)
323  if self.doReturnResults also
324  results (`np.array` with absolute zeropoint offsets)
325  """
326  task = self.TaskClass(butler=butler, config=self.config, log=self.log)
327 
328  exitStatus = 0
329  if self.doRaise:
330  results = task.runDataRef(butler)
331  else:
332  try:
333  results = task.runDataRef(butler)
334  except Exception as e:
335  exitStatus = 1
336  task.log.fatal("Failed: %s" % e)
337  if not isinstance(e, pipeBase.TaskError):
338  traceback.print_exc(file=sys.stderr)
339 
340  task.writeMetadata(butler)
341 
342  if self.doReturnResults:
343  # The results here are the zeropoint offsets for each band
344  return [pipeBase.Struct(exitStatus=exitStatus,
345  results=results)]
346  else:
347  return [pipeBase.Struct(exitStatus=exitStatus)]
348 
349  def run(self, parsedCmd):
350  """
351  Run the task, with no multiprocessing
352 
353  Parameters
354  ----------
355  parsedCmd: `lsst.pipe.base.ArgumentParser` parsed command line
356  """
357 
358  resultList = []
359 
360  if self.precall(parsedCmd):
361  targetList = self.getTargetList(parsedCmd)
362  # make sure that we only get 1
363  resultList = self(targetList[0])
364 
365  return resultList
366 
367 
368 class FgcmOutputProductsTask(pipeBase.PipelineTask, pipeBase.CmdLineTask):
369  """
370  Output products from FGCM global calibration.
371  """
372 
373  ConfigClass = FgcmOutputProductsConfig
374  RunnerClass = FgcmOutputProductsRunner
375  _DefaultName = "fgcmOutputProducts"
376 
377  def __init__(self, butler=None, **kwargs):
378  super().__init__(**kwargs)
379 
380  # no saving of metadata for now
381  def _getMetadataName(self):
382  return None
383 
384  def runQuantum(self, butlerQC, inputRefs, outputRefs):
385  dataRefDict = {}
386  dataRefDict['camera'] = butlerQC.get(inputRefs.camera)
387  dataRefDict['fgcmLookUpTable'] = butlerQC.get(inputRefs.fgcmLookUpTable)
388  dataRefDict['fgcmVisitCatalog'] = butlerQC.get(inputRefs.fgcmVisitCatalog)
389  dataRefDict['fgcmStandardStars'] = butlerQC.get(inputRefs.fgcmStandardStars)
390 
391  if self.config.doZeropointOutput:
392  dataRefDict['fgcmZeropoints'] = butlerQC.get(inputRefs.fgcmZeropoints)
393  photoCalibRefDict = {photoCalibRef.dataId.byName()['visit']:
394  photoCalibRef for photoCalibRef in outputRefs.fgcmPhotoCalib}
395 
396  if self.config.doAtmosphereOutput:
397  dataRefDict['fgcmAtmosphereParameters'] = butlerQC.get(inputRefs.fgcmAtmosphereParameters)
398  atmRefDict = {atmRef.dataId.byName()['visit']: atmRef for
399  atmRef in outputRefs.fgcmTransmissionAtmosphere}
400 
401  if self.config.doReferenceCalibration:
402  refConfig = self.config.refObjLoader
403  self.refObjLoader = ReferenceObjectLoader(dataIds=[ref.datasetRef.dataId
404  for ref in inputRefs.refCat],
405  refCats=butlerQC.get(inputRefs.refCat),
406  config=refConfig,
407  log=self.log)
408  else:
409  self.refObjLoader = None
410 
411  dataRefDict['fgcmBuildStarsTableConfig'] = butlerQC.get(inputRefs.fgcmBuildStarsTableConfig)
412 
413  fgcmBuildStarsConfig = butlerQC.get(inputRefs.fgcmBuildStarsTableConfig)
414  physicalFilterMap = fgcmBuildStarsConfig.physicalFilterMap
415 
416  if self.config.doComposeWcsJacobian and not fgcmBuildStarsConfig.doApplyWcsJacobian:
417  raise RuntimeError("Cannot compose the WCS jacobian if it hasn't been applied "
418  "in fgcmBuildStarsTask.")
419  if not self.config.doComposeWcsJacobian and fgcmBuildStarsConfig.doApplyWcsJacobian:
420  self.log.warn("Jacobian was applied in build-stars but doComposeWcsJacobian is not set.")
421 
422  struct = self.run(dataRefDict, physicalFilterMap, returnCatalogs=True)
423 
424  # Output the photoCalib exposure catalogs
425  if struct.photoCalibCatalogs is not None:
426  self.log.info("Outputting photoCalib catalogs.")
427  for visit, expCatalog in struct.photoCalibCatalogs:
428  butlerQC.put(expCatalog, photoCalibRefDict[visit])
429  self.log.info("Done outputting photoCalib catalogs.")
430 
431  # Output the atmospheres
432  if struct.atmospheres is not None:
433  self.log.info("Outputting atmosphere transmission files.")
434  for visit, atm in struct.atmospheres:
435  butlerQC.put(atm, atmRefDict[visit])
436  self.log.info("Done outputting atmosphere files.")
437 
438  if self.config.doReferenceCalibration:
439  # Turn offset into simple catalog for persistence if necessary
440  schema = afwTable.Schema()
441  schema.addField('offset', type=np.float64,
442  doc="Post-process calibration offset (mag)")
443  offsetCat = afwTable.BaseCatalog(schema)
444  offsetCat.resize(len(struct.offsets))
445  offsetCat['offset'][:] = struct.offsets
446 
447  butlerQC.put(offsetCat, outputRefs.fgcmOffsets)
448 
449  return
450 
451  @pipeBase.timeMethod
452  def runDataRef(self, butler):
453  """
454  Make FGCM output products for use in the stack
455 
456  Parameters
457  ----------
458  butler: `lsst.daf.persistence.Butler`
459  cycleNumber: `int`
460  Final fit cycle number, override config.
461 
462  Returns
463  -------
464  offsets: `lsst.pipe.base.Struct`
465  A structure with array of zeropoint offsets
466 
467  Raises
468  ------
469  RuntimeError:
470  Raised if any one of the following is true:
471 
472  - butler cannot find "fgcmBuildStars_config" or
473  "fgcmBuildStarsTable_config".
474  - butler cannot find "fgcmFitCycle_config".
475  - "fgcmFitCycle_config" does not refer to
476  `self.config.cycleNumber`.
477  - butler cannot find "fgcmAtmosphereParameters" and
478  `self.config.doAtmosphereOutput` is `True`.
479  - butler cannot find "fgcmStandardStars" and
480  `self.config.doReferenceCalibration` is `True` or
481  `self.config.doRefcatOutput` is `True`.
482  - butler cannot find "fgcmZeropoints" and
483  `self.config.doZeropointOutput` is `True`.
484  """
485  if self.config.doReferenceCalibration:
486  # We need the ref obj loader to get the flux field
487  self.makeSubtask("refObjLoader", butler=butler)
488 
489  # Check to make sure that the fgcmBuildStars config exists, to retrieve
490  # the visit and ccd dataset tags
491  if not butler.datasetExists('fgcmBuildStarsTable_config') and \
492  not butler.datasetExists('fgcmBuildStars_config'):
493  raise RuntimeError("Cannot find fgcmBuildStarsTable_config or fgcmBuildStars_config, "
494  "which is prereq for fgcmOutputProducts")
495 
496  if butler.datasetExists('fgcmBuildStarsTable_config'):
497  fgcmBuildStarsConfig = butler.get('fgcmBuildStarsTable_config')
498  else:
499  fgcmBuildStarsConfig = butler.get('fgcmBuildStars_config')
500  visitDataRefName = fgcmBuildStarsConfig.visitDataRefName
501  ccdDataRefName = fgcmBuildStarsConfig.ccdDataRefName
502  physicalFilterMap = fgcmBuildStarsConfig.physicalFilterMap
503 
504  if self.config.doComposeWcsJacobian and not fgcmBuildStarsConfig.doApplyWcsJacobian:
505  raise RuntimeError("Cannot compose the WCS jacobian if it hasn't been applied "
506  "in fgcmBuildStarsTask.")
507 
508  if not self.config.doComposeWcsJacobian and fgcmBuildStarsConfig.doApplyWcsJacobian:
509  self.log.warn("Jacobian was applied in build-stars but doComposeWcsJacobian is not set.")
510 
511  # And make sure that the atmosphere was output properly
512  if (self.config.doAtmosphereOutput
513  and not butler.datasetExists('fgcmAtmosphereParameters', fgcmcycle=self.config.cycleNumber)):
514  raise RuntimeError(f"Atmosphere parameters are missing for cycle {self.config.cycleNumber}.")
515 
516  if not butler.datasetExists('fgcmStandardStars',
517  fgcmcycle=self.config.cycleNumber):
518  raise RuntimeError("Standard stars are missing for cycle %d." %
519  (self.config.cycleNumber))
520 
521  if (self.config.doZeropointOutput
522  and (not butler.datasetExists('fgcmZeropoints', fgcmcycle=self.config.cycleNumber))):
523  raise RuntimeError("Zeropoints are missing for cycle %d." %
524  (self.config.cycleNumber))
525 
526  dataRefDict = {}
527  # This is the _actual_ camera
528  dataRefDict['camera'] = butler.get('camera')
529  dataRefDict['fgcmLookUpTable'] = butler.dataRef('fgcmLookUpTable')
530  dataRefDict['fgcmVisitCatalog'] = butler.dataRef('fgcmVisitCatalog')
531  dataRefDict['fgcmStandardStars'] = butler.dataRef('fgcmStandardStars',
532  fgcmcycle=self.config.cycleNumber)
533 
534  if self.config.doZeropointOutput:
535  dataRefDict['fgcmZeropoints'] = butler.dataRef('fgcmZeropoints',
536  fgcmcycle=self.config.cycleNumber)
537  if self.config.doAtmosphereOutput:
538  dataRefDict['fgcmAtmosphereParameters'] = butler.dataRef('fgcmAtmosphereParameters',
539  fgcmcycle=self.config.cycleNumber)
540 
541  struct = self.run(dataRefDict, physicalFilterMap, butler=butler, returnCatalogs=False)
542 
543  if struct.photoCalibs is not None:
544  self.log.info("Outputting photoCalib files.")
545 
546  for visit, detector, physicalFilter, photoCalib in struct.photoCalibs:
547  butler.put(photoCalib, 'fgcm_photoCalib',
548  dataId={visitDataRefName: visit,
549  ccdDataRefName: detector,
550  'filter': physicalFilter})
551 
552  self.log.info("Done outputting photoCalib files.")
553 
554  if struct.atmospheres is not None:
555  self.log.info("Outputting atmosphere transmission files.")
556  for visit, atm in struct.atmospheres:
557  butler.put(atm, "transmission_atmosphere_fgcm",
558  dataId={visitDataRefName: visit})
559  self.log.info("Done outputting atmosphere transmissions.")
560 
561  return pipeBase.Struct(offsets=struct.offsets)
562 
563  def run(self, dataRefDict, physicalFilterMap, returnCatalogs=True, butler=None):
564  """Run the output products task.
565 
566  Parameters
567  ----------
568  dataRefDict : `dict`
569  All dataRefs are `lsst.daf.persistence.ButlerDataRef` (gen2) or
570  `lsst.daf.butler.DeferredDatasetHandle` (gen3)
571  dataRef dictionary with keys:
572 
573  ``"camera"``
574  Camera object (`lsst.afw.cameraGeom.Camera`)
575  ``"fgcmLookUpTable"``
576  dataRef for the FGCM look-up table.
577  ``"fgcmVisitCatalog"``
578  dataRef for visit summary catalog.
579  ``"fgcmStandardStars"``
580  dataRef for the output standard star catalog.
581  ``"fgcmZeropoints"``
582  dataRef for the zeropoint data catalog.
583  ``"fgcmAtmosphereParameters"``
584  dataRef for the atmosphere parameter catalog.
585  ``"fgcmBuildStarsTableConfig"``
586  Config for `lsst.fgcmcal.fgcmBuildStarsTableTask`.
587  physicalFilterMap : `dict`
588  Dictionary of mappings from physical filter to FGCM band.
589  returnCatalogs : `bool`, optional
590  Return photoCalibs as per-visit exposure catalogs.
591  butler : `lsst.daf.persistence.Butler`, optional
592  Gen2 butler used for reference star outputs
593 
594  Returns
595  -------
596  retStruct : `lsst.pipe.base.Struct`
597  Output structure with keys:
598 
599  offsets : `np.ndarray`
600  Final reference offsets, per band.
601  atmospheres : `generator` [(`int`, `lsst.afw.image.TransmissionCurve`)]
602  Generator that returns (visit, transmissionCurve) tuples.
603  photoCalibs : `generator` [(`int`, `int`, `str`, `lsst.afw.image.PhotoCalib`)]
604  Generator that returns (visit, ccd, filtername, photoCalib) tuples.
605  (returned if returnCatalogs is False).
606  photoCalibCatalogs : `generator` [(`int`, `lsst.afw.table.ExposureCatalog`)]
607  Generator that returns (visit, exposureCatalog) tuples.
608  (returned if returnCatalogs is True).
609  """
610  stdCat = dataRefDict['fgcmStandardStars'].get()
611  md = stdCat.getMetadata()
612  bands = md.getArray('BANDS')
613 
614  if self.config.doReferenceCalibration:
615  offsets = self._computeReferenceOffsets(stdCat, bands)
616  else:
617  offsets = np.zeros(len(bands))
618 
619  # This is Gen2 only, and requires the butler.
620  if self.config.doRefcatOutput and butler is not None:
621  self._outputStandardStars(butler, stdCat, offsets, bands, self.config.datasetConfig)
622 
623  del stdCat
624 
625  if self.config.doZeropointOutput:
626  zptCat = dataRefDict['fgcmZeropoints'].get()
627  visitCat = dataRefDict['fgcmVisitCatalog'].get()
628 
629  pcgen = self._outputZeropoints(dataRefDict['camera'], zptCat, visitCat, offsets, bands,
630  physicalFilterMap, returnCatalogs=returnCatalogs)
631  else:
632  pcgen = None
633 
634  if self.config.doAtmosphereOutput:
635  atmCat = dataRefDict['fgcmAtmosphereParameters'].get()
636  atmgen = self._outputAtmospheres(dataRefDict, atmCat)
637  else:
638  atmgen = None
639 
640  retStruct = pipeBase.Struct(offsets=offsets,
641  atmospheres=atmgen)
642  if returnCatalogs:
643  retStruct.photoCalibCatalogs = pcgen
644  else:
645  retStruct.photoCalibs = pcgen
646 
647  return retStruct
648 
649  def generateTractOutputProducts(self, dataRefDict, tract,
650  visitCat, zptCat, atmCat, stdCat,
651  fgcmBuildStarsConfig,
652  returnCatalogs=True,
653  butler=None):
654  """
655  Generate the output products for a given tract, as specified in the config.
656 
657  This method is here to have an alternate entry-point for
658  FgcmCalibrateTract.
659 
660  Parameters
661  ----------
662  dataRefDict : `dict`
663  All dataRefs are `lsst.daf.persistence.ButlerDataRef` (gen2) or
664  `lsst.daf.butler.DeferredDatasetHandle` (gen3)
665  dataRef dictionary with keys:
666 
667  ``"camera"``
668  Camera object (`lsst.afw.cameraGeom.Camera`)
669  ``"fgcmLookUpTable"``
670  dataRef for the FGCM look-up table.
671  tract : `int`
672  Tract number
673  visitCat : `lsst.afw.table.BaseCatalog`
674  FGCM visitCat from `FgcmBuildStarsTask`
675  zptCat : `lsst.afw.table.BaseCatalog`
676  FGCM zeropoint catalog from `FgcmFitCycleTask`
677  atmCat : `lsst.afw.table.BaseCatalog`
678  FGCM atmosphere parameter catalog from `FgcmFitCycleTask`
679  stdCat : `lsst.afw.table.SimpleCatalog`
680  FGCM standard star catalog from `FgcmFitCycleTask`
681  fgcmBuildStarsConfig : `lsst.fgcmcal.FgcmBuildStarsConfig`
682  Configuration object from `FgcmBuildStarsTask`
683  returnCatalogs : `bool`, optional
684  Return photoCalibs as per-visit exposure catalogs.
685  butler: `lsst.daf.persistence.Butler`, optional
686  Gen2 butler used for reference star outputs
687 
688  Returns
689  -------
690  retStruct : `lsst.pipe.base.Struct`
691  Output structure with keys:
692 
693  offsets : `np.ndarray`
694  Final reference offsets, per band.
695  atmospheres : `generator` [(`int`, `lsst.afw.image.TransmissionCurve`)]
696  Generator that returns (visit, transmissionCurve) tuples.
697  photoCalibs : `generator` [(`int`, `int`, `str`, `lsst.afw.image.PhotoCalib`)]
698  Generator that returns (visit, ccd, filtername, photoCalib) tuples.
699  (returned if returnCatalogs is False).
700  photoCalibCatalogs : `generator` [(`int`, `lsst.afw.table.ExposureCatalog`)]
701  Generator that returns (visit, exposureCatalog) tuples.
702  (returned if returnCatalogs is True).
703  """
704  physicalFilterMap = fgcmBuildStarsConfig.physicalFilterMap
705 
706  md = stdCat.getMetadata()
707  bands = md.getArray('BANDS')
708 
709  if self.config.doComposeWcsJacobian and not fgcmBuildStarsConfig.doApplyWcsJacobian:
710  raise RuntimeError("Cannot compose the WCS jacobian if it hasn't been applied "
711  "in fgcmBuildStarsTask.")
712 
713  if not self.config.doComposeWcsJacobian and fgcmBuildStarsConfig.doApplyWcsJacobian:
714  self.log.warn("Jacobian was applied in build-stars but doComposeWcsJacobian is not set.")
715 
716  if self.config.doReferenceCalibration:
717  offsets = self._computeReferenceOffsets(stdCat, bands)
718  else:
719  offsets = np.zeros(len(bands))
720 
721  if self.config.doRefcatOutput and butler is not None:
722  # Create a special config that has the tract number in it
723  datasetConfig = copy.copy(self.config.datasetConfig)
724  datasetConfig.ref_dataset_name = '%s_%d' % (self.config.datasetConfig.ref_dataset_name,
725  tract)
726  self._outputStandardStars(butler, stdCat, offsets, bands, datasetConfig)
727 
728  if self.config.doZeropointOutput:
729  pcgen = self._outputZeropoints(dataRefDict['camera'], zptCat, visitCat, offsets, bands,
730  physicalFilterMap, returnCatalogs=returnCatalogs)
731  else:
732  pcgen = None
733 
734  if self.config.doAtmosphereOutput:
735  atmgen = self._outputAtmospheres(dataRefDict, atmCat)
736  else:
737  atmgen = None
738 
739  retStruct = pipeBase.Struct(offsets=offsets,
740  atmospheres=atmgen)
741  if returnCatalogs:
742  retStruct.photoCalibCatalogs = pcgen
743  else:
744  retStruct.photoCalibs = pcgen
745 
746  return retStruct
747 
748  def _computeReferenceOffsets(self, stdCat, bands):
749  """
750  Compute offsets relative to a reference catalog.
751 
752  This method splits the star catalog into healpix pixels
753  and computes the calibration transfer for a sample of
754  these pixels to approximate the 'absolute' calibration
755  values (on for each band) to apply to transfer the
756  absolute scale.
757 
758  Parameters
759  ----------
760  stdCat : `lsst.afw.table.SimpleCatalog`
761  FGCM standard stars
762  bands : `list` [`str`]
763  List of band names from FGCM output
764  Returns
765  -------
766  offsets : `numpy.array` of floats
767  Per band zeropoint offsets
768  """
769 
770  # Only use stars that are observed in all the bands that were actually used
771  # This will ensure that we use the same healpix pixels for the absolute
772  # calibration of each band.
773  minObs = stdCat['ngood'].min(axis=1)
774 
775  goodStars = (minObs >= 1)
776  stdCat = stdCat[goodStars]
777 
778  self.log.info("Found %d stars with at least 1 good observation in each band" %
779  (len(stdCat)))
780 
781  # We have to make a table for each pixel with flux/fluxErr
782  # This is a temporary table generated for input to the photoCal task.
783  # These fluxes are not instFlux (they are top-of-the-atmosphere approximate and
784  # have had chromatic corrections applied to get to the standard system
785  # specified by the atmosphere/instrumental parameters), nor are they
786  # in Jansky (since they don't have a proper absolute calibration: the overall
787  # zeropoint is estimated from the telescope size, etc.)
788  sourceMapper = afwTable.SchemaMapper(stdCat.schema)
789  sourceMapper.addMinimalSchema(afwTable.SimpleTable.makeMinimalSchema())
790  sourceMapper.editOutputSchema().addField('instFlux', type=np.float64,
791  doc="instrumental flux (counts)")
792  sourceMapper.editOutputSchema().addField('instFluxErr', type=np.float64,
793  doc="instrumental flux error (counts)")
794  badStarKey = sourceMapper.editOutputSchema().addField('flag_badStar',
795  type='Flag',
796  doc="bad flag")
797 
798  # Split up the stars
799  # Note that there is an assumption here that the ra/dec coords stored
800  # on-disk are in radians, and therefore that starObs['coord_ra'] /
801  # starObs['coord_dec'] return radians when used as an array of numpy float64s.
802  theta = np.pi/2. - stdCat['coord_dec']
803  phi = stdCat['coord_ra']
804 
805  ipring = hp.ang2pix(self.config.referencePixelizationNside, theta, phi)
806  h, rev = esutil.stat.histogram(ipring, rev=True)
807 
808  gdpix, = np.where(h >= self.config.referencePixelizationMinStars)
809 
810  self.log.info("Found %d pixels (nside=%d) with at least %d good stars" %
811  (gdpix.size,
812  self.config.referencePixelizationNside,
813  self.config.referencePixelizationMinStars))
814 
815  if gdpix.size < self.config.referencePixelizationNPixels:
816  self.log.warn("Found fewer good pixels (%d) than preferred in configuration (%d)" %
817  (gdpix.size, self.config.referencePixelizationNPixels))
818  else:
819  # Sample out the pixels we want to use
820  gdpix = np.random.choice(gdpix, size=self.config.referencePixelizationNPixels, replace=False)
821 
822  results = np.zeros(gdpix.size, dtype=[('hpix', 'i4'),
823  ('nstar', 'i4', len(bands)),
824  ('nmatch', 'i4', len(bands)),
825  ('zp', 'f4', len(bands)),
826  ('zpErr', 'f4', len(bands))])
827  results['hpix'] = ipring[rev[rev[gdpix]]]
828 
829  # We need a boolean index to deal with catalogs...
830  selected = np.zeros(len(stdCat), dtype=bool)
831 
832  refFluxFields = [None]*len(bands)
833 
834  for p, pix in enumerate(gdpix):
835  i1a = rev[rev[pix]: rev[pix + 1]]
836 
837  # the stdCat afwTable can only be indexed with boolean arrays,
838  # and not numpy index arrays (see DM-16497). This little trick
839  # converts the index array into a boolean array
840  selected[:] = False
841  selected[i1a] = True
842 
843  for b, band in enumerate(bands):
844 
845  struct = self._computeOffsetOneBand(sourceMapper, badStarKey, b, band, stdCat,
846  selected, refFluxFields)
847  results['nstar'][p, b] = len(i1a)
848  results['nmatch'][p, b] = len(struct.arrays.refMag)
849  results['zp'][p, b] = struct.zp
850  results['zpErr'][p, b] = struct.sigma
851 
852  # And compute the summary statistics
853  offsets = np.zeros(len(bands))
854 
855  for b, band in enumerate(bands):
856  # make configurable
857  ok, = np.where(results['nmatch'][:, b] >= self.config.referenceMinMatch)
858  offsets[b] = np.median(results['zp'][ok, b])
859  # use median absolute deviation to estimate Normal sigma
860  # see https://en.wikipedia.org/wiki/Median_absolute_deviation
861  madSigma = 1.4826*np.median(np.abs(results['zp'][ok, b] - offsets[b]))
862  self.log.info("Reference catalog offset for %s band: %.12f +/- %.12f" %
863  (band, offsets[b], madSigma))
864 
865  return offsets
866 
867  def _computeOffsetOneBand(self, sourceMapper, badStarKey,
868  b, band, stdCat, selected, refFluxFields):
869  """
870  Compute the zeropoint offset between the fgcm stdCat and the reference
871  stars for one pixel in one band
872 
873  Parameters
874  ----------
875  sourceMapper: `lsst.afw.table.SchemaMapper`
876  Mapper to go from stdCat to calibratable catalog
877  badStarKey: `lsst.afw.table.Key`
878  Key for the field with bad stars
879  b: `int`
880  Index of the band in the star catalog
881  band: `str`
882  Name of band for reference catalog
883  stdCat: `lsst.afw.table.SimpleCatalog`
884  FGCM standard stars
885  selected: `numpy.array(dtype=bool)`
886  Boolean array of which stars are in the pixel
887  refFluxFields: `list`
888  List of names of flux fields for reference catalog
889  """
890 
891  sourceCat = afwTable.SimpleCatalog(sourceMapper.getOutputSchema())
892  sourceCat.reserve(selected.sum())
893  sourceCat.extend(stdCat[selected], mapper=sourceMapper)
894  sourceCat['instFlux'] = 10.**(stdCat['mag_std_noabs'][selected, b]/(-2.5))
895  sourceCat['instFluxErr'] = (np.log(10.)/2.5)*(stdCat['magErr_std'][selected, b]
896  * sourceCat['instFlux'])
897  # Make sure we only use stars that have valid measurements
898  # (This is perhaps redundant with requirements above that the
899  # stars be observed in all bands, but it can't hurt)
900  badStar = (stdCat['mag_std_noabs'][selected, b] > 90.0)
901  for rec in sourceCat[badStar]:
902  rec.set(badStarKey, True)
903 
904  exposure = afwImage.ExposureF()
905  exposure.setFilterLabel(afwImage.FilterLabel(band=band))
906 
907  if refFluxFields[b] is None:
908  # Need to find the flux field in the reference catalog
909  # to work around limitations of DirectMatch in PhotoCal
910  ctr = stdCat[0].getCoord()
911  rad = 0.05*lsst.geom.degrees
912  refDataTest = self.refObjLoader.loadSkyCircle(ctr, rad, band)
913  refFluxFields[b] = refDataTest.fluxField
914 
915  # Make a copy of the config so that we can modify it
916  calConfig = copy.copy(self.config.photoCal.value)
917  calConfig.match.referenceSelection.signalToNoise.fluxField = refFluxFields[b]
918  calConfig.match.referenceSelection.signalToNoise.errField = refFluxFields[b] + 'Err'
919  calTask = self.config.photoCal.target(refObjLoader=self.refObjLoader,
920  config=calConfig,
921  schema=sourceCat.getSchema())
922 
923  struct = calTask.run(exposure, sourceCat)
924 
925  return struct
926 
927  def _outputStandardStars(self, butler, stdCat, offsets, bands, datasetConfig):
928  """
929  Output standard stars in indexed reference catalog format.
930  This is not currently supported in Gen3.
931 
932  Parameters
933  ----------
934  butler : `lsst.daf.persistence.Butler`
935  stdCat : `lsst.afw.table.SimpleCatalog`
936  FGCM standard star catalog from fgcmFitCycleTask
937  offsets : `numpy.array` of floats
938  Per band zeropoint offsets
939  bands : `list` [`str`]
940  List of band names from FGCM output
941  datasetConfig : `lsst.meas.algorithms.DatasetConfig`
942  Config for reference dataset
943  """
944 
945  self.log.info("Outputting standard stars to %s" % (datasetConfig.ref_dataset_name))
946 
947  indexer = IndexerRegistry[self.config.datasetConfig.indexer.name](
948  self.config.datasetConfig.indexer.active)
949 
950  # We determine the conversion from the native units (typically radians) to
951  # degrees for the first star. This allows us to treat coord_ra/coord_dec as
952  # numpy arrays rather than Angles, which would we approximately 600x slower.
953  # TODO: Fix this after DM-16524 (HtmIndexer.indexPoints should take coords
954  # (as Angles) for input
955  conv = stdCat[0]['coord_ra'].asDegrees()/float(stdCat[0]['coord_ra'])
956  indices = np.array(indexer.indexPoints(stdCat['coord_ra']*conv,
957  stdCat['coord_dec']*conv))
958 
959  formattedCat = self._formatCatalog(stdCat, offsets, bands)
960 
961  # Write the master schema
962  dataId = indexer.makeDataId('master_schema',
963  datasetConfig.ref_dataset_name)
964  masterCat = afwTable.SimpleCatalog(formattedCat.schema)
965  addRefCatMetadata(masterCat)
966  butler.put(masterCat, 'ref_cat', dataId=dataId)
967 
968  # Break up the pixels using a histogram
969  h, rev = esutil.stat.histogram(indices, rev=True)
970  gd, = np.where(h > 0)
971  selected = np.zeros(len(formattedCat), dtype=bool)
972  for i in gd:
973  i1a = rev[rev[i]: rev[i + 1]]
974 
975  # the formattedCat afwTable can only be indexed with boolean arrays,
976  # and not numpy index arrays (see DM-16497). This little trick
977  # converts the index array into a boolean array
978  selected[:] = False
979  selected[i1a] = True
980 
981  # Write the individual pixel
982  dataId = indexer.makeDataId(indices[i1a[0]],
983  datasetConfig.ref_dataset_name)
984  butler.put(formattedCat[selected], 'ref_cat', dataId=dataId)
985 
986  # And save the dataset configuration
987  dataId = indexer.makeDataId(None, datasetConfig.ref_dataset_name)
988  butler.put(datasetConfig, 'ref_cat_config', dataId=dataId)
989 
990  self.log.info("Done outputting standard stars.")
991 
992  def _formatCatalog(self, fgcmStarCat, offsets, bands):
993  """
994  Turn an FGCM-formatted star catalog, applying zeropoint offsets.
995 
996  Parameters
997  ----------
998  fgcmStarCat : `lsst.afw.Table.SimpleCatalog`
999  SimpleCatalog as output by fgcmcal
1000  offsets : `list` with len(self.bands) entries
1001  Zeropoint offsets to apply
1002  bands : `list` [`str`]
1003  List of band names from FGCM output
1004 
1005  Returns
1006  -------
1007  formattedCat: `lsst.afw.table.SimpleCatalog`
1008  SimpleCatalog suitable for using as a reference catalog
1009  """
1010 
1011  sourceMapper = afwTable.SchemaMapper(fgcmStarCat.schema)
1012  minSchema = LoadIndexedReferenceObjectsTask.makeMinimalSchema(bands,
1013  addCentroid=False,
1014  addIsResolved=True,
1015  coordErrDim=0)
1016  sourceMapper.addMinimalSchema(minSchema)
1017  for band in bands:
1018  sourceMapper.editOutputSchema().addField('%s_nGood' % (band), type=np.int32)
1019  sourceMapper.editOutputSchema().addField('%s_nTotal' % (band), type=np.int32)
1020  sourceMapper.editOutputSchema().addField('%s_nPsfCandidate' % (band), type=np.int32)
1021 
1022  formattedCat = afwTable.SimpleCatalog(sourceMapper.getOutputSchema())
1023  formattedCat.reserve(len(fgcmStarCat))
1024  formattedCat.extend(fgcmStarCat, mapper=sourceMapper)
1025 
1026  # Note that we don't have to set `resolved` because the default is False
1027 
1028  for b, band in enumerate(bands):
1029  mag = fgcmStarCat['mag_std_noabs'][:, b].astype(np.float64) + offsets[b]
1030  # We want fluxes in nJy from calibrated AB magnitudes
1031  # (after applying offset). Updated after RFC-549 and RFC-575.
1032  flux = (mag*units.ABmag).to_value(units.nJy)
1033  fluxErr = (np.log(10.)/2.5)*flux*fgcmStarCat['magErr_std'][:, b].astype(np.float64)
1034 
1035  formattedCat['%s_flux' % (band)][:] = flux
1036  formattedCat['%s_fluxErr' % (band)][:] = fluxErr
1037  formattedCat['%s_nGood' % (band)][:] = fgcmStarCat['ngood'][:, b]
1038  formattedCat['%s_nTotal' % (band)][:] = fgcmStarCat['ntotal'][:, b]
1039  formattedCat['%s_nPsfCandidate' % (band)][:] = fgcmStarCat['npsfcand'][:, b]
1040 
1041  addRefCatMetadata(formattedCat)
1042 
1043  return formattedCat
1044 
1045  def _outputZeropoints(self, camera, zptCat, visitCat, offsets, bands,
1046  physicalFilterMap, returnCatalogs=True,
1047  tract=None):
1048  """Output the zeropoints in fgcm_photoCalib format.
1049 
1050  Parameters
1051  ----------
1052  camera : `lsst.afw.cameraGeom.Camera`
1053  Camera from the butler.
1054  zptCat : `lsst.afw.table.BaseCatalog`
1055  FGCM zeropoint catalog from `FgcmFitCycleTask`.
1056  visitCat : `lsst.afw.table.BaseCatalog`
1057  FGCM visitCat from `FgcmBuildStarsTask`.
1058  offsets : `numpy.array`
1059  Float array of absolute calibration offsets, one for each filter.
1060  bands : `list` [`str`]
1061  List of band names from FGCM output.
1062  physicalFilterMap : `dict`
1063  Dictionary of mappings from physical filter to FGCM band.
1064  returnCatalogs : `bool`, optional
1065  Return photoCalibs as per-visit exposure catalogs.
1066  tract: `int`, optional
1067  Tract number to output. Default is None (global calibration)
1068 
1069  Returns
1070  -------
1071  photoCalibs : `generator` [(`int`, `int`, `str`, `lsst.afw.image.PhotoCalib`)]
1072  Generator that returns (visit, ccd, filtername, photoCalib) tuples.
1073  (returned if returnCatalogs is False).
1074  photoCalibCatalogs : `generator` [(`int`, `lsst.afw.table.ExposureCatalog`)]
1075  Generator that returns (visit, exposureCatalog) tuples.
1076  (returned if returnCatalogs is True).
1077  """
1078  # Select visit/ccds where we have a calibration
1079  # This includes ccds where we were able to interpolate from neighboring
1080  # ccds.
1081  cannot_compute = fgcm.fgcmUtilities.zpFlagDict['CANNOT_COMPUTE_ZEROPOINT']
1082  selected = (((zptCat['fgcmFlag'] & cannot_compute) == 0)
1083  & (zptCat['fgcmZptVar'] > 0.0))
1084 
1085  # Log warnings for any visit which has no valid zeropoints
1086  badVisits = np.unique(zptCat['visit'][~selected])
1087  goodVisits = np.unique(zptCat['visit'][selected])
1088  allBadVisits = badVisits[~np.isin(badVisits, goodVisits)]
1089  for allBadVisit in allBadVisits:
1090  self.log.warn(f'No suitable photoCalib for visit {allBadVisit}')
1091 
1092  # Get a mapping from filtername to the offsets
1093  offsetMapping = {}
1094  for f in physicalFilterMap:
1095  # Not every filter in the map will necesarily have a band.
1096  if physicalFilterMap[f] in bands:
1097  offsetMapping[f] = offsets[bands.index(physicalFilterMap[f])]
1098 
1099  # Get a mapping from "ccd" to the ccd index used for the scaling
1100  ccdMapping = {}
1101  for ccdIndex, detector in enumerate(camera):
1102  ccdMapping[detector.getId()] = ccdIndex
1103 
1104  # And a mapping to get the flat-field scaling values
1105  scalingMapping = {}
1106  for rec in visitCat:
1107  scalingMapping[rec['visit']] = rec['scaling']
1108 
1109  if self.config.doComposeWcsJacobian:
1110  approxPixelAreaFields = computeApproxPixelAreaFields(camera)
1111 
1112  # The zptCat is sorted by visit, which is useful
1113  lastVisit = -1
1114  zptVisitCatalog = None
1115 
1116  metadata = dafBase.PropertyList()
1117  metadata.add("COMMENT", "Catalog id is detector id, sorted.")
1118  metadata.add("COMMENT", "Only detectors with data have entries.")
1119 
1120  for rec in zptCat[selected]:
1121  # Retrieve overall scaling
1122  scaling = scalingMapping[rec['visit']][ccdMapping[rec['detector']]]
1123 
1124  # The postCalibrationOffset describe any zeropoint offsets
1125  # to apply after the fgcm calibration. The first part comes
1126  # from the reference catalog match (used in testing). The
1127  # second part comes from the mean chromatic correction
1128  # (if configured).
1129  postCalibrationOffset = offsetMapping[rec['filtername']]
1130  if self.config.doApplyMeanChromaticCorrection:
1131  postCalibrationOffset += rec['fgcmDeltaChrom']
1132 
1133  fgcmSuperStarField = self._getChebyshevBoundedField(rec['fgcmfZptSstarCheb'],
1134  rec['fgcmfZptChebXyMax'])
1135  # Convert from FGCM AB to nJy
1136  fgcmZptField = self._getChebyshevBoundedField((rec['fgcmfZptCheb']*units.AB).to_value(units.nJy),
1137  rec['fgcmfZptChebXyMax'],
1138  offset=postCalibrationOffset,
1139  scaling=scaling)
1140 
1141  if self.config.doComposeWcsJacobian:
1142 
1143  fgcmField = afwMath.ProductBoundedField([approxPixelAreaFields[rec['detector']],
1144  fgcmSuperStarField,
1145  fgcmZptField])
1146  else:
1147  # The photoCalib is just the product of the fgcmSuperStarField and the
1148  # fgcmZptField
1149  fgcmField = afwMath.ProductBoundedField([fgcmSuperStarField, fgcmZptField])
1150 
1151  # The "mean" calibration will be set to the center of the ccd for reference
1152  calibCenter = fgcmField.evaluate(fgcmField.getBBox().getCenter())
1153  calibErr = (np.log(10.0)/2.5)*calibCenter*np.sqrt(rec['fgcmZptVar'])
1154  photoCalib = afwImage.PhotoCalib(calibrationMean=calibCenter,
1155  calibrationErr=calibErr,
1156  calibration=fgcmField,
1157  isConstant=False)
1158 
1159  if not returnCatalogs:
1160  # Return individual photoCalibs
1161  yield (int(rec['visit']), int(rec['detector']), rec['filtername'], photoCalib)
1162  else:
1163  # Return full per-visit exposure catalogs
1164  if rec['visit'] != lastVisit:
1165  # This is a new visit. If the last visit was not -1, yield
1166  # the ExposureCatalog
1167  if lastVisit > -1:
1168  # ensure that the detectors are in sorted order, for fast lookups
1169  zptVisitCatalog.sort()
1170  yield (int(lastVisit), zptVisitCatalog)
1171  else:
1172  # We need to create a new schema
1173  zptExpCatSchema = afwTable.ExposureTable.makeMinimalSchema()
1174  zptExpCatSchema.addField('visit', type='I', doc='Visit number')
1175 
1176  # And start a new one
1177  zptVisitCatalog = afwTable.ExposureCatalog(zptExpCatSchema)
1178  zptVisitCatalog.setMetadata(metadata)
1179 
1180  lastVisit = int(rec['visit'])
1181 
1182  catRecord = zptVisitCatalog.addNew()
1183  catRecord['id'] = int(rec['detector'])
1184  catRecord['visit'] = rec['visit']
1185  catRecord.setPhotoCalib(photoCalib)
1186 
1187  # Final output of last exposure catalog
1188  if returnCatalogs:
1189  # ensure that the detectors are in sorted order, for fast lookups
1190  zptVisitCatalog.sort()
1191  yield (int(lastVisit), zptVisitCatalog)
1192 
1193  def _getChebyshevBoundedField(self, coefficients, xyMax, offset=0.0, scaling=1.0):
1194  """
1195  Make a ChebyshevBoundedField from fgcm coefficients, with optional offset
1196  and scaling.
1197 
1198  Parameters
1199  ----------
1200  coefficients: `numpy.array`
1201  Flattened array of chebyshev coefficients
1202  xyMax: `list` of length 2
1203  Maximum x and y of the chebyshev bounding box
1204  offset: `float`, optional
1205  Absolute calibration offset. Default is 0.0
1206  scaling: `float`, optional
1207  Flat scaling value from fgcmBuildStars. Default is 1.0
1208 
1209  Returns
1210  -------
1211  boundedField: `lsst.afw.math.ChebyshevBoundedField`
1212  """
1213 
1214  orderPlus1 = int(np.sqrt(coefficients.size))
1215  pars = np.zeros((orderPlus1, orderPlus1))
1216 
1217  bbox = lsst.geom.Box2I(lsst.geom.Point2I(0.0, 0.0),
1218  lsst.geom.Point2I(*xyMax))
1219 
1220  pars[:, :] = (coefficients.reshape(orderPlus1, orderPlus1)
1221  * (10.**(offset/-2.5))*scaling)
1222 
1223  boundedField = afwMath.ChebyshevBoundedField(bbox, pars)
1224 
1225  return boundedField
1226 
1227  def _outputAtmospheres(self, dataRefDict, atmCat):
1228  """
1229  Output the atmospheres.
1230 
1231  Parameters
1232  ----------
1233  dataRefDict : `dict`
1234  All dataRefs are `lsst.daf.persistence.ButlerDataRef` (gen2) or
1235  `lsst.daf.butler.DeferredDatasetHandle` (gen3)
1236  dataRef dictionary with keys:
1237 
1238  ``"fgcmLookUpTable"``
1239  dataRef for the FGCM look-up table.
1240  atmCat : `lsst.afw.table.BaseCatalog`
1241  FGCM atmosphere parameter catalog from fgcmFitCycleTask.
1242 
1243  Returns
1244  -------
1245  atmospheres : `generator` [(`int`, `lsst.afw.image.TransmissionCurve`)]
1246  Generator that returns (visit, transmissionCurve) tuples.
1247  """
1248  # First, we need to grab the look-up table and key info
1249  lutCat = dataRefDict['fgcmLookUpTable'].get()
1250 
1251  atmosphereTableName = lutCat[0]['tablename']
1252  elevation = lutCat[0]['elevation']
1253  atmLambda = lutCat[0]['atmLambda']
1254  lutCat = None
1255 
1256  # Make the atmosphere table if possible
1257  try:
1258  atmTable = fgcm.FgcmAtmosphereTable.initWithTableName(atmosphereTableName)
1259  atmTable.loadTable()
1260  except IOError:
1261  atmTable = None
1262 
1263  if atmTable is None:
1264  # Try to use MODTRAN instead
1265  try:
1266  modGen = fgcm.ModtranGenerator(elevation)
1267  lambdaRange = np.array([atmLambda[0], atmLambda[-1]])/10.
1268  lambdaStep = (atmLambda[1] - atmLambda[0])/10.
1269  except (ValueError, IOError) as e:
1270  raise RuntimeError("FGCM look-up-table generated with modtran, "
1271  "but modtran not configured to run.") from e
1272 
1273  zenith = np.degrees(np.arccos(1./atmCat['secZenith']))
1274 
1275  for i, visit in enumerate(atmCat['visit']):
1276  if atmTable is not None:
1277  # Interpolate the atmosphere table
1278  atmVals = atmTable.interpolateAtmosphere(pmb=atmCat[i]['pmb'],
1279  pwv=atmCat[i]['pwv'],
1280  o3=atmCat[i]['o3'],
1281  tau=atmCat[i]['tau'],
1282  alpha=atmCat[i]['alpha'],
1283  zenith=zenith[i],
1284  ctranslamstd=[atmCat[i]['cTrans'],
1285  atmCat[i]['lamStd']])
1286  else:
1287  # Run modtran
1288  modAtm = modGen(pmb=atmCat[i]['pmb'],
1289  pwv=atmCat[i]['pwv'],
1290  o3=atmCat[i]['o3'],
1291  tau=atmCat[i]['tau'],
1292  alpha=atmCat[i]['alpha'],
1293  zenith=zenith[i],
1294  lambdaRange=lambdaRange,
1295  lambdaStep=lambdaStep,
1296  ctranslamstd=[atmCat[i]['cTrans'],
1297  atmCat[i]['lamStd']])
1298  atmVals = modAtm['COMBINED']
1299 
1300  # Now need to create something to persist...
1301  curve = TransmissionCurve.makeSpatiallyConstant(throughput=atmVals,
1302  wavelengths=atmLambda,
1303  throughputAtMin=atmVals[0],
1304  throughputAtMax=atmVals[-1])
1305 
1306  yield (int(visit), curve)
def computeApproxPixelAreaFields(camera)
Definition: utilities.py:486