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  float dEta = fabs(geometry->getCorners()[0].eta()-
195  geometry->getCorners()[2].eta());
196  float dPhi = fabs(geometry->getCorners()[0].phi() -
197  geometry->getCorners()[2].phi());
198  jetArea += dEta * dPhi;
199  }
200  else {
201  edm::LogWarning("DataNotFound") <<"reco::makeCaloJetSpecific: Geometry for cell "
202  <<tower->id()<<" can not be found. Ignoring cell\n";
203  }
204  }
205  else {
206  edm::LogWarning("DataNotFound")<<"reco::makeCaloJetSpecific: Constituent is not of "
207  <<"CaloTower type\n";
208  }
209  }
210 
211  double towerEnergy = eInHad + eInEm;
212  caloJetSpecific->mHadEnergyInHO = eInHO;
213  caloJetSpecific->mHadEnergyInHB = eInHB;
214  caloJetSpecific->mHadEnergyInHE = eInHE;
215  caloJetSpecific->mHadEnergyInHF = eHadInHF;
216  caloJetSpecific->mEmEnergyInHF = eEmInHF;
217  caloJetSpecific->mEmEnergyInEB = eInEB;
218  caloJetSpecific->mEmEnergyInEE = eInEE;
219  if (towerEnergy > 0) {
220  caloJetSpecific->mEnergyFractionHadronic = eInHad / towerEnergy;
221  caloJetSpecific->mEnergyFractionEm = eInEm / towerEnergy;
222  }
223  else { // HO only jet
224  caloJetSpecific->mEnergyFractionHadronic = 1.;
225  caloJetSpecific->mEnergyFractionEm = 0.;
226  }
227  caloJetSpecific->mTowersArea = jetArea;
228  caloJetSpecific->mMaxEInEmTowers = 0;
229  caloJetSpecific->mMaxEInHadTowers = 0;
230 
231  //Sort the arrays
232  sort(eECal_i.begin(), eECal_i.end(), greater<double>());
233  sort(eHCal_i.begin(), eHCal_i.end(), greater<double>());
234 
235  if (!towers.empty()) {
236  //Highest value in the array is the first element of the array
237  caloJetSpecific->mMaxEInEmTowers = eECal_i.front();
238  caloJetSpecific->mMaxEInHadTowers = eHCal_i.front();
239  }
240 
241  return true;
242 }
243 
244 
245 //______________________________________________________________________________
246 bool reco::makeSpecific(vector<reco::CandidatePtr> const & particles,
247  PFJet::Specific* pfJetSpecific)
248 {
249  if (0==pfJetSpecific) return false;
250 
251  // 1.- Loop over PFCandidates,
252  // 2.- Get the corresponding PFCandidate
253  // 3.- Calculate the different PFJet specific quantities
254 
255  float chargedHadronEnergy=0.;
256  float neutralHadronEnergy=0.;
257  float photonEnergy=0.;
258  float electronEnergy=0.;
259  float muonEnergy=0.;
260  float HFHadronEnergy=0.;
261  float HFEMEnergy=0.;
262  int chargedHadronMultiplicity=0;
263  int neutralHadronMultiplicity=0;
264  int photonMultiplicity=0;
265  int electronMultiplicity=0;
266  int muonMultiplicity=0;
267  int HFHadronMultiplicity=0;
268  int HFEMMultiplicity=0;
269 
270  float chargedEmEnergy=0.;
271  float neutralEmEnergy=0.;
272  float chargedMuEnergy=0.;
273  int chargedMultiplicity=0;
274  int neutralMultiplicity=0;
275 
276  vector<reco::CandidatePtr>::const_iterator itParticle;
277  for (itParticle=particles.begin();itParticle!=particles.end();++itParticle){
278  if ( itParticle->isNull() || !itParticle->isAvailable() ) {
279  edm::LogWarning("DataNotFound") << " JetSpecific: PF Particle is invalid\n";
280  continue;
281  }
282  const PFCandidate* pfCand = dynamic_cast<const PFCandidate*> (itParticle->get());
283  if (pfCand) {
284  switch (PFCandidate::ParticleType(pfCand->particleId())) {
285  case PFCandidate::h: // charged hadron
286  chargedHadronEnergy += pfCand->energy();
287  chargedHadronMultiplicity++;
288  chargedMultiplicity++;
289  break;
290 
291  case PFCandidate::h0 : // neutral hadron
292  neutralHadronEnergy += pfCand->energy();
293  neutralHadronMultiplicity++;
294  neutralMultiplicity++;
295  break;
296 
297  case PFCandidate::gamma: // photon
298  photonEnergy += pfCand->energy();
299  photonMultiplicity++;
300  neutralEmEnergy += pfCand->energy();
301  neutralMultiplicity++;
302  break;
303 
304  case PFCandidate::e: // electron
305  electronEnergy += pfCand->energy();
306  electronMultiplicity++;
307  chargedEmEnergy += pfCand->energy();
308  chargedMultiplicity++;
309  break;
310 
311  case PFCandidate::mu: // muon
312  muonEnergy += pfCand->energy();
313  muonMultiplicity++;
314  chargedMuEnergy += pfCand->energy();
315  chargedMultiplicity++;
316  break;
317 
318  case PFCandidate::h_HF : // hadron in HF
319  HFHadronEnergy += pfCand->energy();
320  HFHadronMultiplicity++;
321  neutralMultiplicity++;
322  break;
323 
324  case PFCandidate::egamma_HF : // electromagnetic in HF
325  HFEMEnergy += pfCand->energy();
326  HFEMMultiplicity++;
327  neutralEmEnergy += pfCand->energy();
328  neutralMultiplicity++;
329  break;
330 
331 
332  default:
333  edm::LogWarning("DataNotFound") <<"reco::makePFJetSpecific: Unknown PFCandidate::ParticleType: "
334  <<pfCand->particleId()<<" is ignored\n";
335  break;
336  }
337  }
338  else {
339  edm::LogWarning("DataNotFound") <<"reco::makePFJetSpecific: Referred constituent is not "
340  <<"a PFCandidate\n";
341  }
342  }
343 
344  pfJetSpecific->mChargedHadronEnergy=chargedHadronEnergy;
345  pfJetSpecific->mNeutralHadronEnergy= neutralHadronEnergy;
346  pfJetSpecific->mPhotonEnergy= photonEnergy;
347  pfJetSpecific->mElectronEnergy= electronEnergy;
348  pfJetSpecific->mMuonEnergy= muonEnergy;
349  pfJetSpecific->mHFHadronEnergy= HFHadronEnergy;
350  pfJetSpecific->mHFEMEnergy= HFEMEnergy;
351 
352  pfJetSpecific->mChargedHadronMultiplicity=chargedHadronMultiplicity;
353  pfJetSpecific->mNeutralHadronMultiplicity= neutralHadronMultiplicity;
354  pfJetSpecific->mPhotonMultiplicity= photonMultiplicity;
355  pfJetSpecific->mElectronMultiplicity= electronMultiplicity;
356  pfJetSpecific->mMuonMultiplicity= muonMultiplicity;
357  pfJetSpecific->mHFHadronMultiplicity= HFHadronMultiplicity;
358  pfJetSpecific->mHFEMMultiplicity= HFEMMultiplicity;
359 
360  pfJetSpecific->mChargedEmEnergy=chargedEmEnergy;
361  pfJetSpecific->mChargedMuEnergy=chargedMuEnergy;
362  pfJetSpecific->mNeutralEmEnergy=neutralEmEnergy;
363  pfJetSpecific->mChargedMultiplicity=chargedMultiplicity;
364  pfJetSpecific->mNeutralMultiplicity=neutralMultiplicity;
365 
366  return true;
367 }
368 
369 
370 //______________________________________________________________________________
371 bool reco::makeSpecific(vector<reco::CandidatePtr> const & mcparticles,
372  GenJet::Specific* genJetSpecific)
373 {
374  if (0==genJetSpecific) return false;
375 
376  vector<reco::CandidatePtr>::const_iterator itMcParticle=mcparticles.begin();
377  for (;itMcParticle!=mcparticles.end();++itMcParticle) {
378  if ( itMcParticle->isNull() || !itMcParticle->isAvailable() ) {
379  edm::LogWarning("DataNotFound") << " JetSpecific: MC Particle is invalid\n";
380  continue;
381  }
382  const Candidate* candidate = itMcParticle->get();
383  if (candidate->hasMasterClone()) candidate = candidate->masterClone().get();
384  const GenParticle* genParticle = GenJet::genParticle(candidate);
385  if (genParticle) {
386  double e = genParticle->energy();
387  switch (abs (genParticle->pdgId ())) {
388  case 22: // photon
389  case 11: // e
390  genJetSpecific->m_EmEnergy += e;
391  break;
392  case 211: // pi
393  case 321: // K
394  case 130: // KL
395  case 2212: // p
396  case 2112: // n
397  genJetSpecific->m_HadEnergy += e;
398  break;
399  case 13: // muon
400  case 12: // nu_e
401  case 14: // nu_mu
402  case 16: // nu_tau
403 
404  genJetSpecific->m_InvisibleEnergy += e;
405  break;
406  default:
407  genJetSpecific->m_AuxiliaryEnergy += e;
408  }
409  }
410  else {
411  edm::LogWarning("DataNotFound") <<"reco::makeGenJetSpecific: Referred GenParticleCandidate "
412  <<"is not available in the event\n";
413  }
414  }
415 
416  return true;
417 }
418 
419 
420 //______________________________________________________________________________
422 {
423  int eta = std::abs(iEta);
424  if (eta <= topology.lastHBRing()) return HcalBarrel;
425  else if (eta <= topology.lastHERing()) return HcalEndcap;
426  else if (eta <= topology.lastHFRing()) return HcalForward;
427  return HcalEmpty;
428 }
429 
430 
431 
float mMaxEInEmTowers
Maximum energy in EM towers.
Definition: CaloJet.h:54
virtual double energy() const GCC11_FINAL
energy
dictionary specific
virtual void setCharge(Charge q) GCC11_FINAL
set electric charge
ParticleType
particle types
Definition: PFCandidate.h:40
Jets made from CaloTowers.
Definition: CaloJet.h:30
float mEmEnergyInHF
Em energy in HF.
Definition: CaloJet.h:70
float mEnergyFractionHadronic
Hadronic energy fraction.
Definition: CaloJet.h:72
int mChargedMultiplicity
Definition: PFJet.h:74
virtual int pdgId() const GCC11_FINAL
PDG identifier.
float mEmEnergyInEB
Em energy in EB.
Definition: CaloJet.h:66
HcalSubdetector hcalSubdetector(int iEta, const HcalTopology &topology)
converts eta to the corresponding HCAL subdetector.
Definition: JetSpecific.cc:421
#define abs(x)
Definition: mlp_lapack.h:159
int lastHBRing() const
Definition: HcalTopology.h:79
int mPhotonMultiplicity
Definition: PFJet.h:63
float mNeutralHadronEnergy
Definition: PFJet.h:54
float mChargedMuEnergy
Definition: PFJet.h:72
float mHadEnergyInHB
Hadronic energy in HB.
Definition: CaloJet.h:60
T eta() const
Jets made from CaloTowers.
Definition: BasicJet.h:21
double charge(const std::vector< uint8_t > &Ampls)
Jets made from PFObjects.
Definition: PFJet.h:22
float mChargedHadronEnergy
Definition: PFJet.h:53
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:71
virtual bool hasMasterClone() const =0
Jets made out of PFClusters.
Definition: PFClusterJet.h:27
double emEnergy() const
Definition: CaloTower.h:79
float mElectronEnergy
Definition: PFJet.h:56
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:61
float mPhotonEnergy
Definition: PFJet.h:55
float mNeutralEmEnergy
Definition: PFJet.h:73
int lastHFRing() const
Definition: HcalTopology.h:83
HcalSubdetector
Definition: HcalAssistant.h:32
math::XYZPoint Point
point in the space
Definition: Particle.h:29
Jets made from MC generator particles.
Definition: GenJet.h:25
float mEnergyFractionEm
Em energy fraction.
Definition: CaloJet.h:74
float mHadEnergyInHF
Hadronic energy in HF.
Definition: CaloJet.h:62
const int mu
Definition: Constants.h:23
float m_InvisibleEnergy
Invisible energy (mu, nu, ...)
Definition: GenJet.h:39
float mMaxEInHadTowers
Maximum energy in HCAL towers.
Definition: CaloJet.h:56
int mElectronMultiplicity
Definition: PFJet.h:64
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:62
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
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:75
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:58
ESHandle< TrackerGeometry > geometry
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:35
float mHFHadronEnergy
Definition: PFJet.h:58
float m_HadEnergy
Energy of Hadrons.
Definition: GenJet.h:37
T get() const
get a component
Definition: Candidate.h:219
float mTowersArea
Area of contributing CaloTowers.
Definition: CaloJet.h:76
int ieta() const
get the tower ieta
virtual ParticleType particleId() const
Definition: PFCandidate.h:347
int mHFHadronMultiplicity
Definition: PFJet.h:66
float m_EmEnergy
Energy of EM particles.
Definition: GenJet.h:35
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:25
value_type const * get() const
Definition: RefToBase.h:212
float mEmEnergyInEE
Em energy in EE.
Definition: CaloJet.h:68
*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:81
virtual const CornersVec & getCorners() const =0
Returns the corner points of this cell&#39;s volume.
float mHadEnergyInHE
Hadronic energy in HE.
Definition: CaloJet.h:64
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