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)
 
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 15 of file HGCalShowerShape.h.

Constructor & Destructor Documentation

HGCalShowerShape::HGCalShowerShape ( )
inline

Definition at line 17 of file HGCalShowerShape.h.

17 {}
HGCalShowerShape::~HGCalShowerShape ( )
inline

Definition at line 19 of file HGCalShowerShape.h.

19 {}

Member Function Documentation

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

Definition at line 85 of file HGCalShowerShape.cc.

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

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

86 {
87  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
89  std::vector<bool> layers(nlayers);
90  for(const auto& id_cluster : clustersPtrs)
91  {
92  int layer = triggerGeometry.triggerLayer(id_cluster.second->detId());
93  layers[layer-1] = true;
94  }
95  int length = 0;
96  int maxlength = 0;
97  for(bool layer : layers)
98  {
99  if(layer) length++;
100  else length = 0;
101  if(length>maxlength) maxlength = length;
102  }
103  return maxlength;
104 }
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:52
float HGCalShowerShape::eMax ( const l1t::HGCalMulticluster c3d) const

Definition at line 371 of file HGCalShowerShape.cc.

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

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

371  {
372 
373  std::unordered_map<int, float> layer_energy;
374 
375  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
376 
377  for(const auto& id_clu : clustersPtrs){
378 
379  unsigned layer = triggerTools_.layerWithOffset(id_clu.second->detId());
380  layer_energy[layer] += id_clu.second->energy();
381 
382  }
383 
384  float EMax=0;
385 
386  for(const auto& layer : layer_energy){
387 
388  if(layer.second>EMax) EMax = layer.second;
389 
390  }
391 
392  return EMax;
393 
394 }
HGCalTriggerTools triggerTools_
unsigned layerWithOffset(const DetId &) const
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:52
void HGCalShowerShape::eventSetup ( const edm::EventSetup es)
inline

Definition at line 21 of file HGCalShowerShape.h.

References HGCalTriggerTools::eventSetup(), firstLayer(), lastLayer(), maxLayer(), and triggerTools_.

Referenced by HGCalMulticlusteringImpl::eventSetup().

void eventSetup(const edm::EventSetup &)
HGCalTriggerTools triggerTools_
int HGCalShowerShape::firstLayer ( const l1t::HGCalMulticluster c3d) const

Definition at line 31 of file HGCalShowerShape.cc.

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

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

31  {
32 
33  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
34 
35  int firstLayer=999;
36 
37  for(const auto& id_clu : clustersPtrs){
38 
39  int layer = triggerTools_.layerWithOffset(id_clu.second->detId());
40  if(layer<firstLayer) firstLayer=layer;
41 
42  }
43 
44  return firstLayer;
45 
46 }
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:52
int HGCalShowerShape::lastLayer ( const l1t::HGCalMulticluster c3d) const

Definition at line 68 of file HGCalShowerShape.cc.

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

Referenced by eventSetup(), and showerLength().

68  {
69 
70  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
71 
72  int lastLayer=-999;
73 
74  for(const auto& id_clu : clustersPtrs){
75 
76  int layer = triggerTools_.layerWithOffset(id_clu.second->detId());
77  if(layer>lastLayer) lastLayer=layer;
78 
79  }
80 
81  return lastLayer;
82 
83 }
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:52
int HGCalShowerShape::maxLayer ( const l1t::HGCalMulticluster c3d) const

Definition at line 49 of file HGCalShowerShape.cc.

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

Referenced by eventSetup(), and HGCalMulticlusteringImpl::finalizeClusters().

49  {
50 
51  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
52  std::unordered_map<int, float> layers_pt;
53  float max_pt = 0.;
54  int max_layer = 0;
55  for(const auto& id_cluster : clustersPtrs){
56  unsigned layer = triggerTools_.layerWithOffset(id_cluster.second->detId());
57  auto itr_insert = layers_pt.emplace(layer, 0.);
58  itr_insert.first->second += id_cluster.second->pt();
59  if(itr_insert.first->second>max_pt){
60  max_pt = itr_insert.first->second;
61  max_layer = layer;
62  }
63  }
64  return max_layer;
65 }
HGCalTriggerTools triggerTools_
unsigned layerWithOffset(const DetId &) const
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:52
float HGCalShowerShape::meanX ( const std::vector< pair< float, float > > &  energy_X_tc) const
private

Definition at line 12 of file HGCalShowerShape.cc.

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

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

Definition at line 26 of file HGCalShowerShape.h.

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

Referenced by HGCalMulticlusteringImpl::finalizeClusters().

26 {return lastLayer(c3d)-firstLayer(c3d)+1; }//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 192 of file HGCalShowerShape.cc.

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

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

192  {
193 
194  std::unordered_map<int, std::vector<std::pair<float,float> > > tc_layer_energy_eta;
195  std::unordered_map<int, LorentzVector> layer_LV;
196 
197  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
198 
199  for(const auto& id_clu : clustersPtrs){
200 
201  unsigned layer = triggerTools_.layerWithOffset(id_clu.second->detId());
202 
203  layer_LV[layer] += id_clu.second->p4();
204 
205  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
206 
207  for(const auto& id_tc : triggerCells){
208 
209  tc_layer_energy_eta[layer].emplace_back( std::make_pair(id_tc.second->energy(),id_tc.second->eta()) );
210 
211  }
212 
213  }
214 
215 
216  float SigmaEtaEtaMax=0;
217 
218  for(auto& tc_iter : tc_layer_energy_eta){
219 
220  const std::vector<std::pair<float, float> >& energy_eta_layer = tc_iter.second;
221  const LorentzVector& LV_layer = layer_LV[tc_iter.first];
222  float SigmaEtaEtaLayer = sigmaXX(energy_eta_layer,LV_layer.eta()); //RMS wrt layer eta, not wrt c3d eta
223  if(SigmaEtaEtaLayer > SigmaEtaEtaMax) SigmaEtaEtaMax = SigmaEtaEtaLayer;
224 
225  }
226 
227 
228  return SigmaEtaEtaMax;
229 
230 
231 }
HGCalTriggerTools triggerTools_
unsigned layerWithOffset(const DetId &) const
float sigmaXX(const std::vector< pair< float, float > > &energy_X_tc, const float X_cluster) const
math::XYZTLorentzVector LorentzVector
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:52
float HGCalShowerShape::sigmaEtaEtaTot ( const l1t::HGCalMulticluster c3d) const

Definition at line 107 of file HGCalShowerShape.cc.

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

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

107  {
108 
109  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
110 
111  std::vector<std::pair<float,float> > tc_energy_eta ;
112 
113  for(const auto& id_clu : clustersPtrs){
114 
115  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
116 
117  for(const auto& id_tc : triggerCells){
118 
119  tc_energy_eta.emplace_back( std::make_pair(id_tc.second->energy(),id_tc.second->eta()) );
120 
121  }
122 
123  }
124 
125  float SeeTot = sigmaXX(tc_energy_eta,c3d.eta());
126 
127  return SeeTot;
128 
129 }
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:52
float HGCalShowerShape::sigmaEtaEtaTot ( const l1t::HGCalCluster c2d) const

Definition at line 429 of file HGCalShowerShape.cc.

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

429  {
430 
431  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& cellsPtrs = c2d.constituents();
432 
433  std::vector<std::pair<float,float> > tc_energy_eta ;
434 
435  for(const auto& id_cell : cellsPtrs){
436 
437  tc_energy_eta.emplace_back( std::make_pair(id_cell.second->energy(),id_cell.second->eta()) );
438 
439  }
440 
441  float See = sigmaXX(tc_energy_eta,c2d.eta());
442 
443  return See;
444 
445 }
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:52
float HGCalShowerShape::sigmaPhiPhi ( const std::vector< pair< float, float > > &  energy_phi_tc,
const float  phi_cluster 
) const
inlineprivate

Definition at line 72 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 236 of file HGCalShowerShape.cc.

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

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

236  {
237 
238  std::unordered_map<int, std::vector<std::pair<float,float> > > tc_layer_energy_phi;
239  std::unordered_map<int, LorentzVector> layer_LV;
240 
241  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
242 
243  for(const auto& id_clu : clustersPtrs){
244 
245  unsigned layer = triggerTools_.layerWithOffset(id_clu.second->detId());
246 
247  layer_LV[layer] += id_clu.second->p4();
248 
249  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
250 
251  for(const auto& id_tc : triggerCells){
252 
253  tc_layer_energy_phi[layer].emplace_back( std::make_pair(id_tc.second->energy(),id_tc.second->phi()) );
254 
255  }
256 
257  }
258 
259 
260  float SigmaPhiPhiMax=0;
261 
262  for(auto& tc_iter : tc_layer_energy_phi){
263 
264  const std::vector<std::pair<float, float> >& energy_phi_layer = tc_iter.second;
265  const LorentzVector& LV_layer = layer_LV[tc_iter.first];
266  float SigmaPhiPhiLayer = sigmaPhiPhi(energy_phi_layer,LV_layer.phi()); //RMS wrt layer phi, not wrt c3d phi
267  if(SigmaPhiPhiLayer > SigmaPhiPhiMax) SigmaPhiPhiMax = SigmaPhiPhiLayer;
268 
269  }
270 
271 
272  return SigmaPhiPhiMax;
273 
274 
275 }
HGCalTriggerTools triggerTools_
unsigned layerWithOffset(const DetId &) const
float sigmaPhiPhi(const std::vector< pair< float, float > > &energy_phi_tc, const float phi_cluster) const
math::XYZTLorentzVector LorentzVector
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:52
float HGCalShowerShape::sigmaPhiPhiTot ( const l1t::HGCalMulticluster c3d) const

Definition at line 134 of file HGCalShowerShape.cc.

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

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

134  {
135 
136  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
137 
138  std::vector<std::pair<float,float> > tc_energy_phi ;
139 
140  for(const auto& id_clu : clustersPtrs){
141 
142  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
143 
144  for(const auto& id_tc : triggerCells){
145 
146  tc_energy_phi.emplace_back( std::make_pair(id_tc.second->energy(),id_tc.second->phi()) );
147 
148  }
149 
150  }
151 
152  float SppTot = sigmaPhiPhi(tc_energy_phi,c3d.phi());
153 
154  return SppTot;
155 
156 }
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:52
float HGCalShowerShape::sigmaPhiPhiTot ( const l1t::HGCalCluster c2d) const

Definition at line 449 of file HGCalShowerShape.cc.

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

449  {
450 
451  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& cellsPtrs = c2d.constituents();
452 
453  std::vector<std::pair<float,float> > tc_energy_phi ;
454 
455  for(const auto& id_cell : cellsPtrs){
456 
457  tc_energy_phi.emplace_back( std::make_pair(id_cell.second->energy(),id_cell.second->phi()) );
458 
459  }
460 
461  float Spp = sigmaPhiPhi(tc_energy_phi,c2d.phi());
462 
463  return Spp;
464 
465 }
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:52
float HGCalShowerShape::sigmaRRMax ( const l1t::HGCalMulticluster c3d) const

Definition at line 280 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 HGCalMulticlusteringImpl::finalizeClusters(), and showerLength().

280  {
281 
282  std::unordered_map<int, std::vector<std::pair<float,float> > > tc_layer_energy_r;
283 
284  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
285 
286  for(const auto& id_clu : clustersPtrs){
287 
288  unsigned layer = triggerTools_.layerWithOffset(id_clu.second->detId());
289 
290  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
291 
292  for(const auto& id_tc : triggerCells){
293 
294  float r = (id_tc.second->position().z()!=0. ? std::sqrt( pow(id_tc.second->position().x(),2) + pow(id_tc.second->position().y(),2) )/std::abs(id_tc.second->position().z()) : 0.);
295  tc_layer_energy_r[layer].emplace_back( std::make_pair(id_tc.second->energy(),r) );
296 
297  }
298 
299  }
300 
301 
302  float SigmaRRMax=0;
303 
304  for(auto& tc_iter : tc_layer_energy_r){
305 
306  const std::vector<std::pair<float, float> >& energy_r_layer = tc_iter.second;
307  float r_mean_layer = meanX(energy_r_layer);
308  float SigmaRRLayer = sigmaXX(energy_r_layer,r_mean_layer);
309  if(SigmaRRLayer > SigmaRRMax) SigmaRRMax = SigmaRRLayer;
310 
311  }
312 
313 
314  return SigmaRRMax;
315 
316 
317 }
HGCalTriggerTools triggerTools_
unsigned layerWithOffset(const DetId &) const
float meanX(const std::vector< pair< float, float > > &energy_X_tc) const
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
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:52
float HGCalShowerShape::sigmaRRMean ( const l1t::HGCalMulticluster c3d,
float  radius = 5. 
) const

Definition at line 320 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 HGCalMulticlusteringImpl::finalizeClusters(), and showerLength().

320  {
321 
322  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
323  // group trigger cells by layer
324  std::unordered_map<int, std::vector<edm::Ptr<l1t::HGCalTriggerCell>> > layers_tcs;
325  for(const auto& id_clu : clustersPtrs){
326  unsigned layer = triggerTools_.layerWithOffset(id_clu.second->detId());
327  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
328  for(const auto& id_tc : triggerCells){
329  layers_tcs[layer].emplace_back(id_tc.second);
330  }
331  }
332 
333  // Select trigger cells within X cm of the max TC in the layer
334  std::unordered_map<int, std::vector<std::pair<float,float> > > tc_layers_energy_r;
335  for(const auto& layer_tcs : layers_tcs){
336  int layer = layer_tcs.first;
337  edm::Ptr<l1t::HGCalTriggerCell> max_tc = layer_tcs.second.front();
338  for(const auto& tc : layer_tcs.second){
339  if(tc->energy()>max_tc->energy()) max_tc = tc;
340  }
341  for(const auto& tc : layer_tcs.second){
342  double dx = tc->position().x() - max_tc->position().x();
343  double dy = tc->position().y() - max_tc->position().y();
344  double distance_to_max = std::sqrt(dx*dx+dy*dy);
345  if(distance_to_max<radius){
346  float r = (tc->position().z()!=0. ? std::sqrt(tc->position().x()*tc->position().x() + tc->position().y()*tc->position().y())/std::abs(tc->position().z()) : 0.);
347  tc_layers_energy_r[layer].emplace_back( std::make_pair(tc->energy(), r));
348  }
349  }
350  }
351 
352  // Compute srr layer by layer
353  std::vector<std::pair<float,float>> layers_energy_srr2;
354  for(const auto& layer_energy_r : tc_layers_energy_r){
355  const auto& energies_r = layer_energy_r.second;
356  float r_mean_layer = meanX(energies_r);
357  float srr = sigmaXX(energies_r,r_mean_layer);
358  double energy_sum = 0.;
359  for(const auto& energy_r : energies_r){
360  energy_sum += energy_r.first;
361  }
362  layers_energy_srr2.emplace_back(std::make_pair(energy_sum, srr*srr));
363  }
364  // Combine all layer srr
365  float srr2_mean = meanX(layers_energy_srr2);
366  return std::sqrt(srr2_mean);
367 }
T y() const
Definition: PV3DBase.h:63
HGCalTriggerTools triggerTools_
unsigned layerWithOffset(const DetId &) const
float meanX(const std::vector< pair< float, float > > &energy_X_tc) const
float sigmaXX(const std::vector< pair< float, float > > &energy_X_tc, const float X_cluster) 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
T x() const
Definition: PV3DBase.h:62
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:52
float HGCalShowerShape::sigmaRRTot ( const l1t::HGCalMulticluster c3d) const

Definition at line 162 of file HGCalShowerShape.cc.

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

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

162  {
163 
164  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
165 
166  std::vector<std::pair<float,float> > tc_energy_r ;
167 
168  for(const auto& id_clu : clustersPtrs){
169 
170  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
171 
172  for(const auto& id_tc : triggerCells){
173 
174  float r = (id_tc.second->position().z()!=0. ? std::sqrt( pow(id_tc.second->position().x(),2) + pow(id_tc.second->position().y(),2) )/std::abs(id_tc.second->position().z()) : 0.);
175  tc_energy_r.emplace_back( std::make_pair(id_tc.second->energy(),r) );
176 
177  }
178 
179  }
180 
181  float r_mean = meanX(tc_energy_r);
182  float Szz = sigmaXX(tc_energy_r,r_mean);
183 
184  return Szz;
185 
186 }
float meanX(const std::vector< pair< float, float > > &energy_X_tc) const
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
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:52
float HGCalShowerShape::sigmaRRTot ( const l1t::HGCalCluster c2d) const

Definition at line 470 of file HGCalShowerShape.cc.

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

470  {
471 
472  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& cellsPtrs = c2d.constituents();
473 
474  std::vector<std::pair<float,float> > tc_energy_r ;
475 
476  for(const auto& id_cell : cellsPtrs){
477 
478  float r = (id_cell.second->position().z()!=0. ? std::sqrt( pow(id_cell.second->position().x(),2) + pow(id_cell.second->position().y(),2) )/std::abs(id_cell.second->position().z()) : 0.);
479  tc_energy_r.emplace_back( std::make_pair(id_cell.second->energy(),r) );
480 
481  }
482 
483  float r_mean = meanX(tc_energy_r);
484  float Srr = sigmaXX(tc_energy_r,r_mean);
485 
486  return Srr;
487 
488 }
float meanX(const std::vector< pair< float, float > > &energy_X_tc) const
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
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:52
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 52 of file HGCalShowerShape.h.

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

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

54  {
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) X_MSE=deltaX2_sum/Etot;
65  float X_RMS=sqrt(X_MSE);
66  return X_RMS;
67  }
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 400 of file HGCalShowerShape.cc.

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

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

400  {
401 
402  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clustersPtrs = c3d.constituents();
403 
404  std::vector<std::pair<float,float> > tc_energy_z ;
405 
406  for(const auto& id_clu : clustersPtrs){
407 
408  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& triggerCells = id_clu.second->constituents();
409 
410  for(const auto& id_tc : triggerCells){
411 
412  tc_energy_z.emplace_back( std::make_pair(id_tc.second->energy(),id_tc.second->position().z()) );
413 
414  }
415 
416  }
417 
418  float z_mean = meanX(tc_energy_z);
419  float Szz = sigmaXX(tc_energy_z,z_mean);
420 
421  return Szz;
422 
423 }
float meanX(const std::vector< pair< float, float > > &energy_X_tc) const
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:52

Member Data Documentation

HGCalTriggerTools HGCalShowerShape::triggerTools_
private