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