CMS 3D CMS Logo

HGCGraph.cc
Go to the documentation of this file.
1 // Author: Felice Pantaleo - felice.pantaleo@cern.ch
2 // Date: 11/2018
6 #include "HGCDoublet.h"
7 #include "HGCGraph.h"
10 
11 template <typename TILES>
13  const std::vector<TICLSeedingRegion> &regions,
14  int nEtaBins,
15  int nPhiBins,
16  const std::vector<reco::CaloCluster> &layerClusters,
17  const std::vector<float> &mask,
18  const edm::ValueMap<std::pair<float, float>> &layerClustersTime,
19  int deltaIEta,
20  int deltaIPhi,
21  float minCosTheta,
22  float minCosPointing,
25  int skip_layers,
26  int maxNumberOfLayers,
27  float maxDeltaTime,
28  int lastLayerEE,
29  int lastLayerFH,
30  const std::vector<double> &siblings_maxRSquared) {
31  isOuterClusterOfDoublets_.clear();
32  isOuterClusterOfDoublets_.resize(layerClusters.size());
33  allDoublets_.clear();
34  theRootDoublets_.clear();
35  bool checkDistanceRootDoubletVsSeed = root_doublet_max_distance_from_seed_squared < 9999;
36  float origin_eta;
37  float origin_phi;
38  float maxRSquared;
39  for (const auto &r : regions) {
40  bool isGlobal = (r.index == -1);
41  auto zSide = r.zSide;
42  int startEtaBin, endEtaBin, startPhiBin, endPhiBin;
43 
44  if (isGlobal) {
45  startEtaBin = 0;
46  startPhiBin = 0;
47  endEtaBin = nEtaBins;
48  endPhiBin = nPhiBins;
49  origin_eta = 0;
50  origin_phi = 0;
51  } else {
52  auto firstLayerOnZSide = maxNumberOfLayers * zSide;
53  const auto &firstLayerHisto = histo[firstLayerOnZSide];
54  origin_eta = r.origin.eta();
55  origin_phi = r.origin.phi();
56  int entryEtaBin = firstLayerHisto.etaBin(origin_eta);
57  int entryPhiBin = firstLayerHisto.phiBin(origin_phi);
58  // For track-seeded iterations, if the impact point is below a certain
59  // eta-threshold, i.e., it has higher eta, make the initial search
60  // window bigger in both eta and phi by one bin, to contain better low
61  // energy showers.
62  auto etaWindow = deltaIEta;
63  auto phiWindow = deltaIPhi;
64  if (std::abs(origin_eta) > etaLimitIncreaseWindow) {
65  etaWindow++;
66  phiWindow++;
67  LogDebug("HGCGraph") << "Limit of Eta for increase: " << etaLimitIncreaseWindow
68  << " reached! Increasing inner search window" << std::endl;
69  }
70  startEtaBin = std::max(entryEtaBin - etaWindow, 0);
71  endEtaBin = std::min(entryEtaBin + etaWindow + 1, nEtaBins);
72  startPhiBin = entryPhiBin - phiWindow;
73  endPhiBin = entryPhiBin + phiWindow + 1;
74  if (verbosity_ > ticl::VerbosityLevel::Guru) {
75  LogDebug("HGCGraph") << " Entrance eta, phi: " << origin_eta << ", " << origin_phi
76  << " entryEtaBin: " << entryEtaBin << " entryPhiBin: " << entryPhiBin
77  << " globalBin: " << firstLayerHisto.globalBin(origin_eta, origin_phi)
78  << " on layer: " << firstLayerOnZSide << " startEtaBin: " << startEtaBin
79  << " endEtaBin: " << endEtaBin << " startPhiBin: " << startPhiBin
80  << " endPhiBin: " << endPhiBin << " phiBin(0): " << firstLayerHisto.phiBin(0.)
81  << " phiBin(" << M_PI / 2. << "): " << firstLayerHisto.phiBin(M_PI / 2.) << " phiBin("
82  << M_PI << "): " << firstLayerHisto.phiBin(M_PI) << " phiBin(" << -M_PI / 2.
83  << "): " << firstLayerHisto.phiBin(-M_PI / 2.) << " phiBin(" << -M_PI
84  << "): " << firstLayerHisto.phiBin(-M_PI) << " phiBin(" << 2. * M_PI
85  << "): " << firstLayerHisto.phiBin(2. * M_PI) << std::endl;
86  }
87  }
88 
89  for (int il = 0; il < maxNumberOfLayers - 1; ++il) {
90  for (int outer_layer = 0; outer_layer < std::min(1 + skip_layers, maxNumberOfLayers - 1 - il); ++outer_layer) {
91  int currentInnerLayerId = il + maxNumberOfLayers * zSide;
92  int currentOuterLayerId = currentInnerLayerId + 1 + outer_layer;
93  auto const &outerLayerHisto = histo[currentOuterLayerId];
94  auto const &innerLayerHisto = histo[currentInnerLayerId];
95  maxRSquared = (il <= lastLayerEE) ? siblings_maxRSquared[0]
96  : (il <= lastLayerFH) ? siblings_maxRSquared[1]
98  const int etaLimitIncreaseWindowBin = innerLayerHisto.etaBin(etaLimitIncreaseWindow);
99  if (verbosity_ > ticl::VerbosityLevel::Advanced) {
100  LogDebug("HGCGraph") << "Limit of Eta for increase: " << etaLimitIncreaseWindow
101  << " at etaBin: " << etaLimitIncreaseWindowBin << std::endl;
102  }
103 
104  for (int ieta = startEtaBin; ieta < endEtaBin; ++ieta) {
105  auto offset = ieta * nPhiBins;
106  for (int iphi_it = startPhiBin; iphi_it < endPhiBin; ++iphi_it) {
107  int iphi = ((iphi_it % nPhiBins + nPhiBins) % nPhiBins);
108  if (verbosity_ > ticl::VerbosityLevel::Advanced) {
109  LogDebug("HGCGraph") << "Inner Global Bin: " << (offset + iphi)
110  << " on layers I/O: " << currentInnerLayerId << "/" << currentOuterLayerId
111  << " with clusters: " << innerLayerHisto[offset + iphi].size() << std::endl;
112  }
113  for (auto innerClusterId : innerLayerHisto[offset + iphi]) {
114  // Skip masked clusters
115  if (mask[innerClusterId] == 0.) {
116  if (verbosity_ > ticl::VerbosityLevel::Advanced)
117  LogDebug("HGCGraph") << "Skipping inner masked cluster " << innerClusterId << std::endl;
118  continue;
119  }
120 
121  // For global-seeded iterations, if the inner cluster is above a certain
122  // eta-threshold, i.e., it has higher eta, make the outer search
123  // window bigger in both eta and phi by one bin, to contain better low
124  // energy showers. Track-Seeded iterations are excluded since, in
125  // that case, the inner search window has already been enlarged.
126  auto etaWindow = deltaIEta;
127  auto phiWindow = deltaIPhi;
128  if (isGlobal && ieta > etaLimitIncreaseWindowBin) {
129  etaWindow++;
130  phiWindow++;
131  if (verbosity_ > ticl::VerbosityLevel::Advanced) {
132  LogDebug("HGCGraph") << "Eta and Phi window increased by one" << std::endl;
133  }
134  }
135  const auto etaRangeMin = std::max(0, ieta - etaWindow);
136  const auto etaRangeMax = std::min(ieta + etaWindow + 1, nEtaBins);
137 
138  for (int oeta = etaRangeMin; oeta < etaRangeMax; ++oeta) {
139  // wrap phi bin
140  for (int phiRange = 0; phiRange < 2 * phiWindow + 1; ++phiRange) {
141  // The first wrapping is to take into account the
142  // cases in which we would have to seach in
143  // negative bins. The second wrap is mandatory to
144  // account for all other cases, since we add in
145  // between a full nPhiBins slot.
146  auto ophi = ((iphi + phiRange - phiWindow) % nPhiBins + nPhiBins) % nPhiBins;
147  if (verbosity_ > ticl::VerbosityLevel::Guru) {
148  LogDebug("HGCGraph") << "Outer Global Bin: " << (oeta * nPhiBins + ophi)
149  << " on layers I/O: " << currentInnerLayerId << "/" << currentOuterLayerId
150  << " with clusters: " << innerLayerHisto[oeta * nPhiBins + ophi].size()
151  << std::endl;
152  }
153  for (auto outerClusterId : outerLayerHisto[oeta * nPhiBins + ophi]) {
154  // Skip masked clusters
155  if (mask[outerClusterId] == 0.) {
156  if (verbosity_ > ticl::VerbosityLevel::Advanced)
157  LogDebug("HGCGraph") << "Skipping outer masked cluster " << outerClusterId << std::endl;
158  continue;
159  }
160  auto doubletId = allDoublets_.size();
161  if (maxDeltaTime != -1 &&
162  !areTimeCompatible(innerClusterId, outerClusterId, layerClustersTime, maxDeltaTime)) {
163  if (verbosity_ > ticl::VerbosityLevel::Advanced)
164  LogDebug("HGCGraph") << "Rejecting doublets due to timing!" << std::endl;
165  continue;
166  }
167  if (currentOuterLayerId - currentInnerLayerId == 1) {
168  if (areOverlappingOnSiblingLayers(innerClusterId, outerClusterId, layerClusters, maxRSquared)) {
169  allDoublets_.emplace_back(
170  innerClusterId, outerClusterId, doubletId, &layerClusters, r.index, true);
171  } else {
172  continue;
173  }
174  } else {
175  allDoublets_.emplace_back(
176  innerClusterId, outerClusterId, doubletId, &layerClusters, r.index, false);
177  }
178  if (verbosity_ > ticl::VerbosityLevel::Advanced) {
179  LogDebug("HGCGraph")
180  << "Creating doubletsId: " << doubletId << " layerLink in-out: [" << currentInnerLayerId
181  << ", " << currentOuterLayerId << "] clusterLink in-out: [" << innerClusterId << ", "
182  << outerClusterId << "]" << std::endl;
183  }
184  isOuterClusterOfDoublets_[outerClusterId].push_back(doubletId);
185  auto &neigDoublets = isOuterClusterOfDoublets_[innerClusterId];
186  auto &thisDoublet = allDoublets_[doubletId];
187  if (verbosity_ > ticl::VerbosityLevel::Expert) {
188  LogDebug("HGCGraph")
189  << "Checking compatibility of doubletId: " << doubletId
190  << " with all possible inners doublets link by the innerClusterId: " << innerClusterId
191  << std::endl;
192  }
193  bool isRootDoublet =
194  thisDoublet.checkCompatibilityAndTag(allDoublets_,
195  neigDoublets,
196  r.directionAtOrigin,
197  minCosTheta,
198  minCosPointing,
199  verbosity_ > ticl::VerbosityLevel::Advanced);
200  if (isRootDoublet and checkDistanceRootDoubletVsSeed) {
201  if (reco::deltaR2(layerClusters[innerClusterId].eta(),
202  layerClusters[innerClusterId].phi(),
203  origin_eta,
205  isRootDoublet = false;
206  }
207  }
208  if (isRootDoublet) {
209  theRootDoublets_.push_back(doubletId);
210  }
211  }
212  }
213  }
214  }
215  }
216  }
217  }
218  }
219  }
220  // #ifdef FP_DEBUG
221  if (verbosity_ > ticl::VerbosityLevel::None) {
222  LogDebug("HGCGraph") << "number of Root doublets " << theRootDoublets_.size() << " over a total number of doublets "
223  << allDoublets_.size() << std::endl;
224  }
225  // #endif
226 }
227 
228 template <typename TILES>
230  int outerIdx,
231  const edm::ValueMap<std::pair<float, float>> &layerClustersTime,
232  float maxDeltaTime) {
233  float timeIn = layerClustersTime.get(innerIdx).first;
234  float timeInE = layerClustersTime.get(innerIdx).second;
235  float timeOut = layerClustersTime.get(outerIdx).first;
236  float timeOutE = layerClustersTime.get(outerIdx).second;
237 
238  return (timeIn == -99. || timeOut == -99. ||
239  std::abs(timeIn - timeOut) < maxDeltaTime * sqrt(timeInE * timeInE + timeOutE * timeOutE));
240 }
241 
242 template <typename TILES>
244  int outerIdx,
245  const std::vector<reco::CaloCluster> &layerClusters,
246  float maxRSquared) {
247  return reco::deltaR2(layerClusters[outerIdx].eta(),
248  layerClusters[outerIdx].phi(),
249  layerClusters[innerIdx].eta(),
250  layerClusters[innerIdx].phi()) < maxRSquared;
251 }
252 
253 //also return a vector of seedIndex for the reconstructed tracksters
254 template <typename TILES>
255 void HGCGraphT<TILES>::findNtuplets(std::vector<HGCDoublet::HGCntuplet> &foundNtuplets,
256  std::vector<int> &seedIndices,
257  const unsigned int minClustersPerNtuplet,
258  const bool outInDFS,
259  unsigned int maxOutInHops) {
260  HGCDoublet::HGCntuplet tmpNtuplet;
261  tmpNtuplet.reserve(minClustersPerNtuplet);
262  std::vector<std::pair<unsigned int, unsigned int>> outInToVisit;
263  for (auto rootDoublet : theRootDoublets_) {
264  tmpNtuplet.clear();
265  outInToVisit.clear();
266  int seedIndex = allDoublets_[rootDoublet].seedIndex();
267  int outInHops = 0;
268  allDoublets_[rootDoublet].findNtuplets(
269  allDoublets_, tmpNtuplet, seedIndex, outInDFS, outInHops, maxOutInHops, outInToVisit);
270  while (!outInToVisit.empty()) {
271  allDoublets_[outInToVisit.back().first].findNtuplets(
272  allDoublets_, tmpNtuplet, seedIndex, outInDFS, outInToVisit.back().second, maxOutInHops, outInToVisit);
273  outInToVisit.pop_back();
274  }
275 
276  if (tmpNtuplet.size() > minClustersPerNtuplet) {
277  foundNtuplets.push_back(tmpNtuplet);
278  seedIndices.push_back(seedIndex);
279  }
280  }
281 }
282 
283 template class HGCGraphT<TICLLayerTiles>;
284 template class HGCGraphT<TICLLayerTilesHFNose>;
size
Write out results.
std::vector< unsigned int > HGCntuplet
Definition: HGCDoublet.h:16
void makeAndConnectDoublets(const TILES &h, const std::vector< TICLSeedingRegion > &regions, int nEtaBins, int nPhiBins, const std::vector< reco::CaloCluster > &layerClusters, const std::vector< float > &mask, const edm::ValueMap< std::pair< float, float >> &layerClustersTime, int deltaIEta, int deltaIPhi, float minCosThetai, float maxCosPointing, float root_doublet_max_distance_from_seed_squared, float etaLimitIncreaseWindow, int skip_layers, int maxNumberOfLayers, float maxDeltaTime, int lastLayerEE, int lastLayerFH, const std::vector< double > &siblings_maxRSquared)
Definition: HGCGraph.cc:12
T sqrt(T t)
Definition: SSEVec.h:19
bool areOverlappingOnSiblingLayers(int innerIdx, int outerIdx, const std::vector< reco::CaloCluster > &layerClusters, float maxRSquared)
Definition: HGCGraph.cc:243
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define M_PI
bool areTimeCompatible(int innerIdx, int outerIdx, const edm::ValueMap< std::pair< float, float >> &layerClustersTime, float maxDeltaTime)
Definition: HGCGraph.cc:229
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
void findNtuplets(std::vector< HGCDoublet::HGCntuplet > &foundNtuplets, std::vector< int > &seedIndices, const unsigned int minClustersPerNtuplet, const bool outInDFS, const unsigned int maxOutInHops)
Definition: HGCGraph.cc:255
#define LogDebug(id)