CMS 3D CMS Logo

List of all members | Classes | Public Types | Public Member Functions | Private Member Functions | Private Attributes
HGCalShowerShape Class Reference

#include <HGCalShowerShape.h>

Classes

struct  DeltaPhi
 

Public Types

typedef math::XYZTLorentzVector LorentzVector
 

Public Member Functions

int coreShowerLength (const l1t::HGCalMulticluster &c3d, const HGCalTriggerGeometryBase &triggerGeometry) const
 
float eMax (const l1t::HGCalMulticluster &c3d) const
 
void eventSetup (const edm::EventSetup &es)
 
void fillShapes (l1t::HGCalMulticluster &, const HGCalTriggerGeometryBase &) const
 
int firstLayer (const l1t::HGCalMulticluster &c3d) const
 
 HGCalShowerShape ()
 
int lastLayer (const l1t::HGCalMulticluster &c3d) const
 
int maxLayer (const l1t::HGCalMulticluster &c3d) const
 
int showerLength (const l1t::HGCalMulticluster &c3d) const
 
float sigmaEtaEtaMax (const l1t::HGCalMulticluster &c3d) const
 
float sigmaEtaEtaTot (const l1t::HGCalMulticluster &c3d) const
 
float sigmaEtaEtaTot (const l1t::HGCalCluster &c2d) const
 
float sigmaPhiPhiMax (const l1t::HGCalMulticluster &c3d) const
 
float sigmaPhiPhiTot (const l1t::HGCalMulticluster &c3d) const
 
float sigmaPhiPhiTot (const l1t::HGCalCluster &c2d) const
 
float sigmaRRMax (const l1t::HGCalMulticluster &c3d) const
 
float sigmaRRMean (const l1t::HGCalMulticluster &c3d, float radius=5.) const
 
float sigmaRRTot (const l1t::HGCalMulticluster &c3d) const
 
float sigmaRRTot (const l1t::HGCalCluster &c2d) const
 
float sigmaZZ (const l1t::HGCalMulticluster &c3d) const
 
 ~HGCalShowerShape ()
 

Private Member Functions

float meanX (const std::vector< pair< float, float >> &energy_X_tc) const
 
float sigmaPhiPhi (const std::vector< pair< float, float >> &energy_phi_tc, const float phi_cluster) const
 
template<typename Delta = std::minus<float>>
float sigmaXX (const std::vector< pair< float, float >> &energy_X_tc, const float X_cluster) const
 

Private Attributes

HGCalTriggerTools triggerTools_
 

Detailed Description

Definition at line 12 of file HGCalShowerShape.h.

Member Typedef Documentation

Definition at line 14 of file HGCalShowerShape.h.

Constructor & Destructor Documentation

HGCalShowerShape::HGCalShowerShape ( )
inline

Definition at line 16 of file HGCalShowerShape.h.

16 {}
HGCalShowerShape::~HGCalShowerShape ( )
inline

Definition at line 18 of file HGCalShowerShape.h.

18 {}

Member Function Documentation

int HGCalShowerShape::coreShowerLength ( const l1t::HGCalMulticluster c3d,
const HGCalTriggerGeometryBase triggerGeometry 
) const

Definition at line 70 of file HGCalShowerShape.cc.

References l1t::HGCalClusterT< C >::constituents(), ForwardEmpty, LayerTriplets::layers(), HGCalTriggerTools::layers(), nlayers, HGCalTriggerGeometryBase::triggerLayer(), and triggerTools_.

Referenced by fillShapes(), and showerLength().

71  {
72  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
74  std::vector<bool> layers(nlayers);
75  for (const auto& id_cluster : clustersPtrs) {
76  unsigned layer = triggerGeometry.triggerLayer(id_cluster.second->detId());
77  if (layer == 0 || layer > nlayers)
78  continue;
79  layers[layer - 1] = true; //layer 0 doesn't exist, so shift by -1
80  }
81  int length = 0;
82  int maxlength = 0;
83  for (bool layer : layers) {
84  if (layer)
85  length++;
86  else
87  length = 0;
88  if (length > maxlength)
89  maxlength = length;
90  }
91  return maxlength;
92 }
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
HGCalTriggerTools triggerTools_
virtual unsigned triggerLayer(const unsigned id) const =0
unsigned layers(ForwardSubdetector type) const
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50
float HGCalShowerShape::eMax ( const l1t::HGCalMulticluster c3d) const

Definition at line 299 of file HGCalShowerShape.cc.

References l1t::HGCalClusterT< C >::constituents(), HGCalTriggerTools::layerWithOffset(), and triggerTools_.

Referenced by fillShapes(), and showerLength().

299  {
300  std::unordered_map<int, float> layer_energy;
301 
302  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
303 
304  for (const auto& id_clu : clustersPtrs) {
305  unsigned layer = triggerTools_.layerWithOffset(id_clu.second->detId());
306  layer_energy[layer] += id_clu.second->energy();
307  }
308 
309  float EMax = 0;
310 
311  for (const auto& layer : layer_energy) {
312  if (layer.second > EMax)
313  EMax = layer.second;
314  }
315 
316  return EMax;
317 }
HGCalTriggerTools triggerTools_
unsigned layerWithOffset(const DetId &) const
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50
void HGCalShowerShape::eventSetup ( const edm::EventSetup es)
inline
void HGCalShowerShape::fillShapes ( l1t::HGCalMulticluster c3d,
const HGCalTriggerGeometryBase triggerGeometry 
) const

Definition at line 385 of file HGCalShowerShape.cc.

References coreShowerLength(), l1t::HGCalClusterT< C >::coreShowerLength(), eMax(), l1t::HGCalClusterT< C >::eMax(), firstLayer(), l1t::HGCalClusterT< C >::firstLayer(), maxLayer(), l1t::HGCalClusterT< C >::maxLayer(), showerLength(), l1t::HGCalClusterT< C >::showerLength(), sigmaEtaEtaMax(), l1t::HGCalClusterT< C >::sigmaEtaEtaMax(), sigmaEtaEtaTot(), l1t::HGCalClusterT< C >::sigmaEtaEtaTot(), sigmaPhiPhiMax(), l1t::HGCalClusterT< C >::sigmaPhiPhiMax(), sigmaPhiPhiTot(), l1t::HGCalClusterT< C >::sigmaPhiPhiTot(), sigmaRRMax(), l1t::HGCalClusterT< C >::sigmaRRMax(), sigmaRRMean(), l1t::HGCalClusterT< C >::sigmaRRMean(), sigmaRRTot(), l1t::HGCalClusterT< C >::sigmaRRTot(), sigmaZZ(), and l1t::HGCalClusterT< C >::sigmaZZ().

Referenced by HGCalMulticlusteringImpl::finalizeClusters(), HGCalHistoClusteringImpl::finalizeClusters(), and showerLength().

385  {
386  c3d.showerLength(showerLength(c3d));
387  c3d.coreShowerLength(coreShowerLength(c3d, triggerGeometry));
388  c3d.firstLayer(firstLayer(c3d));
389  c3d.maxLayer(maxLayer(c3d));
390  c3d.sigmaEtaEtaTot(sigmaEtaEtaTot(c3d));
391  c3d.sigmaEtaEtaMax(sigmaEtaEtaMax(c3d));
392  c3d.sigmaPhiPhiTot(sigmaPhiPhiTot(c3d));
393  c3d.sigmaPhiPhiMax(sigmaPhiPhiMax(c3d));
394  c3d.sigmaZZ(sigmaZZ(c3d));
395  c3d.sigmaRRTot(sigmaRRTot(c3d));
396  c3d.sigmaRRMax(sigmaRRMax(c3d));
397  c3d.sigmaRRMean(sigmaRRMean(c3d));
398  c3d.eMax(eMax(c3d));
399 }
float sigmaRRMax(const l1t::HGCalMulticluster &c3d) const
int maxLayer() const
float sigmaPhiPhiTot() const
float sigmaEtaEtaMax() const
float sigmaRRTot() const
int coreShowerLength(const l1t::HGCalMulticluster &c3d, const HGCalTriggerGeometryBase &triggerGeometry) const
float sigmaRRTot(const l1t::HGCalMulticluster &c3d) const
float sigmaPhiPhiTot(const l1t::HGCalMulticluster &c3d) const
int firstLayer(const l1t::HGCalMulticluster &c3d) const
float sigmaZZ() const
float sigmaRRMean(const l1t::HGCalMulticluster &c3d, float radius=5.) const
int showerLength() const
float sigmaEtaEtaMax(const l1t::HGCalMulticluster &c3d) const
float sigmaEtaEtaTot() const
float eMax() const
int coreShowerLength() const
int firstLayer() const
int showerLength(const l1t::HGCalMulticluster &c3d) const
float sigmaRRMean() const
float sigmaZZ(const l1t::HGCalMulticluster &c3d) const
float sigmaRRMax() const
float sigmaPhiPhiMax(const l1t::HGCalMulticluster &c3d) const
int maxLayer(const l1t::HGCalMulticluster &c3d) const
float sigmaPhiPhiMax() const
float eMax(const l1t::HGCalMulticluster &c3d) const
float sigmaEtaEtaTot(const l1t::HGCalMulticluster &c3d) const
int HGCalShowerShape::firstLayer ( const l1t::HGCalMulticluster c3d) const

Definition at line 25 of file HGCalShowerShape.cc.

References l1t::HGCalClusterT< C >::constituents(), HGCalTriggerTools::layerWithOffset(), and triggerTools_.

Referenced by eventSetup(), fillShapes(), and showerLength().

25  {
26  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
27 
28  int firstLayer = 999;
29 
30  for (const auto& id_clu : clustersPtrs) {
31  int layer = triggerTools_.layerWithOffset(id_clu.second->detId());
32  if (layer < firstLayer)
33  firstLayer = layer;
34  }
35 
36  return firstLayer;
37 }
HGCalTriggerTools triggerTools_
int firstLayer(const l1t::HGCalMulticluster &c3d) const
unsigned layerWithOffset(const DetId &) const
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50
int HGCalShowerShape::lastLayer ( const l1t::HGCalMulticluster c3d) const

Definition at line 56 of file HGCalShowerShape.cc.

References l1t::HGCalClusterT< C >::constituents(), HGCalTriggerTools::layerWithOffset(), and triggerTools_.

Referenced by eventSetup(), and showerLength().

56  {
57  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
58 
59  int lastLayer = -999;
60 
61  for (const auto& id_clu : clustersPtrs) {
62  int layer = triggerTools_.layerWithOffset(id_clu.second->detId());
63  if (layer > lastLayer)
64  lastLayer = layer;
65  }
66 
67  return lastLayer;
68 }
HGCalTriggerTools triggerTools_
unsigned layerWithOffset(const DetId &) const
int lastLayer(const l1t::HGCalMulticluster &c3d) const
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50
int HGCalShowerShape::maxLayer ( const l1t::HGCalMulticluster c3d) const

Definition at line 39 of file HGCalShowerShape.cc.

References l1t::HGCalClusterT< C >::constituents(), HGCalTriggerTools::layerWithOffset(), and triggerTools_.

Referenced by eventSetup(), and fillShapes().

39  {
40  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
41  std::unordered_map<int, float> layers_pt;
42  float max_pt = 0.;
43  int max_layer = 0;
44  for (const auto& id_cluster : clustersPtrs) {
45  unsigned layer = triggerTools_.layerWithOffset(id_cluster.second->detId());
46  auto itr_insert = layers_pt.emplace(layer, 0.);
47  itr_insert.first->second += id_cluster.second->pt();
48  if (itr_insert.first->second > max_pt) {
49  max_pt = itr_insert.first->second;
50  max_layer = layer;
51  }
52  }
53  return max_layer;
54 }
HGCalTriggerTools triggerTools_
unsigned layerWithOffset(const DetId &) const
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50
float HGCalShowerShape::meanX ( const std::vector< pair< float, float >> &  energy_X_tc) const
private

Definition at line 10 of file HGCalShowerShape.cc.

Referenced by showerLength(), sigmaRRMax(), sigmaRRMean(), sigmaRRTot(), and sigmaZZ().

10  {
11  float Etot = 0;
12  float X_sum = 0;
13 
14  for (const auto& energy_X : energy_X_tc) {
15  X_sum += energy_X.first * energy_X.second;
16  Etot += energy_X.first;
17  }
18 
19  float X_mean = 0;
20  if (Etot > 0)
21  X_mean = X_sum / Etot;
22  return X_mean;
23 }
int HGCalShowerShape::showerLength ( const l1t::HGCalMulticluster c3d) const
inline

Definition at line 25 of file HGCalShowerShape.h.

References coreShowerLength(), eMax(), fillShapes(), firstLayer(), lastLayer(), meanX(), TCMET_cfi::radius, sigmaEtaEtaMax(), sigmaEtaEtaTot(), sigmaPhiPhiMax(), sigmaPhiPhiTot(), sigmaRRMax(), sigmaRRMean(), sigmaRRTot(), and sigmaZZ().

Referenced by fillShapes().

25  {
26  return lastLayer(c3d) - firstLayer(c3d) + 1;
27  } //in number of layers
int firstLayer(const l1t::HGCalMulticluster &c3d) const
int lastLayer(const l1t::HGCalMulticluster &c3d) const
float HGCalShowerShape::sigmaEtaEtaMax ( const l1t::HGCalMulticluster c3d) const

Definition at line 153 of file HGCalShowerShape.cc.

References l1t::HGCalClusterT< C >::constituents(), HGCalTriggerTools::layerWithOffset(), sigmaXX(), and triggerTools_.

Referenced by fillShapes(), and showerLength().

153  {
154  std::unordered_map<int, std::vector<std::pair<float, float>>> tc_layer_energy_eta;
155  std::unordered_map<int, LorentzVector> layer_LV;
156 
157  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
158 
159  for (const auto& id_clu : clustersPtrs) {
160  unsigned layer = triggerTools_.layerWithOffset(id_clu.second->detId());
161 
162  layer_LV[layer] += id_clu.second->p4();
163 
164  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
165 
166  for (const auto& id_tc : triggerCells) {
167  tc_layer_energy_eta[layer].emplace_back(std::make_pair(id_tc.second->energy(), id_tc.second->eta()));
168  }
169  }
170 
171  float SigmaEtaEtaMax = 0;
172 
173  for (auto& tc_iter : tc_layer_energy_eta) {
174  const std::vector<std::pair<float, float>>& energy_eta_layer = tc_iter.second;
175  const LorentzVector& LV_layer = layer_LV[tc_iter.first];
176  float SigmaEtaEtaLayer = sigmaXX(energy_eta_layer, LV_layer.eta()); //RMS wrt layer eta, not wrt c3d eta
177  if (SigmaEtaEtaLayer > SigmaEtaEtaMax)
178  SigmaEtaEtaMax = SigmaEtaEtaLayer;
179  }
180 
181  return SigmaEtaEtaMax;
182 }
float sigmaXX(const std::vector< pair< float, float >> &energy_X_tc, const float X_cluster) const
HGCalTriggerTools triggerTools_
unsigned layerWithOffset(const DetId &) const
math::XYZTLorentzVector LorentzVector
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50
float HGCalShowerShape::sigmaEtaEtaTot ( const l1t::HGCalMulticluster c3d) const

Definition at line 94 of file HGCalShowerShape.cc.

References l1t::HGCalClusterT< C >::constituents(), reco::LeafCandidate::eta(), and sigmaXX().

Referenced by fillShapes(), and showerLength().

94  {
95  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
96 
97  std::vector<std::pair<float, float>> tc_energy_eta;
98 
99  for (const auto& id_clu : clustersPtrs) {
100  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
101 
102  for (const auto& id_tc : triggerCells) {
103  tc_energy_eta.emplace_back(std::make_pair(id_tc.second->energy(), id_tc.second->eta()));
104  }
105  }
106 
107  float SeeTot = sigmaXX(tc_energy_eta, c3d.eta());
108 
109  return SeeTot;
110 }
double eta() const final
momentum pseudorapidity
float sigmaXX(const std::vector< pair< float, float >> &energy_X_tc, const float X_cluster) const
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50
float HGCalShowerShape::sigmaEtaEtaTot ( const l1t::HGCalCluster c2d) const

Definition at line 338 of file HGCalShowerShape.cc.

References l1t::HGCalClusterT< C >::constituents(), reco::LeafCandidate::eta(), and sigmaXX().

338  {
339  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& cellsPtrs = c2d.constituents();
340 
341  std::vector<std::pair<float, float>> tc_energy_eta;
342 
343  for (const auto& id_cell : cellsPtrs) {
344  tc_energy_eta.emplace_back(std::make_pair(id_cell.second->energy(), id_cell.second->eta()));
345  }
346 
347  float See = sigmaXX(tc_energy_eta, c2d.eta());
348 
349  return See;
350 }
double eta() const final
momentum pseudorapidity
float sigmaXX(const std::vector< pair< float, float >> &energy_X_tc, const float X_cluster) const
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50
float HGCalShowerShape::sigmaPhiPhi ( const std::vector< pair< float, float >> &  energy_phi_tc,
const float  phi_cluster 
) const
inlineprivate

Definition at line 74 of file HGCalShowerShape.h.

Referenced by sigmaPhiPhiMax(), and sigmaPhiPhiTot().

74  {
75  return sigmaXX<DeltaPhi<float>>(energy_phi_tc, phi_cluster);
76  }
float HGCalShowerShape::sigmaPhiPhiMax ( const l1t::HGCalMulticluster c3d) const

Definition at line 184 of file HGCalShowerShape.cc.

References l1t::HGCalClusterT< C >::constituents(), HGCalTriggerTools::layerWithOffset(), sigmaPhiPhi(), and triggerTools_.

Referenced by fillShapes(), and showerLength().

184  {
185  std::unordered_map<int, std::vector<std::pair<float, float>>> tc_layer_energy_phi;
186  std::unordered_map<int, LorentzVector> layer_LV;
187 
188  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
189 
190  for (const auto& id_clu : clustersPtrs) {
191  unsigned layer = triggerTools_.layerWithOffset(id_clu.second->detId());
192 
193  layer_LV[layer] += id_clu.second->p4();
194 
195  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
196 
197  for (const auto& id_tc : triggerCells) {
198  tc_layer_energy_phi[layer].emplace_back(std::make_pair(id_tc.second->energy(), id_tc.second->phi()));
199  }
200  }
201 
202  float SigmaPhiPhiMax = 0;
203 
204  for (auto& tc_iter : tc_layer_energy_phi) {
205  const std::vector<std::pair<float, float>>& energy_phi_layer = tc_iter.second;
206  const LorentzVector& LV_layer = layer_LV[tc_iter.first];
207  float SigmaPhiPhiLayer = sigmaPhiPhi(energy_phi_layer, LV_layer.phi()); //RMS wrt layer phi, not wrt c3d phi
208  if (SigmaPhiPhiLayer > SigmaPhiPhiMax)
209  SigmaPhiPhiMax = SigmaPhiPhiLayer;
210  }
211 
212  return SigmaPhiPhiMax;
213 }
float sigmaPhiPhi(const std::vector< pair< float, float >> &energy_phi_tc, const float phi_cluster) const
HGCalTriggerTools triggerTools_
unsigned layerWithOffset(const DetId &) const
math::XYZTLorentzVector LorentzVector
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50
float HGCalShowerShape::sigmaPhiPhiTot ( const l1t::HGCalMulticluster c3d) const

Definition at line 112 of file HGCalShowerShape.cc.

References l1t::HGCalClusterT< C >::constituents(), reco::LeafCandidate::phi(), and sigmaPhiPhi().

Referenced by fillShapes(), and showerLength().

112  {
113  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
114 
115  std::vector<std::pair<float, float>> tc_energy_phi;
116 
117  for (const auto& id_clu : clustersPtrs) {
118  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
119 
120  for (const auto& id_tc : triggerCells) {
121  tc_energy_phi.emplace_back(std::make_pair(id_tc.second->energy(), id_tc.second->phi()));
122  }
123  }
124 
125  float SppTot = sigmaPhiPhi(tc_energy_phi, c3d.phi());
126 
127  return SppTot;
128 }
float sigmaPhiPhi(const std::vector< pair< float, float >> &energy_phi_tc, const float phi_cluster) const
double phi() const final
momentum azimuthal angle
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50
float HGCalShowerShape::sigmaPhiPhiTot ( const l1t::HGCalCluster c2d) const

Definition at line 352 of file HGCalShowerShape.cc.

References l1t::HGCalClusterT< C >::constituents(), reco::LeafCandidate::phi(), and sigmaPhiPhi().

352  {
353  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& cellsPtrs = c2d.constituents();
354 
355  std::vector<std::pair<float, float>> tc_energy_phi;
356 
357  for (const auto& id_cell : cellsPtrs) {
358  tc_energy_phi.emplace_back(std::make_pair(id_cell.second->energy(), id_cell.second->phi()));
359  }
360 
361  float Spp = sigmaPhiPhi(tc_energy_phi, c2d.phi());
362 
363  return Spp;
364 }
float sigmaPhiPhi(const std::vector< pair< float, float >> &energy_phi_tc, const float phi_cluster) const
double phi() const final
momentum azimuthal angle
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50
float HGCalShowerShape::sigmaRRMax ( const l1t::HGCalMulticluster c3d) const

Definition at line 215 of file HGCalShowerShape.cc.

References funct::abs(), l1t::HGCalClusterT< C >::constituents(), HGCalTriggerTools::layerWithOffset(), meanX(), funct::pow(), alignCSCRings::r, sigmaXX(), mathSSE::sqrt(), and triggerTools_.

Referenced by fillShapes(), and showerLength().

215  {
216  std::unordered_map<int, std::vector<std::pair<float, float>>> tc_layer_energy_r;
217 
218  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
219 
220  for (const auto& id_clu : clustersPtrs) {
221  unsigned layer = triggerTools_.layerWithOffset(id_clu.second->detId());
222 
223  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
224 
225  for (const auto& id_tc : triggerCells) {
226  float r = (id_tc.second->position().z() != 0.
227  ? std::sqrt(pow(id_tc.second->position().x(), 2) + pow(id_tc.second->position().y(), 2)) /
228  std::abs(id_tc.second->position().z())
229  : 0.);
230  tc_layer_energy_r[layer].emplace_back(std::make_pair(id_tc.second->energy(), r));
231  }
232  }
233 
234  float SigmaRRMax = 0;
235 
236  for (auto& tc_iter : tc_layer_energy_r) {
237  const std::vector<std::pair<float, float>>& energy_r_layer = tc_iter.second;
238  float r_mean_layer = meanX(energy_r_layer);
239  float SigmaRRLayer = sigmaXX(energy_r_layer, r_mean_layer);
240  if (SigmaRRLayer > SigmaRRMax)
241  SigmaRRMax = SigmaRRLayer;
242  }
243 
244  return SigmaRRMax;
245 }
float sigmaXX(const std::vector< pair< float, float >> &energy_X_tc, const float X_cluster) const
HGCalTriggerTools triggerTools_
unsigned layerWithOffset(const DetId &) const
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float meanX(const std::vector< pair< float, float >> &energy_X_tc) const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50
float HGCalShowerShape::sigmaRRMean ( const l1t::HGCalMulticluster c3d,
float  radius = 5. 
) const

Definition at line 247 of file HGCalShowerShape.cc.

References funct::abs(), l1t::HGCalClusterT< C >::constituents(), PVValHelper::dx, PVValHelper::dy, reco::LeafCandidate::energy(), HGCalTriggerTools::layerWithOffset(), meanX(), l1t::HGCalTriggerCell::position(), alignCSCRings::r, sigmaXX(), mathSSE::sqrt(), triggerTools_, PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

Referenced by fillShapes(), and showerLength().

247  {
248  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
249  // group trigger cells by layer
250  std::unordered_map<int, std::vector<edm::Ptr<l1t::HGCalTriggerCell>>> layers_tcs;
251  for (const auto& id_clu : clustersPtrs) {
252  unsigned layer = triggerTools_.layerWithOffset(id_clu.second->detId());
253  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
254  for (const auto& id_tc : triggerCells) {
255  layers_tcs[layer].emplace_back(id_tc.second);
256  }
257  }
258 
259  // Select trigger cells within X cm of the max TC in the layer
260  std::unordered_map<int, std::vector<std::pair<float, float>>> tc_layers_energy_r;
261  for (const auto& layer_tcs : layers_tcs) {
262  int layer = layer_tcs.first;
263  edm::Ptr<l1t::HGCalTriggerCell> max_tc = layer_tcs.second.front();
264  for (const auto& tc : layer_tcs.second) {
265  if (tc->energy() > max_tc->energy())
266  max_tc = tc;
267  }
268  for (const auto& tc : layer_tcs.second) {
269  double dx = tc->position().x() - max_tc->position().x();
270  double dy = tc->position().y() - max_tc->position().y();
271  double distance_to_max = std::sqrt(dx * dx + dy * dy);
272  if (distance_to_max < radius) {
273  float r = (tc->position().z() != 0.
274  ? std::sqrt(tc->position().x() * tc->position().x() + tc->position().y() * tc->position().y()) /
275  std::abs(tc->position().z())
276  : 0.);
277  tc_layers_energy_r[layer].emplace_back(std::make_pair(tc->energy(), r));
278  }
279  }
280  }
281 
282  // Compute srr layer by layer
283  std::vector<std::pair<float, float>> layers_energy_srr2;
284  for (const auto& layer_energy_r : tc_layers_energy_r) {
285  const auto& energies_r = layer_energy_r.second;
286  float r_mean_layer = meanX(energies_r);
287  float srr = sigmaXX(energies_r, r_mean_layer);
288  double energy_sum = 0.;
289  for (const auto& energy_r : energies_r) {
290  energy_sum += energy_r.first;
291  }
292  layers_energy_srr2.emplace_back(std::make_pair(energy_sum, srr * srr));
293  }
294  // Combine all layer srr
295  float srr2_mean = meanX(layers_energy_srr2);
296  return std::sqrt(srr2_mean);
297 }
float sigmaXX(const std::vector< pair< float, float >> &energy_X_tc, const float X_cluster) const
T y() const
Definition: PV3DBase.h:63
HGCalTriggerTools triggerTools_
unsigned layerWithOffset(const DetId &) const
const GlobalPoint & position() const
T sqrt(T t)
Definition: SSEVec.h:18
double energy() const final
energy
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float meanX(const std::vector< pair< float, float >> &energy_X_tc) const
T x() const
Definition: PV3DBase.h:62
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50
float HGCalShowerShape::sigmaRRTot ( const l1t::HGCalMulticluster c3d) const

Definition at line 130 of file HGCalShowerShape.cc.

References funct::abs(), l1t::HGCalClusterT< C >::constituents(), meanX(), funct::pow(), alignCSCRings::r, sigmaXX(), and mathSSE::sqrt().

Referenced by fillShapes(), and showerLength().

130  {
131  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
132 
133  std::vector<std::pair<float, float>> tc_energy_r;
134 
135  for (const auto& id_clu : clustersPtrs) {
136  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
137 
138  for (const auto& id_tc : triggerCells) {
139  float r = (id_tc.second->position().z() != 0.
140  ? std::sqrt(pow(id_tc.second->position().x(), 2) + pow(id_tc.second->position().y(), 2)) /
141  std::abs(id_tc.second->position().z())
142  : 0.);
143  tc_energy_r.emplace_back(std::make_pair(id_tc.second->energy(), r));
144  }
145  }
146 
147  float r_mean = meanX(tc_energy_r);
148  float Szz = sigmaXX(tc_energy_r, r_mean);
149 
150  return Szz;
151 }
float sigmaXX(const std::vector< pair< float, float >> &energy_X_tc, const float X_cluster) const
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float meanX(const std::vector< pair< float, float >> &energy_X_tc) const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50
float HGCalShowerShape::sigmaRRTot ( const l1t::HGCalCluster c2d) const

Definition at line 366 of file HGCalShowerShape.cc.

References funct::abs(), l1t::HGCalClusterT< C >::constituents(), meanX(), funct::pow(), alignCSCRings::r, sigmaXX(), and mathSSE::sqrt().

366  {
367  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& cellsPtrs = c2d.constituents();
368 
369  std::vector<std::pair<float, float>> tc_energy_r;
370 
371  for (const auto& id_cell : cellsPtrs) {
372  float r = (id_cell.second->position().z() != 0.
373  ? std::sqrt(pow(id_cell.second->position().x(), 2) + pow(id_cell.second->position().y(), 2)) /
374  std::abs(id_cell.second->position().z())
375  : 0.);
376  tc_energy_r.emplace_back(std::make_pair(id_cell.second->energy(), r));
377  }
378 
379  float r_mean = meanX(tc_energy_r);
380  float Srr = sigmaXX(tc_energy_r, r_mean);
381 
382  return Srr;
383 }
float sigmaXX(const std::vector< pair< float, float >> &energy_X_tc, const float X_cluster) const
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float meanX(const std::vector< pair< float, float >> &energy_X_tc) const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50
template<typename Delta = std::minus<float>>
float HGCalShowerShape::sigmaXX ( const std::vector< pair< float, float >> &  energy_X_tc,
const float  X_cluster 
) const
inlineprivate

Definition at line 55 of file HGCalShowerShape.h.

References delta, funct::pow(), and mathSSE::sqrt().

Referenced by sigmaEtaEtaMax(), sigmaEtaEtaTot(), sigmaRRMax(), sigmaRRMean(), sigmaRRTot(), and sigmaZZ().

55  {
56  Delta delta;
57  float Etot = 0;
58  float deltaX2_sum = 0;
59  for (const auto& energy_X : energy_X_tc) {
60  deltaX2_sum += energy_X.first * pow(delta(energy_X.second, X_cluster), 2);
61  Etot += energy_X.first;
62  }
63  float X_MSE = 0;
64  if (Etot > 0)
65  X_MSE = deltaX2_sum / Etot;
66  float X_RMS = sqrt(X_MSE);
67  return X_RMS;
68  }
dbl * delta
Definition: mlp_gen.cc:36
T sqrt(T t)
Definition: SSEVec.h:18
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
float HGCalShowerShape::sigmaZZ ( const l1t::HGCalMulticluster c3d) const

Definition at line 319 of file HGCalShowerShape.cc.

References l1t::HGCalClusterT< C >::constituents(), meanX(), and sigmaXX().

Referenced by fillShapes(), and showerLength().

319  {
320  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
321 
322  std::vector<std::pair<float, float>> tc_energy_z;
323 
324  for (const auto& id_clu : clustersPtrs) {
325  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
326 
327  for (const auto& id_tc : triggerCells) {
328  tc_energy_z.emplace_back(std::make_pair(id_tc.second->energy(), id_tc.second->position().z()));
329  }
330  }
331 
332  float z_mean = meanX(tc_energy_z);
333  float Szz = sigmaXX(tc_energy_z, z_mean);
334 
335  return Szz;
336 }
float sigmaXX(const std::vector< pair< float, float >> &energy_X_tc, const float X_cluster) const
float meanX(const std::vector< pair< float, float >> &energy_X_tc) const
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50

Member Data Documentation

HGCalTriggerTools HGCalShowerShape::triggerTools_
private