lsst.fgcmcal 22.0.1-18-g5158719+7737614e3e
fgcmLoadReferenceCatalog.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"""Load reference catalog objects for input to FGCM.
24
25This task will load multi-band reference objects and apply color terms (if
26configured). This wrapper around LoadReferenceObjects task also allows loading
27by healpix pixel (the native pixelization of fgcm), and is self-contained so
28the task can be called by third-party code.
29"""
30
31import numpy as np
32import healpy as hp
33from astropy import units
34
35import lsst.pex.config as pexConfig
36import lsst.pipe.base as pipeBase
37from lsst.meas.algorithms import LoadIndexedReferenceObjectsTask, ReferenceSourceSelectorTask
38from lsst.meas.algorithms import getRefFluxField
39from lsst.pipe.tasks.colorterms import ColortermLibrary
40from lsst.afw.image import abMagErrFromFluxErr
41from lsst.meas.algorithms import ReferenceObjectLoader
42
43import lsst.geom
44
45__all__ = ['FgcmLoadReferenceCatalogConfig', 'FgcmLoadReferenceCatalogTask']
46
47
48class FgcmLoadReferenceCatalogConfig(pexConfig.Config):
49 """Config for FgcmLoadReferenceCatalogTask"""
50
51 refObjLoader = pexConfig.ConfigurableField(
52 target=LoadIndexedReferenceObjectsTask,
53 doc="Reference object loader for photometry",
54 )
55 filterMap = pexConfig.DictField(
56 doc="Mapping from physicalFilter label to reference filter name.",
57 keytype=str,
58 itemtype=str,
59 default={},
60 )
61 applyColorTerms = pexConfig.Field(
62 doc=("Apply photometric color terms to reference stars?"
63 "Requires that colorterms be set to a ColorTermLibrary"),
64 dtype=bool,
65 default=True
66 )
67 colorterms = pexConfig.ConfigField(
68 doc="Library of photometric reference catalog name to color term dict.",
69 dtype=ColortermLibrary,
70 )
71 referenceSelector = pexConfig.ConfigurableField(
72 target=ReferenceSourceSelectorTask,
73 doc="Selection of reference sources",
74 )
75
76 def validate(self):
77 super().validate()
78 if not self.filterMapfilterMap:
79 msg = 'Must set filterMap'
80 raise pexConfig.FieldValidationError(FgcmLoadReferenceCatalogConfig.filterMap, self, msg)
81 if self.applyColorTermsapplyColorTerms and len(self.colortermscolorterms.data) == 0:
82 msg = "applyColorTerms=True requires the `colorterms` field be set to a ColortermLibrary."
83 raise pexConfig.FieldValidationError(FgcmLoadReferenceCatalogConfig.colorterms, self, msg)
84
85
86class FgcmLoadReferenceCatalogTask(pipeBase.Task):
87 """
88 Load multi-band reference objects from a reference catalog.
89
90 Parameters
91 ----------
92 butler: `lsst.daf.persistence.Butler`
93 Data butler for reading catalogs
94 """
95 ConfigClass = FgcmLoadReferenceCatalogConfig
96 _DefaultName = 'fgcmLoadReferenceCatalog'
97
98 def __init__(self, butler=None, refObjLoader=None, **kwargs):
99 """Construct an FgcmLoadReferenceCatalogTask
100
101 Parameters
102 ----------
103 butler: `lsst.daf.persistence.Buter`
104 Data butler for reading catalogs.
105 """
106 pipeBase.Task.__init__(self, **kwargs)
107 if refObjLoader is None and butler is not None:
108 self.makeSubtask('refObjLoader', butler=butler)
109 else:
110 self.refObjLoaderrefObjLoader = refObjLoader
111
112 self.makeSubtask('referenceSelector')
113 self._fluxFilters_fluxFilters = None
114 self._fluxFields_fluxFields = None
115 self._referenceFilter_referenceFilter = None
116
117 def getFgcmReferenceStarsHealpix(self, nside, pixel, filterList, nest=False):
118 """
119 Get a reference catalog that overlaps a healpix pixel, using multiple
120 filters. In addition, apply colorterms if available.
121
122 Return format is a numpy recarray for use with fgcm, with the format:
123
124 dtype = ([('ra', `np.float64`),
125 ('dec', `np.float64`),
126 ('refMag', `np.float32`, len(filterList)),
127 ('refMagErr', `np.float32`, len(filterList)])
128
129 Reference magnitudes (AB) will be 99 for non-detections.
130
131 Parameters
132 ----------
133 nside: `int`
134 Healpix nside of pixel to load
135 pixel: `int`
136 Healpix pixel of pixel to load
137 filterList: `list`
138 list of `str` of camera filter names.
139 nest: `bool`, optional
140 Is the pixel in nest format? Default is False.
141
142 Returns
143 -------
144 fgcmRefCat: `np.recarray`
145 """
146
147 # Determine the size of the sky circle to load
148 theta, phi = hp.pix2ang(nside, pixel, nest=nest)
149 center = lsst.geom.SpherePoint(phi * lsst.geom.radians, (np.pi/2. - theta) * lsst.geom.radians)
150
151 corners = hp.boundaries(nside, pixel, step=1, nest=nest)
152 theta_phi = hp.vec2ang(np.transpose(corners))
153
154 radius = 0.0 * lsst.geom.radians
155 for ctheta, cphi in zip(*theta_phi):
156 rad = center.separation(lsst.geom.SpherePoint(cphi * lsst.geom.radians,
157 (np.pi/2. - ctheta) * lsst.geom.radians))
158 if (rad > radius):
159 radius = rad
160
161 # Load the fgcm-format reference catalog
162 fgcmRefCat = self.getFgcmReferenceStarsSkyCirclegetFgcmReferenceStarsSkyCircle(center.getRa().asDegrees(),
163 center.getDec().asDegrees(),
164 radius.asDegrees(),
165 filterList)
166 catPix = hp.ang2pix(nside, np.radians(90.0 - fgcmRefCat['dec']),
167 np.radians(fgcmRefCat['ra']), nest=nest)
168
169 inPix, = np.where(catPix == pixel)
170
171 return fgcmRefCat[inPix]
172
173 def getFgcmReferenceStarsSkyCircle(self, ra, dec, radius, filterList):
174 """
175 Get a reference catalog that overlaps a circular sky region, using
176 multiple filters. In addition, apply colorterms if available.
177
178 Return format is a numpy recarray for use with fgcm.
179
180 dtype = ([('ra', `np.float64`),
181 ('dec', `np.float64`),
182 ('refMag', `np.float32`, len(filterList)),
183 ('refMagErr', `np.float32`, len(filterList)])
184
185 Reference magnitudes (AB) will be 99 for non-detections.
186
187 Parameters
188 ----------
189 ra: `float`
190 ICRS right ascension, degrees.
191 dec: `float`
192 ICRS declination, degrees.
193 radius: `float`
194 Radius to search, degrees.
195 filterList: `list`
196 list of `str` of camera filter names.
197
198 Returns
199 -------
200 fgcmRefCat: `np.recarray`
201 """
202
203 center = lsst.geom.SpherePoint(ra * lsst.geom.degrees, dec * lsst.geom.degrees)
204
205 # Check if we haev previously cached values for the fluxFields
206 if self._fluxFilters_fluxFilters is None or self._fluxFilters_fluxFilters != filterList:
207 self._determine_flux_fields_determine_flux_fields(center, filterList)
208
209 skyCircle = self.refObjLoaderrefObjLoader.loadSkyCircle(center,
210 radius * lsst.geom.degrees,
211 self._referenceFilter_referenceFilter)
212
213 if not skyCircle.refCat.isContiguous():
214 refCat = skyCircle.refCat.copy(deep=True)
215 else:
216 refCat = skyCircle.refCat
217
218 # Select on raw (uncorrected) catalog, where the errors should make more sense
219 goodSources = self.referenceSelector.selectSources(refCat)
220 selected = goodSources.selected
221
222 fgcmRefCat = np.zeros(np.sum(selected), dtype=[('ra', 'f8'),
223 ('dec', 'f8'),
224 ('refMag', 'f4', len(filterList)),
225 ('refMagErr', 'f4', len(filterList))])
226 if fgcmRefCat.size == 0:
227 # Return an empty catalog if we don't have any selected sources
228 return fgcmRefCat
229
230 # The ra/dec native Angle format is radians
231 # We determine the conversion from the native units (typically
232 # radians) to degrees for the first observation. This allows us
233 # to treate ra/dec as numpy arrays rather than Angles, which would
234 # be approximately 600x slower.
235
236 conv = refCat[0]['coord_ra'].asDegrees() / float(refCat[0]['coord_ra'])
237 fgcmRefCat['ra'] = refCat['coord_ra'][selected] * conv
238 fgcmRefCat['dec'] = refCat['coord_dec'][selected] * conv
239
240 # Default (unset) values are 99.0
241 fgcmRefCat['refMag'][:, :] = 99.0
242 fgcmRefCat['refMagErr'][:, :] = 99.0
243
244 if self.config.applyColorTerms:
245 if isinstance(self.refObjLoaderrefObjLoader, ReferenceObjectLoader):
246 # Gen3
247 refCatName = self.refObjLoaderrefObjLoader.config.value.ref_dataset_name
248 else:
249 # Gen2
250 refCatName = self.refObjLoaderrefObjLoader.ref_dataset_name
251
252 for i, (filterName, fluxField) in enumerate(zip(self._fluxFilters_fluxFilters, self._fluxFields_fluxFields)):
253 if fluxField is None:
254 continue
255
256 self.log.debug("Applying color terms for filtername=%r" % (filterName))
257
258 colorterm = self.config.colorterms.getColorterm(filterName, refCatName, doRaise=True)
259
260 refMag, refMagErr = colorterm.getCorrectedMagnitudes(refCat)
261
262 # nan_to_num replaces nans with zeros, and this ensures that we select
263 # magnitudes that both filter out nans and are not very large (corresponding
264 # to very small fluxes), as "99" is a common sentinel for illegal magnitudes.
265
266 good, = np.where((np.nan_to_num(refMag[selected]) < 90.0)
267 & (np.nan_to_num(refMagErr[selected]) < 90.0)
268 & (np.nan_to_num(refMagErr[selected]) > 0.0))
269
270 fgcmRefCat['refMag'][good, i] = refMag[selected][good]
271 fgcmRefCat['refMagErr'][good, i] = refMagErr[selected][good]
272
273 else:
274 # No colorterms
275
276 for i, (filterName, fluxField) in enumerate(zip(self._fluxFilters_fluxFilters, self._fluxFields_fluxFields)):
277 # nan_to_num replaces nans with zeros, and this ensures that we select
278 # fluxes that both filter out nans and are positive.
279 good, = np.where((np.nan_to_num(refCat[fluxField][selected]) > 0.0)
280 & (np.nan_to_num(refCat[fluxField+'Err'][selected]) > 0.0))
281 refMag = (refCat[fluxField][selected][good] * units.nJy).to_value(units.ABmag)
282 refMagErr = abMagErrFromFluxErr(refCat[fluxField+'Err'][selected][good],
283 refCat[fluxField][selected][good])
284 fgcmRefCat['refMag'][good, i] = refMag
285 fgcmRefCat['refMagErr'][good, i] = refMagErr
286
287 return fgcmRefCat
288
289 def _determine_flux_fields(self, center, filterList):
290 """
291 Determine the flux field names for a reference catalog.
292
293 Will set self._fluxFields_fluxFields, self._referenceFilter_referenceFilter.
294
295 Parameters
296 ----------
297 center: `lsst.geom.SpherePoint`
298 The center around which to load test sources.
299 filterList: `list`
300 list of `str` of camera filter names.
301 """
302
303 # Record self._fluxFilters for checks on subsequent calls
304 self._fluxFilters_fluxFilters = filterList
305
306 # Search for a good filter to use to load the reference catalog
307 # via the refObjLoader task which requires a valid filterName
308 foundReferenceFilter = False
309 for filterName in filterList:
310 refFilterName = self.config.filterMap.get(filterName)
311 if refFilterName is None:
312 continue
313
314 try:
315 results = self.refObjLoaderrefObjLoader.loadSkyCircle(center,
316 0.05 * lsst.geom.degrees,
317 refFilterName)
318 foundReferenceFilter = True
319 self._referenceFilter_referenceFilter = refFilterName
320 break
321 except RuntimeError:
322 # This just means that the filterName wasn't listed
323 # in the reference catalog. This is okay.
324 pass
325
326 if not foundReferenceFilter:
327 raise RuntimeError("Could not find any valid flux field(s) %s" %
328 (", ".join(filterList)))
329
330 # Retrieve all the fluxField names
331 self._fluxFields_fluxFields = []
332 for filterName in filterList:
333 fluxField = None
334
335 refFilterName = self.config.filterMap.get(filterName)
336
337 if refFilterName is not None:
338 try:
339 fluxField = getRefFluxField(results.refCat.schema, filterName=refFilterName)
340 except RuntimeError:
341 # This flux field isn't available. Set to None
342 fluxField = None
343
344 if fluxField is None:
345 self.log.warning(f'No reference flux field for camera filter {filterName}')
346
347 self._fluxFields_fluxFields.append(fluxField)
def __init__(self, butler=None, refObjLoader=None, **kwargs)
def getFgcmReferenceStarsHealpix(self, nside, pixel, filterList, nest=False)