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 
19 #include "CalorimeterDefinitions.h"
20 
22  using namespace particleFlowRecHitProducer;
23 
24  template <typename CAL>
26  public:
28  auto cc = setWhatProduced(this);
29  geomToken_ = cc.consumes();
30  if constexpr (std::is_same_v<CAL, HCAL>)
31  hcalToken_ = cc.consumes();
32  }
33 
34  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
36  descriptions.addWithDefaultLabel(desc);
37  }
38 
39  std::unique_ptr<typename CAL::TopologyTypeHost> produce(const typename CAL::TopologyRecordType& iRecord) {
40  const auto& geom = iRecord.get(geomToken_);
41  auto product = std::make_unique<typename CAL::TopologyTypeHost>(CAL::kSize, cms::alpakatools::host());
42  auto view = product->view();
43 
44  const int calEnums[2] = {CAL::kSubdetectorBarrelId, CAL::kSubdetectorEndcapId};
45  for (const auto subdet : calEnums) {
46  // Construct topology
47  // for HCAL: using dedicated record
48  // for ECAL: from CaloGeometry (separate for barrel and endcap)
49  const CaloSubdetectorGeometry* geo = geom.getSubdetectorGeometry(CAL::kDetectorId, subdet);
50  const CaloSubdetectorTopology* topo;
51  std::variant<EcalBarrelTopology, EcalEndcapTopology> topoVar; // need to store ECAL topology temporarily
52  if constexpr (std::is_same_v<CAL, HCAL>)
53  topo = &iRecord.get(hcalToken_);
54  else if (subdet == EcalSubdetector::EcalBarrel)
55  topo = &topoVar.emplace<EcalBarrelTopology>(geom);
56  else
57  topo = &topoVar.emplace<EcalEndcapTopology>(geom);
58 
59  // Fill product
60  for (auto const detId : geom.getValidDetIds(CAL::kDetectorId, subdet)) {
61  const uint32_t denseId = CAL::detId2denseId(detId);
62  assert(denseId < CAL::kSize);
63 
64  const GlobalPoint pos = geo->getGeometry(detId)->getPosition();
65  view.positionX(denseId) = pos.x();
66  view.positionY(denseId) = pos.y();
67  view.positionZ(denseId) = pos.z();
68 
69  for (uint32_t n = 0; n < 8; n++) {
70  uint32_t neighDetId = getNeighbourDetId(detId, n, *topo);
71  if (CAL::detIdInRange(neighDetId))
72  view.neighbours(denseId)(n) = CAL::detId2denseId(neighDetId);
73  else
74  view.neighbours(denseId)(n) = 0xffffffff;
75  }
76  }
77  }
78 
79  // Remove neighbours that are not backward compatible (only for HCAL)
80  if constexpr (std::is_same_v<CAL, HCAL>) {
81  for (const auto subdet : calEnums)
82  for (auto const detId : geom.getValidDetIds(CAL::kDetectorId, subdet)) {
83  const uint32_t denseId = CAL::detId2denseId(detId);
84  for (uint32_t n = 0; n < 8; n++) {
85  const ::reco::PFRecHitsTopologyNeighbours& neighboursOfNeighbour =
86  view.neighbours(view.neighbours(denseId)[n]);
87  if (std::find(neighboursOfNeighbour.begin(), neighboursOfNeighbour.end(), denseId) ==
88  neighboursOfNeighbour.end())
89  view.neighbours(denseId)[n] = 0xffffffff;
90  }
91  }
92  }
93 
94  // Print results (for debugging)
95  LogDebug("PFRecHitTopologyESProducer").log([&](auto& log) {
96  for (const auto subdet : calEnums)
97  for (const auto detId : geom.getValidDetIds(CAL::kDetectorId, subdet)) {
98  const uint32_t denseId = CAL::detId2denseId(detId);
99  log.format("detId:{} denseId:{} pos:{},{},{} neighbours:{},{},{},{};{},{},{},{}\n",
100  (uint32_t)detId,
101  denseId,
102  view[denseId].positionX(),
103  view[denseId].positionY(),
104  view[denseId].positionZ(),
105  view[denseId].neighbours()(0),
106  view[denseId].neighbours()(1),
107  view[denseId].neighbours()(2),
108  view[denseId].neighbours()(3),
109  view[denseId].neighbours()(4),
110  view[denseId].neighbours()(5),
111  view[denseId].neighbours()(6),
112  view[denseId].neighbours()(7));
113  }
114  });
115 
116  return product;
117  }
118 
119  private:
122 
123  // specialised for HCAL/ECAL, because non-nearest neighbours are defined differently
124  uint32_t getNeighbourDetId(const uint32_t detId, const uint32_t direction, const CaloSubdetectorTopology& topo);
125  };
126 
127  template <>
129  const uint32_t direction,
130  const CaloSubdetectorTopology& topo) {
131  // desired order for PF: NORTH, SOUTH, EAST, WEST, NORTHEAST, SOUTHWEST, SOUTHEAST, NORTHWEST
132  if (detId == 0)
133  return 0;
134 
135  if (direction == 0) // NORTH
136  return topo.goNorth(detId); // larger iphi values (except phi boundary)
137  if (direction == 1) // SOUTH
138  return topo.goSouth(detId); // smaller iphi values (except phi boundary)
139  if (direction == 2) // EAST
140  return topo.goEast(detId); // smaller ieta values
141  if (direction == 3) // WEST
142  return topo.goWest(detId); // larger ieta values
143 
144  if (direction == 4) { // NORTHEAST
145  const uint32_t NE = getNeighbourDetId(getNeighbourDetId(detId, 0, topo), 2, topo);
146  if (NE)
147  return NE;
148  return getNeighbourDetId(getNeighbourDetId(detId, 2, topo), 0, topo);
149  }
150  if (direction == 5) { // SOUTHWEST
151  const uint32_t SW = getNeighbourDetId(getNeighbourDetId(detId, 1, topo), 3, topo);
152  if (SW)
153  return SW;
154  return getNeighbourDetId(getNeighbourDetId(detId, 3, topo), 1, topo);
155  }
156  if (direction == 6) { // SOUTHEAST
157  const uint32_t ES = getNeighbourDetId(getNeighbourDetId(detId, 2, topo), 1, topo);
158  if (ES)
159  return ES;
160  return getNeighbourDetId(getNeighbourDetId(detId, 1, topo), 2, topo);
161  }
162  if (direction == 7) { // NORTHWEST
163  const uint32_t WN = getNeighbourDetId(getNeighbourDetId(detId, 3, topo), 0, topo);
164  if (WN)
165  return WN;
166  return getNeighbourDetId(getNeighbourDetId(detId, 0, topo), 3, topo);
167  }
168  return 0;
169  }
170 
171  template <>
173  const uint32_t direction,
174  const CaloSubdetectorTopology& topo) {
175  // desired order for PF: NORTH, SOUTH, EAST, WEST, NORTHEAST, SOUTHWEST, SOUTHEAST, NORTHWEST
176  if (detId == 0)
177  return 0;
178 
179  if (direction == 0) // NORTH
180  return topo.goNorth(detId); // larger iphi values (except phi boundary)
181  if (direction == 1) // SOUTH
182  return topo.goSouth(detId); // smaller iphi values (except phi boundary)
183  if (direction == 2) // EAST
184  return topo.goEast(detId); // smaller ieta values
185  if (direction == 3) // WEST
186  return topo.goWest(detId); // larger ieta values
187 
188  std::pair<uint32_t, uint32_t> directions;
189  if (direction == 4) { // NORTHEAST
190  if (HCAL::getZside(detId) > 0)
191  directions = {2, 0}; // positive eta: east -> move to smaller |ieta| (finner phi granularity) first
192  else
193  directions = {0, 2}; // negative eta: move in phi first then move to east (coarser phi granularity)
194  } else if (direction == 5) { // SOUTHWEST
195  if (HCAL::getZside(detId) > 0)
196  directions = {1, 3}; // positive eta: move in phi first then move to west (coarser phi granularity)
197  else
198  directions = {3, 1}; // negative eta: west -> move to smaller |ieta| (finner phi granularity) first
199  } else if (direction == 6) { // SOUTHEAST
200  if (HCAL::getZside(detId) > 0)
201  directions = {2, 1}; // positive eta: east -> move to smaller |ieta| (finner phi granularity) first
202  else
203  directions = {1, 2}; // negative eta: move in phi first then move to east (coarser phi granularity)
204  } else if (direction == 7) { // NORTHWEST
205  if (HCAL::getZside(detId) > 0)
206  directions = {0, 3}; // positive eta: move in phi first then move to west (coarser phi granularity)
207  else
208  directions = {3, 0}; // negative eta: west -> move to smaller |ieta| (finner phi granularity) first
209  } else
210  return 0;
211  const uint32_t nn1 = getNeighbourDetId(detId, directions.first, topo); // nearest neighbour in direction 1
212  const uint32_t nn2 = getNeighbourDetId(detId, directions.second, topo); // nearest neighbour in direction 2
213  const uint32_t nnn = getNeighbourDetId(nn1, directions.second, topo); // next-to-nearest neighbour
214  if (nnn == nn1 || nnn == nn2) // avoid duplicates
215  return 0;
216  return nnn;
217  }
218 
221 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
222 
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:26
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
#define LogDebug(id)