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