CMS 3D CMS Logo

HGCalShowerShape.cc
Go to the documentation of this file.
5 
6 #include <unordered_map>
7 #include <numeric>
8 
10  : threshold_(conf.getParameter<double>("shape_threshold")),
11  distance_(conf.getParameter<double>("shape_distance")) {}
12 
13 //Compute energy-weighted mean of any variable X in the cluster
14 
15 float HGCalShowerShape::meanX(const std::vector<pair<float, float>>& energy_X_tc) const {
16  float Etot = 0;
17  float X_sum = 0;
18 
19  for (const auto& energy_X : energy_X_tc) {
20  X_sum += energy_X.first * energy_X.second;
21  Etot += energy_X.first;
22  }
23 
24  float X_mean = 0;
25  if (Etot > 0)
26  X_mean = X_sum / Etot;
27  return X_mean;
28 }
29 
31  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
32 
33  int firstLayer = 999;
34 
35  for (const auto& id_clu : clustersPtrs) {
36  if (!pass(*id_clu.second, c3d))
37  continue;
38  int layer = triggerTools_.layerWithOffset(id_clu.second->detId());
39  if (layer < firstLayer)
40  firstLayer = layer;
41  }
42 
43  return firstLayer;
44 }
45 
47  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
48  std::unordered_map<int, float> layers_pt;
49  float max_pt = 0.;
50  int max_layer = 0;
51  for (const auto& id_cluster : clustersPtrs) {
52  if (!pass(*id_cluster.second, c3d))
53  continue;
54  unsigned layer = triggerTools_.layerWithOffset(id_cluster.second->detId());
55  auto itr_insert = layers_pt.emplace(layer, 0.);
56  itr_insert.first->second += id_cluster.second->pt();
57  if (itr_insert.first->second > max_pt) {
58  max_pt = itr_insert.first->second;
59  max_layer = layer;
60  }
61  }
62  return max_layer;
63 }
64 
66  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
67 
68  int lastLayer = -999;
69 
70  for (const auto& id_clu : clustersPtrs) {
71  if (!pass(*id_clu.second, c3d))
72  continue;
73  int layer = triggerTools_.layerWithOffset(id_clu.second->detId());
74  if (layer > lastLayer)
75  lastLayer = layer;
76  }
77 
78  return lastLayer;
79 }
80 
82  const HGCalTriggerGeometryBase& triggerGeometry) const {
83  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
85  std::vector<bool> layers(nlayers);
86  for (const auto& id_cluster : clustersPtrs) {
87  if (!pass(*id_cluster.second, c3d))
88  continue;
89  unsigned layer = triggerGeometry.triggerLayer(id_cluster.second->detId());
90  if (triggerTools_.isNose(id_cluster.second->detId()))
92  else {
94  }
95  if (layer == 0 || layer > nlayers)
96  continue;
97  layers[layer - 1] = true; //layer 0 doesn't exist, so shift by -1
98  }
99  int length = 0;
100  int maxlength = 0;
101  for (bool layer : layers) {
102  if (layer)
103  length++;
104  else
105  length = 0;
106  if (length > maxlength)
107  maxlength = length;
108  }
109  return maxlength;
110 }
111 
113  const HGCalTriggerGeometryBase& triggerGeometry,
114  float quantile) const {
115  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
117  std::vector<double> layers(nlayers, 0);
118  for (const auto& id_clu : clustersPtrs) {
119  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
120 
121  for (const auto& id_tc : triggerCells) {
122  if (!pass(*id_tc.second, c3d))
123  continue;
124  unsigned layer = triggerGeometry.triggerLayer(id_tc.second->detId());
125  if (triggerTools_.isNose(id_tc.second->detId()))
127  else {
129  }
130  if (layer == 0 || layer > nlayers)
131  continue;
132  layers[layer - 1] += id_tc.second->pt(); //layer 0 doesn't exist, so shift by -1
133  }
134  }
135  std::partial_sum(layers.begin(), layers.end(), layers.begin());
136  double pt_threshold = layers.back() * quantile;
137  unsigned percentile = 0;
138  for (double pt : layers) {
139  if (pt > pt_threshold) {
140  break;
141  }
142  percentile++;
143  }
144  // Linear interpolation of percentile value
145  double pt0 = (percentile > 0 ? layers[percentile - 1] : 0.);
146  double pt1 = (percentile < layers.size() ? layers[percentile] : layers.back());
147  return percentile + (pt1 - pt0 > 0. ? (pt_threshold - pt0) / (pt1 - pt0) : 0.);
148 }
149 
151  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
152  std::set<double> ordered_tcs;
153  double pt_sum = 0.;
154  for (const auto& id_clu : clustersPtrs) {
155  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
156 
157  for (const auto& id_tc : triggerCells) {
158  if (!pass(*id_tc.second, c3d))
159  continue;
160  ordered_tcs.emplace(id_tc.second->pt());
161  pt_sum += id_tc.second->pt();
162  }
163  }
164  double pt_threshold = pt_sum * quantile;
165  double partial_sum = 0.;
166  double partial_sum_prev = 0.;
167  int ntc = 0;
168  for (auto itr = ordered_tcs.rbegin(); itr != ordered_tcs.rend(); ++itr) {
169  partial_sum_prev = partial_sum;
170  partial_sum += *itr;
171  ntc++;
172  if (partial_sum > pt_threshold) {
173  break;
174  }
175  }
176  // Linear interpolation of ntc
177  return ntc - 1 +
178  (partial_sum - partial_sum_prev > 0. ? (pt_threshold - partial_sum_prev) / (partial_sum - partial_sum_prev)
179  : 0.);
180 }
181 
183  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
184 
185  std::vector<std::pair<float, float>> tc_energy_eta;
186 
187  for (const auto& id_clu : clustersPtrs) {
188  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
189 
190  for (const auto& id_tc : triggerCells) {
191  if (!pass(*id_tc.second, c3d))
192  continue;
193  tc_energy_eta.emplace_back(std::make_pair(id_tc.second->energy(), id_tc.second->eta()));
194  }
195  }
196 
197  float SeeTot = sigmaXX(tc_energy_eta, c3d.eta());
198 
199  return SeeTot;
200 }
201 
203  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
204 
205  std::vector<std::pair<float, float>> tc_energy_phi;
206 
207  for (const auto& id_clu : clustersPtrs) {
208  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
209 
210  for (const auto& id_tc : triggerCells) {
211  if (!pass(*id_tc.second, c3d))
212  continue;
213  tc_energy_phi.emplace_back(std::make_pair(id_tc.second->energy(), id_tc.second->phi()));
214  }
215  }
216 
217  float SppTot = sigmaPhiPhi(tc_energy_phi, c3d.phi());
218 
219  return SppTot;
220 }
221 
223  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
224 
225  std::vector<std::pair<float, float>> tc_energy_r;
226 
227  for (const auto& id_clu : clustersPtrs) {
228  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
229 
230  for (const auto& id_tc : triggerCells) {
231  if (!pass(*id_tc.second, c3d))
232  continue;
233  float r = (id_tc.second->position().z() != 0.
234  ? std::sqrt(pow(id_tc.second->position().x(), 2) + pow(id_tc.second->position().y(), 2)) /
235  std::abs(id_tc.second->position().z())
236  : 0.);
237  tc_energy_r.emplace_back(std::make_pair(id_tc.second->energy(), r));
238  }
239  }
240 
241  float r_mean = meanX(tc_energy_r);
242  float Szz = sigmaXX(tc_energy_r, r_mean);
243 
244  return Szz;
245 }
246 
248  std::unordered_map<int, std::vector<std::pair<float, float>>> tc_layer_energy_eta;
249  std::unordered_map<int, LorentzVector> layer_LV;
250 
251  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
252 
253  for (const auto& id_clu : clustersPtrs) {
254  unsigned layer = triggerTools_.layerWithOffset(id_clu.second->detId());
255 
256  layer_LV[layer] += id_clu.second->p4();
257 
258  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
259 
260  for (const auto& id_tc : triggerCells) {
261  if (!pass(*id_tc.second, c3d))
262  continue;
263  tc_layer_energy_eta[layer].emplace_back(std::make_pair(id_tc.second->energy(), id_tc.second->eta()));
264  }
265  }
266 
267  float SigmaEtaEtaMax = 0;
268 
269  for (auto& tc_iter : tc_layer_energy_eta) {
270  const std::vector<std::pair<float, float>>& energy_eta_layer = tc_iter.second;
271  const LorentzVector& LV_layer = layer_LV[tc_iter.first];
272  float SigmaEtaEtaLayer = sigmaXX(energy_eta_layer, LV_layer.eta()); //RMS wrt layer eta, not wrt c3d eta
273  if (SigmaEtaEtaLayer > SigmaEtaEtaMax)
274  SigmaEtaEtaMax = SigmaEtaEtaLayer;
275  }
276 
277  return SigmaEtaEtaMax;
278 }
279 
281  std::unordered_map<int, std::vector<std::pair<float, float>>> tc_layer_energy_phi;
282  std::unordered_map<int, LorentzVector> layer_LV;
283 
284  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
285 
286  for (const auto& id_clu : clustersPtrs) {
287  unsigned layer = triggerTools_.layerWithOffset(id_clu.second->detId());
288 
289  layer_LV[layer] += id_clu.second->p4();
290 
291  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
292 
293  for (const auto& id_tc : triggerCells) {
294  if (!pass(*id_tc.second, c3d))
295  continue;
296  tc_layer_energy_phi[layer].emplace_back(std::make_pair(id_tc.second->energy(), id_tc.second->phi()));
297  }
298  }
299 
300  float SigmaPhiPhiMax = 0;
301 
302  for (auto& tc_iter : tc_layer_energy_phi) {
303  const std::vector<std::pair<float, float>>& energy_phi_layer = tc_iter.second;
304  const LorentzVector& LV_layer = layer_LV[tc_iter.first];
305  float SigmaPhiPhiLayer = sigmaPhiPhi(energy_phi_layer, LV_layer.phi()); //RMS wrt layer phi, not wrt c3d phi
306  if (SigmaPhiPhiLayer > SigmaPhiPhiMax)
307  SigmaPhiPhiMax = SigmaPhiPhiLayer;
308  }
309 
310  return SigmaPhiPhiMax;
311 }
312 
314  std::unordered_map<int, std::vector<std::pair<float, float>>> tc_layer_energy_r;
315 
316  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
317 
318  for (const auto& id_clu : clustersPtrs) {
319  unsigned layer = triggerTools_.layerWithOffset(id_clu.second->detId());
320 
321  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
322 
323  for (const auto& id_tc : triggerCells) {
324  if (!pass(*id_tc.second, c3d))
325  continue;
326  float r = (id_tc.second->position().z() != 0.
327  ? std::sqrt(pow(id_tc.second->position().x(), 2) + pow(id_tc.second->position().y(), 2)) /
328  std::abs(id_tc.second->position().z())
329  : 0.);
330  tc_layer_energy_r[layer].emplace_back(std::make_pair(id_tc.second->energy(), r));
331  }
332  }
333 
334  float SigmaRRMax = 0;
335 
336  for (auto& tc_iter : tc_layer_energy_r) {
337  const std::vector<std::pair<float, float>>& energy_r_layer = tc_iter.second;
338  float r_mean_layer = meanX(energy_r_layer);
339  float SigmaRRLayer = sigmaXX(energy_r_layer, r_mean_layer);
340  if (SigmaRRLayer > SigmaRRMax)
341  SigmaRRMax = SigmaRRLayer;
342  }
343 
344  return SigmaRRMax;
345 }
346 
348  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
349  // group trigger cells by layer
350  std::unordered_map<int, std::vector<edm::Ptr<l1t::HGCalTriggerCell>>> layers_tcs;
351  for (const auto& id_clu : clustersPtrs) {
352  unsigned layer = triggerTools_.layerWithOffset(id_clu.second->detId());
353  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
354  for (const auto& id_tc : triggerCells) {
355  if (!pass(*id_tc.second, c3d))
356  continue;
357  layers_tcs[layer].emplace_back(id_tc.second);
358  }
359  }
360 
361  // Select trigger cells within X cm of the max TC in the layer
362  std::unordered_map<int, std::vector<std::pair<float, float>>> tc_layers_energy_r;
363  for (const auto& layer_tcs : layers_tcs) {
364  int layer = layer_tcs.first;
365  edm::Ptr<l1t::HGCalTriggerCell> max_tc = layer_tcs.second.front();
366  for (const auto& tc : layer_tcs.second) {
367  if (tc->energy() > max_tc->energy())
368  max_tc = tc;
369  }
370  for (const auto& tc : layer_tcs.second) {
371  double dx = tc->position().x() - max_tc->position().x();
372  double dy = tc->position().y() - max_tc->position().y();
373  double distance_to_max = std::sqrt(dx * dx + dy * dy);
374  if (distance_to_max < radius) {
375  float r = (tc->position().z() != 0.
376  ? std::sqrt(tc->position().x() * tc->position().x() + tc->position().y() * tc->position().y()) /
377  std::abs(tc->position().z())
378  : 0.);
379  tc_layers_energy_r[layer].emplace_back(std::make_pair(tc->energy(), r));
380  }
381  }
382  }
383 
384  // Compute srr layer by layer
385  std::vector<std::pair<float, float>> layers_energy_srr2;
386  for (const auto& layer_energy_r : tc_layers_energy_r) {
387  const auto& energies_r = layer_energy_r.second;
388  float r_mean_layer = meanX(energies_r);
389  float srr = sigmaXX(energies_r, r_mean_layer);
390  double energy_sum = 0.;
391  for (const auto& energy_r : energies_r) {
392  energy_sum += energy_r.first;
393  }
394  layers_energy_srr2.emplace_back(std::make_pair(energy_sum, srr * srr));
395  }
396  // Combine all layer srr
397  float srr2_mean = meanX(layers_energy_srr2);
398  return std::sqrt(srr2_mean);
399 }
400 
402  std::unordered_map<int, float> layer_energy;
403 
404  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
405 
406  for (const auto& id_clu : clustersPtrs) {
407  if (!pass(*id_clu.second, c3d))
408  continue;
409  unsigned layer = triggerTools_.layerWithOffset(id_clu.second->detId());
410  layer_energy[layer] += id_clu.second->energy();
411  }
412 
413  float EMax = 0;
414 
415  for (const auto& layer : layer_energy) {
416  if (layer.second > EMax)
417  EMax = layer.second;
418  }
419 
420  return EMax;
421 }
422 
424  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
425 
426  std::vector<std::pair<float, float>> tc_energy_z;
427 
428  for (const auto& id_clu : clustersPtrs) {
429  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
430 
431  for (const auto& id_tc : triggerCells) {
432  if (!pass(*id_tc.second, c3d))
433  continue;
434  tc_energy_z.emplace_back(id_tc.second->energy(), id_tc.second->position().z());
435  }
436  }
437 
438  return meanX(tc_energy_z);
439 }
440 
442  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
443 
444  std::vector<std::pair<float, float>> tc_energy_z;
445 
446  for (const auto& id_clu : clustersPtrs) {
447  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
448 
449  for (const auto& id_tc : triggerCells) {
450  if (!pass(*id_tc.second, c3d))
451  continue;
452  tc_energy_z.emplace_back(std::make_pair(id_tc.second->energy(), id_tc.second->position().z()));
453  }
454  }
455 
456  float z_mean = meanX(tc_energy_z);
457  float Szz = sigmaXX(tc_energy_z, z_mean);
458 
459  return Szz;
460 }
461 
463  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& cellsPtrs = c2d.constituents();
464 
465  std::vector<std::pair<float, float>> tc_energy_eta;
466 
467  for (const auto& id_cell : cellsPtrs) {
468  if (!pass(*id_cell.second, c2d))
469  continue;
470  tc_energy_eta.emplace_back(std::make_pair(id_cell.second->energy(), id_cell.second->eta()));
471  }
472 
473  float See = sigmaXX(tc_energy_eta, c2d.eta());
474 
475  return See;
476 }
477 
479  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& cellsPtrs = c2d.constituents();
480 
481  std::vector<std::pair<float, float>> tc_energy_phi;
482 
483  for (const auto& id_cell : cellsPtrs) {
484  if (!pass(*id_cell.second, c2d))
485  continue;
486  tc_energy_phi.emplace_back(std::make_pair(id_cell.second->energy(), id_cell.second->phi()));
487  }
488 
489  float Spp = sigmaPhiPhi(tc_energy_phi, c2d.phi());
490 
491  return Spp;
492 }
493 
495  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& cellsPtrs = c2d.constituents();
496 
497  std::vector<std::pair<float, float>> tc_energy_r;
498 
499  for (const auto& id_cell : cellsPtrs) {
500  if (!pass(*id_cell.second, c2d))
501  continue;
502  float r = (id_cell.second->position().z() != 0.
503  ? std::sqrt(pow(id_cell.second->position().x(), 2) + pow(id_cell.second->position().y(), 2)) /
504  std::abs(id_cell.second->position().z())
505  : 0.);
506  tc_energy_r.emplace_back(std::make_pair(id_cell.second->energy(), r));
507  }
508 
509  float r_mean = meanX(tc_energy_r);
510  float Srr = sigmaXX(tc_energy_r, r_mean);
511 
512  return Srr;
513 }
514 
516  c3d.showerLength(showerLength(c3d));
517  c3d.coreShowerLength(coreShowerLength(c3d, triggerGeometry));
518  c3d.firstLayer(firstLayer(c3d));
519  c3d.maxLayer(maxLayer(c3d));
520  c3d.sigmaEtaEtaTot(sigmaEtaEtaTot(c3d));
521  c3d.sigmaEtaEtaMax(sigmaEtaEtaMax(c3d));
522  c3d.sigmaPhiPhiTot(sigmaPhiPhiTot(c3d));
523  c3d.sigmaPhiPhiMax(sigmaPhiPhiMax(c3d));
524  c3d.sigmaZZ(sigmaZZ(c3d));
525  c3d.sigmaRRTot(sigmaRRTot(c3d));
526  c3d.sigmaRRMax(sigmaRRMax(c3d));
527  c3d.sigmaRRMean(sigmaRRMean(c3d));
528  c3d.eMax(eMax(c3d));
529  c3d.zBarycenter(meanZ(c3d));
530  c3d.layer10percent(percentileLayer(c3d, triggerGeometry, 0.10));
531  c3d.layer50percent(percentileLayer(c3d, triggerGeometry, 0.50));
532  c3d.layer90percent(percentileLayer(c3d, triggerGeometry, 0.90));
535 }
HGCalTriggerGeometryBase
Definition: HGCalTriggerGeometryBase.h:19
HGCalShowerShape::fillShapes
void fillShapes(l1t::HGCalMulticluster &, const HGCalTriggerGeometryBase &) const
Definition: HGCalShowerShape.cc:515
l1t::HGCalCluster
Definition: HGCalCluster.h:11
l1t::HGCalClusterT::eMax
float eMax() const
Definition: HGCalClusterT.h:149
HGCalShowerShape::pass
bool pass(const T &obj, const Tref &ref) const
Definition: HGCalShowerShape.h:85
ForwardEmpty
Definition: ForwardSubdetector.h:5
l1t::HGCalClusterT::zBarycenter
float zBarycenter() const
Definition: HGCalClusterT.h:158
PV3DBase::x
T x() const
Definition: PV3DBase.h:59
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
deltaPhi.h
l1t::HGCalClusterT::showerLength
int showerLength() const
Definition: HGCalClusterT.h:145
l1t::HGCalClusterT::sigmaEtaEtaTot
float sigmaEtaEtaTot() const
Definition: HGCalClusterT.h:152
HGCalShowerShape::sigmaRRMax
float sigmaRRMax(const l1t::HGCalMulticluster &c3d) const
Definition: HGCalShowerShape.cc:313
l1t::HGCalClusterT::layer10percent
float layer10percent() const
Definition: HGCalClusterT.h:159
l1t::HGCalClusterT::triggerCells90percent
float triggerCells90percent() const
Definition: HGCalClusterT.h:163
HGCalShowerShape::maxLayer
int maxLayer(const l1t::HGCalMulticluster &c3d) const
Definition: HGCalShowerShape.cc:46
HGCalShowerShape::sigmaRRTot
float sigmaRRTot(const l1t::HGCalMulticluster &c3d) const
Definition: HGCalShowerShape.cc:222
HGCalShowerShape.h
l1t::HGCalClusterT::sigmaRRTot
float sigmaRRTot() const
Definition: HGCalClusterT.h:155
l1t::HGCalMulticluster
Definition: HGCalMulticluster.h:13
l1t::HGCalClusterT::sigmaPhiPhiMax
float sigmaPhiPhiMax() const
Definition: HGCalClusterT.h:151
HGCalTriggerGeometryBase::triggerLayer
virtual unsigned triggerLayer(const unsigned id) const =0
HGCalMulticluster.h
HFNose
Definition: ForwardSubdetector.h:11
l1t::HGCalClusterT::sigmaPhiPhiTot
float sigmaPhiPhiTot() const
Definition: HGCalClusterT.h:153
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
HGCalShowerShape::sigmaPhiPhiTot
float sigmaPhiPhiTot(const l1t::HGCalMulticluster &c3d) const
Definition: HGCalShowerShape.cc:202
l1t::HGCalClusterT::layer90percent
float layer90percent() const
Definition: HGCalClusterT.h:161
HGCalShowerShape::sigmaRRMean
float sigmaRRMean(const l1t::HGCalMulticluster &c3d, float radius=5.) const
Definition: HGCalShowerShape.cc:347
HGCalShowerShape::sigmaEtaEtaMax
float sigmaEtaEtaMax(const l1t::HGCalMulticluster &c3d) const
Definition: HGCalShowerShape.cc:247
HLT_FULL_cff.pt1
pt1
Definition: HLT_FULL_cff.py:9870
HGCalTriggerTools::layerWithOffset
unsigned layerWithOffset(const DetId &) const
Definition: HGCalTriggerTools.cc:134
l1t::HGCalClusterT::constituents
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50
HGCalShowerShape::sigmaXX
float sigmaXX(const std::vector< pair< float, float >> &energy_X_tc, const float X_cluster) const
Definition: HGCalShowerShape.h:62
HGCalShowerShape::meanZ
float meanZ(const l1t::HGCalMulticluster &c3d) const
Definition: HGCalShowerShape.cc:423
HGCalShowerShape::firstLayer
int firstLayer(const l1t::HGCalMulticluster &c3d) const
Definition: HGCalShowerShape.cc:30
HGCalTriggerTools::isNose
bool isNose(const DetId &) const
Definition: HGCalTriggerTools.cc:170
SusyPostProcessor_cff.quantile
quantile
Definition: SusyPostProcessor_cff.py:7
distance_
double distance_
Definition: PFRecoTauChargedHadronFromGenericTrackPlugin.cc:198
edm::ParameterSet
Definition: ParameterSet.h:47
HGCalShowerShape::sigmaPhiPhi
float sigmaPhiPhi(const std::vector< pair< float, float >> &energy_phi_tc, const float phi_cluster) const
Definition: HGCalShowerShape.h:81
reco::LeafCandidate::eta
double eta() const final
momentum pseudorapidity
Definition: LeafCandidate.h:152
HGCalShowerShape::meanX
float meanX(const std::vector< pair< float, float >> &energy_X_tc) const
Definition: HGCalShowerShape.cc:15
PV3DBase::y
T y() const
Definition: PV3DBase.h:60
l1t::HGCalClusterT::firstLayer
int firstLayer() const
Definition: HGCalClusterT.h:147
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
l1t::HGCalClusterT::layer50percent
float layer50percent() const
Definition: HGCalClusterT.h:160
PVValHelper::dy
Definition: PVValidationHelpers.h:49
HGCalTriggerTools::layers
unsigned layers(ForwardSubdetector type) const
Definition: HGCalTriggerTools.cc:65
HGCalCluster.h
l1t::HGCalClusterT::coreShowerLength
int coreShowerLength() const
Definition: HGCalClusterT.h:146
itr
std::vector< std::pair< float, float > >::iterator itr
Definition: HGCDigitizer.cc:29
HGCalShowerShape::triggerTools_
HGCalTriggerTools triggerTools_
Definition: HGCalShowerShape.h:92
HGCalShowerShape::percentileTriggerCells
float percentileTriggerCells(const l1t::HGCalMulticluster &c3d, float quantile=0.5) const
Definition: HGCalShowerShape.cc:150
edm::Ptr
Definition: AssociationVector.h:31
alignCSCRings.r
r
Definition: alignCSCRings.py:93
HGCalShowerShape::lastLayer
int lastLayer(const l1t::HGCalMulticluster &c3d) const
Definition: HGCalShowerShape.cc:65
l1t::HGCalClusterT::maxLayer
int maxLayer() const
Definition: HGCalClusterT.h:148
HGCalShowerShape::percentileLayer
float percentileLayer(const l1t::HGCalMulticluster &c3d, const HGCalTriggerGeometryBase &triggerGeometry, float quantile=0.5) const
Definition: HGCalShowerShape.cc:112
l1t::HGCalClusterT::sigmaRRMean
float sigmaRRMean() const
Definition: HGCalClusterT.h:157
reco::LeafCandidate::phi
double phi() const final
momentum azimuthal angle
Definition: LeafCandidate.h:148
l1t::HGCalClusterT::sigmaRRMax
float sigmaRRMax() const
Definition: HGCalClusterT.h:156
l1t::HGCalClusterT::sigmaZZ
float sigmaZZ() const
Definition: HGCalClusterT.h:154
l1t::HGCalTriggerCell::position
const GlobalPoint & position() const
Definition: HGCalTriggerCell.h:26
HGCalShowerShape::eMax
float eMax(const l1t::HGCalMulticluster &c3d) const
Definition: HGCalShowerShape.cc:401
HGCalShowerShape::sigmaZZ
float sigmaZZ(const l1t::HGCalMulticluster &c3d) const
Definition: HGCalShowerShape.cc:441
CosmicsPD_Skims.radius
radius
Definition: CosmicsPD_Skims.py:135
HGCalShowerShape::showerLength
int showerLength(const l1t::HGCalMulticluster &c3d) const
Definition: HGCalShowerShape.h:26
reco::LeafCandidate::energy
double energy() const final
energy
Definition: LeafCandidate.h:125
HGCalShowerShape::sigmaEtaEtaTot
float sigmaEtaEtaTot(const l1t::HGCalMulticluster &c3d) const
Definition: HGCalShowerShape.cc:182
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
PVValHelper::dx
Definition: PVValidationHelpers.h:48
hgcalTopologyTester_cfi.layers
layers
Definition: hgcalTopologyTester_cfi.py:8
l1t::HGCalClusterT::triggerCells67percent
float triggerCells67percent() const
Definition: HGCalClusterT.h:162
HGCalShowerShape::coreShowerLength
int coreShowerLength(const l1t::HGCalMulticluster &c3d, const HGCalTriggerGeometryBase &triggerGeometry) const
Definition: HGCalShowerShape.cc:81
caloTruthCellsProducer_cfi.triggerCells
triggerCells
Definition: caloTruthCellsProducer_cfi.py:7
l1t::HGCalClusterT::sigmaEtaEtaMax
float sigmaEtaEtaMax() const
Definition: HGCalClusterT.h:150
HGCalShowerShape::LorentzVector
math::XYZTLorentzVector LorentzVector
Definition: HGCalShowerShape.h:14
nlayers
Definition: HIMultiTrackSelector.h:48
HGCalShowerShape::HGCalShowerShape
HGCalShowerShape()
Definition: HGCalShowerShape.h:16
HGCalShowerShape::sigmaPhiPhiMax
float sigmaPhiPhiMax(const l1t::HGCalMulticluster &c3d) const
Definition: HGCalShowerShape.cc:280