CMS 3D CMS Logo

PFRecHitTopologyESProducer.cc
Go to the documentation of this file.
1 #include <algorithm>
2 #include <memory>
3 #include <type_traits>
4 #include <utility>
5 #include <variant>
6 
23 
24 #include "CalorimeterDefinitions.h"
25 
27  using namespace particleFlowRecHitProducer;
28 
29  template <typename CAL>
31  public:
33  : ESProducer(iConfig), cutsFromDB_(iConfig.getParameter<bool>("usePFThresholdsFromDB")) {
34  auto cc = setWhatProduced(this);
35  geomToken_ = cc.consumes();
36  if constexpr (std::is_same_v<CAL, HCAL>) {
37  hcalToken_ = cc.consumes();
38  hcalCutsToken_ = cc.consumes();
39  }
40  }
41 
42  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
44  if constexpr (std::is_same_v<CAL, HCAL>)
45  desc.add<bool>("usePFThresholdsFromDB", true);
46  else // only needs to be true for HBHE
47  desc.add<bool>("usePFThresholdsFromDB", false);
48  descriptions.addWithDefaultLabel(desc);
49  }
50 
51  std::unique_ptr<typename CAL::TopologyTypeHost> produce(const typename CAL::TopologyRecordType& iRecord) {
52  const auto& geom = iRecord.get(geomToken_);
53  auto product = std::make_unique<typename CAL::TopologyTypeHost>(CAL::kSize, cms::alpakatools::host());
54  auto view = product->view();
55  const int calEnums[2] = {CAL::kSubdetectorBarrelId, CAL::kSubdetectorEndcapId};
56  for (const auto subdet : calEnums) {
57  // Construct topology
58  // for HCAL: using dedicated record
59  // for ECAL: from CaloGeometry (separate for barrel and endcap)
60  const CaloSubdetectorGeometry* geo = geom.getSubdetectorGeometry(CAL::kDetectorId, subdet);
61  const CaloSubdetectorTopology* topo;
62  std::variant<EcalBarrelTopology, EcalEndcapTopology> topoVar; // need to store ECAL topology temporarily
63  if constexpr (std::is_same_v<CAL, HCAL>)
64  topo = &iRecord.get(hcalToken_);
65  else if (subdet == EcalSubdetector::EcalBarrel)
66  topo = &topoVar.emplace<EcalBarrelTopology>(geom);
67  else
68  topo = &topoVar.emplace<EcalEndcapTopology>(geom);
69 
70  // Fill product
71  for (auto const detId : geom.getValidDetIds(CAL::kDetectorId, subdet)) {
72  const uint32_t denseId = CAL::detId2denseId(detId);
73  assert(denseId < CAL::kSize);
74 
75  // Fill SoA members with HCAL PF Thresholds from GT
76  if constexpr (std::is_same_v<CAL, HCAL>) {
77  view.cutsFromDB() = false;
78  if (cutsFromDB_) {
79  view.cutsFromDB() = true;
80  const HcalPFCuts& pfCuts = iRecord.get(hcalCutsToken_);
81  const HcalTopology& htopo = iRecord.get(hcalToken_);
82  std::unique_ptr<HcalPFCuts> prod = std::make_unique<HcalPFCuts>(pfCuts);
83  prod->setTopo(&htopo);
84  view.noiseThreshold(denseId) = prod->getValues(detId.rawId())->noiseThreshold();
85  view.seedThreshold(denseId) = prod->getValues(detId.rawId())->seedThreshold();
86  }
87  }
88 
89  const GlobalPoint pos = geo->getGeometry(detId)->getPosition();
90  view.positionX(denseId) = pos.x();
91  view.positionY(denseId) = pos.y();
92  view.positionZ(denseId) = pos.z();
93 
94  for (uint32_t n = 0; n < 8; n++) {
95  uint32_t neighDetId = getNeighbourDetId(detId, n, *topo);
96  if (CAL::detIdInRange(neighDetId))
97  view.neighbours(denseId)(n) = CAL::detId2denseId(neighDetId);
98  else
99  view.neighbours(denseId)(n) = 0xffffffff;
100  }
101  }
102  }
103 
104  // Remove neighbours that are not backward compatible (only for HCAL)
105  if constexpr (std::is_same_v<CAL, HCAL>) {
106  for (const auto subdet : calEnums)
107  for (auto const detId : geom.getValidDetIds(CAL::kDetectorId, subdet)) {
108  const uint32_t denseId = CAL::detId2denseId(detId);
109  for (uint32_t n = 0; n < 8; n++) {
110  if (view.neighbours(denseId)[n] == 0xffffffff)
111  continue;
112  const ::reco::PFRecHitsTopologyNeighbours& neighboursOfNeighbour =
113  view.neighbours(view.neighbours(denseId)[n]);
114  if (std::find(neighboursOfNeighbour.begin(), neighboursOfNeighbour.end(), denseId) ==
115  neighboursOfNeighbour.end())
116  view.neighbours(denseId)[n] = 0xffffffff;
117  }
118  }
119  }
120 
121  // Print results (for debugging)
122  LogDebug("PFRecHitTopologyESProducer").log([&](auto& log) {
123  for (const auto subdet : calEnums)
124  for (const auto detId : geom.getValidDetIds(CAL::kDetectorId, subdet)) {
125  const uint32_t denseId = CAL::detId2denseId(detId);
126  log.format("detId:{} denseId:{} pos:{},{},{} neighbours:{},{},{},{};{},{},{},{}\n",
127  (uint32_t)detId,
128  denseId,
129  view[denseId].positionX(),
130  view[denseId].positionY(),
131  view[denseId].positionZ(),
132  view[denseId].neighbours()(0),
133  view[denseId].neighbours()(1),
134  view[denseId].neighbours()(2),
135  view[denseId].neighbours()(3),
136  view[denseId].neighbours()(4),
137  view[denseId].neighbours()(5),
138  view[denseId].neighbours()(6),
139  view[denseId].neighbours()(7));
140  }
141  });
142 
143  return product;
144  }
145 
146  private:
150  const bool cutsFromDB_;
151 
152  // specialised for HCAL/ECAL, because non-nearest neighbours are defined differently
153  uint32_t getNeighbourDetId(const uint32_t detId, const uint32_t direction, const CaloSubdetectorTopology& topo);
154  };
155 
156  template <>
158  const uint32_t direction,
159  const CaloSubdetectorTopology& topo) {
160  // desired order for PF: NORTH, SOUTH, EAST, WEST, NORTHEAST, SOUTHWEST, SOUTHEAST, NORTHWEST
161  if (detId == 0)
162  return 0;
163 
164  if (direction == 0) // NORTH
165  return topo.goNorth(detId); // larger iphi values (except phi boundary)
166  if (direction == 1) // SOUTH
167  return topo.goSouth(detId); // smaller iphi values (except phi boundary)
168  if (direction == 2) // EAST
169  return topo.goEast(detId); // smaller ieta values
170  if (direction == 3) // WEST
171  return topo.goWest(detId); // larger ieta values
172 
173  if (direction == 4) { // NORTHEAST
174  const uint32_t NE = getNeighbourDetId(getNeighbourDetId(detId, 0, topo), 2, topo);
175  if (NE)
176  return NE;
177  return getNeighbourDetId(getNeighbourDetId(detId, 2, topo), 0, topo);
178  }
179  if (direction == 5) { // SOUTHWEST
180  const uint32_t SW = getNeighbourDetId(getNeighbourDetId(detId, 1, topo), 3, topo);
181  if (SW)
182  return SW;
183  return getNeighbourDetId(getNeighbourDetId(detId, 3, topo), 1, topo);
184  }
185  if (direction == 6) { // SOUTHEAST
186  const uint32_t ES = getNeighbourDetId(getNeighbourDetId(detId, 2, topo), 1, topo);
187  if (ES)
188  return ES;
189  return getNeighbourDetId(getNeighbourDetId(detId, 1, topo), 2, topo);
190  }
191  if (direction == 7) { // NORTHWEST
192  const uint32_t WN = getNeighbourDetId(getNeighbourDetId(detId, 3, topo), 0, topo);
193  if (WN)
194  return WN;
195  return getNeighbourDetId(getNeighbourDetId(detId, 0, topo), 3, topo);
196  }
197  return 0;
198  }
199 
200  template <>
202  const uint32_t direction,
203  const CaloSubdetectorTopology& topo) {
204  // desired order for PF: NORTH, SOUTH, EAST, WEST, NORTHEAST, SOUTHWEST, SOUTHEAST, NORTHWEST
205  if (detId == 0)
206  return 0;
207 
208  if (direction == 0) // NORTH
209  return topo.goNorth(detId); // larger iphi values (except phi boundary)
210  if (direction == 1) // SOUTH
211  return topo.goSouth(detId); // smaller iphi values (except phi boundary)
212  if (direction == 2) // EAST
213  return topo.goEast(detId); // smaller ieta values
214  if (direction == 3) // WEST
215  return topo.goWest(detId); // larger ieta values
216 
217  std::pair<uint32_t, uint32_t> directions;
218  if (direction == 4) { // NORTHEAST
219  if (HCAL::getZside(detId) > 0)
220  directions = {2, 0}; // positive eta: east -> move to smaller |ieta| (finner phi granularity) first
221  else
222  directions = {0, 2}; // negative eta: move in phi first then move to east (coarser phi granularity)
223  } else if (direction == 5) { // SOUTHWEST
224  if (HCAL::getZside(detId) > 0)
225  directions = {1, 3}; // positive eta: move in phi first then move to west (coarser phi granularity)
226  else
227  directions = {3, 1}; // negative eta: west -> move to smaller |ieta| (finner phi granularity) first
228  } else if (direction == 6) { // SOUTHEAST
229  if (HCAL::getZside(detId) > 0)
230  directions = {2, 1}; // positive eta: east -> move to smaller |ieta| (finner phi granularity) first
231  else
232  directions = {1, 2}; // negative eta: move in phi first then move to east (coarser phi granularity)
233  } else if (direction == 7) { // NORTHWEST
234  if (HCAL::getZside(detId) > 0)
235  directions = {0, 3}; // positive eta: move in phi first then move to west (coarser phi granularity)
236  else
237  directions = {3, 0}; // negative eta: west -> move to smaller |ieta| (finner phi granularity) first
238  } else
239  return 0;
240  const uint32_t nn1 = getNeighbourDetId(detId, directions.first, topo); // nearest neighbour in direction 1
241  const uint32_t nn2 = getNeighbourDetId(detId, directions.second, topo); // nearest neighbour in direction 2
242  const uint32_t nnn = getNeighbourDetId(nn1, directions.second, topo); // next-to-nearest neighbour
243  if (nnn == nn1 || nnn == nn2) // avoid duplicates
244  return 0;
245  return nnn;
246  }
247 
250 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
251 
virtual DetId goNorth(const DetId &id) const
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
PFRecHitTopologyESProducer< ECAL > PFRecHitECALTopologyESProducer
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
assert(be >=bs)
virtual DetId goSouth(const DetId &id) const
std::unique_ptr< typename CAL::TopologyTypeHost > produce(const typename CAL::TopologyRecordType &iRecord)
Eigen::Matrix< uint32_t, 8, 1 > PFRecHitsTopologyNeighbours
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > hcalToken_
virtual DetId goEast(const DetId &id) const
PFRecHitTopologyESProducer< HCAL > PFRecHitHCALTopologyESProducer
virtual std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
alpaka::DevCpu const & host()
Definition: host.h:14
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
uint32_t getNeighbourDetId(const uint32_t detId, const uint32_t direction, const CaloSubdetectorTopology &topo)
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geomToken_
#define DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(type)
Definition: ModuleFactory.h:17
virtual DetId goWest(const DetId &id) const
edm::ESGetToken< HcalPFCuts, HcalPFCutsRcd > hcalCutsToken_
#define LogDebug(id)