CMS 3D CMS Logo

ClusterTools.cc
Go to the documentation of this file.
2 
8 
10 
14 
15 #include "vdt/vdtMath.h"
16 
17 #include <iostream>
18 
19 using namespace hgcal;
21 }
22 
24  edm::ConsumesCollector& sumes) :
25  eetok( sumes.consumes<HGCRecHitCollection>(conf.getParameter<edm::InputTag>("HGCEEInput")) ),
26  fhtok( sumes.consumes<HGCRecHitCollection>(conf.getParameter<edm::InputTag>("HGCFHInput")) ),
27  bhtok( sumes.consumes<HGCRecHitCollection>(conf.getParameter<edm::InputTag>("HGCBHInput")) ) {
28  }
29 
31  rhtools_.getEvent(ev);
33  ev.getByToken(eetok, temp);
34  eerh_ = temp.product();
35  ev.getByToken(fhtok, temp);
36  fhrh_ = temp.product();
37  ev.getByToken(bhtok, temp);
38  bhrh_ = temp.product();
39 }
40 
43 }
44 
46  float energy=0.f, energyHad=0.f;
47  const auto& hits = clus.hitsAndFractions();
48  for( const auto& hit : hits ) {
49  const auto& id = hit.first;
50  const float fraction = hit.second;
51  if( id.det() == DetId::HGCalEE ) {
52  energy += eerh_->find(id)->energy()*fraction;
53  } else if( id.det() == DetId::HGCalHSi ) {
54  const float temp = fhrh_->find(id)->energy();
55  energy += temp*fraction;
56  energyHad += temp*fraction;
57  } else if( id.det() == DetId::HGCalHSc ) {
58  const float temp = bhrh_->find(id)->energy();
59  energy += temp*fraction;
60  energyHad += temp*fraction;
61  } else if( id.det() == DetId::Forward ) {
62  switch( id.subdetId() ) {
63  case HGCEE:
64  energy += eerh_->find(id)->energy()*fraction;
65  break;
66  case HGCHEF:
67  {
68  const float temp = fhrh_->find(id)->energy();
69  energy += temp*fraction;
70  energyHad += temp*fraction;
71  }
72  break;
73  default:
74  throw cms::Exception("HGCalClusterTools")
75  << " Cluster contains hits that are not from HGCal! " << std::endl;
76  }
77  } else if ( id.det() == DetId::Hcal && id.subdetId() == HcalEndcap ) {
78  const float temp = bhrh_->find(id)->energy();
79  energy += temp*fraction;
80  energyHad += temp*fraction;
81  } else {
82  throw cms::Exception("HGCalClusterTools")
83  << " Cluster contains hits that are not from HGCal! " << std::endl;
84  }
85  }
86  float fraction = -1.f;
87  if( energy > 0.f ) {
88  fraction = energyHad/energy;
89  }
90  return fraction;
91 }
92 
93 
95  if( clu.clusters().empty() ) return math::XYZPoint();
96 
97  double acc_x = 0.0;
98  double acc_y = 0.0;
99  double acc_z = 0.0;
100  double totweight = 0.;
101  double mcenergy = getMultiClusterEnergy(clu);
102  for( const auto& ptr : clu.clusters() ) {
103  if (mcenergy != 0) {
104  if(ptr->energy()<.01*mcenergy) continue; //cutoff < 1% layer contribution
105  }
106  const double weight = ptr->energy(); // weigth each corrdinate only by the total energy of the layer cluster
107  acc_x += ptr->x()*weight;
108  acc_y += ptr->y()*weight;
109  acc_z += ptr->z()*weight;
110  totweight += weight;
111  }
112  if (totweight != 0) {
113  acc_x /= totweight;
114  acc_y /= totweight;
115  acc_z /= totweight;
116  }
117  // return x/y/z in absolute coordinates
118  return math::XYZPoint(acc_x,acc_y,acc_z);
119 }
120 
121 int ClusterTools::getLayer(const DetId detid) const {
122  return rhtools_.getLayerWithOffset(detid);
123 }
124 
126  double acc = 0.0;
127  for(const auto& ptr : clu.clusters() ) {
128  acc += ptr->energy();
129  }
130  return acc;
131 }
132 
133 bool ClusterTools::getWidths(const reco::CaloCluster & clus,double & sigmaetaeta, double & sigmaphiphi, double & sigmaetaetal, double & sigmaphiphil ) const {
134 
135  if (getLayer(clus.hitsAndFractions()[0].first) > lastLayerEE) return false;
136  const math::XYZPoint & position(clus.position());
137  unsigned nhit=clus.hitsAndFractions().size();
138 
139  sigmaetaeta=0.;
140  sigmaphiphi=0.;
141  sigmaetaetal=0.;
142  sigmaphiphil=0.;
143 
144  double sumw=0.;
145  double sumlogw=0.;
146 
147  for (unsigned int ih=0; ih<nhit; ++ih) {
148  const DetId & id = (clus.hitsAndFractions())[ih].first;
149  if ((clus.hitsAndFractions())[ih].second==0.) continue;
150 
151  if ((id.det() == DetId::HGCalEE) ||
152  (id.det()==DetId::Forward && id.subdetId()==HGCEE)) {
153  const HGCRecHit * theHit = &(*eerh_->find(id));
154 
155  GlobalPoint cellPos = rhtools_.getPosition(id);
156  double weight = theHit->energy();
157  // take w0=2 To be optimized
158  double logweight = 0;
159  if (clus.energy() != 0) {
160  logweight = std::max(0.,2 + log(theHit->energy()/clus.energy()));
161  }
162  double deltaetaeta2 = (cellPos.eta()-position.eta())*(cellPos.eta()-position.eta());
163  double deltaphiphi2 = (cellPos.phi()-position.phi())*(cellPos.phi()-position.phi());
164  sigmaetaeta += deltaetaeta2* weight;
165  sigmaphiphi += deltaphiphi2 * weight;
166  sigmaetaetal += deltaetaeta2* logweight;
167  sigmaphiphil += deltaphiphi2 * logweight;
168  sumw += weight;
169  sumlogw += logweight;
170  }
171  }
172 
173 
174  if (sumw<=0.) return false;
175 
176  sigmaetaeta /= sumw;
177  sigmaetaeta = std::sqrt(sigmaetaeta);
178  sigmaphiphi /= sumw;
179  sigmaphiphi = std::sqrt(sigmaphiphi);
180 
181  if (sumlogw != 0) {
182  sigmaetaetal /= sumlogw;
183  sigmaetaetal = std::sqrt(sigmaetaetal);
184  sigmaphiphil /= sumlogw;
185  sigmaphiphil = std::sqrt(sigmaphiphil);
186  }
187 
188  return true;
189 }
constexpr float energy() const
Definition: CaloRecHit.h:31
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:131
void getEvent(const edm::Event &)
Definition: RecHitTools.cc:70
void getEvent(const edm::Event &)
Definition: ClusterTools.cc:30
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
Definition: weight.py:1
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:197
bool ev
math::XYZPoint getMultiClusterPosition(const reco::HGCalMultiCluster &) const
Definition: ClusterTools.cc:94
void getEventSetup(const edm::EventSetup &)
Definition: RecHitTools.cc:73
U second(std::pair< T, U > const &p)
RecHitTools rhtools_
Definition: ClusterTools.h:56
const edm::PtrVector< reco::BasicCluster > & clusters() const
T sqrt(T t)
Definition: SSEVec.h:18
double f[11][100]
double energy() const
cluster energy
Definition: CaloCluster.h:126
int getLayer(const DetId) const
const HGCRecHitCollection * eerh_
Definition: ClusterTools.h:58
Definition: DetId.h:18
double getMultiClusterEnergy(const reco::HGCalMultiCluster &) const
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:283
T const * product() const
Definition: Handle.h:74
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:112
float getClusterHadronFraction(const reco::CaloCluster &) const
Definition: ClusterTools.cc:45
const edm::EDGetTokenT< HGCRecHitCollection > fhtok
Definition: ClusterTools.h:57
T eta() const
Definition: PV3DBase.h:76
const edm::EDGetTokenT< HGCRecHitCollection > eetok
Definition: ClusterTools.h:57
iterator find(key_type k)
const HGCRecHitCollection * bhrh_
Definition: ClusterTools.h:58
HLT enums.
static const int lastLayerEE
Definition: ClusterTools.h:59
static int position[264][3]
Definition: ReadPGInfo.cc:509
const edm::EDGetTokenT< HGCRecHitCollection > bhtok
Definition: ClusterTools.h:57
bool getWidths(const reco::CaloCluster &clus, double &sigmaetaeta, double &sigmaphiphi, double &sigmaetaetalog, double &sigmaphiphilog) const
const HGCRecHitCollection * fhrh_
Definition: ClusterTools.h:58
void getEventSetup(const edm::EventSetup &)
Definition: ClusterTools.cc:41