lsst.pipe.tasks  21.0.0-142-gef555c1e+bcbe3dafff
drpAssociationPipe.py
Go to the documentation of this file.
1 # This file is part of pipe_tasks.
2 
3 # Developed for the LSST Data Management System.
4 # This product includes software developed by the LSST Project
5 # (https://www.lsst.org).
6 # See the COPYRIGHT file at the top-level directory of this distribution
7 # for details of code ownership.
8 
9 # This program is free software: you can redistribute it and/or modify
10 # it under the terms of the GNU General Public License as published by
11 # the Free Software Foundation, either version 3 of the License, or
12 # (at your option) any later version.
13 
14 # This program is distributed in the hope that it will be useful,
15 # but WITHOUT ANY WARRANTY; without even the implied warranty of
16 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 # GNU General Public License for more details.
18 
19 # You should have received a copy of the GNU General Public License
20 # along with this program. If not, see <https://www.gnu.org/licenses/>.
21 #
22 
23 """Pipeline for running DiaSource association in a DRP context.
24 """
25 
26 import numpy as np
27 import pandas as pd
28 
29 import lsst.geom as geom
30 import lsst.pex.config as pexConfig
31 import lsst.pipe.base as pipeBase
32 from lsst.skymap import BaseSkyMap
33 
34 from .coaddBase import makeSkyInfo
35 from .simpleAssociation import SimpleAssociationTask
36 
37 __all__ = ["DrpAssociationPipeTask",
38  "DrpAssociationPipeConfig",
39  "DrpAssociationPipeConnections"]
40 
41 
42 class DrpAssociationPipeConnections(pipeBase.PipelineTaskConnections,
43  dimensions=("tract", "patch", "skymap"),
44  defaultTemplates={"coaddName": "deep",
45  "warpTypeSuffix": "",
46  "fakesType": ""}):
47  diaSourceTables = pipeBase.connectionTypes.Input(
48  doc="Set of catalogs of calibrated DiaSources.",
49  name="{fakesType}{coaddName}Diff_diaSrcTable",
50  storageClass="DataFrame",
51  dimensions=("instrument", "visit", "detector"),
52  deferLoad=True,
53  multiple=True
54  )
55  skyMap = pipeBase.connectionTypes.Input(
56  doc="Input definition of geometry/bbox and projection/wcs for coadded "
57  "exposures",
58  name=BaseSkyMap.SKYMAP_DATASET_TYPE_NAME,
59  storageClass="SkyMap",
60  dimensions=("skymap", ),
61  )
62  assocDiaSourceTable = pipeBase.connectionTypes.Output(
63  doc="Catalog of DiaSources covering the patch and associated with a "
64  "DiaObject.",
65  name="{fakesType}{coaddName}Diff_assocDiaSrcTable",
66  storageClass="DataFrame",
67  dimensions=("tract", "patch"),
68  )
69  diaObjectTable = pipeBase.connectionTypes.Output(
70  doc="Catalog of DiaObjects created from spatially associating "
71  "DiaSources.",
72  name="{fakesType}{coaddName}Diff_diaObjTable",
73  storageClass="DataFrame",
74  dimensions=("tract", "patch"),
75  )
76 
77 
78 class DrpAssociationPipeConfig(
79  pipeBase.PipelineTaskConfig,
80  pipelineConnections=DrpAssociationPipeConnections):
81  associator = pexConfig.ConfigurableField(
82  target=SimpleAssociationTask,
83  doc="Task used to associate DiaSources with DiaObjects.",
84  )
85 
86 
87 class DrpAssociationPipeTask(pipeBase.PipelineTask):
88  """Driver pipeline for loading DiaSource catalogs in a patch/tract
89  region and associating them.
90  """
91  ConfigClass = DrpAssociationPipeConfig
92  _DefaultName = "drpAssociation"
93 
94  def __init__(self, **kwargs):
95  super().__init__(**kwargs)
96  self.makeSubtask('associator')
97 
98  def runQuantum(self, butlerQC, inputRefs, outputRefs):
99  inputs = butlerQC.get(inputRefs)
100 
101  inputs["tractId"] = butlerQC.quantum.dataId["tract"]
102  inputs["patchId"] = butlerQC.quantum.dataId["patch"]
103  tractPatchId, skymapBits = butlerQC.quantum.dataId.pack(
104  "tract_patch",
105  returnMaxBits=True)
106  inputs["tractPatchId"] = tractPatchId
107  inputs["skymapBits"] = skymapBits
108 
109  outputs = self.run(**inputs)
110  butlerQC.put(outputs, outputRefs)
111 
112  def run(self,
113  diaSourceTables,
114  skyMap,
115  tractId,
116  patchId,
117  tractPatchId,
118  skymapBits):
119  """Trim DiaSources to the current Patch and run association.
120 
121  Takes in the set of DiaSource catalogs that covers the current patch,
122  trims them to the dimensions of the patch, and [TODO: eventually]
123  runs association on the concatenated DiaSource Catalog.
124 
125  Parameters
126  ----------
127  diaSourceTables : `list` of `lst.daf.butler.DeferredDatasetHandle`
128  Set of DiaSource catalogs potentially covering this patch/tract.
129  skyMap : `lsst.skymap.BaseSkyMap`
130  SkyMap defining the patch/tract
131  tractId : `int`
132  Id of current tract being processed.
133  patchId : `int`
134  Id of current patch being processed
135 
136  Returns
137  -------
138  output : `lsst.pipe.base.Struct`
139  Results struct with attributes:
140 
141  ``assocDiaSourceTable``
142  Table of DiaSources with updated value for diaObjectId.
143  (`pandas.DataFrame`)
144  ``diaObjectTable``
145  Table of DiaObjects from matching DiaSources
146  (`pandas.DataFrame`).
147  """
148  self.log.info("Running DPR Association on patch %i, tract %i...",
149  patchId, tractId)
150 
151  skyInfo = makeSkyInfo(skyMap, tractId, patchId)
152 
153  # Get the patch bounding box.
154  innerPatchBox = geom.Box2D(skyInfo.patchInfo.getInnerBBox())
155 
156  diaSourceHistory = []
157  for catRef in diaSourceTables:
158  cat = catRef.get(
159  datasetType=self.config.connections.diaSourceTables,
160  immediate=True)
161 
162  isInTractPatch = self._trimToPatch(cat,
163  innerPatchBox,
164  skyInfo.wcs)
165 
166  nDiaSrc = isInTractPatch.sum()
167  self.log.info(
168  "Read DiaSource catalog of length %i from visit %i, "
169  "detector %i. Found %i sources within the patch/tract "
170  "footprint.",
171  len(cat), catRef.dataId["visit"],
172  catRef.dataId["detector"], nDiaSrc)
173 
174  if nDiaSrc <= 0:
175  diaSourceHistory.append(pd.DataFrame(columns=cat.columns))
176  continue
177 
178  cutCat = cat[isInTractPatch]
179  diaSourceHistory.append(cutCat)
180 
181  diaSourceHistoryCat = pd.concat(diaSourceHistory)
182  self.log.info("Found %i DiaSources overlapping patch %i, tract %i",
183  len(diaSourceHistoryCat), patchId, tractId)
184 
185  assocResult = self.associator.run(diaSourceHistoryCat,
186  tractPatchId,
187  skymapBits)
188 
189  self.log.info("Associated DiaSources into %i DiaObjects",
190  len(assocResult.diaObjects))
191 
192  return pipeBase.Struct(
193  diaObjectTable=assocResult.diaObjects,
194  assocDiaSourceTable=assocResult.assocDiaSources)
195 
196  def _trimToPatch(self, cat, innerPatchBox, wcs):
197  """Create generator testing if a set of DiaSources are in the
198  patch/tract.
199 
200  Parameters
201  ----------
202  cat : `pandas.DataFrame`
203  Catalog of DiaSources to test within patch/tract.
204  innerPatchBox : `lsst.geom.Box2D`
205  Bounding box of the patch.
206  wcs : `lsst.geom.SkyWcs`
207  Wcs of the tract.
208 
209  Returns
210  ------
211  isInPatch : `numpy.ndarray`, (N,)
212  Booleans representing if the DiaSources are contained within the
213  current patch and tract.
214  """
215  isInPatch = np.array([
216  innerPatchBox.contains(
217  wcs.skyToPixel(
218  geom.SpherePoint(row["ra"], row["decl"], geom.degrees)))
219  for idx, row in cat.iterrows()])
220  return isInPatch
def makeSkyInfo(skyMap, tractId, patchId)
Definition: coaddBase.py:289