LSST Applications g180d380827+78227d2bc4,g2079a07aa2+86d27d4dc4,g2305ad1205+bdd7851fe3,g2bbee38e9b+c6a8a0fb72,g337abbeb29+c6a8a0fb72,g33d1c0ed96+c6a8a0fb72,g3a166c0a6a+c6a8a0fb72,g3d1719c13e+260d7c3927,g3ddfee87b4+723a6db5f3,g487adcacf7+29e55ea757,g50ff169b8f+96c6868917,g52b1c1532d+585e252eca,g591dd9f2cf+9443c4b912,g62aa8f1a4b+7e2ea9cd42,g858d7b2824+260d7c3927,g864b0138d7+8498d97249,g95921f966b+dffe86973d,g991b906543+260d7c3927,g99cad8db69+4809d78dd9,g9c22b2923f+e2510deafe,g9ddcbc5298+9a081db1e4,ga1e77700b3+03d07e1c1f,gb0e22166c9+60f28cb32d,gb23b769143+260d7c3927,gba4ed39666+c2a2e4ac27,gbb8dafda3b+e22341fd87,gbd998247f1+585e252eca,gc120e1dc64+713f94b854,gc28159a63d+c6a8a0fb72,gc3e9b769f7+385ea95214,gcf0d15dbbd+723a6db5f3,gdaeeff99f8+f9a426f77a,ge6526c86ff+fde82a80b9,ge79ae78c31+c6a8a0fb72,gee10cc3b42+585e252eca,w.2024.18
LSST Data Management Base Package
Loading...
Searching...
No Matches
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 hpgeom as hpg
33from astropy import units
34
35import lsst.pex.config as pexConfig
36import lsst.pipe.base as pipeBase
37from lsst.meas.algorithms import ReferenceSourceSelectorTask
38from lsst.meas.algorithms import getRefFluxField
39from lsst.pipe.tasks.colorterms import ColortermLibrary
40from lsst.afw.image import abMagErrFromFluxErr
41
42import lsst.geom
43
44__all__ = ['FgcmLoadReferenceCatalogConfig', 'FgcmLoadReferenceCatalogTask']
45
46
47class FgcmLoadReferenceCatalogConfig(pexConfig.Config):
48 """Config for FgcmLoadReferenceCatalogTask"""
49
50 filterMap = pexConfig.DictField(
51 doc="Mapping from physicalFilter label to reference filter name.",
52 keytype=str,
53 itemtype=str,
54 default={},
55 )
56 applyColorTerms = pexConfig.Field(
57 doc=("Apply photometric color terms to reference stars?"
58 "Requires that colorterms be set to a ColorTermLibrary"),
59 dtype=bool,
60 default=True
61 )
62 colorterms = pexConfig.ConfigField(
63 doc="Library of photometric reference catalog name to color term dict.",
64 dtype=ColortermLibrary,
65 )
66 referenceSelector = pexConfig.ConfigurableField(
67 target=ReferenceSourceSelectorTask,
68 doc="Selection of reference sources",
69 )
70
71 def validate(self):
72 super().validate()
73 if not self.filterMap:
74 msg = 'Must set filterMap'
75 raise pexConfig.FieldValidationError(FgcmLoadReferenceCatalogConfig.filterMap, self, msg)
76 if self.applyColorTerms and len(self.colorterms.data) == 0:
77 msg = "applyColorTerms=True requires the `colorterms` field be set to a ColortermLibrary."
78 raise pexConfig.FieldValidationError(FgcmLoadReferenceCatalogConfig.colorterms, self, msg)
79
80
81class FgcmLoadReferenceCatalogTask(pipeBase.Task):
82 """
83 Load multi-band reference objects from a reference catalog.
84
85 Parameters
86 ----------
87 refObjLoader : `lsst.meas.algorithms.ReferenceObjectLoader`
88 Reference object loader.
89 refCatName : `str`
90 Name of reference catalog (for color term lookups).
91 """
92 ConfigClass = FgcmLoadReferenceCatalogConfig
93 _DefaultName = 'fgcmLoadReferenceCatalog'
94
95 def __init__(self, refObjLoader=None, refCatName=None, **kwargs):
96 """Construct an FgcmLoadReferenceCatalogTask
97 """
98 pipeBase.Task.__init__(self, **kwargs)
99 self.refObjLoader = refObjLoader
100 self.refCatName = refCatName
101
102 if refObjLoader is None or refCatName is None:
103 raise RuntimeError("FgcmLoadReferenceCatalogTask requires a refObjLoader and refCatName.")
104
105 self.makeSubtask('referenceSelector')
106 self._fluxFilters = None
107 self._fluxFields = None
109
110 def getFgcmReferenceStarsHealpix(self, nside, pixel, filterList, nest=False):
111 """
112 Get a reference catalog that overlaps a healpix pixel, using multiple
113 filters. In addition, apply colorterms if available.
114
115 Return format is a numpy recarray for use with fgcm, with the format:
116
117 dtype = ([('ra', `np.float64`),
118 ('dec', `np.float64`),
119 ('refMag', `np.float32`, len(filterList)),
120 ('refMagErr', `np.float32`, len(filterList)])
121
122 Reference magnitudes (AB) will be 99 for non-detections.
123
124 Parameters
125 ----------
126 nside: `int`
127 Healpix nside of pixel to load
128 pixel: `int`
129 Healpix pixel of pixel to load
130 filterList: `list`
131 list of `str` of camera filter names.
132 nest: `bool`, optional
133 Is the pixel in nest format? Default is False.
134
135 Returns
136 -------
137 fgcmRefCat: `np.recarray`
138 """
139
140 # Determine the size of the sky circle to load
141 lon, lat = hpg.pixel_to_angle(nside, pixel, nest=nest, degrees=False)
142 center = lsst.geom.SpherePoint(lon * lsst.geom.radians, lat * lsst.geom.radians)
143
144 theta_phi = hpg.boundaries(nside, pixel, step=1, nest=nest, lonlat=False)
145
146 radius = 0.0 * lsst.geom.radians
147 for ctheta, cphi in zip(*theta_phi):
148 rad = center.separation(lsst.geom.SpherePoint(cphi * lsst.geom.radians,
149 (np.pi/2. - ctheta) * lsst.geom.radians))
150 if (rad > radius):
151 radius = rad
152
153 # Load the fgcm-format reference catalog
154 fgcmRefCat = self.getFgcmReferenceStarsSkyCircle(center.getRa().asDegrees(),
155 center.getDec().asDegrees(),
156 radius.asDegrees(),
157 filterList)
158 catPix = hpg.angle_to_pixel(nside, fgcmRefCat['ra'], fgcmRefCat['dec'], nest=nest)
159
160 inPix, = np.where(catPix == pixel)
161
162 return fgcmRefCat[inPix]
163
164 def getFgcmReferenceStarsSkyCircle(self, ra, dec, radius, filterList):
165 """
166 Get a reference catalog that overlaps a circular sky region, using
167 multiple filters. In addition, apply colorterms if available.
168
169 Return format is a numpy recarray for use with fgcm.
170
171 dtype = ([('ra', `np.float64`),
172 ('dec', `np.float64`),
173 ('refMag', `np.float32`, len(filterList)),
174 ('refMagErr', `np.float32`, len(filterList)])
175
176 Reference magnitudes (AB) will be 99 for non-detections.
177
178 Parameters
179 ----------
180 ra: `float`
181 ICRS right ascension, degrees.
182 dec: `float`
183 ICRS declination, degrees.
184 radius: `float`
185 Radius to search, degrees.
186 filterList: `list`
187 list of `str` of camera filter names.
188
189 Returns
190 -------
191 fgcmRefCat: `np.recarray`
192 """
193
194 center = lsst.geom.SpherePoint(ra * lsst.geom.degrees, dec * lsst.geom.degrees)
195
196 # Check if we haev previously cached values for the fluxFields
197 if self._fluxFilters is None or self._fluxFilters != filterList:
198 self._determine_flux_fields(center, filterList)
199
200 skyCircle = self.refObjLoader.loadSkyCircle(center,
201 radius * lsst.geom.degrees,
202 self._referenceFilter)
203
204 if not skyCircle.refCat.isContiguous():
205 refCat = skyCircle.refCat.copy(deep=True)
206 else:
207 refCat = skyCircle.refCat
208
209 # Select on raw (uncorrected) catalog, where the errors should make more sense
210 goodSources = self.referenceSelector.selectSources(refCat)
211 selected = goodSources.selected
212
213 fgcmRefCat = np.zeros(np.sum(selected), dtype=[('ra', 'f8'),
214 ('dec', 'f8'),
215 ('refMag', 'f4', len(filterList)),
216 ('refMagErr', 'f4', len(filterList))])
217 if fgcmRefCat.size == 0:
218 # Return an empty catalog if we don't have any selected sources
219 return fgcmRefCat
220
221 # The ra/dec native Angle format is radians
222 # We determine the conversion from the native units (typically
223 # radians) to degrees for the first observation. This allows us
224 # to treate ra/dec as numpy arrays rather than Angles, which would
225 # be approximately 600x slower.
226
227 conv = refCat[0]['coord_ra'].asDegrees() / float(refCat[0]['coord_ra'])
228 fgcmRefCat['ra'] = refCat['coord_ra'][selected] * conv
229 fgcmRefCat['dec'] = refCat['coord_dec'][selected] * conv
230
231 # Default (unset) values are 99.0
232 fgcmRefCat['refMag'][:, :] = 99.0
233 fgcmRefCat['refMagErr'][:, :] = 99.0
234
235 if self.config.applyColorTerms:
236 for i, (filterName, fluxField) in enumerate(zip(self._fluxFilters, self._fluxFields)):
237 if fluxField is None:
238 continue
239
240 self.log.debug("Applying color terms for filtername=%r" % (filterName))
241
242 colorterm = self.config.colorterms.getColorterm(filterName, self.refCatName, doRaise=True)
243
244 refMag, refMagErr = colorterm.getCorrectedMagnitudes(refCat)
245
246 # nan_to_num replaces nans with zeros, and this ensures that we select
247 # magnitudes that both filter out nans and are not very large (corresponding
248 # to very small fluxes), as "99" is a common sentinel for illegal magnitudes.
249
250 good, = np.where((np.nan_to_num(refMag[selected]) < 90.0)
251 & (np.nan_to_num(refMagErr[selected]) < 90.0)
252 & (np.nan_to_num(refMagErr[selected]) > 0.0))
253
254 fgcmRefCat['refMag'][good, i] = refMag[selected][good]
255 fgcmRefCat['refMagErr'][good, i] = refMagErr[selected][good]
256
257 else:
258 # No colorterms
259
260 for i, (filterName, fluxField) in enumerate(zip(self._fluxFilters, self._fluxFields)):
261 # nan_to_num replaces nans with zeros, and this ensures that we select
262 # fluxes that both filter out nans and are positive.
263 good, = np.where((np.nan_to_num(refCat[fluxField][selected]) > 0.0)
264 & (np.nan_to_num(refCat[fluxField+'Err'][selected]) > 0.0))
265 refMag = (refCat[fluxField][selected][good] * units.nJy).to_value(units.ABmag)
266 refMagErr = abMagErrFromFluxErr(refCat[fluxField+'Err'][selected][good],
267 refCat[fluxField][selected][good])
268 fgcmRefCat['refMag'][good, i] = refMag
269 fgcmRefCat['refMagErr'][good, i] = refMagErr
270
271 return fgcmRefCat
272
273 def _determine_flux_fields(self, center, filterList):
274 """
275 Determine the flux field names for a reference catalog.
276
277 Will set self._fluxFields, self._referenceFilter.
278
279 Parameters
280 ----------
281 center: `lsst.geom.SpherePoint`
282 The center around which to load test sources.
283 filterList: `list`
284 list of `str` of camera filter names.
285 """
286
287 # Record self._fluxFilters for checks on subsequent calls
288 self._fluxFilters = filterList
289
290 # Search for a good filter to use to load the reference catalog
291 # via the refObjLoader task which requires a valid filterName
292 foundReferenceFilter = False
293 for filterName in filterList:
294 refFilterName = self.config.filterMap.get(filterName)
295 if refFilterName is None:
296 continue
297
298 try:
299 results = self.refObjLoader.loadSkyCircle(center,
300 0.05 * lsst.geom.degrees,
301 refFilterName)
302 foundReferenceFilter = True
303 self._referenceFilter = refFilterName
304 break
305 except RuntimeError:
306 # This just means that the filterName wasn't listed
307 # in the reference catalog. This is okay.
308 pass
309
310 if not foundReferenceFilter:
311 raise RuntimeError("Could not find any valid flux field(s) %s" %
312 (", ".join(filterList)))
313
314 # Retrieve all the fluxField names
315 self._fluxFields = []
316 for filterName in filterList:
317 fluxField = None
318
319 refFilterName = self.config.filterMap.get(filterName)
320
321 if refFilterName is not None:
322 try:
323 fluxField = getRefFluxField(results.refCat.schema, filterName=refFilterName)
324 except RuntimeError:
325 # This flux field isn't available. Set to None
326 fluxField = None
327
328 if fluxField is None:
329 self.log.warning(f'No reference flux field for camera filter {filterName}')
330
331 self._fluxFields.append(fluxField)
__init__(self, refObjLoader=None, refCatName=None, **kwargs)
Point in an unspecified spherical coordinate system.
Definition SpherePoint.h:57