LSSTApplications  20.0.0
LSSTDataManagementBasePackage
Functions
lsst.fgcmcal.utilities Namespace Reference

Functions

def makeConfigDict (config, log, camera, maxIter, resetFitParameters, outputZeropoints, tract=None)
 
def translateFgcmLut (lutCat, filterMap)
 
def translateVisitCatalog (visitCat)
 
def computeCcdOffsets (camera, defaultOrientation)
 
def computeReferencePixelScale (camera)
 
def computeApproxPixelAreaFields (camera)
 
def makeZptSchema (superStarChebyshevSize, zptChebyshevSize)
 
def makeZptCat (zptSchema, zpStruct)
 
def makeAtmSchema ()
 
def makeAtmCat (atmSchema, atmStruct)
 
def makeStdSchema (nBands)
 
def makeStdCat (stdSchema, stdStruct, goodBands)
 
def computeApertureRadius (schema, fluxField)
 
def extractReferenceMags (refStars, bands, filterMap)
 

Function Documentation

◆ computeApertureRadius()

def lsst.fgcmcal.utilities.computeApertureRadius (   schema,
  fluxField 
)
Compute the radius associated with a CircularApertureFlux field or
associated slot.

Parameters
----------
schema : `lsst.afw.table.schema`
fluxField : `str`
   CircularApertureFlux field or associated slot.

Returns
-------
apertureRadius: `float`
   Radius of the aperture field, in pixels.

Raises
------
RuntimeError: Raised if flux field is not a CircularApertureFlux
   or associated slot.

Definition at line 783 of file utilities.py.

783 def computeApertureRadius(schema, fluxField):
784  """
785  Compute the radius associated with a CircularApertureFlux field or
786  associated slot.
787 
788  Parameters
789  ----------
790  schema : `lsst.afw.table.schema`
791  fluxField : `str`
792  CircularApertureFlux field or associated slot.
793 
794  Returns
795  -------
796  apertureRadius: `float`
797  Radius of the aperture field, in pixels.
798 
799  Raises
800  ------
801  RuntimeError: Raised if flux field is not a CircularApertureFlux
802  or associated slot.
803  """
804  fluxFieldName = schema[fluxField].asField().getName()
805 
806  m = re.search(r'CircularApertureFlux_(\d+)_(\d+)_', fluxFieldName)
807 
808  if m is None:
809  raise RuntimeError("Flux field %s does not correspond to a circular aperture"
810  % (fluxField))
811 
812  apertureRadius = float(m.groups()[0]) + float(m.groups()[1])/10.
813 
814  return apertureRadius
815 
816 

◆ computeApproxPixelAreaFields()

def lsst.fgcmcal.utilities.computeApproxPixelAreaFields (   camera)
Compute the approximate pixel area bounded fields from the camera
geometry.

Parameters
----------
camera: `lsst.afw.cameraGeom.Camera`

Returns
-------
approxPixelAreaFields: `dict`
   Dictionary of approximate area fields, keyed with detector ID

Definition at line 474 of file utilities.py.

474 def computeApproxPixelAreaFields(camera):
475  """
476  Compute the approximate pixel area bounded fields from the camera
477  geometry.
478 
479  Parameters
480  ----------
481  camera: `lsst.afw.cameraGeom.Camera`
482 
483  Returns
484  -------
485  approxPixelAreaFields: `dict`
486  Dictionary of approximate area fields, keyed with detector ID
487  """
488 
489  areaScaling = 1. / computeReferencePixelScale(camera)**2.
490 
491  # Generate fake WCSs centered at 180/0 to avoid the RA=0/360 problem,
492  # since we are looking for relative scales
493  boresight = geom.SpherePoint(180.0*geom.degrees, 0.0*geom.degrees)
494 
495  flipX = False
496  # Create a temporary visitInfo for input to createInitialSkyWcs
497  # The orientation does not matter for the area computation
498  visitInfo = afwImage.VisitInfo(boresightRaDec=boresight,
499  boresightRotAngle=0.0*geom.degrees,
500  rotType=afwImage.visitInfo.RotType.SKY)
501 
502  approxPixelAreaFields = {}
503 
504  for i, detector in enumerate(camera):
505  key = detector.getId()
506 
507  wcs = createInitialSkyWcs(visitInfo, detector, flipX)
508  bbox = detector.getBBox()
509 
510  areaField = afwMath.PixelAreaBoundedField(bbox, wcs,
511  unit=geom.arcseconds, scaling=areaScaling)
512  approxAreaField = afwMath.ChebyshevBoundedField.approximate(areaField)
513 
514  approxPixelAreaFields[key] = approxAreaField
515 
516  return approxPixelAreaFields
517 
518 

◆ computeCcdOffsets()

def lsst.fgcmcal.utilities.computeCcdOffsets (   camera,
  defaultOrientation 
)
Compute the CCD offsets in ra/dec and x/y space

Parameters
----------
camera: `lsst.afw.cameraGeom.Camera`
defaultOrientation: `float`
   Default camera orientation (degrees)

Returns
-------
ccdOffsets: `numpy.ndarray`
   Numpy array with ccd offset information for input to FGCM.
   Angular units are degrees, and x/y units are pixels.

Definition at line 378 of file utilities.py.

378 def computeCcdOffsets(camera, defaultOrientation):
379  """
380  Compute the CCD offsets in ra/dec and x/y space
381 
382  Parameters
383  ----------
384  camera: `lsst.afw.cameraGeom.Camera`
385  defaultOrientation: `float`
386  Default camera orientation (degrees)
387 
388  Returns
389  -------
390  ccdOffsets: `numpy.ndarray`
391  Numpy array with ccd offset information for input to FGCM.
392  Angular units are degrees, and x/y units are pixels.
393  """
394  # TODO: DM-21215 will fully generalize to arbitrary camera orientations
395 
396  # and we need to know the ccd offsets from the camera geometry
397  ccdOffsets = np.zeros(len(camera), dtype=[('CCDNUM', 'i4'),
398  ('DELTA_RA', 'f8'),
399  ('DELTA_DEC', 'f8'),
400  ('RA_SIZE', 'f8'),
401  ('DEC_SIZE', 'f8'),
402  ('X_SIZE', 'i4'),
403  ('Y_SIZE', 'i4')])
404 
405  # Generate fake WCSs centered at 180/0 to avoid the RA=0/360 problem,
406  # since we are looking for relative positions
407  boresight = geom.SpherePoint(180.0*geom.degrees, 0.0*geom.degrees)
408 
409  # TODO: DM-17597 will update testdata_jointcal so that the test data
410  # does not have nan as the boresight angle for HSC data. For the
411  # time being, there is this ungainly hack.
412  if camera.getName() == 'HSC' and np.isnan(defaultOrientation):
413  orientation = 270*geom.degrees
414  else:
415  orientation = defaultOrientation*geom.degrees
416  flipX = False
417 
418  # Create a temporary visitInfo for input to createInitialSkyWcs
419  visitInfo = afwImage.VisitInfo(boresightRaDec=boresight,
420  boresightRotAngle=orientation,
421  rotType=afwImage.visitInfo.RotType.SKY)
422 
423  for i, detector in enumerate(camera):
424  ccdOffsets['CCDNUM'][i] = detector.getId()
425 
426  wcs = createInitialSkyWcs(visitInfo, detector, flipX)
427 
428  detCenter = wcs.pixelToSky(detector.getCenter(afwCameraGeom.PIXELS))
429  ccdOffsets['DELTA_RA'][i] = (detCenter.getRa() - boresight.getRa()).asDegrees()
430  ccdOffsets['DELTA_DEC'][i] = (detCenter.getDec() - boresight.getDec()).asDegrees()
431 
432  bbox = detector.getBBox()
433 
434  detCorner1 = wcs.pixelToSky(geom.Point2D(bbox.getMin()))
435  detCorner2 = wcs.pixelToSky(geom.Point2D(bbox.getMax()))
436 
437  ccdOffsets['RA_SIZE'][i] = np.abs((detCorner2.getRa() - detCorner1.getRa()).asDegrees())
438  ccdOffsets['DEC_SIZE'][i] = np.abs((detCorner2.getDec() - detCorner1.getDec()).asDegrees())
439 
440  ccdOffsets['X_SIZE'][i] = bbox.getMaxX()
441  ccdOffsets['Y_SIZE'][i] = bbox.getMaxY()
442 
443  return ccdOffsets
444 
445 

◆ computeReferencePixelScale()

def lsst.fgcmcal.utilities.computeReferencePixelScale (   camera)
Compute the median pixel scale in the camera

Returns
-------
pixelScale: `float`
   Average pixel scale (arcsecond) over the camera

Definition at line 446 of file utilities.py.

446 def computeReferencePixelScale(camera):
447  """
448  Compute the median pixel scale in the camera
449 
450  Returns
451  -------
452  pixelScale: `float`
453  Average pixel scale (arcsecond) over the camera
454  """
455 
456  boresight = geom.SpherePoint(180.0*geom.degrees, 0.0*geom.degrees)
457  orientation = 0.0*geom.degrees
458  flipX = False
459 
460  # Create a temporary visitInfo for input to createInitialSkyWcs
461  visitInfo = afwImage.VisitInfo(boresightRaDec=boresight,
462  boresightRotAngle=orientation,
463  rotType=afwImage.visitInfo.RotType.SKY)
464 
465  pixelScales = np.zeros(len(camera))
466  for i, detector in enumerate(camera):
467  wcs = createInitialSkyWcs(visitInfo, detector, flipX)
468  pixelScales[i] = wcs.getPixelScale().asArcseconds()
469 
470  ok, = np.where(pixelScales > 0.0)
471  return np.median(pixelScales[ok])
472 
473 

◆ extractReferenceMags()

def lsst.fgcmcal.utilities.extractReferenceMags (   refStars,
  bands,
  filterMap 
)
Extract reference magnitudes from refStars for given bands and
associated filterMap.

Parameters
----------
refStars : `lsst.afw.table.BaseCatalog`
   FGCM reference star catalog
bands : `list`
   List of bands for calibration
filterMap: `dict`
   FGCM mapping of filter to band

Returns
-------
refMag : `np.ndarray`
   nstar x nband array of reference magnitudes
refMagErr : `np.ndarray`
   nstar x nband array of reference magnitude errors

Definition at line 817 of file utilities.py.

817 def extractReferenceMags(refStars, bands, filterMap):
818  """
819  Extract reference magnitudes from refStars for given bands and
820  associated filterMap.
821 
822  Parameters
823  ----------
824  refStars : `lsst.afw.table.BaseCatalog`
825  FGCM reference star catalog
826  bands : `list`
827  List of bands for calibration
828  filterMap: `dict`
829  FGCM mapping of filter to band
830 
831  Returns
832  -------
833  refMag : `np.ndarray`
834  nstar x nband array of reference magnitudes
835  refMagErr : `np.ndarray`
836  nstar x nband array of reference magnitude errors
837  """
838  # After DM-23331 fgcm reference catalogs have FILTERNAMES to prevent
839  # against index errors and allow more flexibility in fitting after
840  # the build stars step.
841 
842  md = refStars.getMetadata()
843  if 'FILTERNAMES' in md:
844  filternames = md.getArray('FILTERNAMES')
845 
846  # The reference catalog that fgcm wants has one entry per band
847  # in the config file
848  refMag = np.zeros((len(refStars), len(bands)),
849  dtype=refStars['refMag'].dtype) + 99.0
850  refMagErr = np.zeros_like(refMag) + 99.0
851  for i, filtername in enumerate(filternames):
852  # We are allowed to run the fit configured so that we do not
853  # use every column in the reference catalog.
854  try:
855  band = filterMap[filtername]
856  except KeyError:
857  continue
858  try:
859  ind = bands.index(band)
860  except ValueError:
861  continue
862 
863  refMag[:, ind] = refStars['refMag'][:, i]
864  refMagErr[:, ind] = refStars['refMagErr'][:, i]
865 
866  else:
867  # Continue to use old catalogs as before.
868  refMag = refStars['refMag'][:, :]
869  refMagErr = refStars['refMagErr'][:, :]
870 
871  return refMag, refMagErr

◆ makeAtmCat()

def lsst.fgcmcal.utilities.makeAtmCat (   atmSchema,
  atmStruct 
)
Make the atmosphere catalog for persistence

Parameters
----------
atmSchema: `lsst.afw.table.Schema`
   Atmosphere catalog schema
atmStruct: `numpy.ndarray`
   Atmosphere structure from fgcm

Returns
-------
atmCat: `lsst.afw.table.BaseCatalog`
   Atmosphere catalog for persistence

Definition at line 675 of file utilities.py.

675 def makeAtmCat(atmSchema, atmStruct):
676  """
677  Make the atmosphere catalog for persistence
678 
679  Parameters
680  ----------
681  atmSchema: `lsst.afw.table.Schema`
682  Atmosphere catalog schema
683  atmStruct: `numpy.ndarray`
684  Atmosphere structure from fgcm
685 
686  Returns
687  -------
688  atmCat: `lsst.afw.table.BaseCatalog`
689  Atmosphere catalog for persistence
690  """
691 
692  atmCat = afwTable.BaseCatalog(atmSchema)
693  atmCat.reserve(atmStruct.size)
694  for i in range(atmStruct.size):
695  atmCat.addNew()
696 
697  atmCat['visit'][:] = atmStruct['VISIT']
698  atmCat['pmb'][:] = atmStruct['PMB']
699  atmCat['pwv'][:] = atmStruct['PWV']
700  atmCat['tau'][:] = atmStruct['TAU']
701  atmCat['alpha'][:] = atmStruct['ALPHA']
702  atmCat['o3'][:] = atmStruct['O3']
703  atmCat['secZenith'][:] = atmStruct['SECZENITH']
704  atmCat['cTrans'][:] = atmStruct['CTRANS']
705  atmCat['lamStd'][:] = atmStruct['LAMSTD']
706 
707  return atmCat
708 
709 

◆ makeAtmSchema()

def lsst.fgcmcal.utilities.makeAtmSchema ( )
Make the atmosphere schema

Returns
-------
atmSchema: `lsst.afw.table.Schema`

Definition at line 651 of file utilities.py.

651 def makeAtmSchema():
652  """
653  Make the atmosphere schema
654 
655  Returns
656  -------
657  atmSchema: `lsst.afw.table.Schema`
658  """
659 
660  atmSchema = afwTable.Schema()
661 
662  atmSchema.addField('visit', type=np.int32, doc='Visit number')
663  atmSchema.addField('pmb', type=np.float64, doc='Barometric pressure (mb)')
664  atmSchema.addField('pwv', type=np.float64, doc='Water vapor (mm)')
665  atmSchema.addField('tau', type=np.float64, doc='Aerosol optical depth')
666  atmSchema.addField('alpha', type=np.float64, doc='Aerosol slope')
667  atmSchema.addField('o3', type=np.float64, doc='Ozone (dobson)')
668  atmSchema.addField('secZenith', type=np.float64, doc='Secant(zenith) (~ airmass)')
669  atmSchema.addField('cTrans', type=np.float64, doc='Transmission correction factor')
670  atmSchema.addField('lamStd', type=np.float64, doc='Wavelength for transmission correction')
671 
672  return atmSchema
673 
674 

◆ makeConfigDict()

def lsst.fgcmcal.utilities.makeConfigDict (   config,
  log,
  camera,
  maxIter,
  resetFitParameters,
  outputZeropoints,
  tract = None 
)
Make the FGCM fit cycle configuration dict

Parameters
----------
config: `lsst.fgcmcal.FgcmFitCycleConfig`
    Configuration object
log: `lsst.log.Log`
    LSST log object
camera: `lsst.afw.cameraGeom.Camera`
    Camera from the butler
maxIter: `int`
    Maximum number of iterations
resetFitParameters: `bool`
    Reset fit parameters before fitting?
outputZeropoints: `bool`
    Compute zeropoints for output?
tract: `int`, optional
    Tract number for extending the output file name for debugging.
    Default is None.

Returns
-------
configDict: `dict`
    Configuration dictionary for fgcm

Definition at line 43 of file utilities.py.

43 def makeConfigDict(config, log, camera, maxIter,
44  resetFitParameters, outputZeropoints, tract=None):
45  """
46  Make the FGCM fit cycle configuration dict
47 
48  Parameters
49  ----------
50  config: `lsst.fgcmcal.FgcmFitCycleConfig`
51  Configuration object
52  log: `lsst.log.Log`
53  LSST log object
54  camera: `lsst.afw.cameraGeom.Camera`
55  Camera from the butler
56  maxIter: `int`
57  Maximum number of iterations
58  resetFitParameters: `bool`
59  Reset fit parameters before fitting?
60  outputZeropoints: `bool`
61  Compute zeropoints for output?
62  tract: `int`, optional
63  Tract number for extending the output file name for debugging.
64  Default is None.
65 
66  Returns
67  -------
68  configDict: `dict`
69  Configuration dictionary for fgcm
70  """
71  # Extract the bands that are _not_ being fit for fgcm configuration
72  notFitBands = [b for b in config.bands if b not in config.fitBands]
73 
74  # process the starColorCuts
75  starColorCutList = []
76  for ccut in config.starColorCuts:
77  parts = ccut.split(',')
78  starColorCutList.append([parts[0], parts[1], float(parts[2]), float(parts[3])])
79 
80  # TODO: Having direct access to the mirror area from the camera would be
81  # useful. See DM-16489.
82  # Mirror area in cm**2
83  mirrorArea = np.pi*(camera.telescopeDiameter*100./2.)**2.
84 
85  # Get approximate average camera gain:
86  gains = [amp.getGain() for detector in camera for amp in detector.getAmplifiers()]
87  cameraGain = float(np.median(gains))
88 
89  if tract is None:
90  outfileBase = config.outfileBase
91  else:
92  outfileBase = '%s-%06d' % (config.outfileBase, tract)
93 
94  # create a configuration dictionary for fgcmFitCycle
95  configDict = {'outfileBase': outfileBase,
96  'logger': log,
97  'exposureFile': None,
98  'obsFile': None,
99  'indexFile': None,
100  'lutFile': None,
101  'mirrorArea': mirrorArea,
102  'cameraGain': cameraGain,
103  'ccdStartIndex': camera[0].getId(),
104  'expField': 'VISIT',
105  'ccdField': 'CCD',
106  'seeingField': 'DELTA_APER',
107  'fwhmField': 'PSFSIGMA',
108  'skyBrightnessField': 'SKYBACKGROUND',
109  'deepFlag': 'DEEPFLAG', # unused
110  'bands': list(config.bands),
111  'fitBands': list(config.fitBands),
112  'notFitBands': notFitBands,
113  'requiredBands': list(config.requiredBands),
114  'filterToBand': dict(config.filterMap),
115  'logLevel': 'INFO', # FIXME
116  'nCore': config.nCore,
117  'nStarPerRun': config.nStarPerRun,
118  'nExpPerRun': config.nExpPerRun,
119  'reserveFraction': config.reserveFraction,
120  'freezeStdAtmosphere': config.freezeStdAtmosphere,
121  'precomputeSuperStarInitialCycle': config.precomputeSuperStarInitialCycle,
122  'superStarSubCCDDict': dict(config.superStarSubCcdDict),
123  'superStarSubCCDChebyshevOrder': config.superStarSubCcdChebyshevOrder,
124  'superStarSubCCDTriangular': config.superStarSubCcdTriangular,
125  'superStarSigmaClip': config.superStarSigmaClip,
126  'ccdGraySubCCDDict': dict(config.ccdGraySubCcdDict),
127  'ccdGraySubCCDChebyshevOrder': config.ccdGraySubCcdChebyshevOrder,
128  'ccdGraySubCCDTriangular': config.ccdGraySubCcdTriangular,
129  'cycleNumber': config.cycleNumber,
130  'maxIter': maxIter,
131  'UTBoundary': config.utBoundary,
132  'washMJDs': config.washMjds,
133  'epochMJDs': config.epochMjds,
134  'coatingMJDs': config.coatingMjds,
135  'minObsPerBand': config.minObsPerBand,
136  'latitude': config.latitude,
137  'brightObsGrayMax': config.brightObsGrayMax,
138  'minStarPerCCD': config.minStarPerCcd,
139  'minCCDPerExp': config.minCcdPerExp,
140  'maxCCDGrayErr': config.maxCcdGrayErr,
141  'minStarPerExp': config.minStarPerExp,
142  'minExpPerNight': config.minExpPerNight,
143  'expGrayInitialCut': config.expGrayInitialCut,
144  'expGrayPhotometricCutDict': dict(config.expGrayPhotometricCutDict),
145  'expGrayHighCutDict': dict(config.expGrayHighCutDict),
146  'expGrayRecoverCut': config.expGrayRecoverCut,
147  'expVarGrayPhotometricCutDict': dict(config.expVarGrayPhotometricCutDict),
148  'expGrayErrRecoverCut': config.expGrayErrRecoverCut,
149  'refStarSnMin': config.refStarSnMin,
150  'refStarOutlierNSig': config.refStarOutlierNSig,
151  'applyRefStarColorCuts': config.applyRefStarColorCuts,
152  'illegalValue': -9999.0, # internally used by fgcm.
153  'starColorCuts': starColorCutList,
154  'aperCorrFitNBins': config.aperCorrFitNBins,
155  'aperCorrInputSlopeDict': dict(config.aperCorrInputSlopeDict),
156  'sedBoundaryTermDict': config.sedboundaryterms.toDict()['data'],
157  'sedTermDict': config.sedterms.toDict()['data'],
158  'colorSplitBands': list(config.colorSplitBands),
159  'sigFgcmMaxErr': config.sigFgcmMaxErr,
160  'sigFgcmMaxEGrayDict': dict(config.sigFgcmMaxEGrayDict),
161  'ccdGrayMaxStarErr': config.ccdGrayMaxStarErr,
162  'approxThroughputDict': dict(config.approxThroughputDict),
163  'sigmaCalRange': list(config.sigmaCalRange),
164  'sigmaCalFitPercentile': list(config.sigmaCalFitPercentile),
165  'sigmaCalPlotPercentile': list(config.sigmaCalPlotPercentile),
166  'sigma0Phot': config.sigma0Phot,
167  'mapLongitudeRef': config.mapLongitudeRef,
168  'mapNSide': config.mapNSide,
169  'varNSig': 100.0, # Turn off 'variable star selection' which doesn't work yet
170  'varMinBand': 2,
171  'useRetrievedPwv': False,
172  'useNightlyRetrievedPwv': False,
173  'pwvRetrievalSmoothBlock': 25,
174  'useQuadraticPwv': config.useQuadraticPwv,
175  'useRetrievedTauInit': False,
176  'tauRetrievalMinCCDPerNight': 500,
177  'modelMagErrors': config.modelMagErrors,
178  'instrumentParsPerBand': config.instrumentParsPerBand,
179  'instrumentSlopeMinDeltaT': config.instrumentSlopeMinDeltaT,
180  'fitMirrorChromaticity': config.fitMirrorChromaticity,
181  'useRepeatabilityForExpGrayCutsDict': dict(config.useRepeatabilityForExpGrayCutsDict),
182  'autoPhotometricCutNSig': config.autoPhotometricCutNSig,
183  'autoHighCutNSig': config.autoHighCutNSig,
184  'printOnly': False,
185  'quietMode': config.quietMode,
186  'outputStars': False,
187  'clobber': True,
188  'useSedLUT': False,
189  'resetParameters': resetFitParameters,
190  'outputFgcmcalZpts': True, # when outputting zpts, use fgcmcal format
191  'outputZeropoints': outputZeropoints}
192 
193  return configDict
194 
195 

◆ makeStdCat()

def lsst.fgcmcal.utilities.makeStdCat (   stdSchema,
  stdStruct,
  goodBands 
)
Make the standard star catalog for persistence

Parameters
----------
stdSchema: `lsst.afw.table.Schema`
   Standard star catalog schema
stdStruct: `numpy.ndarray`
   Standard star structure in FGCM format
goodBands: `list`
   List of good band names used in stdStruct

Returns
-------
stdCat: `lsst.afw.table.BaseCatalog`
   Standard star catalog for persistence

Definition at line 742 of file utilities.py.

742 def makeStdCat(stdSchema, stdStruct, goodBands):
743  """
744  Make the standard star catalog for persistence
745 
746  Parameters
747  ----------
748  stdSchema: `lsst.afw.table.Schema`
749  Standard star catalog schema
750  stdStruct: `numpy.ndarray`
751  Standard star structure in FGCM format
752  goodBands: `list`
753  List of good band names used in stdStruct
754 
755  Returns
756  -------
757  stdCat: `lsst.afw.table.BaseCatalog`
758  Standard star catalog for persistence
759  """
760 
761  stdCat = afwTable.SimpleCatalog(stdSchema)
762 
763  stdCat.reserve(stdStruct.size)
764  for i in range(stdStruct.size):
765  stdCat.addNew()
766 
767  stdCat['id'][:] = stdStruct['FGCM_ID']
768  stdCat['coord_ra'][:] = stdStruct['RA'] * geom.degrees
769  stdCat['coord_dec'][:] = stdStruct['DEC'] * geom.degrees
770  stdCat['ngood'][:, :] = stdStruct['NGOOD'][:, :]
771  stdCat['ntotal'][:, :] = stdStruct['NTOTAL'][:, :]
772  stdCat['mag_std_noabs'][:, :] = stdStruct['MAG_STD'][:, :]
773  stdCat['magErr_std'][:, :] = stdStruct['MAGERR_STD'][:, :]
774  stdCat['npsfcand'][:, :] = stdStruct['NPSFCAND'][:, :]
775 
776  md = PropertyList()
777  md.set("BANDS", list(goodBands))
778  stdCat.setMetadata(md)
779 
780  return stdCat
781 
782 

◆ makeStdSchema()

def lsst.fgcmcal.utilities.makeStdSchema (   nBands)
Make the standard star schema

Parameters
----------
nBands: `int`
   Number of bands in standard star catalog

Returns
-------
stdSchema: `lsst.afw.table.Schema`

Definition at line 710 of file utilities.py.

710 def makeStdSchema(nBands):
711  """
712  Make the standard star schema
713 
714  Parameters
715  ----------
716  nBands: `int`
717  Number of bands in standard star catalog
718 
719  Returns
720  -------
721  stdSchema: `lsst.afw.table.Schema`
722  """
723 
724  stdSchema = afwTable.SimpleTable.makeMinimalSchema()
725  stdSchema.addField('ngood', type='ArrayI', doc='Number of good observations',
726  size=nBands)
727  stdSchema.addField('ntotal', type='ArrayI', doc='Number of total observations',
728  size=nBands)
729  stdSchema.addField('mag_std_noabs', type='ArrayF',
730  doc='Standard magnitude (no absolute calibration)',
731  size=nBands)
732  stdSchema.addField('magErr_std', type='ArrayF',
733  doc='Standard magnitude error',
734  size=nBands)
735  stdSchema.addField('npsfcand', type='ArrayI',
736  doc='Number of observations flagged as psf candidates',
737  size=nBands)
738 
739  return stdSchema
740 
741 

◆ makeZptCat()

def lsst.fgcmcal.utilities.makeZptCat (   zptSchema,
  zpStruct 
)
Make the zeropoint catalog for persistence

Parameters
----------
zptSchema: `lsst.afw.table.Schema`
   Zeropoint catalog schema
zpStruct: `numpy.ndarray`
   Zeropoint structure from fgcm

Returns
-------
zptCat: `afwTable.BaseCatalog`
   Zeropoint catalog for persistence

Definition at line 598 of file utilities.py.

598 def makeZptCat(zptSchema, zpStruct):
599  """
600  Make the zeropoint catalog for persistence
601 
602  Parameters
603  ----------
604  zptSchema: `lsst.afw.table.Schema`
605  Zeropoint catalog schema
606  zpStruct: `numpy.ndarray`
607  Zeropoint structure from fgcm
608 
609  Returns
610  -------
611  zptCat: `afwTable.BaseCatalog`
612  Zeropoint catalog for persistence
613  """
614 
615  zptCat = afwTable.BaseCatalog(zptSchema)
616  zptCat.reserve(zpStruct.size)
617 
618  for filterName in zpStruct['FILTERNAME']:
619  rec = zptCat.addNew()
620  rec['filtername'] = filterName.decode('utf-8')
621 
622  zptCat['visit'][:] = zpStruct['VISIT']
623  zptCat['ccd'][:] = zpStruct['CCD']
624  zptCat['fgcmFlag'][:] = zpStruct['FGCM_FLAG']
625  zptCat['fgcmZpt'][:] = zpStruct['FGCM_ZPT']
626  zptCat['fgcmZptErr'][:] = zpStruct['FGCM_ZPTERR']
627  zptCat['fgcmfZptChebXyMax'][:, :] = zpStruct['FGCM_FZPT_XYMAX']
628  zptCat['fgcmfZptCheb'][:, :] = zpStruct['FGCM_FZPT_CHEB']
629  zptCat['fgcmfZptSstarCheb'][:, :] = zpStruct['FGCM_FZPT_SSTAR_CHEB']
630  zptCat['fgcmI0'][:] = zpStruct['FGCM_I0']
631  zptCat['fgcmI10'][:] = zpStruct['FGCM_I10']
632  zptCat['fgcmR0'][:] = zpStruct['FGCM_R0']
633  zptCat['fgcmR10'][:] = zpStruct['FGCM_R10']
634  zptCat['fgcmGry'][:] = zpStruct['FGCM_GRY']
635  zptCat['fgcmZptVar'][:] = zpStruct['FGCM_ZPTVAR']
636  zptCat['fgcmTilings'][:] = zpStruct['FGCM_TILINGS']
637  zptCat['fgcmFpGry'][:] = zpStruct['FGCM_FPGRY']
638  zptCat['fgcmFpGryBlue'][:] = zpStruct['FGCM_FPGRY_CSPLIT'][:, 0]
639  zptCat['fgcmFpGryBlueErr'][:] = zpStruct['FGCM_FPGRY_CSPLITERR'][:, 0]
640  zptCat['fgcmFpGryRed'][:] = zpStruct['FGCM_FPGRY_CSPLIT'][:, 2]
641  zptCat['fgcmFpGryRedErr'][:] = zpStruct['FGCM_FPGRY_CSPLITERR'][:, 2]
642  zptCat['fgcmFpVar'][:] = zpStruct['FGCM_FPVAR']
643  zptCat['fgcmDust'][:] = zpStruct['FGCM_DUST']
644  zptCat['fgcmFlat'][:] = zpStruct['FGCM_FLAT']
645  zptCat['fgcmAperCorr'][:] = zpStruct['FGCM_APERCORR']
646  zptCat['exptime'][:] = zpStruct['EXPTIME']
647 
648  return zptCat
649 
650 

◆ makeZptSchema()

def lsst.fgcmcal.utilities.makeZptSchema (   superStarChebyshevSize,
  zptChebyshevSize 
)
Make the zeropoint schema

Parameters
----------
superStarChebyshevSize: `int`
   Length of the superstar chebyshev array
zptChebyshevSize: `int`
   Length of the zeropoint chebyshev array

Returns
-------
zptSchema: `lsst.afw.table.schema`

Definition at line 519 of file utilities.py.

519 def makeZptSchema(superStarChebyshevSize, zptChebyshevSize):
520  """
521  Make the zeropoint schema
522 
523  Parameters
524  ----------
525  superStarChebyshevSize: `int`
526  Length of the superstar chebyshev array
527  zptChebyshevSize: `int`
528  Length of the zeropoint chebyshev array
529 
530  Returns
531  -------
532  zptSchema: `lsst.afw.table.schema`
533  """
534 
535  zptSchema = afwTable.Schema()
536 
537  zptSchema.addField('visit', type=np.int32, doc='Visit number')
538  zptSchema.addField('ccd', type=np.int32, doc='CCD number')
539  zptSchema.addField('fgcmFlag', type=np.int32, doc=('FGCM flag value: '
540  '1: Photometric, used in fit; '
541  '2: Photometric, not used in fit; '
542  '4: Non-photometric, on partly photometric night; '
543  '8: Non-photometric, on non-photometric night; '
544  '16: No zeropoint could be determined; '
545  '32: Too few stars for reliable gray computation'))
546  zptSchema.addField('fgcmZpt', type=np.float64, doc='FGCM zeropoint (center of CCD)')
547  zptSchema.addField('fgcmZptErr', type=np.float64,
548  doc='Error on zeropoint, estimated from repeatability + number of obs')
549  zptSchema.addField('fgcmfZptChebXyMax', type='ArrayD', size=2,
550  doc='maximum x/maximum y to scale to apply chebyshev parameters')
551  zptSchema.addField('fgcmfZptCheb', type='ArrayD',
552  size=zptChebyshevSize,
553  doc='Chebyshev parameters (flattened) for zeropoint')
554  zptSchema.addField('fgcmfZptSstarCheb', type='ArrayD',
555  size=superStarChebyshevSize,
556  doc='Chebyshev parameters (flattened) for superStarFlat')
557  zptSchema.addField('fgcmI0', type=np.float64, doc='Integral of the passband')
558  zptSchema.addField('fgcmI10', type=np.float64, doc='Normalized chromatic integral')
559  zptSchema.addField('fgcmR0', type=np.float64,
560  doc='Retrieved i0 integral, estimated from stars (only for flag 1)')
561  zptSchema.addField('fgcmR10', type=np.float64,
562  doc='Retrieved i10 integral, estimated from stars (only for flag 1)')
563  zptSchema.addField('fgcmGry', type=np.float64,
564  doc='Estimated gray extinction relative to atmospheric solution; '
565  'only for flag <= 4')
566  zptSchema.addField('fgcmZptVar', type=np.float64, doc='Variance of zeropoint over ccd')
567  zptSchema.addField('fgcmTilings', type=np.float64,
568  doc='Number of photometric tilings used for solution for ccd')
569  zptSchema.addField('fgcmFpGry', type=np.float64,
570  doc='Average gray extinction over the full focal plane '
571  '(same for all ccds in a visit)')
572  zptSchema.addField('fgcmFpGryBlue', type=np.float64,
573  doc='Average gray extinction over the full focal plane '
574  'for 25% bluest stars')
575  zptSchema.addField('fgcmFpGryBlueErr', type=np.float64,
576  doc='Error on Average gray extinction over the full focal plane '
577  'for 25% bluest stars')
578  zptSchema.addField('fgcmFpGryRed', type=np.float64,
579  doc='Average gray extinction over the full focal plane '
580  'for 25% reddest stars')
581  zptSchema.addField('fgcmFpGryRedErr', type=np.float64,
582  doc='Error on Average gray extinction over the full focal plane '
583  'for 25% reddest stars')
584  zptSchema.addField('fgcmFpVar', type=np.float64,
585  doc='Variance of gray extinction over the full focal plane '
586  '(same for all ccds in a visit)')
587  zptSchema.addField('fgcmDust', type=np.float64,
588  doc='Gray dust extinction from the primary/corrector'
589  'at the time of the exposure')
590  zptSchema.addField('fgcmFlat', type=np.float64, doc='Superstarflat illumination correction')
591  zptSchema.addField('fgcmAperCorr', type=np.float64, doc='Aperture correction estimated by fgcm')
592  zptSchema.addField('exptime', type=np.float32, doc='Exposure time')
593  zptSchema.addField('filtername', type=str, size=10, doc='Filter name')
594 
595  return zptSchema
596 
597 

◆ translateFgcmLut()

def lsst.fgcmcal.utilities.translateFgcmLut (   lutCat,
  filterMap 
)
Translate the FGCM look-up-table into an fgcm-compatible object

Parameters
----------
lutCat: `lsst.afw.table.BaseCatalog`
   Catalog describing the FGCM look-up table
filterMap: `dict`
   Filter to band mapping

Returns
-------
fgcmLut: `lsst.fgcm.FgcmLut`
   Lookup table for FGCM
lutIndexVals: `numpy.ndarray`
   Numpy array with LUT index information for FGCM
lutStd: `numpy.ndarray`
   Numpy array with LUT standard throughput values for FGCM

Notes
-----
After running this code, it is wise to `del lutCat` to clear the memory.

Definition at line 196 of file utilities.py.

196 def translateFgcmLut(lutCat, filterMap):
197  """
198  Translate the FGCM look-up-table into an fgcm-compatible object
199 
200  Parameters
201  ----------
202  lutCat: `lsst.afw.table.BaseCatalog`
203  Catalog describing the FGCM look-up table
204  filterMap: `dict`
205  Filter to band mapping
206 
207  Returns
208  -------
209  fgcmLut: `lsst.fgcm.FgcmLut`
210  Lookup table for FGCM
211  lutIndexVals: `numpy.ndarray`
212  Numpy array with LUT index information for FGCM
213  lutStd: `numpy.ndarray`
214  Numpy array with LUT standard throughput values for FGCM
215 
216  Notes
217  -----
218  After running this code, it is wise to `del lutCat` to clear the memory.
219  """
220 
221  # first we need the lutIndexVals
222  # dtype is set for py2/py3/fits/fgcm compatibility
223  lutFilterNames = np.array(lutCat[0]['filterNames'].split(','), dtype='a')
224  lutStdFilterNames = np.array(lutCat[0]['stdFilterNames'].split(','), dtype='a')
225 
226  # Note that any discrepancies between config values will raise relevant
227  # exceptions in the FGCM code.
228 
229  lutIndexVals = np.zeros(1, dtype=[('FILTERNAMES', lutFilterNames.dtype.str,
230  lutFilterNames.size),
231  ('STDFILTERNAMES', lutStdFilterNames.dtype.str,
232  lutStdFilterNames.size),
233  ('PMB', 'f8', lutCat[0]['pmb'].size),
234  ('PMBFACTOR', 'f8', lutCat[0]['pmbFactor'].size),
235  ('PMBELEVATION', 'f8'),
236  ('LAMBDANORM', 'f8'),
237  ('PWV', 'f8', lutCat[0]['pwv'].size),
238  ('O3', 'f8', lutCat[0]['o3'].size),
239  ('TAU', 'f8', lutCat[0]['tau'].size),
240  ('ALPHA', 'f8', lutCat[0]['alpha'].size),
241  ('ZENITH', 'f8', lutCat[0]['zenith'].size),
242  ('NCCD', 'i4')])
243 
244  lutIndexVals['FILTERNAMES'][:] = lutFilterNames
245  lutIndexVals['STDFILTERNAMES'][:] = lutStdFilterNames
246  lutIndexVals['PMB'][:] = lutCat[0]['pmb']
247  lutIndexVals['PMBFACTOR'][:] = lutCat[0]['pmbFactor']
248  lutIndexVals['PMBELEVATION'] = lutCat[0]['pmbElevation']
249  lutIndexVals['LAMBDANORM'] = lutCat[0]['lambdaNorm']
250  lutIndexVals['PWV'][:] = lutCat[0]['pwv']
251  lutIndexVals['O3'][:] = lutCat[0]['o3']
252  lutIndexVals['TAU'][:] = lutCat[0]['tau']
253  lutIndexVals['ALPHA'][:] = lutCat[0]['alpha']
254  lutIndexVals['ZENITH'][:] = lutCat[0]['zenith']
255  lutIndexVals['NCCD'] = lutCat[0]['nCcd']
256 
257  # now we need the Standard Values
258  lutStd = np.zeros(1, dtype=[('PMBSTD', 'f8'),
259  ('PWVSTD', 'f8'),
260  ('O3STD', 'f8'),
261  ('TAUSTD', 'f8'),
262  ('ALPHASTD', 'f8'),
263  ('ZENITHSTD', 'f8'),
264  ('LAMBDARANGE', 'f8', 2),
265  ('LAMBDASTEP', 'f8'),
266  ('LAMBDASTD', 'f8', lutFilterNames.size),
267  ('LAMBDASTDFILTER', 'f8', lutStdFilterNames.size),
268  ('I0STD', 'f8', lutFilterNames.size),
269  ('I1STD', 'f8', lutFilterNames.size),
270  ('I10STD', 'f8', lutFilterNames.size),
271  ('I2STD', 'f8', lutFilterNames.size),
272  ('LAMBDAB', 'f8', lutFilterNames.size),
273  ('ATMLAMBDA', 'f8', lutCat[0]['atmLambda'].size),
274  ('ATMSTDTRANS', 'f8', lutCat[0]['atmStdTrans'].size)])
275  lutStd['PMBSTD'] = lutCat[0]['pmbStd']
276  lutStd['PWVSTD'] = lutCat[0]['pwvStd']
277  lutStd['O3STD'] = lutCat[0]['o3Std']
278  lutStd['TAUSTD'] = lutCat[0]['tauStd']
279  lutStd['ALPHASTD'] = lutCat[0]['alphaStd']
280  lutStd['ZENITHSTD'] = lutCat[0]['zenithStd']
281  lutStd['LAMBDARANGE'][:] = lutCat[0]['lambdaRange'][:]
282  lutStd['LAMBDASTEP'] = lutCat[0]['lambdaStep']
283  lutStd['LAMBDASTD'][:] = lutCat[0]['lambdaStd']
284  lutStd['LAMBDASTDFILTER'][:] = lutCat[0]['lambdaStdFilter']
285  lutStd['I0STD'][:] = lutCat[0]['i0Std']
286  lutStd['I1STD'][:] = lutCat[0]['i1Std']
287  lutStd['I10STD'][:] = lutCat[0]['i10Std']
288  lutStd['I2STD'][:] = lutCat[0]['i2Std']
289  lutStd['LAMBDAB'][:] = lutCat[0]['lambdaB']
290  lutStd['ATMLAMBDA'][:] = lutCat[0]['atmLambda'][:]
291  lutStd['ATMSTDTRANS'][:] = lutCat[0]['atmStdTrans'][:]
292 
293  lutTypes = [row['luttype'] for row in lutCat]
294 
295  # And the flattened look-up-table
296  lutFlat = np.zeros(lutCat[0]['lut'].size, dtype=[('I0', 'f4'),
297  ('I1', 'f4')])
298 
299  lutFlat['I0'][:] = lutCat[lutTypes.index('I0')]['lut'][:]
300  lutFlat['I1'][:] = lutCat[lutTypes.index('I1')]['lut'][:]
301 
302  lutDerivFlat = np.zeros(lutCat[0]['lut'].size, dtype=[('D_LNPWV', 'f4'),
303  ('D_O3', 'f4'),
304  ('D_LNTAU', 'f4'),
305  ('D_ALPHA', 'f4'),
306  ('D_SECZENITH', 'f4'),
307  ('D_LNPWV_I1', 'f4'),
308  ('D_O3_I1', 'f4'),
309  ('D_LNTAU_I1', 'f4'),
310  ('D_ALPHA_I1', 'f4'),
311  ('D_SECZENITH_I1', 'f4')])
312 
313  for name in lutDerivFlat.dtype.names:
314  lutDerivFlat[name][:] = lutCat[lutTypes.index(name)]['lut'][:]
315 
316  # The fgcm.FgcmLUT() class copies all the LUT information into special
317  # shared memory objects that will not blow up the memory usage when used
318  # with python multiprocessing. Once all the numbers are copied, the
319  # references to the temporary objects (lutCat, lutFlat, lutDerivFlat)
320  # will fall out of scope and can be cleaned up by the garbage collector.
321  fgcmLut = fgcm.FgcmLUT(lutIndexVals, lutFlat, lutDerivFlat, lutStd,
322  filterToBand=filterMap)
323 
324  return fgcmLut, lutIndexVals, lutStd
325 
326 

◆ translateVisitCatalog()

def lsst.fgcmcal.utilities.translateVisitCatalog (   visitCat)
Translate the FGCM visit catalog to an fgcm-compatible object

Parameters
----------
visitCat: `lsst.afw.table.BaseCatalog`
   FGCM visitCat from `lsst.fgcmcal.FgcmBuildStarsTask`

Returns
-------
fgcmExpInfo: `numpy.ndarray`
   Numpy array for visit information for FGCM

Notes
-----
After running this code, it is wise to `del visitCat` to clear the memory.

Definition at line 327 of file utilities.py.

327 def translateVisitCatalog(visitCat):
328  """
329  Translate the FGCM visit catalog to an fgcm-compatible object
330 
331  Parameters
332  ----------
333  visitCat: `lsst.afw.table.BaseCatalog`
334  FGCM visitCat from `lsst.fgcmcal.FgcmBuildStarsTask`
335 
336  Returns
337  -------
338  fgcmExpInfo: `numpy.ndarray`
339  Numpy array for visit information for FGCM
340 
341  Notes
342  -----
343  After running this code, it is wise to `del visitCat` to clear the memory.
344  """
345 
346  fgcmExpInfo = np.zeros(len(visitCat), dtype=[('VISIT', 'i8'),
347  ('MJD', 'f8'),
348  ('EXPTIME', 'f8'),
349  ('PSFSIGMA', 'f8'),
350  ('DELTA_APER', 'f8'),
351  ('SKYBACKGROUND', 'f8'),
352  ('DEEPFLAG', 'i2'),
353  ('TELHA', 'f8'),
354  ('TELRA', 'f8'),
355  ('TELDEC', 'f8'),
356  ('TELROT', 'f8'),
357  ('PMB', 'f8'),
358  ('FILTERNAME', 'a10')])
359  fgcmExpInfo['VISIT'][:] = visitCat['visit']
360  fgcmExpInfo['MJD'][:] = visitCat['mjd']
361  fgcmExpInfo['EXPTIME'][:] = visitCat['exptime']
362  fgcmExpInfo['DEEPFLAG'][:] = visitCat['deepFlag']
363  fgcmExpInfo['TELHA'][:] = visitCat['telha']
364  fgcmExpInfo['TELRA'][:] = visitCat['telra']
365  fgcmExpInfo['TELDEC'][:] = visitCat['teldec']
366  fgcmExpInfo['TELROT'][:] = visitCat['telrot']
367  fgcmExpInfo['PMB'][:] = visitCat['pmb']
368  fgcmExpInfo['PSFSIGMA'][:] = visitCat['psfSigma']
369  fgcmExpInfo['DELTA_APER'][:] = visitCat['deltaAper']
370  fgcmExpInfo['SKYBACKGROUND'][:] = visitCat['skyBackground']
371  # Note that we have to go through asAstropy() to get a string
372  # array out of an afwTable. This is faster than a row-by-row loop.
373  fgcmExpInfo['FILTERNAME'][:] = visitCat.asAstropy()['filtername']
374 
375  return fgcmExpInfo
376 
377 
lsst::afw::math::PixelAreaBoundedField
A BoundedField that evaluate the pixel area of a SkyWcs in angular units.
Definition: PixelAreaBoundedField.h:47
lsst.fgcmcal.utilities.computeCcdOffsets
def computeCcdOffsets(camera, defaultOrientation)
Definition: utilities.py:378
lsst.fgcmcal.utilities.makeZptCat
def makeZptCat(zptSchema, zpStruct)
Definition: utilities.py:598
lsst.fgcmcal.utilities.computeReferencePixelScale
def computeReferencePixelScale(camera)
Definition: utilities.py:446
lsst.fgcmcal.utilities.makeStdCat
def makeStdCat(stdSchema, stdStruct, goodBands)
Definition: utilities.py:742
lsst.fgcmcal.utilities.makeConfigDict
def makeConfigDict(config, log, camera, maxIter, resetFitParameters, outputZeropoints, tract=None)
Definition: utilities.py:43
lsst::afw::table::Schema
Defines the fields and offsets for a table.
Definition: Schema.h:50
lsst.fgcmcal.utilities.computeApproxPixelAreaFields
def computeApproxPixelAreaFields(camera)
Definition: utilities.py:474
lsst.fgcmcal.utilities.computeApertureRadius
def computeApertureRadius(schema, fluxField)
Definition: utilities.py:783
lsst.fgcmcal.utilities.translateFgcmLut
def translateFgcmLut(lutCat, filterMap)
Definition: utilities.py:196
lsst::afw::image::VisitInfo
Information about a single exposure of an imaging camera.
Definition: VisitInfo.h:68
lsst.fgcmcal.utilities.translateVisitCatalog
def translateVisitCatalog(visitCat)
Definition: utilities.py:327
lsst.fgcmcal.utilities.makeStdSchema
def makeStdSchema(nBands)
Definition: utilities.py:710
lsst.fgcmcal.utilities.makeAtmSchema
def makeAtmSchema()
Definition: utilities.py:651
list
daf::base::PropertyList * list
Definition: fits.cc:913
lsst::geom::Point< double, 2 >
lsst.fgcmcal.utilities.makeZptSchema
def makeZptSchema(superStarChebyshevSize, zptChebyshevSize)
Definition: utilities.py:519
lsst.fgcmcal.utilities.makeAtmCat
def makeAtmCat(atmSchema, atmStruct)
Definition: utilities.py:675
lsst.obs.base.utils.createInitialSkyWcs
def createInitialSkyWcs(visitInfo, detector, flipX=False)
Definition: utils.py:44
lsst::geom::SpherePoint
Point in an unspecified spherical coordinate system.
Definition: SpherePoint.h:57
lsst::afw::table::CatalogT< BaseRecord >
lsst::afw::table::SortedCatalogT
Custom catalog class for record/table subclasses that are guaranteed to have an ID,...
Definition: fwd.h:63
lsst.fgcmcal.utilities.extractReferenceMags
def extractReferenceMags(refStars, bands, filterMap)
Definition: utilities.py:817