CMS 3D CMS Logo

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