CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 =
51  geometry->getSubdetectorGeometry(DetId::Calo, CaloTowerDetId::SubdetId);
52 
54  c.get<IdealGeometryRecord>().get(topology);
55 
56  // Make the 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
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
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  vector<reco::CandidatePtr>::const_iterator itParticle;
273  for (itParticle=particles.begin();itParticle!=particles.end();++itParticle){
274  if ( itParticle->isNull() || !itParticle->isAvailable() ) {
275  edm::LogWarning("DataNotFound") << " JetSpecific: PF Particle is invalid\n";
276  continue;
277  }
278  const Candidate* pfCand = itParticle->get();
279  if (pfCand) {
280  switch (std::abs(pfCand->pdgId())) {
281  case 211: //PFCandidate::h: // charged hadron
282  chargedHadronEnergy += pfCand->energy();
283  chargedHadronMultiplicity++;
284  chargedMultiplicity++;
285  break;
286 
287  case 130: //PFCandidate::h0 : // neutral hadron
288  neutralHadronEnergy += pfCand->energy();
289  neutralHadronMultiplicity++;
290  neutralMultiplicity++;
291  break;
292 
293  case 22: //PFCandidate::gamma: // photon
294  photonEnergy += pfCand->energy();
295  photonMultiplicity++;
296  neutralEmEnergy += pfCand->energy();
297  neutralMultiplicity++;
298  break;
299 
300  case 11: // PFCandidate::e: // electron
301  electronEnergy += pfCand->energy();
302  electronMultiplicity++;
303  chargedEmEnergy += pfCand->energy();
304  chargedMultiplicity++;
305  break;
306 
307  case 13: //PFCandidate::mu: // muon
308  muonEnergy += pfCand->energy();
309  muonMultiplicity++;
310  chargedMuEnergy += pfCand->energy();
311  chargedMultiplicity++;
312  break;
313 
314  case 1: // PFCandidate::h_HF : // hadron in HF
315  HFHadronEnergy += pfCand->energy();
316  HFHadronMultiplicity++;
317  neutralMultiplicity++;
318  break;
319 
320  case 2: //PFCandidate::egamma_HF : // electromagnetic in HF
321  HFEMEnergy += pfCand->energy();
322  HFEMMultiplicity++;
323  neutralEmEnergy += pfCand->energy();
324  neutralMultiplicity++;
325  break;
326 
327 
328  default:
329  edm::LogWarning("DataNotFound") <<"reco::makePFJetSpecific: Unknown PFCandidate::ParticleType: "
330  <<pfCand->pdgId()<<" is ignored\n";
331  break;
332  }
333  }
334  else {
335  edm::LogWarning("DataNotFound") <<"reco::makePFJetSpecific: Referred constituent is not "
336  <<"a PFCandidate\n";
337  }
338  }
339 
340  pfJetSpecific->mChargedHadronEnergy=chargedHadronEnergy;
341  pfJetSpecific->mNeutralHadronEnergy= neutralHadronEnergy;
342  pfJetSpecific->mPhotonEnergy= photonEnergy;
343  pfJetSpecific->mElectronEnergy= electronEnergy;
344  pfJetSpecific->mMuonEnergy= muonEnergy;
345  pfJetSpecific->mHFHadronEnergy= HFHadronEnergy;
346  pfJetSpecific->mHFEMEnergy= HFEMEnergy;
347 
350  pfJetSpecific->mPhotonMultiplicity= photonMultiplicity;
352  pfJetSpecific->mMuonMultiplicity= muonMultiplicity;
353  pfJetSpecific->mHFHadronMultiplicity= HFHadronMultiplicity;
354  pfJetSpecific->mHFEMMultiplicity= HFEMMultiplicity;
355 
356  pfJetSpecific->mChargedEmEnergy=chargedEmEnergy;
357  pfJetSpecific->mChargedMuEnergy=chargedMuEnergy;
358  pfJetSpecific->mNeutralEmEnergy=neutralEmEnergy;
359  pfJetSpecific->mChargedMultiplicity=chargedMultiplicity;
360  pfJetSpecific->mNeutralMultiplicity=neutralMultiplicity;
361 
362  return true;
363 }
364 
365 
366 //______________________________________________________________________________
367 bool reco::makeSpecific(vector<reco::CandidatePtr> const & mcparticles,
368  GenJet::Specific* genJetSpecific)
369 {
370  if (0==genJetSpecific) return false;
371 
372  vector<reco::CandidatePtr>::const_iterator itMcParticle=mcparticles.begin();
373  for (;itMcParticle!=mcparticles.end();++itMcParticle) {
374  if ( itMcParticle->isNull() || !itMcParticle->isAvailable() ) {
375  edm::LogWarning("DataNotFound") << " JetSpecific: MC Particle is invalid\n";
376  continue;
377  }
378  const Candidate* candidate = itMcParticle->get();
379  if (candidate->hasMasterClone()) candidate = candidate->masterClone().get();
380  //const GenParticle* genParticle = GenJet::genParticle(candidate);
381  if (candidate) {
382  double e = candidate->energy();
383  switch (abs (candidate->pdgId ())) {
384  case 22: // photon
385  case 11: // e
386  genJetSpecific->m_EmEnergy += e;
387  break;
388  case 211: // pi
389  case 321: // K
390  case 130: // KL
391  case 2212: // p
392  case 2112: // n
393  genJetSpecific->m_HadEnergy += e;
394  break;
395  case 13: // muon
396  case 12: // nu_e
397  case 14: // nu_mu
398  case 16: // nu_tau
399 
400  genJetSpecific->m_InvisibleEnergy += e;
401  break;
402  default:
403  genJetSpecific->m_AuxiliaryEnergy += e;
404  }
405  }
406  else {
407  edm::LogWarning("DataNotFound") <<"reco::makeGenJetSpecific: Referred GenParticleCandidate "
408  <<"is not available in the event\n";
409  }
410  }
411 
412  return true;
413 }
414 
415 
416 //______________________________________________________________________________
418 {
419  int eta = std::abs(iEta);
420  if (eta <= topology.lastHBRing()) return HcalBarrel;
421  else if (eta <= topology.lastHERing()) return HcalEndcap;
422  else if (eta <= topology.lastHFRing()) return HcalForward;
423  return HcalEmpty;
424 }
425 
426 
427 
float mMaxEInEmTowers
Maximum energy in EM towers.
Definition: CaloJet.h:53
dictionary specific
virtual void setCharge(Charge q) GCC11_FINAL
set electric charge
virtual double energy() const =0
energy
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:73
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:417
int lastHBRing() const
Definition: HcalTopology.h:80
int mPhotonMultiplicity
Definition: PFJet.h:62
float mNeutralHadronEnergy
Definition: PFJet.h:53
float mChargedMuEnergy
Definition: PFJet.h:71
float mHadEnergyInHB
Hadronic energy in HB.
Definition: CaloJet.h:59
T eta() const
Jets made from CaloTowers.
Definition: BasicJet.h:20
double charge(const std::vector< uint8_t > &Ampls)
Jets made from PFObjects.
Definition: PFJet.h:21
float mChargedHadronEnergy
Definition: PFJet.h:52
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:70
virtual bool hasMasterClone() const =0
Jets made out of PFClusters.
Definition: PFClusterJet.h:26
double emEnergy() const
Definition: CaloTower.h:94
float mElectronEnergy
Definition: PFJet.h:55
static const int SubdetId
double p4[4]
Definition: TauolaWrapper.h:92
int mChargedHadronMultiplicity
Definition: PFJet.h:60
float mPhotonEnergy
Definition: PFJet.h:54
float mNeutralEmEnergy
Definition: PFJet.h:72
int lastHFRing() const
Definition: HcalTopology.h:84
HcalSubdetector
Definition: HcalAssistant.h:31
math::XYZPoint Point
point in the space
Definition: Particle.h:31
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:63
Jets made out of tracks.
Definition: TrackJet.h:27
float phiSpan() const
double hadEnergy() const
Definition: CaloTower.h:95
virtual int pdgId() const =0
PDG identifier.
CaloTowerDetId id() const
Definition: CaloTower.h:87
int mNeutralHadronMultiplicity
Definition: PFJet.h:61
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:74
const T & get() const
Definition: EventSetup.h:55
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
float mHFHadronEnergy
Definition: PFJet.h:57
float m_HadEnergy
Energy of Hadrons.
Definition: GenJet.h:36
T get() const
get a component
Definition: Candidate.h:219
float mTowersArea
Area of contributing CaloTowers.
Definition: CaloJet.h:75
int ieta() const
get the tower ieta
int mHFHadronMultiplicity
Definition: PFJet.h:65
float m_EmEnergy
Energy of EM particles.
Definition: GenJet.h:34
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:27
value_type const * get() const
Definition: RefToBase.h:212
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:96
int lastHERing() const
Definition: HcalTopology.h:82
float mHadEnergyInHE
Hadronic energy in HE.
Definition: CaloJet.h:63
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
virtual const CandidateBaseRef & masterClone() const =0