LSST Applications g0265f82a02+0e5473021a,g02d81e74bb+bd2ed33bd6,g1470d8bcf6+de7501a2e0,g14a832a312+ff425fae3c,g2079a07aa2+86d27d4dc4,g2305ad1205+91a32aca49,g295015adf3+762506a1ad,g2bbee38e9b+0e5473021a,g337abbeb29+0e5473021a,g3ddfee87b4+c34e8be1fa,g487adcacf7+5fae3daba8,g50ff169b8f+96c6868917,g52b1c1532d+585e252eca,g591dd9f2cf+ea1711114f,g5a732f18d5+53520f316c,g64a986408d+bd2ed33bd6,g858d7b2824+bd2ed33bd6,g8a8a8dda67+585e252eca,g99cad8db69+016a06b37a,g9ddcbc5298+9a081db1e4,ga1e77700b3+15fc3df1f7,ga8c6da7877+ef4e3a5875,gb0e22166c9+60f28cb32d,gb6a65358fc+0e5473021a,gba4ed39666+c2a2e4ac27,gbb8dafda3b+09e12c87ab,gc120e1dc64+bc2e06c061,gc28159a63d+0e5473021a,gcf0d15dbbd+c34e8be1fa,gdaeeff99f8+f9a426f77a,ge6526c86ff+508d0e0a30,ge79ae78c31+0e5473021a,gee10cc3b42+585e252eca,gf18bd8381d+8d59551888,gf1cff7945b+bd2ed33bd6,w.2024.16
LSST Data Management Base Package
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | Static Public Attributes | Protected Member Functions | Protected Attributes | Static Protected Attributes | List of all members
lsst.fgcmcal.fgcmLoadReferenceCatalog.FgcmLoadReferenceCatalogTask Class Reference
Inheritance diagram for lsst.fgcmcal.fgcmLoadReferenceCatalog.FgcmLoadReferenceCatalogTask:

Public Member Functions

 __init__ (self, refObjLoader=None, refCatName=None, **kwargs)
 
 getFgcmReferenceStarsHealpix (self, nside, pixel, filterList, nest=False)
 
 getFgcmReferenceStarsSkyCircle (self, ra, dec, radius, filterList)
 

Public Attributes

 refObjLoader
 
 refCatName
 

Static Public Attributes

 ConfigClass = FgcmLoadReferenceCatalogConfig
 

Protected Member Functions

 _determine_flux_fields (self, center, filterList)
 

Protected Attributes

 _fluxFilters
 
 _fluxFields
 
 _referenceFilter
 

Static Protected Attributes

str _DefaultName = 'fgcmLoadReferenceCatalog'
 

Detailed Description

Load multi-band reference objects from a reference catalog.

Parameters
----------
refObjLoader : `lsst.meas.algorithms.ReferenceObjectLoader`
    Reference object loader.
refCatName : `str`
    Name of reference catalog (for color term lookups).

Definition at line 81 of file fgcmLoadReferenceCatalog.py.

Constructor & Destructor Documentation

◆ __init__()

lsst.fgcmcal.fgcmLoadReferenceCatalog.FgcmLoadReferenceCatalogTask.__init__ ( self,
refObjLoader = None,
refCatName = None,
** kwargs )
Construct an FgcmLoadReferenceCatalogTask

Definition at line 95 of file fgcmLoadReferenceCatalog.py.

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
108 self._referenceFilter = None
109

Member Function Documentation

◆ _determine_flux_fields()

lsst.fgcmcal.fgcmLoadReferenceCatalog.FgcmLoadReferenceCatalogTask._determine_flux_fields ( self,
center,
filterList )
protected
Determine the flux field names for a reference catalog.

Will set self._fluxFields, self._referenceFilter.

Parameters
----------
center: `lsst.geom.SpherePoint`
   The center around which to load test sources.
filterList: `list`
   list of `str` of camera filter names.

Definition at line 273 of file fgcmLoadReferenceCatalog.py.

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)

◆ getFgcmReferenceStarsHealpix()

lsst.fgcmcal.fgcmLoadReferenceCatalog.FgcmLoadReferenceCatalogTask.getFgcmReferenceStarsHealpix ( self,
nside,
pixel,
filterList,
nest = False )
Get a reference catalog that overlaps a healpix pixel, using multiple
filters.  In addition, apply colorterms if available.

Return format is a numpy recarray for use with fgcm, with the format:

dtype = ([('ra', `np.float64`),
          ('dec', `np.float64`),
          ('refMag', `np.float32`, len(filterList)),
          ('refMagErr', `np.float32`, len(filterList)])

Reference magnitudes (AB) will be 99 for non-detections.

Parameters
----------
nside: `int`
   Healpix nside of pixel to load
pixel: `int`
   Healpix pixel of pixel to load
filterList: `list`
   list of `str` of camera filter names.
nest: `bool`, optional
   Is the pixel in nest format?  Default is False.

Returns
-------
fgcmRefCat: `np.recarray`

Definition at line 110 of file fgcmLoadReferenceCatalog.py.

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
Point in an unspecified spherical coordinate system.
Definition SpherePoint.h:57

◆ getFgcmReferenceStarsSkyCircle()

lsst.fgcmcal.fgcmLoadReferenceCatalog.FgcmLoadReferenceCatalogTask.getFgcmReferenceStarsSkyCircle ( self,
ra,
dec,
radius,
filterList )
Get a reference catalog that overlaps a circular sky region, using
multiple filters.  In addition, apply colorterms if available.

Return format is a numpy recarray for use with fgcm.

dtype = ([('ra', `np.float64`),
          ('dec', `np.float64`),
          ('refMag', `np.float32`, len(filterList)),
          ('refMagErr', `np.float32`, len(filterList)])

Reference magnitudes (AB) will be 99 for non-detections.

Parameters
----------
ra: `float`
   ICRS right ascension, degrees.
dec: `float`
   ICRS declination, degrees.
radius: `float`
   Radius to search, degrees.
filterList: `list`
   list of `str` of camera filter names.

Returns
-------
fgcmRefCat: `np.recarray`

Definition at line 164 of file fgcmLoadReferenceCatalog.py.

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

Member Data Documentation

◆ _DefaultName

str lsst.fgcmcal.fgcmLoadReferenceCatalog.FgcmLoadReferenceCatalogTask._DefaultName = 'fgcmLoadReferenceCatalog'
staticprotected

Definition at line 93 of file fgcmLoadReferenceCatalog.py.

◆ _fluxFields

lsst.fgcmcal.fgcmLoadReferenceCatalog.FgcmLoadReferenceCatalogTask._fluxFields
protected

Definition at line 107 of file fgcmLoadReferenceCatalog.py.

◆ _fluxFilters

lsst.fgcmcal.fgcmLoadReferenceCatalog.FgcmLoadReferenceCatalogTask._fluxFilters
protected

Definition at line 106 of file fgcmLoadReferenceCatalog.py.

◆ _referenceFilter

lsst.fgcmcal.fgcmLoadReferenceCatalog.FgcmLoadReferenceCatalogTask._referenceFilter
protected

Definition at line 108 of file fgcmLoadReferenceCatalog.py.

◆ ConfigClass

lsst.fgcmcal.fgcmLoadReferenceCatalog.FgcmLoadReferenceCatalogTask.ConfigClass = FgcmLoadReferenceCatalogConfig
static

Definition at line 92 of file fgcmLoadReferenceCatalog.py.

◆ refCatName

lsst.fgcmcal.fgcmLoadReferenceCatalog.FgcmLoadReferenceCatalogTask.refCatName

Definition at line 100 of file fgcmLoadReferenceCatalog.py.

◆ refObjLoader

lsst.fgcmcal.fgcmLoadReferenceCatalog.FgcmLoadReferenceCatalogTask.refObjLoader

Definition at line 99 of file fgcmLoadReferenceCatalog.py.


The documentation for this class was generated from the following file: