CMS 3D CMS Logo

JetSpecific.cc
Go to the documentation of this file.
1 //
3 // JetSpecific
4 // -----------
5 //
7 
8 
10 
11 
25 
26 #include <cmath>
27 
28 
29 using namespace std;
30 
31 
33 // implementation of global functions
35 
36 //______________________________________________________________________________
37 // Overloaded methods to write out specific types
38 
39 
40 // CaloJet
44  std::vector<reco::CandidatePtr> const & constituents,
45  edm::EventSetup const & c )
46 {
47  // Get geometry
49  c.get<CaloGeometryRecord>().get(geometry);
50  const CaloSubdetectorGeometry* towerGeometry =
52 
54  c.get<HcalRecNumberingRecord>().get(topology);
55 
56  // Make the specific
57  reco::CaloJet::Specific specific;
58  makeSpecific (constituents, *towerGeometry, &specific, *topology);
59  // Set the calo jet
60  jet = reco::CaloJet( p4, point, specific, constituents);
61 }
62 
63 
64 
65 // BasicJet
69  std::vector<reco::CandidatePtr> const & constituents,
70  edm::EventSetup const & c )
71 {
72  jet = reco::BasicJet( p4, point, constituents);
73 }
74 
75 // GenJet
79  std::vector<reco::CandidatePtr> const & constituents,
80  edm::EventSetup const & c )
81 {
82 
83  // Make the specific
84  reco::GenJet::Specific specific;
85  makeSpecific (constituents, &specific);
86  // Set to the jet
87  jet = reco::GenJet( p4, point, specific, constituents);
88 }
89 
90 // PFJet
94  std::vector<reco::CandidatePtr> const & constituents,
95  edm::EventSetup const & c )
96 {
97  // Make the specific
98  reco::PFJet::Specific specific;
99  makeSpecific (constituents, &specific);
100  // now make jet charge
101  int charge = 0.;
102  for ( std::vector<reco::CandidatePtr>::const_iterator ic = constituents.begin(),
103  icend = constituents.end();
104  ic != icend; ++ic ) {
105  charge += (*ic)->charge();
106  }
107  jet = reco::PFJet( p4, point, specific, constituents);
108  jet.setCharge( charge );
109 }
110 
111 
112 // TrackJet
115  reco::Particle::Point const & point,
116  std::vector<reco::CandidatePtr> const & constituents,
117  edm::EventSetup const & c )
118 {
119  jet = reco::TrackJet(p4, point, constituents);
120 }
121 
122 // PFClusterJet
125  reco::Particle::Point const & point,
126  std::vector<reco::CandidatePtr> const & constituents,
127  edm::EventSetup const & c )
128 {
129  jet = reco::PFClusterJet( p4, point, constituents);
130 }
131 
132 
133 
134 //______________________________________________________________________________
135 bool reco::makeSpecific(vector<reco::CandidatePtr> const & towers,
136  const CaloSubdetectorGeometry& towerGeometry,
137  CaloJet::Specific* caloJetSpecific,
138  const HcalTopology &topology)
139 {
140  if (0==caloJetSpecific) return false;
141 
142  // 1.- Loop over the tower Ids,
143  // 2.- Get the corresponding CaloTower
144  // 3.- Calculate the different CaloJet specific quantities
145  vector<double> eECal_i;
146  vector<double> eHCal_i;
147  double eInHad = 0.;
148  double eInEm = 0.;
149  double eInHO = 0.;
150  double eInHB = 0.;
151  double eInHE = 0.;
152  double eHadInHF = 0.;
153  double eEmInHF = 0.;
154  double eInEB = 0.;
155  double eInEE = 0.;
156  double jetArea = 0.;
157 
158  vector<reco::CandidatePtr>::const_iterator itTower;
159  for (itTower=towers.begin();itTower!=towers.end();++itTower) {
160  if ( itTower->isNull() || !itTower->isAvailable() ) {
161  edm::LogWarning("DataNotFound") << " JetSpecific: Tower is invalid\n";
162  continue;
163  }
164  const CaloTower* tower = dynamic_cast<const CaloTower*>(itTower->get());
165  if (tower) {
166  //Array of energy in EM Towers:
167  eECal_i.push_back(tower->emEnergy());
168  eInEm += tower->emEnergy();
169  //Array of energy in HCAL Towers:
170  eHCal_i.push_back(tower->hadEnergy());
171  eInHad += tower->hadEnergy();
172 
173  // figure out contributions
174  switch (reco::hcalSubdetector(tower->id().ieta(),topology)) {
175  case HcalBarrel:
176  eInHB += tower->hadEnergy();
177  eInHO += tower->outerEnergy();
178  eInEB += tower->emEnergy();
179  break;
180  case HcalEndcap:
181  eInHE += tower->hadEnergy();
182  eInEE += tower->emEnergy();
183  break;
184  case HcalForward:
185  eHadInHF += tower->hadEnergy();
186  eEmInHF += tower->emEnergy();
187  break;
188  default:
189  break;
190  }
191  // get area of the tower (++ minus --)
192  const CaloCellGeometry* geometry = towerGeometry.getGeometry(tower->id());
193  if (geometry) {
194  jetArea += geometry->etaSpan() * geometry->phiSpan();
195  }
196  else {
197  edm::LogWarning("DataNotFound") <<"reco::makeCaloJetSpecific: Geometry for cell "
198  <<tower->id()<<" can not be found. Ignoring cell\n";
199  }
200  }
201  else {
202  edm::LogWarning("DataNotFound")<<"reco::makeCaloJetSpecific: Constituent is not of "
203  <<"CaloTower type\n";
204  }
205  }
206 
207  double towerEnergy = eInHad + eInEm;
208  caloJetSpecific->mHadEnergyInHO = eInHO;
209  caloJetSpecific->mHadEnergyInHB = eInHB;
210  caloJetSpecific->mHadEnergyInHE = eInHE;
211  caloJetSpecific->mHadEnergyInHF = eHadInHF;
212  caloJetSpecific->mEmEnergyInHF = eEmInHF;
213  caloJetSpecific->mEmEnergyInEB = eInEB;
214  caloJetSpecific->mEmEnergyInEE = eInEE;
215  if (towerEnergy > 0) {
216  caloJetSpecific->mEnergyFractionHadronic = eInHad / towerEnergy;
217  caloJetSpecific->mEnergyFractionEm = eInEm / towerEnergy;
218  }
219  else { // HO only jet
220  caloJetSpecific->mEnergyFractionHadronic = 1.;
221  caloJetSpecific->mEnergyFractionEm = 0.;
222  }
223  caloJetSpecific->mTowersArea = jetArea;
224  caloJetSpecific->mMaxEInEmTowers = 0;
225  caloJetSpecific->mMaxEInHadTowers = 0;
226 
227  //Sort the arrays
228  sort(eECal_i.begin(), eECal_i.end(), greater<double>());
229  sort(eHCal_i.begin(), eHCal_i.end(), greater<double>());
230 
231  if (!towers.empty()) {
232  //Highest value in the array is the first element of the array
233  caloJetSpecific->mMaxEInEmTowers = eECal_i.front();
234  caloJetSpecific->mMaxEInHadTowers = eHCal_i.front();
235  }
236 
237  return true;
238 }
239 
240 
241 //______________________________________________________________________________
242 bool reco::makeSpecific(vector<reco::CandidatePtr> const & particles,
243  PFJet::Specific* pfJetSpecific)
244 {
245  if (0==pfJetSpecific) return false;
246 
247  // 1.- Loop over PFCandidates,
248  // 2.- Get the corresponding PFCandidate
249  // 3.- Calculate the different PFJet specific quantities
250 
251  float chargedHadronEnergy=0.;
252  float neutralHadronEnergy=0.;
253  float photonEnergy=0.;
254  float electronEnergy=0.;
255  float muonEnergy=0.;
256  float HFHadronEnergy=0.;
257  float HFEMEnergy=0.;
260  int photonMultiplicity=0;
261  int electronMultiplicity=0;
262  int muonMultiplicity=0;
263  int HFHadronMultiplicity=0;
264  int HFEMMultiplicity=0;
265 
266  float chargedEmEnergy=0.;
267  float neutralEmEnergy=0.;
268  float chargedMuEnergy=0.;
269  int chargedMultiplicity=0;
270  int neutralMultiplicity=0;
271 
272  float HOEnergy=0.;
273 
274  vector<reco::CandidatePtr>::const_iterator itParticle;
275  for (itParticle=particles.begin();itParticle!=particles.end();++itParticle){
276  if ( itParticle->isNull() || !itParticle->isAvailable() ) {
277  edm::LogWarning("DataNotFound") << " JetSpecific: PF Particle is invalid\n";
278  continue;
279  }
280  const Candidate* pfCand = itParticle->get();
281  if (pfCand) {
282 
283  const PFCandidate* pfCandCast = dynamic_cast<const PFCandidate*>(pfCand);
284  if (pfCandCast)
285  HOEnergy += pfCandCast->hoEnergy();
286 
287  switch (std::abs(pfCand->pdgId())) {
288  case 211: //PFCandidate::h: // charged hadron
289  chargedHadronEnergy += pfCand->energy();
290  chargedHadronMultiplicity++;
291  chargedMultiplicity++;
292  break;
293 
294  case 130: //PFCandidate::h0 : // neutral hadron
295  neutralHadronEnergy += pfCand->energy();
296  neutralHadronMultiplicity++;
297  neutralMultiplicity++;
298  break;
299 
300  case 22: //PFCandidate::gamma: // photon
301  photonEnergy += pfCand->energy();
302  photonMultiplicity++;
303  neutralEmEnergy += pfCand->energy();
304  neutralMultiplicity++;
305  break;
306 
307  case 11: // PFCandidate::e: // electron
308  electronEnergy += pfCand->energy();
309  electronMultiplicity++;
310  chargedEmEnergy += pfCand->energy();
311  chargedMultiplicity++;
312  break;
313 
314  case 13: //PFCandidate::mu: // muon
315  muonEnergy += pfCand->energy();
316  muonMultiplicity++;
317  chargedMuEnergy += pfCand->energy();
318  chargedMultiplicity++;
319  break;
320 
321  case 1: // PFCandidate::h_HF : // hadron in HF
322  HFHadronEnergy += pfCand->energy();
323  HFHadronMultiplicity++;
324  neutralHadronEnergy += pfCand->energy();
325  neutralMultiplicity++;
326  break;
327 
328  case 2: //PFCandidate::egamma_HF : // electromagnetic in HF
329  HFEMEnergy += pfCand->energy();
330  HFEMMultiplicity++;
331  neutralEmEnergy += pfCand->energy();
332  neutralMultiplicity++;
333  break;
334 
335 
336  default:
337  edm::LogWarning("DataNotFound") <<"reco::makePFJetSpecific: Unknown PFCandidate::ParticleType: "
338  <<pfCand->pdgId()<<" is ignored\n";
339  break;
340  }
341  }
342  else {
343  edm::LogWarning("DataNotFound") <<"reco::makePFJetSpecific: Referred constituent is not "
344  <<"a PFCandidate\n";
345  }
346  }
347 
348  pfJetSpecific->mChargedHadronEnergy=chargedHadronEnergy;
349  pfJetSpecific->mNeutralHadronEnergy= neutralHadronEnergy;
350  pfJetSpecific->mPhotonEnergy= photonEnergy;
351  pfJetSpecific->mElectronEnergy= electronEnergy;
352  pfJetSpecific->mMuonEnergy= muonEnergy;
353  pfJetSpecific->mHFHadronEnergy= HFHadronEnergy;
354  pfJetSpecific->mHFEMEnergy= HFEMEnergy;
355 
358  pfJetSpecific->mPhotonMultiplicity= photonMultiplicity;
360  pfJetSpecific->mMuonMultiplicity= muonMultiplicity;
361  pfJetSpecific->mHFHadronMultiplicity= HFHadronMultiplicity;
362  pfJetSpecific->mHFEMMultiplicity= HFEMMultiplicity;
363 
364  pfJetSpecific->mChargedEmEnergy=chargedEmEnergy;
365  pfJetSpecific->mChargedMuEnergy=chargedMuEnergy;
366  pfJetSpecific->mNeutralEmEnergy=neutralEmEnergy;
367  pfJetSpecific->mChargedMultiplicity=chargedMultiplicity;
368  pfJetSpecific->mNeutralMultiplicity=neutralMultiplicity;
369 
370  pfJetSpecific->mHOEnergy= HOEnergy;
371 
372  return true;
373 }
374 
375 
376 //______________________________________________________________________________
377 bool reco::makeSpecific(vector<reco::CandidatePtr> const & mcparticles,
378  GenJet::Specific* genJetSpecific)
379 {
380  if (0==genJetSpecific) return false;
381 
382  vector<reco::CandidatePtr>::const_iterator itMcParticle=mcparticles.begin();
383  for (;itMcParticle!=mcparticles.end();++itMcParticle) {
384  if ( itMcParticle->isNull() || !itMcParticle->isAvailable() ) {
385  edm::LogWarning("DataNotFound") << " JetSpecific: MC Particle is invalid\n";
386  continue;
387  }
388  const Candidate* candidate = itMcParticle->get();
389  if (candidate->hasMasterClone()) candidate = candidate->masterClone().get();
390  //const GenParticle* genParticle = GenJet::genParticle(candidate);
391  if (candidate) {
392  double e = candidate->energy();
393  switch (abs (candidate->pdgId ())) {
394  case 22: // photon
395  case 11: // e
396  genJetSpecific->m_EmEnergy += e;
397  break;
398  case 211: // pi
399  case 321: // K
400  case 130: // KL
401  case 2212: // p
402  case 2112: // n
403  genJetSpecific->m_HadEnergy += e;
404  break;
405  case 13: // muon
406  case 12: // nu_e
407  case 14: // nu_mu
408  case 16: // nu_tau
409 
410  genJetSpecific->m_InvisibleEnergy += e;
411  break;
412  default:
413  genJetSpecific->m_AuxiliaryEnergy += e;
414  }
415  }
416  else {
417  edm::LogWarning("DataNotFound") <<"reco::makeGenJetSpecific: Referred GenParticleCandidate "
418  <<"is not available in the event\n";
419  }
420  }
421 
422  return true;
423 }
424 
425 
426 //______________________________________________________________________________
428 {
429  int eta = std::abs(iEta);
430  if (eta <= topology.lastHBRing()) return HcalBarrel;
431  else if (eta <= topology.lastHERing()) return HcalEndcap;
432  else if (eta <= topology.lastHFRing()) return HcalForward;
433  return HcalEmpty;
434 }
435 
436 
437 
float mMaxEInEmTowers
Maximum energy in EM towers.
Definition: CaloJet.h:53
value_type const * get() const
Definition: RefToBase.h:234
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:45
Jets made from CaloTowers.
Definition: CaloJet.h:29
float etaSpan() const
float mEmEnergyInHF
Em energy in HF.
Definition: CaloJet.h:69
CaloTopology const * topology(0)
float mEnergyFractionHadronic
Hadronic energy fraction.
Definition: CaloJet.h:71
int mChargedMultiplicity
Definition: PFJet.h:75
float mEmEnergyInEB
Em energy in EB.
Definition: CaloJet.h:65
HcalSubdetector hcalSubdetector(int iEta, const HcalTopology &topology)
converts eta to the corresponding HCAL subdetector.
Definition: JetSpecific.cc:427
int lastHBRing() const
Definition: HcalTopology.h:88
int mPhotonMultiplicity
Definition: PFJet.h:64
float mNeutralHadronEnergy
Definition: PFJet.h:55
float mChargedMuEnergy
Definition: PFJet.h:73
float mHadEnergyInHB
Hadronic energy in HB.
Definition: CaloJet.h:59
Jets made from CaloTowers.
Definition: BasicJet.h:20
Jets made from PFObjects.
Definition: PFJet.h:21
float mChargedHadronEnergy
Definition: PFJet.h:54
virtual const CaloCellGeometry * getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
float mChargedEmEnergy
Definition: PFJet.h:72
virtual double energy() const =0
energy
Jets made out of PFClusters.
Definition: PFClusterJet.h:26
double emEnergy() const
Definition: CaloTower.h:110
float mElectronEnergy
Definition: PFJet.h:57
static const int SubdetId
virtual int pdgId() const =0
PDG identifier.
double p4[4]
Definition: TauolaWrapper.h:92
int mChargedHadronMultiplicity
Definition: PFJet.h:62
float mPhotonEnergy
Definition: PFJet.h:56
float mNeutralEmEnergy
Definition: PFJet.h:74
virtual void setCharge(Charge q) final
set electric charge
Definition: LeafCandidate.h:93
int lastHFRing() const
Definition: HcalTopology.h:92
HcalSubdetector
Definition: HcalAssistant.h:31
math::XYZPoint Point
point in the space
Definition: Particle.h:25
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Jets made from MC generator particles.
Definition: GenJet.h:24
float mEnergyFractionEm
Em energy fraction.
Definition: CaloJet.h:73
float mHadEnergyInHF
Hadronic energy in HF.
Definition: CaloJet.h:61
float m_InvisibleEnergy
Invisible energy (mu, nu, ...)
Definition: GenJet.h:38
float mMaxEInHadTowers
Maximum energy in HCAL towers.
Definition: CaloJet.h:55
int mElectronMultiplicity
Definition: PFJet.h:65
Jets made out of tracks.
Definition: TrackJet.h:27
float phiSpan() const
double hadEnergy() const
Definition: CaloTower.h:111
double hoEnergy() const
return corrected Hcal energy
Definition: PFCandidate.h:241
CaloTowerDetId id() const
Definition: CaloTower.h:103
virtual const CandidateBaseRef & masterClone() const =0
int mNeutralHadronMultiplicity
Definition: PFJet.h:63
bool makeSpecific(std::vector< reco::CandidatePtr > const &towers, const CaloSubdetectorGeometry &towerGeometry, reco::CaloJet::Specific *caloJetSpecific, const HcalTopology &topology)
Make CaloJet specifics. Assumes PseudoJet is made from CaloTowerCandidates.
Definition: JetSpecific.cc:135
int mNeutralMultiplicity
Definition: PFJet.h:76
const T & get() const
Definition: EventSetup.h:56
float m_AuxiliaryEnergy
Anything else (undecayed Sigmas etc.)
Definition: GenJet.h:40
float mHadEnergyInHO
Hadronic nergy fraction in HO.
Definition: CaloJet.h:57
ESHandle< TrackerGeometry > geometry
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
float mHFHadronEnergy
Definition: PFJet.h:59
float m_HadEnergy
Energy of Hadrons.
Definition: GenJet.h:36
T get() const
get a component
Definition: Candidate.h:217
float mTowersArea
Area of contributing CaloTowers.
Definition: CaloJet.h:75
int ieta() const
get the tower ieta
int mHFHadronMultiplicity
Definition: PFJet.h:67
float m_EmEnergy
Energy of EM particles.
Definition: GenJet.h:34
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
float mEmEnergyInEE
Em energy in EE.
Definition: CaloJet.h:67
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
double outerEnergy() const
Definition: CaloTower.h:112
int lastHERing() const
Definition: HcalTopology.h:90
float mHadEnergyInHE
Hadronic energy in HE.
Definition: CaloJet.h:63
virtual bool hasMasterClone() const =0
void writeSpecific(reco::CaloJet &jet, reco::Particle::LorentzVector const &p4, reco::Particle::Point const &point, std::vector< reco::CandidatePtr > const &constituents, edm::EventSetup const &c)
Definition: JetSpecific.cc:41