LSST Applications g02d81e74bb+86cf3d8bc9,g180d380827+7a4e862ed4,g2079a07aa2+86d27d4dc4,g2305ad1205+e1ca1c66fa,g29320951ab+012e1474a1,g295015adf3+341ea1ce94,g2bbee38e9b+0e5473021a,g337abbeb29+0e5473021a,g33d1c0ed96+0e5473021a,g3a166c0a6a+0e5473021a,g3ddfee87b4+c429d67c83,g48712c4677+f88676dd22,g487adcacf7+27e1e21933,g50ff169b8f+96c6868917,g52b1c1532d+585e252eca,g591dd9f2cf+b41db86c35,g5a732f18d5+53520f316c,g64a986408d+86cf3d8bc9,g858d7b2824+86cf3d8bc9,g8a8a8dda67+585e252eca,g99cad8db69+84912a7fdc,g9ddcbc5298+9a081db1e4,ga1e77700b3+15fc3df1f7,ga8c6da7877+a2b54eae19,gb0e22166c9+60f28cb32d,gba4ed39666+c2a2e4ac27,gbb8dafda3b+6681f309db,gc120e1dc64+f0fcc2f6d8,gc28159a63d+0e5473021a,gcf0d15dbbd+c429d67c83,gdaeeff99f8+f9a426f77a,ge6526c86ff+0433e6603d,ge79ae78c31+0e5473021a,gee10cc3b42+585e252eca,gff1a9f87cc+86cf3d8bc9,w.2024.17
LSST Data Management Base Package
Loading...
Searching...
No Matches
Functions | Variables
lsst.fgcmcal.utilities Namespace Reference

Functions

 makeConfigDict (config, log, camera, maxIter, resetFitParameters, outputZeropoints, lutFilterNames, tract=None, nCore=1)
 
 translateFgcmLut (lutCat, physicalFilterMap)
 
 translateVisitCatalog (visitCat)
 
 computeReferencePixelScale (camera)
 
 computeApproxPixelAreaFields (camera)
 
 makeZptSchema (superStarChebyshevSize, zptChebyshevSize)
 
 makeZptCat (zptSchema, zpStruct)
 
 makeAtmSchema ()
 
 makeAtmCat (atmSchema, atmStruct)
 
 makeStdSchema (nBands)
 
 makeStdCat (stdSchema, stdStruct, goodBands)
 
 computeApertureRadiusFromName (fluxField)
 
 extractReferenceMags (refStars, bands, filterMap)
 
 lookupStaticCalibrations (datasetType, registry, quantumDataId, collections)
 

Variables

str FGCM_EXP_FIELD = 'VISIT'
 
str FGCM_CCD_FIELD = 'DETECTOR'
 
float FGCM_ILLEGAL_VALUE = -9999.0
 

Detailed Description

Utility functions for fgcmcal.

This file contains utility functions that are used by more than one task,
and do not need to be part of a task.

Function Documentation

◆ computeApertureRadiusFromName()

lsst.fgcmcal.utilities.computeApertureRadiusFromName ( fluxField)
Compute the radius associated with a CircularApertureFlux or ApFlux field.

Parameters
----------
fluxField : `str`
   CircularApertureFlux or ApFlux

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

Raises
------
 RuntimeError: Raised if flux field is not a CircularApertureFlux,
   ApFlux, or apFlux.

Definition at line 777 of file utilities.py.

777def computeApertureRadiusFromName(fluxField):
778 """
779 Compute the radius associated with a CircularApertureFlux or ApFlux field.
780
781 Parameters
782 ----------
783 fluxField : `str`
784 CircularApertureFlux or ApFlux
785
786 Returns
787 -------
788 apertureRadius : `float`
789 Radius of the aperture field, in pixels.
790
791 Raises
792 ------
793 RuntimeError: Raised if flux field is not a CircularApertureFlux,
794 ApFlux, or apFlux.
795 """
796 # TODO: Move this method to more general stack method in DM-25775
797 m = re.search(r'(CircularApertureFlux|ApFlux|apFlux)_(\d+)_(\d+)_', fluxField)
798
799 if m is None:
800 raise RuntimeError(f"Flux field {fluxField} does not correspond to a CircularApertureFlux or ApFlux")
801
802 apertureRadius = float(m.groups()[1]) + float(m.groups()[2])/10.
803
804 return apertureRadius
805
806

◆ computeApproxPixelAreaFields()

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 459 of file utilities.py.

459def computeApproxPixelAreaFields(camera):
460 """
461 Compute the approximate pixel area bounded fields from the camera
462 geometry.
463
464 Parameters
465 ----------
466 camera: `lsst.afw.cameraGeom.Camera`
467
468 Returns
469 -------
470 approxPixelAreaFields: `dict`
471 Dictionary of approximate area fields, keyed with detector ID
472 """
473
474 areaScaling = 1. / computeReferencePixelScale(camera)**2.
475
476 # Generate fake WCSs centered at 180/0 to avoid the RA=0/360 problem,
477 # since we are looking for relative scales
478 boresight = geom.SpherePoint(180.0*geom.degrees, 0.0*geom.degrees)
479
480 flipX = False
481 # Create a temporary visitInfo for input to createInitialSkyWcs
482 # The orientation does not matter for the area computation
483 visitInfo = afwImage.VisitInfo(boresightRaDec=boresight,
484 boresightRotAngle=0.0*geom.degrees,
485 rotType=afwImage.RotType.SKY)
486
487 approxPixelAreaFields = {}
488
489 for i, detector in enumerate(camera):
490 key = detector.getId()
491
492 wcs = createInitialSkyWcs(visitInfo, detector, flipX)
493 bbox = detector.getBBox()
494
495 areaField = afwMath.PixelAreaBoundedField(bbox, wcs,
496 unit=geom.arcseconds, scaling=areaScaling)
497 approxAreaField = afwMath.ChebyshevBoundedField.approximate(areaField)
498
499 approxPixelAreaFields[key] = approxAreaField
500
501 return approxPixelAreaFields
502
503
Information about a single exposure of an imaging camera.
Definition VisitInfo.h:68
A BoundedField that evaluate the pixel area of a SkyWcs in angular units.
Point in an unspecified spherical coordinate system.
Definition SpherePoint.h:57

◆ computeReferencePixelScale()

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 431 of file utilities.py.

431def computeReferencePixelScale(camera):
432 """
433 Compute the median pixel scale in the camera
434
435 Returns
436 -------
437 pixelScale: `float`
438 Average pixel scale (arcsecond) over the camera
439 """
440
441 boresight = geom.SpherePoint(180.0*geom.degrees, 0.0*geom.degrees)
442 orientation = 0.0*geom.degrees
443 flipX = False
444
445 # Create a temporary visitInfo for input to createInitialSkyWcs
446 visitInfo = afwImage.VisitInfo(boresightRaDec=boresight,
447 boresightRotAngle=orientation,
448 rotType=afwImage.RotType.SKY)
449
450 pixelScales = np.zeros(len(camera))
451 for i, detector in enumerate(camera):
452 wcs = createInitialSkyWcs(visitInfo, detector, flipX)
453 pixelScales[i] = wcs.getPixelScale().asArcseconds()
454
455 ok, = np.where(pixelScales > 0.0)
456 return np.median(pixelScales[ok])
457
458

◆ extractReferenceMags()

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

Parameters
----------
refStars : `astropy.table.Table` or `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 807 of file utilities.py.

807def extractReferenceMags(refStars, bands, filterMap):
808 """
809 Extract reference magnitudes from refStars for given bands and
810 associated filterMap.
811
812 Parameters
813 ----------
814 refStars : `astropy.table.Table` or `lsst.afw.table.BaseCatalog`
815 FGCM reference star catalog.
816 bands : `list`
817 List of bands for calibration.
818 filterMap: `dict`
819 FGCM mapping of filter to band.
820
821 Returns
822 -------
823 refMag : `np.ndarray`
824 nstar x nband array of reference magnitudes.
825 refMagErr : `np.ndarray`
826 nstar x nband array of reference magnitude errors.
827 """
828 hasAstropyMeta = False
829 try:
830 meta = refStars.meta
831 hasAstropyMeta = True
832 except AttributeError:
833 meta = refStars.getMetadata()
834
835 if 'FILTERNAMES' in meta:
836 if hasAstropyMeta:
837 filternames = meta['FILTERNAMES']
838 else:
839 filternames = meta.getArray('FILTERNAMES')
840
841 # The reference catalog that fgcm wants has one entry per band
842 # in the config file
843 refMag = np.zeros((len(refStars), len(bands)),
844 dtype=refStars['refMag'].dtype) + 99.0
845 refMagErr = np.zeros_like(refMag) + 99.0
846 for i, filtername in enumerate(filternames):
847 # We are allowed to run the fit configured so that we do not
848 # use every column in the reference catalog.
849 try:
850 band = filterMap[filtername]
851 except KeyError:
852 continue
853 try:
854 ind = bands.index(band)
855 except ValueError:
856 continue
857
858 refMag[:, ind] = refStars['refMag'][:, i]
859 refMagErr[:, ind] = refStars['refMagErr'][:, i]
860 else:
861 raise RuntimeError("FGCM reference stars missing FILTERNAMES metadata.")
862
863 return refMag, refMagErr
864
865

◆ lookupStaticCalibrations()

lsst.fgcmcal.utilities.lookupStaticCalibrations ( datasetType,
registry,
quantumDataId,
collections )

Definition at line 866 of file utilities.py.

866def lookupStaticCalibrations(datasetType, registry, quantumDataId, collections):
867 # For static calibrations, we search with a timespan that has unbounded
868 # begin and end; we'll get an error if there's more than one match (because
869 # then it's not static).
870 timespan = Timespan(begin=None, end=None)
871 result = []
872 # First iterate over all of the data IDs for this dataset type that are
873 # consistent with the quantum data ID.
874 for dataId in registry.queryDataIds(datasetType.dimensions, dataId=quantumDataId):
875 # Find the dataset with this data ID using the unbounded timespan.
876 if ref := registry.findDataset(datasetType, dataId, collections=collections, timespan=timespan):
877 result.append(ref)
878 return result

◆ makeAtmCat()

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 669 of file utilities.py.

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

◆ makeAtmSchema()

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

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

Definition at line 645 of file utilities.py.

645def makeAtmSchema():
646 """
647 Make the atmosphere schema
648
649 Returns
650 -------
651 atmSchema: `lsst.afw.table.Schema`
652 """
653
654 atmSchema = afwTable.Schema()
655
656 atmSchema.addField('visit', type=np.int64, doc='Visit number')
657 atmSchema.addField('pmb', type=np.float64, doc='Barometric pressure (mb)')
658 atmSchema.addField('pwv', type=np.float64, doc='Water vapor (mm)')
659 atmSchema.addField('tau', type=np.float64, doc='Aerosol optical depth')
660 atmSchema.addField('alpha', type=np.float64, doc='Aerosol slope')
661 atmSchema.addField('o3', type=np.float64, doc='Ozone (dobson)')
662 atmSchema.addField('secZenith', type=np.float64, doc='Secant(zenith) (~ airmass)')
663 atmSchema.addField('cTrans', type=np.float64, doc='Transmission correction factor')
664 atmSchema.addField('lamStd', type=np.float64, doc='Wavelength for transmission correction')
665
666 return atmSchema
667
668
Defines the fields and offsets for a table.
Definition Schema.h:51

◆ makeConfigDict()

lsst.fgcmcal.utilities.makeConfigDict ( config,
log,
camera,
maxIter,
resetFitParameters,
outputZeropoints,
lutFilterNames,
tract = None,
nCore = 1 )
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?
lutFilterNames : array-like, `str`
    Array of physical filter names in the LUT.
tract : `int`, optional
    Tract number for extending the output file name for debugging.
    Default is None.
nCore : `int`, optional
    Number of cores to use.

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

Definition at line 47 of file utilities.py.

49 lutFilterNames, tract=None, nCore=1):
50 """
51 Make the FGCM fit cycle configuration dict
52
53 Parameters
54 ----------
55 config : `lsst.fgcmcal.FgcmFitCycleConfig`
56 Configuration object
57 log : `lsst.log.Log`
58 LSST log object
59 camera : `lsst.afw.cameraGeom.Camera`
60 Camera from the butler
61 maxIter : `int`
62 Maximum number of iterations
63 resetFitParameters: `bool`
64 Reset fit parameters before fitting?
65 outputZeropoints : `bool`
66 Compute zeropoints for output?
67 lutFilterNames : array-like, `str`
68 Array of physical filter names in the LUT.
69 tract : `int`, optional
70 Tract number for extending the output file name for debugging.
71 Default is None.
72 nCore : `int`, optional
73 Number of cores to use.
74
75 Returns
76 -------
77 configDict : `dict`
78 Configuration dictionary for fgcm
79 """
80 # Extract the bands that are _not_ being fit for fgcm configuration
81 notFitBands = [b for b in config.bands if b not in config.fitBands]
82
83 # process the starColorCuts
84 starColorCutList = []
85 for ccut in config.starColorCuts:
86 if ccut == 'NO_DATA':
87 # No color cuts to apply.
88 break
89 parts = ccut.split(',')
90 starColorCutList.append([parts[0], parts[1], float(parts[2]), float(parts[3])])
91
92 # process the refStarColorCuts
93 refStarColorCutList = []
94 for ccut in config.refStarColorCuts:
95 if ccut == 'NO_DATA':
96 # No color cuts to apply.
97 break
98 parts = ccut.split(',')
99 refStarColorCutList.append([parts[0], parts[1], float(parts[2]), float(parts[3])])
100
101 # TODO: Having direct access to the mirror area from the camera would be
102 # useful. See DM-16489.
103 # Mirror area in cm**2
104 if config.mirrorArea is None:
105 mirrorArea = np.pi*(camera.telescopeDiameter*100./2.)**2.
106 else:
107 # Convert to square cm.
108 mirrorArea = config.mirrorArea * 100.**2.
109
110 # Get approximate average camera gain:
111 gains = [amp.getGain() for detector in camera for amp in detector.getAmplifiers()]
112 cameraGain = float(np.median(gains))
113
114 # Cut down the filter map to those that are in the LUT
115 filterToBand = {filterName: config.physicalFilterMap[filterName] for
116 filterName in lutFilterNames}
117
118 if tract is None:
119 outfileBase = config.outfileBase
120 else:
121 outfileBase = '%s-%06d' % (config.outfileBase, tract)
122
123 # create a configuration dictionary for fgcmFitCycle
124 configDict = {'outfileBase': outfileBase,
125 'logger': log,
126 'exposureFile': None,
127 'obsFile': None,
128 'indexFile': None,
129 'lutFile': None,
130 'mirrorArea': mirrorArea,
131 'cameraGain': cameraGain,
132 'ccdStartIndex': camera[0].getId(),
133 'expField': FGCM_EXP_FIELD,
134 'ccdField': FGCM_CCD_FIELD,
135 'seeingField': 'DELTA_APER',
136 'fwhmField': 'PSFSIGMA',
137 'skyBrightnessField': 'SKYBACKGROUND',
138 'deepFlag': 'DEEPFLAG', # unused
139 'bands': list(config.bands),
140 'fitBands': list(config.fitBands),
141 'notFitBands': notFitBands,
142 'requiredBands': list(config.requiredBands),
143 'filterToBand': filterToBand,
144 'logLevel': 'INFO',
145 'nCore': nCore,
146 'nStarPerRun': config.nStarPerRun,
147 'nExpPerRun': config.nExpPerRun,
148 'reserveFraction': config.reserveFraction,
149 'freezeStdAtmosphere': config.freezeStdAtmosphere,
150 'precomputeSuperStarInitialCycle': config.precomputeSuperStarInitialCycle,
151 'superStarSubCCDDict': dict(config.superStarSubCcdDict),
152 'superStarSubCCDChebyshevOrder': config.superStarSubCcdChebyshevOrder,
153 'superStarSubCCDTriangular': config.superStarSubCcdTriangular,
154 'superStarSigmaClip': config.superStarSigmaClip,
155 'superStarPlotCCDResiduals': config.superStarPlotCcdResiduals,
156 'focalPlaneSigmaClip': config.focalPlaneSigmaClip,
157 'ccdGraySubCCDDict': dict(config.ccdGraySubCcdDict),
158 'ccdGraySubCCDChebyshevOrder': config.ccdGraySubCcdChebyshevOrder,
159 'ccdGraySubCCDTriangular': config.ccdGraySubCcdTriangular,
160 'ccdGrayFocalPlaneDict': dict(config.ccdGrayFocalPlaneDict),
161 'ccdGrayFocalPlaneChebyshevOrder': config.ccdGrayFocalPlaneChebyshevOrder,
162 'ccdGrayFocalPlaneFitMinCcd': config.ccdGrayFocalPlaneFitMinCcd,
163 'cycleNumber': config.cycleNumber,
164 'maxIter': maxIter,
165 'deltaMagBkgOffsetPercentile': config.deltaMagBkgOffsetPercentile,
166 'deltaMagBkgPerCcd': config.deltaMagBkgPerCcd,
167 'UTBoundary': config.utBoundary,
168 'washMJDs': config.washMjds,
169 'epochMJDs': config.epochMjds,
170 'coatingMJDs': config.coatingMjds,
171 'minObsPerBand': config.minObsPerBand,
172 'latitude': config.latitude,
173 'defaultCameraOrientation': config.defaultCameraOrientation,
174 'brightObsGrayMax': config.brightObsGrayMax,
175 'minStarPerCCD': config.minStarPerCcd,
176 'minCCDPerExp': config.minCcdPerExp,
177 'maxCCDGrayErr': config.maxCcdGrayErr,
178 'minStarPerExp': config.minStarPerExp,
179 'minExpPerNight': config.minExpPerNight,
180 'expGrayInitialCut': config.expGrayInitialCut,
181 'expGrayPhotometricCutDict': dict(config.expGrayPhotometricCutDict),
182 'expGrayHighCutDict': dict(config.expGrayHighCutDict),
183 'expGrayRecoverCut': config.expGrayRecoverCut,
184 'expVarGrayPhotometricCutDict': dict(config.expVarGrayPhotometricCutDict),
185 'expGrayErrRecoverCut': config.expGrayErrRecoverCut,
186 'refStarSnMin': config.refStarSnMin,
187 'refStarOutlierNSig': config.refStarOutlierNSig,
188 'applyRefStarColorCuts': config.applyRefStarColorCuts,
189 'refStarMaxFracUse': config.refStarMaxFracUse,
190 'useExposureReferenceOffset': config.useExposureReferenceOffset,
191 'illegalValue': FGCM_ILLEGAL_VALUE, # internally used by fgcm.
192 'starColorCuts': starColorCutList,
193 'refStarColorCuts': refStarColorCutList,
194 'aperCorrFitNBins': config.aperCorrFitNBins,
195 'aperCorrInputSlopeDict': dict(config.aperCorrInputSlopeDict),
196 'sedBoundaryTermDict': config.sedboundaryterms.toDict()['data'],
197 'sedTermDict': config.sedterms.toDict()['data'],
198 'colorSplitBands': list(config.colorSplitBands),
199 'sigFgcmMaxErr': config.sigFgcmMaxErr,
200 'sigFgcmMaxEGrayDict': dict(config.sigFgcmMaxEGrayDict),
201 'ccdGrayMaxStarErr': config.ccdGrayMaxStarErr,
202 'approxThroughputDict': dict(config.approxThroughputDict),
203 'sigmaCalRange': list(config.sigmaCalRange),
204 'sigmaCalFitPercentile': list(config.sigmaCalFitPercentile),
205 'sigmaCalPlotPercentile': list(config.sigmaCalPlotPercentile),
206 'sigma0Phot': config.sigma0Phot,
207 'mapLongitudeRef': config.mapLongitudeRef,
208 'mapNSide': config.mapNSide,
209 'varNSig': 100.0, # Turn off 'variable star selection' which doesn't work yet
210 'varMinBand': 2,
211 'useRetrievedPwv': False,
212 'useNightlyRetrievedPwv': False,
213 'pwvRetrievalSmoothBlock': 25,
214 'useQuadraticPwv': config.useQuadraticPwv,
215 'useRetrievedTauInit': False,
216 'tauRetrievalMinCCDPerNight': 500,
217 'modelMagErrors': config.modelMagErrors,
218 'instrumentParsPerBand': config.instrumentParsPerBand,
219 'instrumentSlopeMinDeltaT': config.instrumentSlopeMinDeltaT,
220 'fitMirrorChromaticity': config.fitMirrorChromaticity,
221 'fitCCDChromaticityDict': dict(config.fitCcdChromaticityDict),
222 'useRepeatabilityForExpGrayCutsDict': dict(config.useRepeatabilityForExpGrayCutsDict),
223 'autoPhotometricCutNSig': config.autoPhotometricCutNSig,
224 'autoHighCutNSig': config.autoHighCutNSig,
225 'deltaAperInnerRadiusArcsec': config.deltaAperInnerRadiusArcsec,
226 'deltaAperOuterRadiusArcsec': config.deltaAperOuterRadiusArcsec,
227 'deltaAperFitMinNgoodObs': config.deltaAperFitMinNgoodObs,
228 'deltaAperFitPerCcdNx': config.deltaAperFitPerCcdNx,
229 'deltaAperFitPerCcdNy': config.deltaAperFitPerCcdNy,
230 'deltaAperFitSpatialNside': config.deltaAperFitSpatialNside,
231 'doComputeDeltaAperExposures': config.doComputeDeltaAperPerVisit,
232 'doComputeDeltaAperStars': config.doComputeDeltaAperPerStar,
233 'doComputeDeltaAperMap': config.doComputeDeltaAperMap,
234 'doComputeDeltaAperPerCcd': config.doComputeDeltaAperPerCcd,
235 'printOnly': False,
236 'quietMode': config.quietMode,
237 'randomSeed': config.randomSeed,
238 'outputStars': False,
239 'outputPath': os.path.abspath('.'),
240 'clobber': True,
241 'useSedLUT': False,
242 'resetParameters': resetFitParameters,
243 'doPlots': config.doPlots,
244 'outputFgcmcalZpts': True, # when outputting zpts, use fgcmcal format
245 'outputZeropoints': outputZeropoints}
246
247 return configDict
248
249

◆ makeStdCat()

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 737 of file utilities.py.

737def makeStdCat(stdSchema, stdStruct, goodBands):
738 """
739 Make the standard star catalog for persistence
740
741 Parameters
742 ----------
743 stdSchema: `lsst.afw.table.Schema`
744 Standard star catalog schema
745 stdStruct: `numpy.ndarray`
746 Standard star structure in FGCM format
747 goodBands: `list`
748 List of good band names used in stdStruct
749
750 Returns
751 -------
752 stdCat: `lsst.afw.table.BaseCatalog`
753 Standard star catalog for persistence
754 """
755
756 stdCat = afwTable.SimpleCatalog(stdSchema)
757 stdCat.resize(stdStruct.size)
758
759 stdCat['id'][:] = stdStruct['FGCM_ID']
760 stdCat['coord_ra'][:] = stdStruct['RA'] * geom.degrees
761 stdCat['coord_dec'][:] = stdStruct['DEC'] * geom.degrees
762 stdCat['ngood'][:, :] = stdStruct['NGOOD'][:, :]
763 stdCat['ntotal'][:, :] = stdStruct['NTOTAL'][:, :]
764 stdCat['mag_std_noabs'][:, :] = stdStruct['MAG_STD'][:, :]
765 stdCat['magErr_std'][:, :] = stdStruct['MAGERR_STD'][:, :]
766 if 'NPSFCAND' in stdStruct.dtype.names:
767 stdCat['npsfcand'][:, :] = stdStruct['NPSFCAND'][:, :]
768 stdCat['delta_aper'][:, :] = stdStruct['DELTA_APER'][:, :]
769
770 md = PropertyList()
771 md.set("BANDS", list(goodBands))
772 stdCat.setMetadata(md)
773
774 return stdCat
775
776
Custom catalog class for record/table subclasses that are guaranteed to have an ID,...

◆ makeStdSchema()

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 702 of file utilities.py.

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

◆ makeZptCat()

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 590 of file utilities.py.

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

◆ makeZptSchema()

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 504 of file utilities.py.

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

◆ translateFgcmLut()

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

Parameters
----------
lutCat: `lsst.afw.table.BaseCatalog`
   Catalog describing the FGCM look-up table
physicalFilterMap: `dict`
   Physical 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 250 of file utilities.py.

250def translateFgcmLut(lutCat, physicalFilterMap):
251 """
252 Translate the FGCM look-up-table into an fgcm-compatible object
253
254 Parameters
255 ----------
256 lutCat: `lsst.afw.table.BaseCatalog`
257 Catalog describing the FGCM look-up table
258 physicalFilterMap: `dict`
259 Physical filter to band mapping
260
261 Returns
262 -------
263 fgcmLut: `lsst.fgcm.FgcmLut`
264 Lookup table for FGCM
265 lutIndexVals: `numpy.ndarray`
266 Numpy array with LUT index information for FGCM
267 lutStd: `numpy.ndarray`
268 Numpy array with LUT standard throughput values for FGCM
269
270 Notes
271 -----
272 After running this code, it is wise to `del lutCat` to clear the memory.
273 """
274
275 # first we need the lutIndexVals
276 lutFilterNames = np.array(lutCat[0]['physicalFilters'].split(','), dtype='U')
277 lutStdFilterNames = np.array(lutCat[0]['stdPhysicalFilters'].split(','), dtype='U')
278
279 # Note that any discrepancies between config values will raise relevant
280 # exceptions in the FGCM code.
281
282 lutIndexVals = np.zeros(1, dtype=[('FILTERNAMES', lutFilterNames.dtype.str,
283 lutFilterNames.size),
284 ('STDFILTERNAMES', lutStdFilterNames.dtype.str,
285 lutStdFilterNames.size),
286 ('PMB', 'f8', lutCat[0]['pmb'].size),
287 ('PMBFACTOR', 'f8', lutCat[0]['pmbFactor'].size),
288 ('PMBELEVATION', 'f8'),
289 ('LAMBDANORM', 'f8'),
290 ('PWV', 'f8', lutCat[0]['pwv'].size),
291 ('O3', 'f8', lutCat[0]['o3'].size),
292 ('TAU', 'f8', lutCat[0]['tau'].size),
293 ('ALPHA', 'f8', lutCat[0]['alpha'].size),
294 ('ZENITH', 'f8', lutCat[0]['zenith'].size),
295 ('NCCD', 'i4')])
296
297 lutIndexVals['FILTERNAMES'][:] = lutFilterNames
298 lutIndexVals['STDFILTERNAMES'][:] = lutStdFilterNames
299 lutIndexVals['PMB'][:] = lutCat[0]['pmb']
300 lutIndexVals['PMBFACTOR'][:] = lutCat[0]['pmbFactor']
301 lutIndexVals['PMBELEVATION'] = lutCat[0]['pmbElevation']
302 lutIndexVals['LAMBDANORM'] = lutCat[0]['lambdaNorm']
303 lutIndexVals['PWV'][:] = lutCat[0]['pwv']
304 lutIndexVals['O3'][:] = lutCat[0]['o3']
305 lutIndexVals['TAU'][:] = lutCat[0]['tau']
306 lutIndexVals['ALPHA'][:] = lutCat[0]['alpha']
307 lutIndexVals['ZENITH'][:] = lutCat[0]['zenith']
308 lutIndexVals['NCCD'] = lutCat[0]['nCcd']
309
310 # now we need the Standard Values
311 lutStd = np.zeros(1, dtype=[('PMBSTD', 'f8'),
312 ('PWVSTD', 'f8'),
313 ('O3STD', 'f8'),
314 ('TAUSTD', 'f8'),
315 ('ALPHASTD', 'f8'),
316 ('ZENITHSTD', 'f8'),
317 ('LAMBDARANGE', 'f8', 2),
318 ('LAMBDASTEP', 'f8'),
319 ('LAMBDASTD', 'f8', lutFilterNames.size),
320 ('LAMBDASTDFILTER', 'f8', lutStdFilterNames.size),
321 ('I0STD', 'f8', lutFilterNames.size),
322 ('I1STD', 'f8', lutFilterNames.size),
323 ('I10STD', 'f8', lutFilterNames.size),
324 ('I2STD', 'f8', lutFilterNames.size),
325 ('LAMBDAB', 'f8', lutFilterNames.size),
326 ('ATMLAMBDA', 'f8', lutCat[0]['atmLambda'].size),
327 ('ATMSTDTRANS', 'f8', lutCat[0]['atmStdTrans'].size)])
328 lutStd['PMBSTD'] = lutCat[0]['pmbStd']
329 lutStd['PWVSTD'] = lutCat[0]['pwvStd']
330 lutStd['O3STD'] = lutCat[0]['o3Std']
331 lutStd['TAUSTD'] = lutCat[0]['tauStd']
332 lutStd['ALPHASTD'] = lutCat[0]['alphaStd']
333 lutStd['ZENITHSTD'] = lutCat[0]['zenithStd']
334 lutStd['LAMBDARANGE'][:] = lutCat[0]['lambdaRange'][:]
335 lutStd['LAMBDASTEP'] = lutCat[0]['lambdaStep']
336 lutStd['LAMBDASTD'][:] = lutCat[0]['lambdaStd']
337 lutStd['LAMBDASTDFILTER'][:] = lutCat[0]['lambdaStdFilter']
338 lutStd['I0STD'][:] = lutCat[0]['i0Std']
339 lutStd['I1STD'][:] = lutCat[0]['i1Std']
340 lutStd['I10STD'][:] = lutCat[0]['i10Std']
341 lutStd['I2STD'][:] = lutCat[0]['i2Std']
342 lutStd['LAMBDAB'][:] = lutCat[0]['lambdaB']
343 lutStd['ATMLAMBDA'][:] = lutCat[0]['atmLambda'][:]
344 lutStd['ATMSTDTRANS'][:] = lutCat[0]['atmStdTrans'][:]
345
346 lutTypes = [row['luttype'] for row in lutCat]
347
348 # And the flattened look-up-table
349 lutFlat = np.zeros(lutCat[0]['lut'].size, dtype=[('I0', 'f4'),
350 ('I1', 'f4')])
351
352 lutFlat['I0'][:] = lutCat[lutTypes.index('I0')]['lut'][:]
353 lutFlat['I1'][:] = lutCat[lutTypes.index('I1')]['lut'][:]
354
355 lutDerivFlat = np.zeros(lutCat[0]['lut'].size, dtype=[('D_LNPWV', 'f4'),
356 ('D_O3', 'f4'),
357 ('D_LNTAU', 'f4'),
358 ('D_ALPHA', 'f4'),
359 ('D_SECZENITH', 'f4'),
360 ('D_LNPWV_I1', 'f4'),
361 ('D_O3_I1', 'f4'),
362 ('D_LNTAU_I1', 'f4'),
363 ('D_ALPHA_I1', 'f4'),
364 ('D_SECZENITH_I1', 'f4')])
365
366 for name in lutDerivFlat.dtype.names:
367 lutDerivFlat[name][:] = lutCat[lutTypes.index(name)]['lut'][:]
368
369 # The fgcm.FgcmLUT() class copies all the LUT information into special
370 # shared memory objects that will not blow up the memory usage when used
371 # with python multiprocessing. Once all the numbers are copied, the
372 # references to the temporary objects (lutCat, lutFlat, lutDerivFlat)
373 # will fall out of scope and can be cleaned up by the garbage collector.
374 fgcmLut = fgcm.FgcmLUT(lutIndexVals, lutFlat, lutDerivFlat, lutStd,
375 filterToBand=physicalFilterMap)
376
377 return fgcmLut, lutIndexVals, lutStd
378
379

◆ translateVisitCatalog()

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 380 of file utilities.py.

380def translateVisitCatalog(visitCat):
381 """
382 Translate the FGCM visit catalog to an fgcm-compatible object
383
384 Parameters
385 ----------
386 visitCat: `lsst.afw.table.BaseCatalog`
387 FGCM visitCat from `lsst.fgcmcal.FgcmBuildStarsTask`
388
389 Returns
390 -------
391 fgcmExpInfo: `numpy.ndarray`
392 Numpy array for visit information for FGCM
393
394 Notes
395 -----
396 After running this code, it is wise to `del visitCat` to clear the memory.
397 """
398
399 fgcmExpInfo = np.zeros(len(visitCat), dtype=[('VISIT', 'i8'),
400 ('MJD', 'f8'),
401 ('EXPTIME', 'f8'),
402 ('PSFSIGMA', 'f8'),
403 ('DELTA_APER', 'f8'),
404 ('SKYBACKGROUND', 'f8'),
405 ('DEEPFLAG', 'i2'),
406 ('TELHA', 'f8'),
407 ('TELRA', 'f8'),
408 ('TELDEC', 'f8'),
409 ('TELROT', 'f8'),
410 ('PMB', 'f8'),
411 ('FILTERNAME', 'a50')])
412 fgcmExpInfo['VISIT'][:] = visitCat['visit']
413 fgcmExpInfo['MJD'][:] = visitCat['mjd']
414 fgcmExpInfo['EXPTIME'][:] = visitCat['exptime']
415 fgcmExpInfo['DEEPFLAG'][:] = visitCat['deepFlag']
416 fgcmExpInfo['TELHA'][:] = visitCat['telha']
417 fgcmExpInfo['TELRA'][:] = visitCat['telra']
418 fgcmExpInfo['TELDEC'][:] = visitCat['teldec']
419 fgcmExpInfo['TELROT'][:] = visitCat['telrot']
420 fgcmExpInfo['PMB'][:] = visitCat['pmb']
421 fgcmExpInfo['PSFSIGMA'][:] = visitCat['psfSigma']
422 fgcmExpInfo['DELTA_APER'][:] = visitCat['deltaAper']
423 fgcmExpInfo['SKYBACKGROUND'][:] = visitCat['skyBackground']
424 # Note that we have to go through asAstropy() to get a string
425 # array out of an afwTable. This is faster than a row-by-row loop.
426 fgcmExpInfo['FILTERNAME'][:] = visitCat.asAstropy()['physicalFilter']
427
428 return fgcmExpInfo
429
430

Variable Documentation

◆ FGCM_CCD_FIELD

str lsst.fgcmcal.utilities.FGCM_CCD_FIELD = 'DETECTOR'

Definition at line 43 of file utilities.py.

◆ FGCM_EXP_FIELD

str lsst.fgcmcal.utilities.FGCM_EXP_FIELD = 'VISIT'

Definition at line 42 of file utilities.py.

◆ FGCM_ILLEGAL_VALUE

float lsst.fgcmcal.utilities.FGCM_ILLEGAL_VALUE = -9999.0

Definition at line 44 of file utilities.py.