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