lsst.fgcmcal g1b140448ed+8dc7ba3ca7
Loading...
Searching...
No Matches
fgcmBuildFromIsolatedStars.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"""Build star observations for input to FGCM using sourceTable_visit.
24
25This task finds all the visits and sourceTable_visits in a repository (or a
26subset based on command line parameters) and extracts all the potential
27calibration stars for input into fgcm. This task additionally uses fgcm to
28match star observations into unique stars, and performs as much cleaning of the
29input catalog as possible.
30"""
31import numpy as np
32import esutil
33import hpgeom as hpg
34from smatch.matcher import Matcher
35from astropy.table import Table, vstack
36
37from fgcm.fgcmUtilities import objFlagDict
38
39import lsst.pex.config as pexConfig
40import lsst.pipe.base as pipeBase
41from lsst.pipe.base import connectionTypes
42from lsst.meas.algorithms import ReferenceObjectLoader, LoadReferenceObjectsConfig
43from lsst.pipe.tasks.reserveIsolatedStars import ReserveIsolatedStarsTask
44
45from .fgcmBuildStarsBase import FgcmBuildStarsConfigBase, FgcmBuildStarsBaseTask
46from .utilities import computeApproxPixelAreaFields, computeApertureRadiusFromName
47from .utilities import lookupStaticCalibrations
48
49__all__ = ["FgcmBuildFromIsolatedStarsConfig", "FgcmBuildFromIsolatedStarsTask"]
50
51
52class FgcmBuildFromIsolatedStarsConnections(pipeBase.PipelineTaskConnections,
53 dimensions=("instrument",),
54 defaultTemplates={}):
55 camera = connectionTypes.PrerequisiteInput(
56 doc="Camera instrument",
57 name="camera",
58 storageClass="Camera",
59 dimensions=("instrument",),
60 lookupFunction=lookupStaticCalibrations,
61 isCalibration=True,
62 )
63 fgcm_lookup_table = connectionTypes.PrerequisiteInput(
64 doc=("Atmosphere + instrument look-up-table for FGCM throughput and "
65 "chromatic corrections."),
66 name="fgcmLookUpTable",
67 storageClass="Catalog",
68 dimensions=("instrument",),
69 deferLoad=True,
70 )
71 ref_cat = connectionTypes.PrerequisiteInput(
72 doc="Reference catalog to use for photometric calibration.",
73 name="cal_ref_cat",
74 storageClass="SimpleCatalog",
75 dimensions=("skypix",),
76 deferLoad=True,
77 multiple=True,
78 )
79 isolated_star_cats = pipeBase.connectionTypes.Input(
80 doc=("Catalog of isolated stars with average positions, number of associated "
81 "sources, and indexes to the isolated_star_sources catalogs."),
82 name="isolated_star_cat",
83 storageClass="ArrowAstropy",
84 dimensions=("instrument", "tract", "skymap"),
85 deferLoad=True,
86 multiple=True,
87 )
88 isolated_star_sources = pipeBase.connectionTypes.Input(
89 doc=("Catalog of isolated star sources with sourceIds, and indexes to the "
90 "isolated_star_cats catalogs."),
91 name="isolated_star_sources",
92 storageClass="ArrowAstropy",
93 dimensions=("instrument", "tract", "skymap"),
94 deferLoad=True,
95 multiple=True,
96 )
97 visit_summaries = connectionTypes.Input(
98 doc=("Per-visit consolidated exposure metadata. These catalogs use "
99 "detector id for the id and must be sorted for fast lookups of a "
100 "detector."),
101 name="visitSummary",
102 storageClass="ExposureCatalog",
103 dimensions=("instrument", "visit"),
104 deferLoad=True,
105 multiple=True,
106 )
107 fgcm_visit_catalog = connectionTypes.Output(
108 doc="Catalog of visit information for fgcm.",
109 name="fgcmVisitCatalog",
110 storageClass="Catalog",
111 dimensions=("instrument",),
112 )
113 fgcm_star_observations = connectionTypes.Output(
114 doc="Catalog of star observations for fgcm.",
115 name="fgcm_star_observations",
116 storageClass="ArrowAstropy",
117 dimensions=("instrument",),
118 )
119 fgcm_star_ids = connectionTypes.Output(
120 doc="Catalog of fgcm calibration star IDs.",
121 name="fgcm_star_ids",
122 storageClass="ArrowAstropy",
123 dimensions=("instrument",),
124 )
125 fgcm_reference_stars = connectionTypes.Output(
126 doc="Catalog of fgcm-matched reference stars.",
127 name="fgcm_reference_stars",
128 storageClass="ArrowAstropy",
129 dimensions=("instrument",),
130 )
131
132 def __init__(self, *, config=None):
133 super().__init__(config=config)
134
135 if not config.doReferenceMatches:
136 self.prerequisiteInputs.remove("ref_cat")
137 self.prerequisiteInputs.remove("fgcm_lookup_table")
138 self.outputs.remove("fgcm_reference_stars")
139
140
142 pipelineConnections=FgcmBuildFromIsolatedStarsConnections):
143 """Config for FgcmBuildFromIsolatedStarsTask."""
144 referenceCCD = pexConfig.Field(
145 doc="Reference detector for checking PSF and background.",
146 dtype=int,
147 default=40,
148 )
149 reserve_selection = pexConfig.ConfigurableField(
150 target=ReserveIsolatedStarsTask,
151 doc="Task to select reserved stars.",
152 )
153
154 def setDefaults(self):
155 super().setDefaults()
156
157 self.reserve_selection.reserve_name = "fgcmcal"
158 self.reserve_selection.reserve_fraction = 0.1
159
160 # The names here correspond to the isolated_star_sources.
161 self.instFluxFieldinstFluxField = 'apFlux_12_0_instFlux'
162 self.localBackgroundFluxFieldlocalBackgroundFluxField = 'localBackground_instFlux'
165
166 source_selector = self.sourceSelector["science"]
167 source_selector.setDefaults()
168
169 source_selector.doFlags = True
170 source_selector.doSignalToNoise = True
171 source_selector.doUnresolved = False
172 source_selector.doIsolated = False
173 source_selector.doRequireFiniteRaDec = False
174
175 source_selector.flags.bad = []
176
178 local_background_flag_name = self.localBackgroundFluxFieldlocalBackgroundFluxField[0: -len('instFlux')] + 'flag'
179 source_selector.flags.bad.append(local_background_flag_name)
180
181 source_selector.signalToNoise.minimum = 11.0
182 source_selector.signalToNoise.maximum = 1000.0
183 source_selector.signalToNoise.fluxField = self.instFluxFieldinstFluxField
184 source_selector.signalToNoise.errField = self.instFluxFieldinstFluxField + 'Err'
185
186
188 """Build star catalog for FGCM global calibration, using the isolated star catalogs.
189 """
190 ConfigClass = FgcmBuildFromIsolatedStarsConfig
191 _DefaultName = "fgcmBuildFromIsolatedStars"
192
193 canMultiprocess = False
194
195 def __init__(self, initInputs=None, **kwargs):
196 super().__init__(**kwargs)
197 self.makeSubtask('reserve_selection')
198
199 def runQuantum(self, butlerQC, inputRefs, outputRefs):
200 input_ref_dict = butlerQC.get(inputRefs)
201
202 isolated_star_cat_handles = input_ref_dict["isolated_star_cats"]
203 isolated_star_source_handles = input_ref_dict["isolated_star_sources"]
204
205 isolated_star_cat_handle_dict = {
206 handle.dataId["tract"]: handle for handle in isolated_star_cat_handles
207 }
208 isolated_star_source_handle_dict = {
209 handle.dataId["tract"]: handle for handle in isolated_star_source_handles
210 }
211
212 if len(isolated_star_cat_handle_dict) != len(isolated_star_source_handle_dict):
213 raise RuntimeError("isolated_star_cats and isolate_star_sources must have same length.")
214
215 for tract in isolated_star_cat_handle_dict:
216 if tract not in isolated_star_source_handle_dict:
217 raise RuntimeError(f"tract {tract} in isolated_star_cats but not isolated_star_sources")
218
219 if self.config.doReferenceMatches:
220 lookup_table_handle = input_ref_dict["fgcm_lookup_table"]
221
222 # Prepare the reference catalog loader
223 ref_config = LoadReferenceObjectsConfig()
224 ref_config.filterMap = self.config.fgcmLoadReferenceCatalog.filterMap
225 ref_obj_loader = ReferenceObjectLoader(dataIds=[ref.datasetRef.dataId
226 for ref in inputRefs.ref_cat],
227 refCats=butlerQC.get(inputRefs.ref_cat),
228 name=self.config.connections.ref_cat,
229 log=self.log,
230 config=ref_config)
231 self.makeSubtask('fgcmLoadReferenceCatalog',
232 refObjLoader=ref_obj_loader,
233 refCatName=self.config.connections.ref_cat)
234 else:
235 lookup_table_handle = None
236
237 # The visit summary handles for use with fgcmMakeVisitCatalog must be keyed with
238 # visit, and values are a list with the first value being the visit_summary_handle,
239 # the second is the source catalog (which is not used, but a place-holder is
240 # left for compatibility reasons.)
241 visit_summary_handle_dict = {handle.dataId['visit']: [handle, None] for
242 handle in input_ref_dict['visit_summaries']}
243
244 camera = input_ref_dict["camera"]
245
246 struct = self.run(
247 camera=camera,
248 visit_summary_handle_dict=visit_summary_handle_dict,
249 isolated_star_cat_handle_dict=isolated_star_cat_handle_dict,
250 isolated_star_source_handle_dict=isolated_star_source_handle_dict,
251 lookup_table_handle=lookup_table_handle,
252 )
253
254 butlerQC.put(struct.fgcm_visit_catalog, outputRefs.fgcm_visit_catalog)
255 butlerQC.put(struct.fgcm_star_observations, outputRefs.fgcm_star_observations)
256 butlerQC.put(struct.fgcm_star_ids, outputRefs.fgcm_star_ids)
257 if self.config.doReferenceMatches:
258 butlerQC.put(struct.fgcm_reference_stars, outputRefs.fgcm_reference_stars)
259
260 def run(self, *, camera, visit_summary_handle_dict, isolated_star_cat_handle_dict,
261 isolated_star_source_handle_dict, lookup_table_handle=None):
262 """Run the fgcmBuildFromIsolatedStarsTask.
263
264 Parameters
265 ----------
266 camera : `lsst.afw.cameraGeom.Camera`
267 Camera object.
268 visit_summary_handle_dict : `dict` [`int`, [`lsst.daf.butler.DeferredDatasetHandle`]]
269 Visit summary dataset handles, with the visit as key.
270 isolated_star_cat_handle_dict : `dict` [`int`, `lsst.daf.butler.DeferredDatasetHandle`]
271 Isolated star catalog dataset handles, with the tract as key.
272 isolated_star_source_handle_dict : `dict` [`int`, `lsst.daf.butler.DeferredDatasetHandle`]
273 Isolated star source dataset handles, with the tract as key.
274 lookup_table_handle : `lsst.daf.butler.DeferredDatasetHandle`, optional
275 Data reference to fgcm look-up table (used if matching reference stars).
276
277 Returns
278 -------
279 struct : `lsst.pipe.base.struct`
280 Catalogs for persistence, with attributes:
281
282 ``fgcm_visit_catalog``
283 Catalog of per-visit data (`lsst.afw.table.ExposureCatalog`).
284 ``fgcm_star_observations``
285 Catalog of individual star observations (`astropy.table.Table`).
286 ``fgcm_star_ids``
287 Catalog of per-star ids and positions (`astropy.table.Table`).
288 ``fgcm_reference_stars``
289 Catalog of reference stars matched to fgcm stars (`astropy.table.Table`).
290 """
291 # Compute aperture radius if necessary. This is useful to do now before
292 # any heave lifting has happened (fail early).
293 calib_flux_aperture_radius = None
294 if self.config.doSubtractLocalBackground:
295 try:
296 calib_flux_aperture_radius = computeApertureRadiusFromName(self.config.instFluxField)
297 except RuntimeError as e:
298 raise RuntimeError("Could not determine aperture radius from %s. "
299 "Cannot use doSubtractLocalBackground." %
300 (self.config.instFluxField)) from e
301
302 # Check that we have the lookup_table_handle if we are doing reference matches.
303 if self.config.doReferenceMatches:
304 if lookup_table_handle is None:
305 raise RuntimeError("Must supply lookup_table_handle if config.doReferenceMatches is True.")
306
307 visit_cat = self.fgcmMakeVisitCatalog(camera, visit_summary_handle_dict)
308
309 # Select and concatenate the isolated stars and sources.
310 fgcm_obj, star_obs = self._make_all_star_obs_from_isolated_stars(
311 isolated_star_cat_handle_dict,
312 isolated_star_source_handle_dict,
313 visit_cat,
314 camera,
315 calib_flux_aperture_radius=calib_flux_aperture_radius,
316 )
317
319 visit_cat,
320 star_obs,
321 )
322
323 # Do density down-sampling.
324 self._density_downsample(fgcm_obj, star_obs)
325
326 # Mark reserve stars
327 self._mark_reserve_stars(fgcm_obj)
328
329 # Do reference association.
330 if self.config.doReferenceMatches:
331 fgcm_ref = self._associate_reference_stars(lookup_table_handle, visit_cat, fgcm_obj)
332 else:
333 fgcm_ref = None
334
335 return pipeBase.Struct(
336 fgcm_visit_catalog=visit_cat,
337 fgcm_star_observations=star_obs,
338 fgcm_star_ids=fgcm_obj,
339 fgcm_reference_stars=fgcm_ref,
340 )
341
342 def _make_all_star_obs_from_isolated_stars(
343 self,
344 isolated_star_cat_handle_dict,
345 isolated_star_source_handle_dict,
346 visit_cat,
347 camera,
348 calib_flux_aperture_radius=None,
349 ):
350 """Make all star observations from isolated star catalogs.
351
352 Parameters
353 ----------
354 isolated_star_cat_handle_dict : `dict` [`int`, `lsst.daf.butler.DeferredDatasetHandle`]
355 Isolated star catalog dataset handles, with the tract as key.
356 isolated_star_source_handle_dict : `dict` [`int`, `lsst.daf.butler.DeferredDatasetHandle`]
357 Isolated star source dataset handles, with the tract as key.
358 visit_cat : `lsst.afw.table.ExposureCatalog`
359 Catalog of per-visit data.
360 camera : `lsst.afw.cameraGeom.Camera`
361 The camera object.
362 calib_flux_aperture_radius : `float`, optional
363 Radius for the calibration flux aperture.
364
365 Returns
366 -------
367 fgcm_obj : `astropy.table.Table`
368 Catalog of ids and positions for each star.
369 star_obs : `astropy.table.Table`
370 Catalog of individual star observations.
371 """
372 source_columns = [
373 "sourceId",
374 "visit",
375 "detector",
376 "ra",
377 "decl",
378 "x",
379 "y",
380 "physical_filter",
381 "band",
382 "obj_index",
383 self.config.instFluxField,
384 self.config.instFluxField + "Err",
385 self.config.apertureInnerInstFluxField,
386 self.config.apertureInnerInstFluxField + "Err",
387 self.config.apertureOuterInstFluxField,
388 self.config.apertureOuterInstFluxField + "Err",
389 ]
390
391 if self.config.doSubtractLocalBackground:
392 source_columns.append(self.config.localBackgroundFluxField)
393
394 if self.sourceSelector.config.doFlags:
395 source_columns.extend(self.sourceSelector.config.flags.bad)
396
397 local_background_area = np.pi*calib_flux_aperture_radius**2.
398
399 # Compute the approximate pixel area over the full focal plane
400 # from the WCS jacobian using the camera model.
401 approx_pixel_area_fields = computeApproxPixelAreaFields(camera)
402
403 # Construct mapping from detector number to index.
404 detector_mapping = {}
405 for detector_index, detector in enumerate(camera):
406 detector_mapping[detector.getId()] = detector_index
407
408 star_obs_dtype = [
409 ("ra", "f8"),
410 ("dec", "f8"),
411 ("x", "f8"),
412 ("y", "f8"),
413 ("visit", "i8"),
414 ("detector", "i4"),
415 ("inst_mag", "f4"),
416 ("inst_mag_err", "f4"),
417 ("jacobian", "f4"),
418 ("delta_mag_bkg", "f4"),
419 ("delta_mag_aper", "f4"),
420 ("delta_mag_err_aper", "f4"),
421 ]
422
423 fgcm_obj_dtype = [
424 ("fgcm_id", "i8"),
425 ("isolated_star_id", "i8"),
426 ("ra", "f8"),
427 ("dec", "f8"),
428 ("obs_arr_index", "i4"),
429 ("n_obs", "i4"),
430 ("obj_flag", "i4"),
431 ]
432
433 fgcm_objs = []
434 star_obs_cats = []
435 merge_source_counter = 0
436
437 k = 2.5/np.log(10.)
438
439 visit_cat_table = visit_cat.asAstropy()
440
441 for tract in sorted(isolated_star_cat_handle_dict):
442 stars = isolated_star_cat_handle_dict[tract].get()
443 sources = isolated_star_source_handle_dict[tract].get(parameters={"columns": source_columns})
444
445 # Down-select sources.
446 good_sources = self.sourceSelector.selectSources(sources).selected
447 if self.config.doSubtractLocalBackground:
448 good_sources &= (~sources["localBackground_flag"])
449 local_background = local_background_area*sources[self.config.localBackgroundFluxField]
450 good_sources &= ((sources[self.config.instFluxField] - local_background) > 0)
451
452 if good_sources.sum() == 0:
453 self.log.info("No good sources found in tract %d", tract)
454 continue
455
456 # Need to re-sum numbers of observations of each star in the required bands.
457 if len(self.config.requiredBands) > 0:
458 n_req = np.zeros((len(self.config.requiredBands), len(stars)), dtype=np.int32)
459 for i, band in enumerate(self.config.requiredBands):
460 (band_use,) = (sources[good_sources]["band"] == band).nonzero()
461 np.add.at(
462 n_req,
463 (i, sources[good_sources]["obj_index"][band_use]),
464 1,
465 )
466 (good_stars,) = (n_req.min(axis=0) >= self.config.minPerBand).nonzero()
467 else:
468 good_stars = np.arange(len(stars))
469
470 # With the following matching:
471 # sources[good_sources][b] <-> stars[good_stars[a]]
472 obj_index = sources["obj_index"][good_sources]
473 a, b = esutil.numpy_util.match(good_stars, obj_index)
474
475 # Update indexes and cut to selected stars/sources
476 _, index_new = np.unique(a, return_index=True)
477 stars["source_cat_index"][good_stars] = index_new
478 sources = sources[good_sources][b]
479 sources["obj_index"][:] = a
480 stars = stars[good_stars]
481
482 nsource = np.zeros(len(stars), dtype=np.int32)
483 np.add.at(
484 nsource,
485 sources["obj_index"],
486 1,
487 )
488 stars["nsource"][:] = nsource
489
490 # After these cuts, the catalogs have the following properties:
491 # - ``stars`` only contains isolated stars that have the minimum number of good
492 # sources in the required bands.
493 # - ``sources`` has been cut to the good sources.
494 # - The slice [stars["source_cat_index"]: stars["source_cat_index"]
495 # + stars["nsource"]]
496 # applied to ``sources`` will give all the sources associated with the star.
497 # - For each source, sources["obj_index"] points to the index of the associated
498 # isolated star.
499
500 # We now reformat the sources and compute the ``observations`` that fgcm expects.
501 star_obs = Table(data=np.zeros(len(sources), dtype=star_obs_dtype))
502 star_obs["ra"] = sources["ra"]
503 star_obs["dec"] = sources["decl"]
504 star_obs["x"] = sources["x"]
505 star_obs["y"] = sources["y"]
506 star_obs["visit"] = sources["visit"]
507 star_obs["detector"] = sources["detector"]
508
509 visit_match, obs_match = esutil.numpy_util.match(visit_cat_table["visit"], sources["visit"])
510
511 exp_time = np.zeros(len(star_obs))
512 exp_time[obs_match] = visit_cat_table["exptime"][visit_match]
513
514 with np.warnings.catch_warnings():
515 # Ignore warnings, we will filter infinities and nans below.
516 np.warnings.simplefilter("ignore")
517
518 inst_mag_inner = -2.5*np.log10(sources[self.config.apertureInnerInstFluxField])
519 inst_mag_err_inner = k*(sources[self.config.apertureInnerInstFluxField + "Err"]
520 / sources[self.config.apertureInnerInstFluxField])
521 inst_mag_outer = -2.5*np.log10(sources[self.config.apertureOuterInstFluxField])
522 inst_mag_err_outer = k*(sources[self.config.apertureOuterInstFluxField + "Err"]
523 / sources[self.config.apertureOuterInstFluxField])
524 star_obs["delta_mag_aper"] = inst_mag_inner - inst_mag_outer
525 star_obs["delta_mag_err_aper"] = np.sqrt(inst_mag_err_inner**2. + inst_mag_err_outer**2.)
526 # Set bad values to sentinel values for fgcm.
527 bad = ~np.isfinite(star_obs["delta_mag_aper"])
528 star_obs["delta_mag_aper"][bad] = 99.0
529 star_obs["delta_mag_err_aper"][bad] = 99.0
530
531 if self.config.doSubtractLocalBackground:
532 # At the moment we only adjust the flux and not the flux
533 # error by the background because the error on
534 # base_LocalBackground_instFlux is the rms error in the
535 # background annulus, not the error on the mean in the
536 # background estimate (which is much smaller, by sqrt(n)
537 # pixels used to estimate the background, which we do not
538 # have access to in this task). In the default settings,
539 # the annulus is sufficiently large such that these
540 # additional errors are negligibly small (much less
541 # than a mmag in quadrature).
542
543 # This is the difference between the mag with local background correction
544 # and the mag without local background correction.
545 local_background = local_background_area*sources[self.config.localBackgroundFluxField]
546 star_obs["delta_mag_bkg"] = (-2.5*np.log10(sources[self.config.instFluxField]
547 - local_background) -
548 -2.5*np.log10(sources[self.config.instFluxField]))
549
550 # Need to loop over detectors here.
551 for detector in camera:
552 detector_id = detector.getId()
553 # used index for all observations with a given detector
554 (use,) = (star_obs["detector"][obs_match] == detector_id).nonzero()
555 # Prior to running the calibration, we want to remove the effect
556 # of the jacobian of the WCS because that is a known quantity.
557 # Ideally, this would be done for each individual WCS, but this
558 # is extremely slow and makes small differences that are much
559 # smaller than the variation in the throughput due to other issues.
560 # Therefore, we use the approximate jacobian estimated from the
561 # camera model.
562 jac = approx_pixel_area_fields[detector_id].evaluate(
563 star_obs["x"][obs_match][use],
564 star_obs["y"][obs_match][use],
565 )
566 star_obs["jacobian"][obs_match[use]] = jac
567 scaled_inst_flux = (sources[self.config.instFluxField][obs_match[use]]
568 * visit_cat_table["scaling"][visit_match[use],
569 detector_mapping[detector_id]])
570 star_obs["inst_mag"][obs_match[use]] = (-2.5 * np.log10(scaled_inst_flux
571 / exp_time[use]))
572
573 # Compute instMagErr from inst_flux_err/inst_flux; scaling will cancel out.
574 star_obs["inst_mag_err"] = k*(sources[self.config.instFluxField + "Err"]
575 / sources[self.config.instFluxField])
576
577 # Apply the jacobian if configured to do so.
578 if self.config.doApplyWcsJacobian:
579 star_obs["inst_mag"] -= 2.5*np.log10(star_obs["jacobian"])
580
581 # We now reformat the stars and compute the ''objects'' that fgcm expects.
582 fgcm_obj = Table(data=np.zeros(len(stars), dtype=fgcm_obj_dtype))
583 fgcm_obj["isolated_star_id"] = stars["isolated_star_id"]
584 fgcm_obj["ra"] = stars["ra"]
585 fgcm_obj["dec"] = stars["decl"]
586 fgcm_obj["obs_arr_index"] = stars["source_cat_index"]
587 fgcm_obj["n_obs"] = stars["nsource"]
588
589 # Offset indexes to account for tract merging
590 fgcm_obj["obs_arr_index"] += merge_source_counter
591
592 fgcm_objs.append(fgcm_obj)
593 star_obs_cats.append(star_obs)
594
595 merge_source_counter += len(star_obs)
596
597 # Set the fgcm_id to a unique 64-bit integer for easier sorting.
598 fgcm_obj["fgcm_id"][:] = np.arange(len(fgcm_obj)) + 1
599
600 return vstack(fgcm_objs), vstack(star_obs_cats)
601
602 def _density_downsample(self, fgcm_obj, star_obs):
603 """Downsample stars according to density.
604
605 Catalogs are modified in-place.
606
607 Parameters
608 ----------
609 fgcm_obj : `astropy.table.Table`
610 Catalog of per-star ids and positions.
611 star_obs : `astropy.table.Table`
612 Catalog of individual star observations.
613 """
614 if self.config.randomSeed is not None:
615 np.random.seed(seed=self.config.randomSeed)
616
617 ipnest = hpg.angle_to_pixel(self.config.densityCutNside, fgcm_obj["ra"], fgcm_obj["dec"])
618 # Use the esutil.stat.histogram function to pull out the histogram of
619 # grouped pixels, and the rev_indices which describes which inputs
620 # are grouped together.
621 hist, rev_indices = esutil.stat.histogram(ipnest, rev=True)
622
623 obj_use = np.ones(len(fgcm_obj), dtype=bool)
624
625 (high,) = (hist > self.config.densityCutMaxPerPixel).nonzero()
626 (ok,) = (hist > 0).nonzero()
627 self.log.info("There are %d/%d pixels with high stellar density.", high.size, ok.size)
628 for i in range(high.size):
629 # The pix_indices are the indices of every star in the pixel.
630 pix_indices = rev_indices[rev_indices[high[i]]: rev_indices[high[i] + 1]]
631 # Cut down to the maximum number of stars in the pixel.
632 cut = np.random.choice(
633 pix_indices,
634 size=pix_indices.size - self.config.densityCutMaxPerPixel,
635 replace=False,
636 )
637 obj_use[cut] = False
638
639 fgcm_obj = fgcm_obj[obj_use]
640
641 obs_index = np.zeros(np.sum(fgcm_obj["n_obs"]), dtype=np.int32)
642 ctr = 0
643 for i in range(len(fgcm_obj)):
644 n_obs = fgcm_obj["n_obs"][i]
645 obs_index[ctr: ctr + n_obs] = (
646 np.arange(fgcm_obj["obs_arr_index"][i], fgcm_obj["obs_arr_index"][i] + n_obs)
647 )
648 fgcm_obj["obs_arr_index"][i] = ctr
649 ctr += n_obs
650
651 star_obs = star_obs[obs_index]
652
653 def _mark_reserve_stars(self, fgcm_obj):
654 """Run the star reservation task to flag reserved stars.
655
656 Parameters
657 ----------
658 fgcm_obj : `astropy.table.Table`
659 Catalog of per-star ids and positions.
660 """
661 reserved = self.reserve_selection.run(
662 len(fgcm_obj),
663 extra=','.join(self.config.requiredBands),
664 )
665 fgcm_obj["obj_flag"][reserved] |= objFlagDict["RESERVED"]
666
667 def _associate_reference_stars(self, lookup_table_handle, visit_cat, fgcm_obj):
668 """Associate fgcm isolated stars with reference stars.
669
670 Parameters
671 ----------
672 lookup_table_handle : `lsst.daf.butler.DeferredDatasetHandle`, optional
673 Data reference to fgcm look-up table (used if matching reference stars).
674 visit_cat : `lsst.afw.table.ExposureCatalog`
675 Catalog of per-visit data.
676 fgcm_obj : `astropy.table.Table`
677 Catalog of per-star ids and positions
678
679 Returns
680 -------
681 ref_cat : `astropy.table.Table`
682 Catalog of reference stars matched to fgcm stars.
683 """
684 # Figure out the correct bands/filters that we need based on the data.
685 lut_cat = lookup_table_handle.get()
686
687 std_filter_dict = {filter_name: std_filter for (filter_name, std_filter) in
688 zip(lut_cat[0]["physicalFilters"].split(","),
689 lut_cat[0]["stdPhysicalFilters"].split(","))}
690 std_lambda_dict = {std_filter: std_lambda for (std_filter, std_lambda) in
691 zip(lut_cat[0]["stdPhysicalFilters"].split(","),
692 lut_cat[0]["lambdaStdFilter"])}
693 del lut_cat
694
695 reference_filter_names = self._getReferenceFilterNames(
696 visit_cat,
697 std_filter_dict,
698 std_lambda_dict,
699 )
700 self.log.info("Using the following reference filters: %s", (", ".join(reference_filter_names)))
701
702 # Split into healpix pixels for more efficient matching.
703 ipnest = hpg.angle_to_pixel(self.config.coarseNside, fgcm_obj["ra"], fgcm_obj["dec"])
704 hpix, revpix = esutil.stat.histogram(ipnest, rev=True)
705
706 pixel_cats = []
707
708 # Compute the dtype from the filter names.
709 dtype = [("fgcm_id", "i4"),
710 ("refMag", "f4", (len(reference_filter_names), )),
711 ("refMagErr", "f4", (len(reference_filter_names), ))]
712
713 (gdpix,) = (hpix > 0).nonzero()
714 for ii, gpix in enumerate(gdpix):
715 p1a = revpix[revpix[gpix]: revpix[gpix + 1]]
716
717 # We do a simple wrapping of RA if we need to.
718 ra_wrap = fgcm_obj["ra"][p1a]
719 if (ra_wrap.min() < 10.0) and (ra_wrap.max() > 350.0):
720 ra_wrap[ra_wrap > 180.0] -= 360.0
721 mean_ra = np.mean(ra_wrap) % 360.0
722 else:
723 mean_ra = np.mean(ra_wrap)
724 mean_dec = np.mean(fgcm_obj["dec"][p1a])
725
726 dist = esutil.coords.sphdist(mean_ra, mean_dec, fgcm_obj["ra"][p1a], fgcm_obj["dec"][p1a])
727 rad = dist.max()
728
729 if rad < hpg.nside_to_resolution(self.config.coarseNside)/2.:
730 # Small radius, read just the circle.
731 ref_cat = self.fgcmLoadReferenceCatalog.getFgcmReferenceStarsSkyCircle(
732 mean_ra,
733 mean_dec,
734 rad,
735 reference_filter_names,
736 )
737 else:
738 # Otherwise the less efficient but full coverage.
739 ref_cat = self.fgcmLoadReferenceCatalog.getFgcmReferenceStarsHealpix(
740 self.config.coarseNside,
741 hpg.nest_to_ring(self.config.coarseNside, ipnest[p1a[0]]),
742 reference_filter_names,
743 )
744 if ref_cat.size == 0:
745 # No stars in this pixel; that's okay.
746 continue
747
748 with Matcher(fgcm_obj["ra"][p1a], fgcm_obj["dec"][p1a]) as matcher:
749 idx, i1, i2, d = matcher.query_radius(
750 ref_cat["ra"],
751 ref_cat["dec"],
752 self.config.matchRadius/3600.,
753 return_indices=True,
754 )
755
756 if len(i1) == 0:
757 # No matched stars in this pixel; that's okay.
758 continue
759
760 pixel_cat = Table(data=np.zeros(i1.size, dtype=dtype))
761 pixel_cat["fgcm_id"] = fgcm_obj["fgcm_id"][p1a[i1]]
762 pixel_cat["refMag"][:, :] = ref_cat["refMag"][i2, :]
763 pixel_cat["refMagErr"][:, :] = ref_cat["refMagErr"][i2, :]
764
765 pixel_cats.append(pixel_cat)
766
767 self.log.info(
768 "Found %d reference matches in pixel %d (%d of %d).",
769 len(pixel_cat),
770 ipnest[p1a[0]],
771 ii,
772 gdpix.size - 1,
773 )
774
775 ref_cat_full = vstack(pixel_cats)
776 ref_cat_full.meta.update({'FILTERNAMES': reference_filter_names})
777
778 return ref_cat_full
779
780 def _compute_delta_aper_summary_statistics(self, visit_cat, star_obs):
781 """Compute delta aperture summary statistics (per visit).
782
783 Parameters
784 ----------
785 visit_cat : `lsst.afw.table.ExposureCatalog`
786 Catalog of per-visit data.
787 star_obs : `astropy.table.Table`
788 Catalog of individual star observations.
789 """
790 (ok,) = ((star_obs["delta_mag_aper"] < 99.0)
791 & (star_obs["delta_mag_err_aper"] < 99.0)).nonzero()
792
793 visit_index = np.zeros(len(star_obs[ok]), dtype=np.int32)
794 a, b = esutil.numpy_util.match(visit_cat["visit"], star_obs["visit"][ok])
795 visit_index[b] = a
796
797 h, rev = esutil.stat.histogram(visit_index, rev=True)
798
799 visit_cat["deltaAper"] = -9999.0
800 h_use, = np.where(h >= 3)
801 for index in h_use:
802 i1a = rev[rev[index]: rev[index + 1]]
803 visit_cat["deltaAper"][index] = np.median(star_obs["delta_mag_aper"][ok[i1a]])
def _associate_reference_stars(self, lookup_table_handle, visit_cat, fgcm_obj)
def run(self, *camera, visit_summary_handle_dict, isolated_star_cat_handle_dict, isolated_star_source_handle_dict, lookup_table_handle=None)
def _make_all_star_obs_from_isolated_stars(self, isolated_star_cat_handle_dict, isolated_star_source_handle_dict, visit_cat, camera, calib_flux_aperture_radius=None)
def _getReferenceFilterNames(self, visitCat, stdFilterDict, stdLambdaDict)