CMS 3D CMS Logo

PFHCALDenseIdNavigator.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_PFClusterProducer_PFHCALDenseIdNavigator_h
2 #define RecoParticleFlow_PFClusterProducer_PFHCALDenseIdNavigator_h
3 
5 
10 
16 
21 
24 
25 template <typename DET, typename TOPO, bool ownsTopo = true>
27 public:
29  if (!ownsTopo) {
30  topology_.release();
31  }
32  }
33 
35  : vhcalEnum_(iConfig.getParameter<std::vector<int>>("hcalEnums")),
36  hcalToken_(cc.esConsumes<edm::Transition::BeginRun>()),
37  geomToken_(cc.esConsumes<edm::Transition::BeginRun>()) {}
38 
39  // Initialization
40  void init(const edm::EventSetup& iSetup) override {
41  bool check = theRecNumberWatcher_.check(iSetup);
42  if (!check)
43  return;
44 
45  edm::ESHandle<HcalTopology> hcalTopology = iSetup.getHandle(hcalToken_);
46  topology_.release();
47  topology_.reset(hcalTopology.product());
48 
49  // Fill a vector of valid denseid's
51  const CaloGeometry& caloGeom = *hGeom;
52 
53  std::vector<DetId> vecHcal;
54  std::vector<unsigned int> vDenseIdHcal;
55  neighboursHcal_.clear();
56  for (auto hcalSubdet : vhcalEnum_) {
57  std::vector<DetId> vecDetIds(caloGeom.getValidDetIds(DetId::Hcal, hcalSubdet));
58  vecHcal.insert(vecHcal.end(), vecDetIds.begin(), vecDetIds.end());
59  }
60  vDenseIdHcal.reserve(vecHcal.size());
61  for (auto hDetId : vecHcal) {
62  vDenseIdHcal.push_back(topology_.get()->detId2denseId(hDetId));
63  }
64  std::sort(vDenseIdHcal.begin(), vDenseIdHcal.end());
65 
66  // Fill a vector of cell neighbours
67  denseIdHcalMax_ = *max_element(vDenseIdHcal.begin(), vDenseIdHcal.end());
68  denseIdHcalMin_ = *min_element(vDenseIdHcal.begin(), vDenseIdHcal.end());
70 
71  for (auto denseid : vDenseIdHcal) {
72  std::vector<DetId> neighbours(9, DetId(0));
73 
74  // the centre
75  unsigned denseid_c = denseid;
76  DetId detid_c = topology_.get()->denseId2detId(denseid_c);
77  CaloNavigator<DET> navigator(detid_c, topology_.get());
78 
79  // Using enum in Geometry/CaloTopology/interface/CaloDirection.h
80  // Order: CENTER(NONE),SOUTH,SOUTHEAST,SOUTHWEST,EAST,WEST,NORTHEAST,NORTHWEST,NORTH
81  neighbours.at(NONE) = detid_c;
82 
83  neighbours.at(NORTH) = getNeighbourDetId(detid_c, 0);
84  neighbours.at(SOUTH) = getNeighbourDetId(detid_c, 1);
85  neighbours.at(EAST) = getNeighbourDetId(detid_c, 2);
86  neighbours.at(WEST) = getNeighbourDetId(detid_c, 3);
87 
88  neighbours.at(NORTHEAST) = getNeighbourDetId(detid_c, 4);
89  neighbours.at(SOUTHWEST) = getNeighbourDetId(detid_c, 5);
90  neighbours.at(SOUTHEAST) = getNeighbourDetId(detid_c, 6);
91  neighbours.at(NORTHWEST) = getNeighbourDetId(detid_c, 7);
92 
93  unsigned index = getIdx(denseid_c);
94  neighboursHcal_[index] = neighbours;
95 
96  } // loop over vDenseIdHcal
97 
98  if (debug) {
99  backwardCompatibilityCheck(vDenseIdHcal);
100  printNeighbourInfo(vDenseIdHcal);
101  }
102  }
103 
104  // Check neighbour
105  void backwardCompatibilityCheck(const std::vector<unsigned int> vDenseIdHcal) {
106  //
107  // Check backward compatibility (does a neighbour of a channel have the channel as a neighbour?)
108  //
109  for (auto denseid : vDenseIdHcal) {
110  DetId detid = topology_.get()->denseId2detId(denseid);
111  HcalDetId hid = HcalDetId(detid);
112  if (detid == DetId(0)) {
113  edm::LogWarning("PFHCALDenseIdNavigator") << "Found an invalid DetId";
114  continue;
115  }
116  if (!validNeighbours(denseid)) {
117  edm::LogWarning("PFHCALDenseIdNavigator") << "The vector for neighbour information has an invalid length";
118  continue;
119  }
120  std::vector<DetId> neighbours(9, DetId(0));
121  unsigned index = getIdx(denseid);
122  if (index >= neighboursHcal_.size()) {
123  edm::LogWarning("PFHCALDenseIdNavigator") << "The vector for neighbour information is not found";
124  continue; // Skip if not found
125  }
126  neighbours = neighboursHcal_.at(index);
127 
128  //
129  // Loop over neighbours
130  int ineighbour = -1;
131  for (auto neighbour : neighbours) {
132  ineighbour++;
133  if (neighbour == DetId(0))
134  continue;
135  HcalDetId hidn = HcalDetId(neighbour);
136  std::vector<DetId> neighboursOfNeighbour(9, DetId(0));
137  std::unordered_set<unsigned int> listOfNeighboursOfNeighbour; // list of neighbours of neighbour
138  unsigned denseidNeighbour = topology_.get()->detId2denseId(neighbour);
139 
140  if (!validNeighbours(denseidNeighbour)) {
141  edm::LogWarning("PFHCALDenseIdNavigator")
142  << "This neighbour does not have a valid neighbour information (another subdetector?). Ignore: "
143  << detid.det() << " " << hid.ieta() << " " << hid.iphi() << " " << hid.depth() << " " << neighbour.det()
144  << " " << hidn.ieta() << " " << hidn.iphi() << " " << hidn.depth();
145  neighboursHcal_[index][ineighbour] = DetId(0); // Not a valid neighbour. Set to DetId(0)
146  continue;
147  }
148  if (getIdx(denseidNeighbour) >= neighboursHcal_.size())
149  continue;
150  neighboursOfNeighbour = neighboursHcal_.at(getIdx(denseidNeighbour));
151 
152  //
153  // Loop over neighbours of neighbours
154  for (auto neighbourOfNeighbour : neighboursOfNeighbour) {
155  if (neighbourOfNeighbour == DetId(0))
156  continue;
157  unsigned denseidNeighbourOfNeighbour = topology_.get()->detId2denseId(neighbourOfNeighbour);
158  if (!validNeighbours(denseidNeighbourOfNeighbour))
159  continue;
160  listOfNeighboursOfNeighbour.insert(denseidNeighbourOfNeighbour);
161  }
162 
163  //
164  if (listOfNeighboursOfNeighbour.find(denseid) == listOfNeighboursOfNeighbour.end()) {
165  // This neighbour is not backward compatible.
166  edm::LogWarning("PFHCALDenseIdNavigator")
167  << "This neighbour does not have the original channel as its neighbour. Ignore: " << detid.det() << " "
168  << hid.ieta() << " " << hid.iphi() << " " << hid.depth() << " " << neighbour.det() << " " << hidn.ieta()
169  << " " << hidn.iphi() << " " << hidn.depth();
170  neighboursHcal_[index][ineighbour] = DetId(0);
171  }
172 
173  } // loop over neighbours
174  } // loop over denseId_
175  }
176 
177  // Print out neighbour DetId's
178  void printNeighbourInfo(const std::vector<unsigned int> vDenseIdHcal) {
179  // Print neighbour definitions
180  for (auto denseid : vDenseIdHcal) {
181  std::vector<DetId> neighbours(9, DetId(0));
182 
183  // the centre
184  unsigned denseid_c = denseid;
185  DetId detid_c = topology_.get()->denseId2detId(denseid_c);
186  CaloNavigator<DET> navigator(detid_c, topology_.get());
187 
188  unsigned index = getIdx(denseid_c);
189 
190  neighbours = neighboursHcal_[index];
191 
192  const DetId N = neighbours.at(NORTH);
193  const DetId S = neighbours.at(SOUTH);
194  const DetId E = neighbours.at(EAST);
195  const DetId W = neighbours.at(WEST);
196 
197  const DetId NE = neighbours.at(NORTHEAST);
198  const DetId SW = neighbours.at(SOUTHWEST);
199  const DetId SE = neighbours.at(SOUTHEAST);
200  const DetId NW = neighbours.at(NORTHWEST);
201 
202  edm::LogPrint("PFHCALDenseIdNavigator")
203  << "PFHCALDenseIdNavigator: " << HcalDetId(detid_c).ieta() << " " << HcalDetId(detid_c).iphi() << " "
204  << HcalDetId(detid_c).depth() << " " << HcalDetId(detid_c).subdetId() << ": " << HcalDetId(N).ieta() << " "
205  << HcalDetId(N).iphi() << " " << HcalDetId(N).depth() << " " << HcalDetId(N).subdetId() << ", "
206  << HcalDetId(E).ieta() << " " << HcalDetId(E).iphi() << " " << HcalDetId(E).depth() << " "
207  << HcalDetId(E).subdetId() << ", " << HcalDetId(S).ieta() << " " << HcalDetId(S).iphi() << " "
208  << HcalDetId(S).depth() << " " << HcalDetId(S).subdetId() << ", " << HcalDetId(W).ieta() << " "
209  << HcalDetId(W).iphi() << " " << HcalDetId(W).depth() << " " << HcalDetId(W).subdetId() << ", "
210  << HcalDetId(NE).ieta() << " " << HcalDetId(NE).iphi() << " " << HcalDetId(NE).depth() << " "
211  << HcalDetId(NE).subdetId() << ", " << HcalDetId(SW).ieta() << " " << HcalDetId(SW).iphi() << " "
212  << HcalDetId(SW).depth() << " " << HcalDetId(SW).subdetId() << ", " << HcalDetId(SE).ieta() << " "
213  << HcalDetId(SE).iphi() << " " << HcalDetId(SE).depth() << " " << HcalDetId(SE).subdetId() << ", "
214  << HcalDetId(NW).ieta() << " " << HcalDetId(NW).iphi() << " " << HcalDetId(NW).depth() << " "
215  << HcalDetId(NW).subdetId();
216 
217  } // print ends
218  }
219 
220  // Check if two DetId's have the same subdetId
221  const bool checkSameSubDet(const DetId detId, const DetId detId2) {
222  HcalDetId hid = HcalDetId(detId);
223  HcalDetId hid2 = HcalDetId(detId2);
224  return hid.subdetId() == hid2.subdetId();
225  }
226 
227  //https://cmssdt.cern.ch/lxr/source/DataFormats/HcalDetId/interface/HcalDetId.h#0141
228  static constexpr int getZside(const DetId detId) { return ((detId & HcalDetId::kHcalZsideMask2) ? (1) : (-1)); }
229 
230  // Obtain the neighbour's DetId
231  const uint32_t getNeighbourDetId(const DetId detId, const uint32_t direction) {
232  // desired order for PF: NORTH, SOUTH, EAST, WEST, NORTHEAST, SOUTHWEST, SOUTHEAST, NORTHWEST
233  if (detId == DetId(0))
234  return 0;
235 
236  if (direction == 0) // NORTH
237  return topology_.get()->goNorth(detId); // larger iphi values (except phi boundary)
238 
239  if (direction == 1) // SOUTH
240  return topology_.get()->goSouth(detId); // smaller iphi values (except phi boundary)
241 
242  // In the case of East/West, make sure we are not moving to another subdetector
243  if (direction == 2 && checkSameSubDet(detId, topology_.get()->goEast(detId))) // EAST
244  return topology_.get()->goEast(detId); // smaller ieta values
245 
246  if (direction == 3 && checkSameSubDet(detId, topology_.get()->goWest(detId))) // WEST
247  return topology_.get()->goWest(detId); // larger ieta values
248 
249  // In the case of cornor cells, ensure backward compatibility.
250  // Also make sure it's not already defined as E(ast) or W(est)
251  if (direction == 4) { // NORTHEAST
252  if (getZside(detId) > 0) { // positive eta: east -> move to smaller |ieta| (finner phi granularity) first
253  const DetId E = getNeighbourDetId(detId, 2);
254  const DetId NE = getNeighbourDetId(E, 0);
255  if (getNeighbourDetId(NE, 1) == E && NE != E)
256  return NE; //
257  } else { // negative eta: move in phi first then move to east (coarser phi granularity)
258  const DetId N = getNeighbourDetId(detId, 0);
259  const DetId NE = getNeighbourDetId(N, 2);
260  const DetId E = getNeighbourDetId(detId, 2);
261  if (getNeighbourDetId(NE, 3) == N && NE != E)
262  return NE;
263  }
264  }
265  if (direction == 5) { // SOUTHWEST
266  if (getZside(detId) > 0) { // positive eta: move in phi first then move to west (coarser phi granularity)
267  const DetId S = getNeighbourDetId(detId, 1);
268  const DetId SW = getNeighbourDetId(S, 3);
269  const DetId W = getNeighbourDetId(detId, 3);
270  if (getNeighbourDetId(SW, 2) == S && SW != W)
271  return SW;
272  } else { // negative eta: west -> move to smaller |ieta| (finner phi granularity) first
273  const DetId W = getNeighbourDetId(detId, 3);
274  const DetId SW = getNeighbourDetId(W, 1);
275  if (getNeighbourDetId(SW, 0) == W && SW != W)
276  return SW;
277  }
278  }
279  if (direction == 6) { // SOUTHEAST
280  if (getZside(detId) > 0) { // positive eta: east -> move to smaller |ieta| (finner phi granularity) first
281  const DetId E = getNeighbourDetId(detId, 2);
282  const DetId SE = getNeighbourDetId(E, 1);
283  if (getNeighbourDetId(SE, 0) == E && SE != E)
284  return SE;
285  } else { // negative eta: move in phi first then move to east (coarser phi granularity)
286  const DetId S = getNeighbourDetId(detId, 1);
287  const DetId SE = getNeighbourDetId(S, 2);
288  const DetId E = getNeighbourDetId(detId, 2);
289  if (getNeighbourDetId(SE, 3) == S && SE != E)
290  return SE;
291  }
292  }
293  if (direction == 7) { // NORTHWEST
294  if (getZside(detId) > 0) { // positive eta: move in phi first then move to west (coarser phi granularity)
295  const DetId N = getNeighbourDetId(detId, 0);
296  const DetId NW = getNeighbourDetId(N, 3);
297  const DetId W = getNeighbourDetId(detId, 3);
298  if (getNeighbourDetId(NW, 2) == N && NW != W)
299  return NW;
300  } else { // negative eta: west -> move to smaller |ieta| (finner phi granularity) first
301  const DetId W = getNeighbourDetId(detId, 3);
302  const DetId NW = getNeighbourDetId(W, 0);
303  if (getNeighbourDetId(NW, 1) == W && NW != W)
304  return NW;
305  }
306  }
307  return 0;
308  }
309 
310  // Associate neighbour PFRecHits
312  std::unique_ptr<reco::PFRecHitCollection>& hits,
313  edm::RefProd<reco::PFRecHitCollection>& refProd) override {
314  DetId detid(hit.detId());
315  HcalDetId hid(detid);
316  unsigned denseid = topology_.get()->detId2denseId(detid);
317 
318  std::vector<DetId> neighbours(9, DetId(0));
319 
320  if (denseid < denseIdHcalMin_ || denseid > denseIdHcalMax_) {
321  edm::LogWarning("PFRecHitHCALCachedNavigator")
322  << " DenseId for this cell is out of the expected range." << std::endl;
323  } else if (!validNeighbours(denseid)) {
324  edm::LogWarning("PFRecHitHCALCachedNavigator")
325  << " DenseId for this cell does not have the neighbour information. " << hid.ieta() << " " << hid.iphi()
326  << " " << hid.depth() << " " << hid.subdetId();
327  } else {
328  unsigned index = getIdx(denseid);
329  neighbours = neighboursHcal_.at(index);
330  }
331 
332  associateNeighbour(neighbours.at(NORTH), hit, hits, refProd, 0, 1, 0); // N
333  associateNeighbour(neighbours.at(NORTHEAST), hit, hits, refProd, 1, 1, 0); // NE
334  associateNeighbour(neighbours.at(SOUTH), hit, hits, refProd, 0, -1, 0); // S
335  associateNeighbour(neighbours.at(SOUTHWEST), hit, hits, refProd, -1, -1, 0); // SW
336  associateNeighbour(neighbours.at(EAST), hit, hits, refProd, 1, 0, 0); // E
337  associateNeighbour(neighbours.at(SOUTHEAST), hit, hits, refProd, 1, -1, 0); // SE
338  associateNeighbour(neighbours.at(WEST), hit, hits, refProd, -1, 0, 0); // W
339  associateNeighbour(neighbours.at(NORTHWEST), hit, hits, refProd, -1, 1, 0); // NW
340  }
341 
342  // Check if we can get the valid neighbour vector for a given denseid
343  const bool validNeighbours(const unsigned int denseid) const {
344  if (denseid < denseIdHcalMin_ || denseid > denseIdHcalMax_)
345  return false; // neighbour denseid is out of the expected range
346  unsigned index = getIdx(denseid);
347  if (neighboursHcal_.at(index).size() != 9)
348  return false; // the neighbour vector size should be 3x3
349  return true;
350  }
351 
352  // Get the index of neighboursHcal_ for a given denseid
353  const unsigned int getIdx(const unsigned int denseid) const {
354  if (denseid < denseIdHcalMin_ || denseid > denseIdHcalMax_)
355  return unsigned(denseIdHcalMax_ - denseIdHcalMin_ +
356  1); // out-of-bounce (different subdetector groups. give a dummy number.)
357  unsigned index = denseid - denseIdHcalMin_;
358  return index;
359  }
360 
361 protected:
363  std::unique_ptr<const TOPO> topology_;
364  std::vector<int> vhcalEnum_;
365  std::vector<std::vector<DetId>> neighboursHcal_;
366  unsigned int denseIdHcalMax_;
367  unsigned int denseIdHcalMin_;
368 
369 private:
372  const bool debug = false;
373 };
374 
375 #endif
void associateNeighbour(const DetId &id, reco::PFRecHit &hit, std::unique_ptr< reco::PFRecHitCollection > &hits, edm::RefProd< reco::PFRecHitCollection > &refProd, short eta, short phi, short depth)
const bool checkSameSubDet(const DetId detId, const DetId detId2)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
void associateNeighbours(reco::PFRecHit &hit, std::unique_ptr< reco::PFRecHitCollection > &hits, edm::RefProd< reco::PFRecHitCollection > &refProd) override
void backwardCompatibilityCheck(const std::vector< unsigned int > vDenseIdHcal)
static constexpr uint32_t kHcalZsideMask2
Definition: HcalDetId.h:21
std::unique_ptr< const TOPO > topology_
edm::ESWatcher< HcalRecNumberingRecord > theRecNumberWatcher_
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geomToken_
T const * product() const
Definition: ESHandle.h:86
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
constexpr int ieta() const
get the cell ieta
Definition: HcalDetId.h:155
static constexpr int getZside(const DetId detId)
Transition
Definition: Transition.h:12
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
void printNeighbourInfo(const std::vector< unsigned int > vDenseIdHcal)
std::vector< std::vector< DetId > > neighboursHcal_
Log< level::Warning, true > LogPrint
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > hcalToken_
Definition: DetId.h:17
#define N
Definition: blowfish.cc:9
PFHCALDenseIdNavigator(const edm::ParameterSet &iConfig, edm::ConsumesCollector &cc)
std::vector< DetId > getValidDetIds() const
Get the list of all valid detector ids.
Definition: CaloGeometry.cc:75
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
const uint32_t getNeighbourDetId(const DetId detId, const uint32_t direction)
HLT enums.
void init(const edm::EventSetup &iSetup) override
const bool validNeighbours(const unsigned int denseid) const
Log< level::Warning, false > LogWarning
Definition: TkAlStyle.h:43
const unsigned int getIdx(const unsigned int denseid) const
constexpr int iphi() const
get the cell iphi
Definition: HcalDetId.h:157
constexpr int depth() const
get the tower depth
Definition: HcalDetId.h:164