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