CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetAnalyzer_HeavyIons_matching.cc
Go to the documentation of this file.
1 //
2 // Jet Analyzer class for heavy ion jets. for DQM jet analysis monitoring
3 // For CMSSW_7_5_X, especially reading background subtracted jets - jet matching 1 - 1
4 // author: Raghav Kunnawalkam Elayavalli,
5 // April 6th 2015
6 // Rutgers University, email: raghav.k.e at CERN dot CH
7 //
8 
9 
11 
12 using namespace edm;
13 using namespace reco;
14 using namespace std;
15 
16 // declare the constructors:
17 
18 
20  mInputJet1Collection (iConfig.getParameter<edm::InputTag> ("src_Jet1")),
21  mInputJet2Collection (iConfig.getParameter<edm::InputTag> ("src_Jet2")),
22  JetType1 (iConfig.getUntrackedParameter<std::string>("Jet1")),
23  JetType2 (iConfig.getUntrackedParameter<std::string>("Jet2")),
24  mRecoJetPtThreshold (iConfig.getParameter<double> ("recoJetPtThreshold")),
25  mRecoDelRMatch (iConfig.getParameter<double> ("recoDelRMatch")),
26  mRecoJetEtaCut (iConfig.getParameter<double> ("recoJetEtaCut"))
27 {
28  std::string inputCollectionLabelJet1(mInputJet1Collection.label());
29  std::string inputCollectionLabelJet2(mInputJet2Collection.label());
30 
31  //consumes
32 
33  if(std::string("VsCalo") == JetType1) caloJet1Token_ = consumes<reco::CaloJetCollection>(mInputJet1Collection);
34  if(std::string("VsPF") == JetType1) pfJetsToken_ = consumes<reco::PFJetCollection>(mInputJet1Collection);
35  if(std::string("PuCalo") == JetType1) caloJet2Token_ = consumes<reco::CaloJetCollection>(mInputJet1Collection);
36  if(std::string("PuPF") == JetType1) basicJetsToken_ = consumes<reco::BasicJetCollection>(mInputJet1Collection);
37 
38  if(std::string("VsCalo") == JetType2) caloJet1Token_ = consumes<reco::CaloJetCollection>(mInputJet2Collection);
39  if(std::string("VsPF") == JetType2) pfJetsToken_ = consumes<reco::PFJetCollection>(mInputJet2Collection);
40  if(std::string("PuCalo") == JetType2) caloJet2Token_ = consumes<reco::CaloJetCollection>(mInputJet2Collection);
41  if(std::string("PuPF") == JetType2) basicJetsToken_ = consumes<reco::BasicJetCollection>(mInputJet2Collection);
42 
43  // initialize the Jet matching histograms
44 
46  mpT_Jet1_matched = 0;
47  mpT_Jet2_matched = 0;
50 
51  // we need to add histograms which will hold the hadronic and electromagnetic energy content for the unmatched Jets.
52  if(std::string("VsCalo") == JetType1 || std::string("PuCalo") == JetType1){
55  }
56  if(std::string("VsCalo") == JetType2 || std::string("PuCalo") == JetType2) {
59  }
60 
61  if(std::string("VsPF") == JetType1){
72  }
73 
74  if(std::string("VsPF") == JetType2){
85  }
86 
87 
88 }
89 
91  {
92 
93  ibooker.setCurrentFolder("JetMET/HIJetValidation/"+mInputJet1Collection.label()+"_DeltaRMatched_"+mInputJet2Collection.label());
94 
95  mpT_ratio_Jet1Jet2 = ibooker.book1D("Ratio_Jet1pT_vs_Jet2pT",Form(";Matched %s Jet pT / %s Jet pT;Counts", "JetType1", "JetType2"),100, 0, 10);
96  mpT_Jet1_matched = ibooker.book1D("Jet1_matched_jet_Spectra",Form(";Matched %s Spectra; counts", "JetType1"),100, 0, 1000);
97  mpT_Jet2_matched = ibooker.book1D("Jet2_matched_jet_Spectra",Form(";Matched %s Spectra; counts", "JetType2"),100, 0, 1000);
98  mpT_Jet1_unmatched = ibooker.book1D("Jet1_unmatched_jet_Spectra",Form(";Unmatched %s spectra;counts","JetType1"),100, 0, 1000);
99  mpT_Jet2_unmatched = ibooker.book1D("Jet2_unmatched_jet_Spectra",Form(";Unmatched %s spectra;counts","JetType2"),100, 0, 1000);
100 
101  if(std::string("VsCalo") == JetType1 || std::string("PuCalo") == JetType1){
102  mHadEnergy_Jet1_unmatched = ibooker.book1D("HadEnergy_Jet1_unmatched", Form("HadEnergy_Jet1_unmatched;HadEnergy unmatched %s;counts", "JetType1"), 50, 0, 200);
103  mEmEnergy_Jet1_unmatched = ibooker.book1D("EmEnergy_Jet1_unmatched", Form("EmEnergy_Jet1_unmatched;EMEnergy unmatched %s;counts", "JetType1"), 50, 0, 200);
104  }
105 
106  if(std::string("VsCalo") == JetType2 || std::string("PuCalo") == JetType2){
107  mHadEnergy_Jet2_unmatched = ibooker.book1D("HadEnergy_Jet2_unmatched", Form("HadEnergy_Jet2_unmatched;HadEnergy unmatched %s;counts", "JetType2"), 50, 0, 200);
108  mEmEnergy_Jet2_unmatched = ibooker.book1D("EmEnergy_Jet2_unmatched", Form("EmEnergy_Jet2_unmatched;EMEnergy unmatched %s;counts", "JetType2"), 50, 0, 200);
109  }
110 
111  if(std::string("VsPF") == JetType1){
112  mChargedHadronEnergy_Jet1_unmatched = ibooker.book1D("ChargedHadronEnergy_Jet1_unmatched", Form(";charged HAD energy unmatched %s;counts","JetType1"), 100, 0, 300);
113  mNeutralHadronEnergy_Jet1_unmatched = ibooker.book1D("neutralHadronEnergy_Jet1_unmatched", Form(";neutral HAD energy unmatched %s;counts", "JetType1"), 100, 0, 300);
114  mChargedEmEnergy_Jet1_unmatched = ibooker.book1D("ChargedEmEnergy_Jet1_unmatched", Form(";charged EM energy unmatched %s;counts", "JetType1"), 100, 0, 300);
115  mNeutralEmEnergy_Jet1_unmatched = ibooker.book1D("neutralEmEnergy_Jet1_unmatched", Form(";neutral EM energy unmatched %s;counts", "JetType1"), 100, 0, 300);
116  mChargedMuEnergy_Jet1_unmatched = ibooker.book1D("ChargedMuEnergy_Jet1_unmatched", Form(";charged Mu energy unmatched %s;counts", "JetType1"), 100, 0, 300);
117 
118  mChargedHadEnergyFraction_Jet1_unmatched = ibooker.book1D("ChargedHadEnergyFraction_Jet1_unmatched",Form(";h^{+/-} Energy Fraction %s;counts", "JetType1"),50, 0, 1);
119  mNeutralHadEnergyFraction_Jet1_unmatched = ibooker.book1D("NeutralHadEnergyFraction_Jet1_unmatched",Form(";h^{0} Energy Fraction %s;counts", "JetType1"),50, 0, 1);
120  mPhotonEnergyFraction_Jet1_unmatched = ibooker.book1D("PhotonEnergyFraction_Jet1_unmatched",Form(";#gamma Energy Fraction %s;counts", "JetType1"),50, 0, 1);
121  mElectronEnergyFraction_Jet1_unmatched = ibooker.book1D("ElectronEnergyFraction_Jet1_unmatched",Form(";e Energy Fraction %s;counts", "JetType1"),50, 0, 1);
122  mMuonEnergyFraction_Jet1_unmatched = ibooker.book1D("MuonoEnergyFraction_Jet1_unmatched",Form(";#mu Energy Fraction %s;counts", "JetType1"),50, 0, 1);
123 
124  }
125 
126  if(std::string("VsPF") == JetType2){
127  mChargedHadronEnergy_Jet2_unmatched = ibooker.book1D("ChargedHadronEnergy_Jet2_unmatched", Form(";charged HAD energy unmatched %s;counts","JetType2"), 100, 0, 300);
128  mNeutralHadronEnergy_Jet2_unmatched = ibooker.book1D("neutralHadronEnergy_Jet2_unmatched", Form(";neutral HAD energy unmatched %s;counts", "JetType2"), 100, 0, 300);
129  mChargedEmEnergy_Jet2_unmatched = ibooker.book1D("ChargedEmEnergy_Jet2_unmatched", Form(";charged EM energy unmatched %s;counts", "JetType2"), 100, 0, 300);
130  mNeutralEmEnergy_Jet2_unmatched = ibooker.book1D("neutralEmEnergy_Jet2_unmatched", Form(";neutral EM energy unmatched %s;counts", "JetType2"), 100, 0, 300);
131  mChargedMuEnergy_Jet2_unmatched = ibooker.book1D("ChargedMuEnergy_Jet2_unmatched", Form(";charged Mu energy unmatched %s;counts", "JetType2"), 100, 0, 300);
132 
133  mChargedHadEnergyFraction_Jet2_unmatched = ibooker.book1D("ChargedHadEnergyFraction_Jet2_unmatched",Form(";h^{+/-} Energy Fraction %s;counts", "JetType2"),50, 0, 1);
134  mNeutralHadEnergyFraction_Jet2_unmatched = ibooker.book1D("NeutralHadEnergyFraction_Jet2_unmatched",Form(";h^{0} Energy Fraction %s;counts", "JetType2"),50, 0, 1);
135  mPhotonEnergyFraction_Jet2_unmatched = ibooker.book1D("PhotonEnergyFraction_Jet2_unmatched",Form(";#gamma Energy Fraction %s;counts", "JetType2"),50, 0, 1);
136  mElectronEnergyFraction_Jet2_unmatched = ibooker.book1D("ElectronEnergyFraction_Jet2_unmatched",Form(";e Energy Fraction %s;counts", "JetType2"),50, 0, 1);
137  mMuonEnergyFraction_Jet2_unmatched = ibooker.book1D("MuonoEnergyFraction_Jet2_unmatched",Form(";#mu Energy Fraction %s;counts", "JetType2"),50, 0, 1);
138 
139  }
140 
141  if (mOutputFile.empty ())
142  LogInfo("OutputInfo") << " Histograms will NOT be saved";
143  else
144  LogInfo("OutputInfo") << " Histograms will be saved to file:" << mOutputFile;
145 
146 
147  }
148 
149 
150 
151 //------------------------------------------------------------------------------
152 // ~JetAnalyzer_HeavyIons
153 //------------------------------------------------------------------------------
155 
156 
157 //------------------------------------------------------------------------------
158 // beginJob
159 //------------------------------------------------------------------------------
160 //void JetAnalyzer_HeavyIons_matching::beginJob() {
161 // std::cout<<"inside the begin job function"<<endl;
162 //}
163 
164 
165 //------------------------------------------------------------------------------
166 // endJob
167 //------------------------------------------------------------------------------
168 //void JetAnalyzer_HeavyIons_matching::endJob()
169 //{
170 // if (!mOutputFile.empty() && &*edm::Service<DQMStore>())
171 // {
172 // edm::Service<DQMStore>()->save(mOutputFile);
173 // }
174 //}
175 
176 
177 //------------------------------------------------------------------------------
178 // analyze
179 //------------------------------------------------------------------------------
181 {
182 
183  // Get the Jet collection
184  //----------------------------------------------------------------------------
185 
186  std::vector<const Jet*> recoJet1;
187  recoJet1.clear();
188  std::vector<const Jet*> recoJet2;
189  recoJet2.clear();
190 
196 
197  if(std::string("VsCalo") == JetType1) {
198  mEvent.getByToken(caloJet1Token_, caloJet1);
199  for (unsigned ijet=0; ijet<caloJet1->size(); ++ijet) recoJet1.push_back(&(*caloJet1)[ijet]);
200  }
201  if(std::string("PuCalo") == JetType1) {
202  mEvent.getByToken(caloJet2Token_, caloJet1);
203  for (unsigned ijet=0; ijet<caloJet1->size(); ++ijet) recoJet1.push_back(&(*caloJet1)[ijet]);
204  }
205  if(std::string("VsPF") == JetType1) {
206  mEvent.getByToken(pfJetsToken_, pfJets);
207  for (unsigned ijet=0; ijet<pfJets->size(); ++ijet) recoJet1.push_back(&(*pfJets)[ijet]);
208  }
209  if(std::string("PuPF") == JetType1) {
210  mEvent.getByToken(basicJetsToken_, basicJets);
211  for (unsigned ijet=0; ijet<basicJets->size(); ++ijet) recoJet1.push_back(&(*basicJets)[ijet]);
212  }
213 
214  if(std::string("VsCalo") == JetType2) {
215  mEvent.getByToken(caloJet1Token_, caloJet2);
216  for (unsigned ijet=0; ijet<caloJet2->size(); ++ijet) recoJet2.push_back(&(*caloJet2)[ijet]);
217  }
218  if(std::string("PuCalo") == JetType2) {
219  mEvent.getByToken(caloJet2Token_, caloJet2);
220  for (unsigned ijet=0; ijet<caloJet2->size(); ++ijet) recoJet2.push_back(&(*caloJet2)[ijet]);
221  }
222  if(std::string("VsPF") == JetType2) {
223  mEvent.getByToken(pfJetsToken_, pfJets);
224  for (unsigned ijet=0; ijet<pfJets->size(); ++ijet) recoJet2.push_back(&(*pfJets)[ijet]);
225  }
226  if(std::string("PuPF") == JetType2) {
227  mEvent.getByToken(basicJetsToken_, basicJets);
228  for (unsigned ijet=0; ijet<basicJets->size(); ++ijet) recoJet2.push_back(&(*basicJets)[ijet]);
229  }
230 
231  // start to perform the matching - between recoJet1 and recoJet2.
232 
233  Int_t Jet1_nref = recoJet1.size();
234  Int_t Jet2_nref = recoJet2.size();
235 
236  int jet1 = 0;
237  int jet2 = 0;
238 
239  std::vector <MyJet> vJet1, vJet2;
240  std::vector <int> Jet1_ID(Jet1_nref), Jet2_ID(Jet2_nref);
241 
242  if(Jet1_nref == 0 || Jet2_nref == 0) return;
243 
244  for(unsigned ijet1 = 0; ijet1 < recoJet1.size(); ++ijet1){
245 
246 
247  if(recoJet1[ijet1]->pt() < mRecoJetPtThreshold) continue;
248  if(fabs(recoJet1[ijet1]->eta()) < mRecoJetEtaCut) continue;
249 
250  MyJet JET1;
251  JET1.eta = recoJet1[ijet1]->eta();
252  JET1.phi = recoJet1[ijet1]->phi();
253  JET1.pt = recoJet1[ijet1]->pt();
254  JET1.id = ijet1;
255 
256  vJet1.push_back(JET1);
257  jet1++;
258 
259  }// first jet loop
260 
261  for(unsigned ijet2 = 0; ijet2 < recoJet2.size(); ++ijet2){
262 
263  if(recoJet2[ijet2]->pt() < mRecoJetPtThreshold) continue;
264  if(fabs(recoJet2[ijet2]->eta()) < mRecoJetEtaCut) continue;
265 
266  MyJet JET2;
267  JET2.eta = recoJet2[ijet2]->eta();
268  JET2.phi = recoJet2[ijet2]->phi();
269  JET2.pt = recoJet2[ijet2]->pt();
270  JET2.id = ijet2;
271 
272  vJet2.push_back(JET2);
273  jet2++;
274 
275  }// second jet loop
276 
277  bool onlyJet1 = (jet1>0 && jet2==0) ? true : false;
278  bool onlyJet2 = (jet1==0 && jet2 >0) ? true : false;
279  bool bothJet1Jet2 = (jet1>0 && jet2 >0) ? true : false;
280 
281  int matchedJets = 0;
282  int unmatchedJet1 = 0;
283  int unmatchedJet2 = 0;
284 
285  std::vector < MyJet >::const_iterator iJet, jJet;
286 
287  if(onlyJet1) {
288 
289  for(iJet = vJet1.begin(); iJet != vJet1.end(); ++iJet){
290 
291  int pj = (*iJet).id;
292 
293  mpT_Jet1_unmatched->Fill(recoJet1[pj]->pt());
294 
295  }
296 
297  }else if(onlyJet2) {
298 
299  for(iJet = vJet2.begin(); iJet != vJet2.end(); ++iJet){
300 
301  int cj = (*iJet).id;
302 
303  mpT_Jet2_unmatched->Fill(recoJet2[cj]->pt());
304  }
305 
306  }else if (bothJet1Jet2){
307 
308  ABMatchedJets mABMatchedJets;
309 
310  for(iJet = vJet1.begin(); iJet != vJet1.end(); ++iJet){
311  for(jJet = vJet2.begin(); jJet != vJet2.end(); ++jJet){
312  mABMatchedJets.insert(std::make_pair(*iJet, *jJet));
313  }
314  }
315 
316 
317  ABItr itr;
318  // matched Jets matching Jet 1 to Jet 2
319  for(itr = mABMatchedJets.begin(); itr != mABMatchedJets.end(); ++itr){
320 
321  ABJetPair jetpair = (*itr);
322  MyJet Aj = jetpair.first;
323  MyJet Bj = jetpair.second;
324 
325  float delr = JetAnalyzer_HeavyIons_matching::deltaRR(Bj.eta, Bj.phi, Aj.eta, Aj.phi);
326 
327  if( delr < mRecoDelRMatch && Jet1_ID[Aj.id] == 0){
328 
329  mpT_ratio_Jet1Jet2->Fill((Float_t) recoJet2[Bj.id]->pt()/recoJet1[Aj.id]->pt());
330 
331  mpT_Jet1_matched->Fill(recoJet1[Aj.id]->pt());
332  mpT_Jet2_matched->Fill(recoJet2[Bj.id]->pt());
333 
334  Jet1_ID[Aj.id] = 1;
335  Jet2_ID[Bj.id] = 1;
336 
337  matchedJets++;
338 
339 
340  }
341 
342  }
343 
344  // for unmatched Jets
345  for(itr = mABMatchedJets.begin(); itr != mABMatchedJets.end(); ++itr){
346 
347  ABJetPair jetpair = (*itr);
348 
349  MyJet Aj = jetpair.first;
350  MyJet Bj = jetpair.second;
351 
352  if(Jet1_ID[Aj.id] == 0) {
353 
354  mpT_Jet1_unmatched->Fill(recoJet1[Aj.id]->pt());
355  unmatchedJet1++;
356  Jet1_ID[Aj.id] = 1;
357 
358  if(std::string("VsCalo") == JetType1 || std::string("PuCalo") == JetType1){
359  mHadEnergy_Jet1_unmatched->Fill((*caloJet1)[Aj.id].hadEnergyInHO()
360  + (*caloJet1)[Aj.id].hadEnergyInHB()
361  + (*caloJet1)[Aj.id].hadEnergyInHF()
362  + (*caloJet1)[Aj.id].hadEnergyInHE()
363  );
364  mEmEnergy_Jet1_unmatched->Fill((*caloJet1)[Aj.id].emEnergyInEB()
365  + (*caloJet1)[Aj.id].emEnergyInEE()
366  + (*caloJet1)[Aj.id].emEnergyInHF()
367  );
368  }
369 
370  if(std::string("VsPF") == JetType1){
371  mChargedHadronEnergy_Jet1_unmatched->Fill((*pfJets)[Aj.id].chargedHadronEnergy());
372  mNeutralHadronEnergy_Jet1_unmatched->Fill((*pfJets)[Aj.id].neutralHadronEnergy());
373  mChargedEmEnergy_Jet1_unmatched->Fill((*pfJets)[Aj.id].chargedEmEnergy());
374  mNeutralEmEnergy_Jet1_unmatched->Fill((*pfJets)[Aj.id].neutralEmEnergy());
375  mChargedMuEnergy_Jet1_unmatched->Fill((*pfJets)[Aj.id].chargedMuEnergy());
376 
377  mChargedHadEnergyFraction_Jet1_unmatched->Fill((*pfJets)[Aj.id].chargedHadronEnergyFraction());
378  mNeutralHadEnergyFraction_Jet1_unmatched->Fill((*pfJets)[Aj.id].neutralHadronEnergyFraction());
379  mPhotonEnergyFraction_Jet1_unmatched->Fill((*pfJets)[Aj.id].photonEnergyFraction());
380  mElectronEnergyFraction_Jet1_unmatched->Fill((*pfJets)[Aj.id].electronEnergyFraction());
381  mMuonEnergyFraction_Jet1_unmatched->Fill((*pfJets)[Aj.id].muonEnergyFraction());
382  }
383 
384  }
385 
386  if(Jet2_ID[Bj.id] == 0) {
387 
388  mpT_Jet2_unmatched->Fill(recoJet2[Bj.id]->pt());
389  unmatchedJet2++;
390  Jet2_ID[Bj.id] = 2;
391  if(std::string("VsCalo") == JetType2 || std::string("PuCalo") == JetType2){
392  mHadEnergy_Jet2_unmatched->Fill((*caloJet2)[Bj.id].hadEnergyInHO()
393  + (*caloJet2)[Bj.id].hadEnergyInHB()
394  + (*caloJet2)[Bj.id].hadEnergyInHF()
395  + (*caloJet2)[Bj.id].hadEnergyInHE()
396  );
397  mEmEnergy_Jet2_unmatched->Fill((*caloJet2)[Bj.id].emEnergyInEB()
398  + (*caloJet2)[Bj.id].emEnergyInEE()
399  + (*caloJet2)[Bj.id].emEnergyInHF()
400  );
401  }
402 
403 
404  if(std::string("VsPF") == JetType2){
405  mChargedHadronEnergy_Jet2_unmatched->Fill((*pfJets)[Bj.id].chargedHadronEnergy());
406  mNeutralHadronEnergy_Jet2_unmatched->Fill((*pfJets)[Bj.id].neutralHadronEnergy());
407  mChargedEmEnergy_Jet2_unmatched->Fill((*pfJets)[Bj.id].chargedEmEnergy());
408  mNeutralEmEnergy_Jet2_unmatched->Fill((*pfJets)[Bj.id].neutralEmEnergy());
409  mChargedMuEnergy_Jet2_unmatched->Fill((*pfJets)[Bj.id].chargedMuEnergy());
410 
411  mChargedHadEnergyFraction_Jet2_unmatched->Fill((*pfJets)[Bj.id].chargedHadronEnergyFraction());
412  mNeutralHadEnergyFraction_Jet2_unmatched->Fill((*pfJets)[Bj.id].neutralHadronEnergyFraction());
413  mPhotonEnergyFraction_Jet2_unmatched->Fill((*pfJets)[Bj.id].photonEnergyFraction());
414  mElectronEnergyFraction_Jet2_unmatched->Fill((*pfJets)[Bj.id].electronEnergyFraction());
415  mMuonEnergyFraction_Jet2_unmatched->Fill((*pfJets)[Bj.id].muonEnergyFraction());
416  }
417 
418  }
419 
420  }
421 
422  }// both Jet1 and Jet2 in the event.
423 
424 
425 }
426 
427 
428 
edm::EDGetTokenT< reco::PFJetCollection > pfJetsToken_
std::multiset< ABJetPair, CompareMatchedJets > ABMatchedJets
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
void analyze(const edm::Event &, const edm::EventSetup &) override
void Fill(long long x)
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:276
edm::EDGetTokenT< reco::CaloJetCollection > caloJet1Token_
static float deltaRR(float eta1, float phi1, float eta2, float phi2)
std::string const & label() const
Definition: InputTag.h:36
JetAnalyzer_HeavyIons_matching(const edm::ParameterSet &)
std::multiset< ABJetPair >::iterator ABItr
edm::EDGetTokenT< reco::CaloJetCollection > caloJet2Token_
tuple pfJets
Definition: pfJets_cff.py:8
Definition: Run.h:43
edm::EDGetTokenT< reco::BasicJetCollection > basicJetsToken_