ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/CmsHi/JetAnalysis/src/HiInclusiveJetAnalyzer.cc
Revision: 1.8
Committed: Tue Mar 6 15:03:10 2012 UTC (13 years, 1 month ago) by mnguyen
Content type: text/plain
Branch: MAIN
Changes since 1.7: +212 -14 lines
Log Message:
Added option to do btagging

File Contents

# User Rev Content
1 yilmaz 1.1 /*
2     Based on the jet response analyzer
3     Modified by Matt Nguyen, November 2010
4    
5     */
6    
7     #include "CmsHi/JetAnalysis/interface/HiInclusiveJetAnalyzer.h"
8    
9     #include "FWCore/Framework/interface/EventSetup.h"
10     #include "FWCore/Framework/interface/ESHandle.h"
11     #include <Math/DistFunc.h>
12     #include "TMath.h"
13    
14     #include "FWCore/Framework/interface/ESHandle.h"
15    
16     #include "FWCore/MessageLogger/interface/MessageLogger.h"
17     #include "FWCore/Utilities/interface/Exception.h"
18     #include "FWCore/Framework/interface/EventSetup.h"
19    
20     #include "DataFormats/HeavyIonEvent/interface/CentralityBins.h"
21     #include "DataFormats/HeavyIonEvent/interface/Centrality.h"
22    
23     #include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
24     #include "DataFormats/JetReco/interface/CaloJetCollection.h"
25     #include "SimDataFormats/HiGenData/interface/GenHIEvent.h"
26     #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
27     #include "DataFormats/JetReco/interface/GenJetCollection.h"
28     #include "DataFormats/VertexReco/interface/Vertex.h"
29     #include "DataFormats/Math/interface/deltaPhi.h"
30     #include "DataFormats/Math/interface/deltaR.h"
31    
32     #include "DataFormats/Common/interface/TriggerResults.h"
33     #include "FWCore/Common/interface/TriggerNames.h"
34     #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutSetup.h"
35     #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
36     #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerObjectMapRecord.h"
37     #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerObjectMap.h"
38     #include "L1Trigger/GlobalTrigger/interface/L1GlobalTrigger.h"
39    
40     using namespace std;
41     using namespace edm;
42     using namespace reco;
43    
44     HiInclusiveJetAnalyzer::HiInclusiveJetAnalyzer(const edm::ParameterSet& iConfig) {
45    
46    
47     jetTag_ = iConfig.getParameter<InputTag>("jetTag");
48 mnguyen 1.8 vtxTag_ = iConfig.getUntrackedParameter<edm::InputTag>("vtxTag",edm::InputTag("hiSelectedVertex"));
49 yilmaz 1.1
50     isMC_ = iConfig.getUntrackedParameter<bool>("isMC",false);
51 mnguyen 1.8
52 yilmaz 1.1 if(isMC_){
53     genjetTag_ = iConfig.getParameter<InputTag>("genjetTag");
54     eventInfoTag_ = iConfig.getParameter<InputTag>("eventInfoTag");
55     }
56     verbose_ = iConfig.getUntrackedParameter<bool>("verbose",false);
57    
58     useCentrality_ = iConfig.getUntrackedParameter<bool>("useCentrality",false);
59 yjlee 1.4 useVtx_ = iConfig.getUntrackedParameter<bool>("useVtx",false);
60 yilmaz 1.1 useJEC_ = iConfig.getUntrackedParameter<bool>("useJEC",true);
61 yilmaz 1.2 usePat_ = iConfig.getUntrackedParameter<bool>("usePAT",true);
62 yilmaz 1.1
63 mnguyen 1.8 doLifeTimeTagging_ = iConfig.getUntrackedParameter<bool>("doLifeTimeTagging",false);
64     pfCandidateLabel_ = iConfig.getUntrackedParameter<edm::InputTag>("pfCandidateLabel",edm::InputTag("particleFlowTmp"));
65    
66 yilmaz 1.1 if(!isMC_){
67     L1gtReadout_ = iConfig.getParameter<edm::InputTag>("L1gtReadout");
68     hltResName_ = iConfig.getUntrackedParameter<string>("hltTrgResults","TriggerResults::HLT");
69    
70    
71     if (iConfig.exists("hltTrgNames"))
72     hltTrgNames_ = iConfig.getUntrackedParameter<vector<string> >("hltTrgNames");
73    
74     if (iConfig.exists("hltProcNames"))
75     hltProcNames_ = iConfig.getUntrackedParameter<vector<string> >("hltProcNames");
76     else {
77     hltProcNames_.push_back("FU");
78     hltProcNames_.push_back("HLT");
79     }
80     }
81    
82     cout<<" jet collection : "<<jetTag_<<endl;
83     if(isMC_)cout<<" genjet collection : "<<genjetTag_<<endl;
84    
85    
86    
87     }
88    
89    
90    
91     HiInclusiveJetAnalyzer::~HiInclusiveJetAnalyzer() { }
92    
93    
94    
95     void
96     HiInclusiveJetAnalyzer::beginRun(const edm::Run& run,
97     const edm::EventSetup & es) {}
98    
99     void
100     HiInclusiveJetAnalyzer::beginJob() {
101    
102     centrality_ = 0;
103    
104     //string jetTagName = jetTag_.label()+"_tree";
105     string jetTagTitle = jetTag_.label()+" Jet Analysis Tree";
106     t = fs1->make<TTree>("t",jetTagTitle.c_str());
107    
108 yjlee 1.6 // TTree* t= new TTree("t","Jet Response Analyzer");
109 yjlee 1.4 //t->Branch("run",&jets_.run,"run/I");
110 yilmaz 1.1 t->Branch("evt",&jets_.evt,"evt/I");
111 yjlee 1.4 //t->Branch("lumi",&jets_.lumi,"lumi/I");
112 yilmaz 1.1 t->Branch("b",&jets_.b,"b/F");
113 yjlee 1.4 if (useVtx_) {
114     t->Branch("vx",&jets_.vx,"vx/F");
115     t->Branch("vy",&jets_.vy,"vy/F");
116     t->Branch("vz",&jets_.vz,"vz/F");
117     }
118     if (useCentrality_) {
119     t->Branch("hf",&jets_.hf,"hf/F");
120     t->Branch("bin",&jets_.bin,"bin/I");
121     }
122 yjlee 1.7
123     t->Branch("nref",&jets_.nref,"nref/I");
124 yilmaz 1.1 t->Branch("rawpt",jets_.rawpt,"rawpt[nref]/F");
125     t->Branch("jtpt",jets_.jtpt,"jtpt[nref]/F");
126     t->Branch("jteta",jets_.jteta,"jteta[nref]/F");
127     t->Branch("jty",jets_.jty,"jty[nref]/F");
128     t->Branch("jtphi",jets_.jtphi,"jtphi[nref]/F");
129 yilmaz 1.2 t->Branch("jtpu",jets_.jtpu,"jtpu[nref]/F");
130 yilmaz 1.1
131 mnguyen 1.8 // b-jet discriminators
132     if (doLifeTimeTagging_) {
133    
134     t->Branch("discr_csvMva",jets_.discr_csvMva,"discr_csvMva[nref]/F");
135     t->Branch("discr_csvSimple",jets_.discr_csvSimple,"discr_csvSimple[nref]/F");
136     t->Branch("discr_muByIp3",jets_.discr_muByIp3,"discr_muByIp3[nref]/F");
137     t->Branch("discr_muByPt",jets_.discr_muByPt,"discr_muByPt[nref]/F");
138     t->Branch("discr_prob",jets_.discr_prob,"discr_prob[nref]/F");
139     t->Branch("discr_probb",jets_.discr_probb,"discr_probb[nref]/F");
140     t->Branch("discr_tcHighEff",jets_.discr_tcHighEff,"discr_tcHighEff[nref]/F");
141     t->Branch("discr_tcHighPur",jets_.discr_tcHighPur,"discr_tcHighPur[nref]/F");
142    
143     t->Branch("nsvtx", jets_.nsvtx, "nsvtx[nref]/b");
144     t->Branch("svtxntrk", jets_.svtxntrk, "svtxntrk[nref]/b");
145     t->Branch("svtxdl", jets_.svtxdl, "svtxdl[nref]/F");
146     t->Branch("svtxdls", jets_.svtxdls, "svtxdls[nref]/F");
147     t->Branch("svtxm", jets_.svtxm, "svtxm[nref]/F");
148     t->Branch("svtxpt", jets_.svtxpt, "svtxpt[nref]/F");
149    
150     t->Branch("nIPtracks",jets_.nIPtracks,"nIPTracks[nref]/b");
151     t->Branch("nselIPtracks",jets_.nselIPtracks,"nselIPTracks[nref]/b");
152    
153     t->Branch("mue", jets_.mue, "mue[nref]/F");
154     t->Branch("mupt", jets_.mupt, "mupt[nref]/F");
155     t->Branch("mueta", jets_.mueta, "mueta[nref]/F");
156     t->Branch("muphi", jets_.muphi, "muphi[nref]/F");
157     t->Branch("mudr", jets_.mudr, "mudr[nref]/F");
158     t->Branch("muptrel", jets_.muptrel, "muptrel[nref]/F");
159     t->Branch("muchg", jets_.muchg, "muchg[nref]/I");
160     }
161    
162 yilmaz 1.1 if(isMC_){
163     t->Branch("pthat",&jets_.pthat,"pthat/F");
164    
165     // Only matched gen jets
166     t->Branch("refpt",jets_.refpt,"refpt[nref]/F");
167     t->Branch("refeta",jets_.refeta,"refeta[nref]/F");
168     t->Branch("refy",jets_.refy,"refy[nref]/F");
169     t->Branch("refphi",jets_.refphi,"refphi[nref]/F");
170     t->Branch("refdphijt",jets_.refdphijt,"refdphijt[nref]/F");
171     t->Branch("refdrjt",jets_.refdrjt,"refdrjt[nref]/F");
172     // matched parton
173     t->Branch("refparton_pt",jets_.refparton_pt,"refparton_pt[nref]/F");
174 mnguyen 1.8 t->Branch("refparton_flavor",jets_.refparton_flavor,"refparton_flavor[nref]/I");
175     t->Branch("refparton_flavorForB",jets_.refparton_flavorForB,"refparton_flavorForB[nref]/I");
176 yilmaz 1.1
177     // For all gen jets, matched or unmatched
178     t->Branch("ngen",&jets_.ngen,"ngen/I");
179     t->Branch("genmatchindex",jets_.genmatchindex,"genmatchindex[ngen]/I");
180     t->Branch("genpt",jets_.genpt,"genpt[ngen]/F");
181     t->Branch("geneta",jets_.geneta,"geneta[ngen]/F");
182     t->Branch("geny",jets_.geny,"geny[ngen]/F");
183     t->Branch("genphi",jets_.genphi,"genphi[ngen]/F");
184     t->Branch("gendphijt",jets_.gendphijt,"gendphijt[ngen]/F");
185     t->Branch("gendrjt",jets_.gendrjt,"gendrjt[ngen]/F");
186     }
187 yjlee 1.7 /*
188 yilmaz 1.1 if(!isMC_){
189     t->Branch("nL1TBit",&jets_.nL1TBit,"nL1TBit/I");
190     t->Branch("l1TBit",jets_.l1TBit,"l1TBit[nL1TBit]/O");
191    
192     t->Branch("nL1ABit",&jets_.nL1ABit,"nL1ABit/I");
193     t->Branch("l1ABit",jets_.l1ABit,"l1ABit[nL1ABit]/O");
194    
195     t->Branch("nHLTBit",&jets_.nHLTBit,"nHLTBit/I");
196     t->Branch("hltBit",jets_.hltBit,"hltBit[nHLTBit]/O");
197    
198     }
199 yjlee 1.7 */
200 yilmaz 1.1 TH1D::SetDefaultSumw2();
201    
202    
203     }
204    
205    
206     void
207     HiInclusiveJetAnalyzer::analyze(const Event& iEvent,
208     const EventSetup& iSetup) {
209    
210     int event = iEvent.id().event();
211     int run = iEvent.id().run();
212     int lumi = iEvent.id().luminosityBlock();
213    
214     jets_.run = run;
215     jets_.evt = event;
216     jets_.lumi = lumi;
217    
218     LogDebug("HiInclusiveJetAnalyzer")<<"START event: "<<event<<" in run "<<run<<endl;
219    
220     int bin = -1;
221     double hf = 0.;
222     double b = 999.;
223    
224     if(useCentrality_){
225     if(!centrality_) centrality_ = new CentralityProvider(iSetup);
226     centrality_->newEvent(iEvent,iSetup); // make sure you do this first in every event
227     //double c = centrality_->centralityValue();
228     const reco::Centrality *cent = centrality_->raw();
229    
230     hf = cent->EtHFhitSum();
231    
232     bin = centrality_->getBin();
233     b = centrality_->bMean();
234     }
235    
236    
237    
238    
239     //double jetPtMin = 35.;
240    
241    
242     // loop the events
243    
244     jets_.bin = bin;
245     jets_.hf = hf;
246    
247    
248     if (useVtx_) {
249     edm::Handle<vector<reco::Vertex> >vertex;
250     iEvent.getByLabel(vtxTag_, vertex);
251    
252     if(vertex->size()>0) {
253     jets_.vx=vertex->begin()->x();
254     jets_.vy=vertex->begin()->y();
255     jets_.vz=vertex->begin()->z();
256     }
257     }
258    
259 yilmaz 1.2 edm::Handle<pat::JetCollection> patjets;
260     if(usePat_)iEvent.getByLabel(jetTag_, patjets);
261    
262     edm::Handle<reco::JetView> jets;
263 yilmaz 1.1 iEvent.getByLabel(jetTag_, jets);
264    
265 mnguyen 1.8 edm::Handle<reco::PFCandidateCollection> pfCandidates;
266     if(doLifeTimeTagging_){
267     iEvent.getByLabel(pfCandidateLabel_,pfCandidates);
268     }
269 yilmaz 1.1 // FILL JRA TREE
270    
271     jets_.b = b;
272     jets_.nref = 0;
273    
274     if(!isMC_){
275     fillL1Bits(iEvent);
276     fillHLTBits(iEvent);
277     }
278    
279     for(unsigned int j = 0 ; j < jets->size(); ++j){
280 yilmaz 1.2 const reco::Jet& jet = (*jets)[j];
281 yilmaz 1.1
282     //cout<<" jet pt "<<jet.pt()<<endl;
283     //if(jet.pt() < jetPtMin) continue;
284 yilmaz 1.2 if (useJEC_ && usePat_){
285     jets_.rawpt[jets_.nref]=(*patjets)[j].correctedJet("Uncorrected").pt();
286     }
287 mnguyen 1.8
288     if(doLifeTimeTagging_){
289    
290     if(jetTag_.label()=="icPu5patJets"){
291     jets_.discr_csvMva[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5CombinedSecondaryVertexMVABJetTags");
292     jets_.discr_csvSimple[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5CombinedSecondaryVertexBJetTags");
293     jets_.discr_muByIp3[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5SoftMuonByIP3dBJetTags");
294     jets_.discr_muByPt[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5SoftMuonByPtBJetTags");
295     jets_.discr_prob[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5JetBProbabilityBJetTags");
296     jets_.discr_probb[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5JetProbabilityBJetTags");
297     jets_.discr_tcHighEff[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5TrackCountingHighEffBJetTags");
298     jets_.discr_tcHighPur[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5TrackCountingHighPurBJetTags");
299     }
300     else if(jetTag_.label()=="akPu3PFpatJets"){
301     jets_.discr_csvMva[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFCombinedSecondaryVertexMVABJetTags");
302     jets_.discr_csvSimple[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFCombinedSecondaryVertexBJetTags");
303     jets_.discr_muByIp3[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFSoftMuonByIP3dBJetTags");
304     jets_.discr_muByPt[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFSoftMuonByPtBJetTags");
305     jets_.discr_prob[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFJetBProbabilityBJetTags");
306     jets_.discr_probb[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFJetProbabilityBJetTags");
307     jets_.discr_tcHighEff[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFTrackCountingHighEffBJetTags");
308     jets_.discr_tcHighPur[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFTrackCountingHighPurBJetTags");
309     }
310     else{
311     cout<<" you fail at btagging "<<endl;
312     }
313    
314     const reco::SecondaryVertexTagInfo& tagInfoSV=*(*patjets)[j].tagInfoSecondaryVertex();
315     if (tagInfoSV.nVertices()>0) {
316     Measurement1D m1D = tagInfoSV.flightDistance(0);
317     jets_.nsvtx[jets_.nref] = tagInfoSV.nVertices();
318     jets_.svtxntrk[jets_.nref] = tagInfoSV.nVertexTracks(0);
319     jets_.svtxdl[jets_.nref] = m1D.value();
320     jets_.svtxdls[jets_.nref] = m1D.significance();
321    
322     const reco::Vertex& svtx = tagInfoSV.secondaryVertex(0);
323    
324     jets_.svtxm[jets_.nref] = svtx.p4().mass();
325     jets_.svtxpt[jets_.nref] = svtx.p4().pt();
326    
327     }
328     else {
329     jets_.nsvtx[jets_.nref] = 0;
330     jets_.svtxntrk[jets_.nref] = 0;
331     jets_.svtxdl[jets_.nref] = 0.0;
332     jets_.svtxdls[jets_.nref] = 0.0;
333     jets_.svtxm[jets_.nref] = 0.0;
334     jets_.svtxpt[jets_.nref] = 0.0;
335     }
336     const reco::TrackIPTagInfo& tagInfoIP=*(*patjets)[j].tagInfoTrackIP();
337    
338     jets_.nIPtracks[jets_.nref] = tagInfoIP.tracks().size();
339     jets_.nselIPtracks[jets_.nref] = tagInfoIP.selectedTracks().size();
340    
341     // would be nice to save the info for the tracks, but requires a more sophisticated tree
342     /*
343     //cout << "Jet pt: " << tagInfoIP.jet()->pt() << endl;
344     cout << "Tot tracks: " << tagInfoIP.tracks().size() << endl;
345     TrackRefVector selTracks=tagInfoIP.selectedTracks();
346     int n=selTracks.size();
347     cout << "Sel tracks: " << n << endl;
348     // false cout << " Pt \t d len \t jet dist \t p3d \t p2d\t ip3d \t ip2d " << endl;
349     GlobalPoint pv(tagInfoIP.primaryVertex()->position().x(),tagInfoIP.primaryVertex()->position().y(),tagInfoIP.primaryVertex()->position().z());
350     cout << pv << " vs " << tagInfoIP.primaryVertex()->position() << endl;
351     for(int j=0;j< n;j++)
352     {
353     TrackIPTagInfo::TrackIPData data = tagInfoIP.impactParameterData()[j];
354     cout << selTracks[j]->pt() << "\t";
355     cout << tagInfoIP.probabilities(0)[j] << "\t";
356     cout << tagInfoIP.probabilities(1)[j] << "\t";
357     cout << data.ip3d.value() << "\t";
358     cout << data.ip3d.significance() << "\t";
359     cout << data.distanceToJetAxis.value() << "\t";
360     cout << data.distanceToJetAxis.significance() << "\t";
361     cout << data.distanceToGhostTrack.value() << "\t";
362     cout << data.distanceToGhostTrack.significance() << "\t";
363     cout << data.closestToJetAxis << "\t";
364     cout << (data.closestToJetAxis - pv).mag() << "\t";
365     cout << data.closestToGhostTrack << "\t";
366     cout << (data.closestToGhostTrack - pv).mag() << "\t";
367     cout << data.ip2d.value() << "\t";
368     cout << data.ip2d.significance() << endl;
369     }
370     */
371    
372     const reco::PFCandidateCollection *pfCandidateColl = &(*pfCandidates);
373     int pfMuonIndex = getPFJetMuon(jet, pfCandidateColl);
374     if(pfMuonIndex >=0){
375     const reco::PFCandidate muon = pfCandidateColl->at(pfMuonIndex);
376     jets_.mupt[jets_.nref] = muon.pt();
377     jets_.mueta[jets_.nref] = muon.eta();
378     jets_.muphi[jets_.nref] = muon.phi();
379     jets_.mue[jets_.nref] = muon.energy();
380     jets_.mudr[jets_.nref] = reco::deltaR(jet,muon);
381     jets_.muptrel[jets_.nref] = getPtRel(muon, jet);
382     jets_.muchg[jets_.nref] = muon.charge();
383     }
384     else{
385     jets_.mupt[jets_.nref] = 0.0;
386     jets_.mueta[jets_.nref] = 0.0;
387     jets_.muphi[jets_.nref] = 0.0;
388     jets_.mue[jets_.nref] = 0.0;
389     jets_.mudr[jets_.nref] = 9.9;
390     jets_.muptrel[jets_.nref] = 0.0;
391     jets_.muchg[jets_.nref] = 0;
392     }
393     }
394    
395 yilmaz 1.1 jets_.jtpt[jets_.nref] = jet.pt();
396     jets_.jteta[jets_.nref] = jet.eta();
397     jets_.jtphi[jets_.nref] = jet.phi();
398     jets_.jty[jets_.nref] = jet.eta();
399 yilmaz 1.2 jets_.jtpu[jets_.nref] = jet.pileup();
400 yilmaz 1.1
401 yilmaz 1.2 if(isMC_ && usePat_){
402 mnguyen 1.8 jets_.refparton_flavorForB[jets_.nref] = (*patjets)[j].partonFlavour();
403 yilmaz 1.2 const reco::GenJet* genjet = (*patjets)[j].genJet();
404 mnguyen 1.8
405 yilmaz 1.2 if(genjet){
406     jets_.refpt[jets_.nref] = genjet->pt();
407     jets_.refeta[jets_.nref] = genjet->eta();
408     jets_.refphi[jets_.nref] = genjet->phi();
409     jets_.refy[jets_.nref] = genjet->eta();
410     jets_.refdphijt[jets_.nref] = reco::deltaPhi(jet.phi(), genjet->phi());
411     jets_.refdrjt[jets_.nref] = reco::deltaR(jet.eta(),jet.phi(),genjet->eta(),genjet->phi());
412     }else{
413     jets_.refpt[jets_.nref] = -999.;
414     jets_.refeta[jets_.nref] = -999.;
415     jets_.refphi[jets_.nref] = -999.;
416     jets_.refy[jets_.nref] = -999.;
417     jets_.refdphijt[jets_.nref] = -999.;
418     jets_.refdrjt[jets_.nref] = -999.;
419     }
420 yilmaz 1.1
421 yilmaz 1.2 // matched partons
422     const reco::GenParticle * parton = (*patjets)[j].genParton();
423     if(parton){
424 mnguyen 1.8 jets_.refparton_pt[jets_.nref] = parton->pt();
425     jets_.refparton_flavor[jets_.nref] = parton->pdgId();
426     } else {
427 yilmaz 1.1 jets_.refparton_pt[jets_.nref] = -999;
428     jets_.refparton_flavor[jets_.nref] = -999;
429 mnguyen 1.8 }
430 yilmaz 1.2 }
431 yilmaz 1.1
432     jets_.nref++;
433    
434    
435     }
436    
437    
438     if(isMC_){
439    
440     edm::Handle<GenEventInfoProduct> hEventInfo;
441     iEvent.getByLabel(eventInfoTag_,hEventInfo);
442     //jets_.pthat = hEventInfo->binningValues()[0];
443    
444     // binning values and qscale appear to be equivalent, but binning values not always present
445     jets_.pthat = hEventInfo->qScale();
446    
447     edm::Handle<vector<reco::GenJet> >genjets;
448     iEvent.getByLabel(genjetTag_, genjets);
449    
450     jets_.ngen = 0;
451     for(unsigned int igen = 0 ; igen < genjets->size(); ++igen){
452     const reco::GenJet & genjet = (*genjets)[igen];
453    
454     float genjet_pt = genjet.pt();
455    
456     // threshold to reduce size of output in minbias PbPb
457     if(genjet_pt>20.){
458    
459     jets_.genpt[jets_.ngen] = genjet_pt;
460     jets_.geneta[jets_.ngen] = genjet.eta();
461     jets_.genphi[jets_.ngen] = genjet.phi();
462     jets_.geny[jets_.ngen] = genjet.eta();
463    
464    
465     // find matching patJet if there is one
466    
467     jets_.gendrjt[jets_.ngen] = -1.0;
468     jets_.genmatchindex[jets_.ngen] = -1;
469    
470     for(unsigned int ijet = 0 ; ijet < jets->size(); ++ijet){
471     const pat::Jet& jet = (*jets)[ijet];
472 mnguyen 1.8 const reco::GenJet* matchedGenJet = (*patjets)[ijet].genJet();
473     if(matchedGenJet){
474     // poor man's matching, someone fix please
475     if(fabs(genjet.pt()-matchedGenJet->pt())<0.0001 &&
476     fabs(genjet.eta()-matchedGenJet->eta())<0.0001 &&
477     (fabs(genjet.phi()-matchedGenJet->phi())<0.0001 || fabs(genjet.phi()-matchedGenJet->phi() - 2.0*TMath::Pi()) < 0.0001 )){
478 yilmaz 1.1
479     jets_.genmatchindex[jets_.ngen] = (int)ijet;
480     jets_.gendphijt[jets_.ngen] = reco::deltaPhi(jet.phi(),genjet.phi());
481     jets_.gendrjt[jets_.ngen] = reco::deltaR(jet,genjet);
482    
483     break;
484     }
485     }
486     }
487     jets_.ngen++;
488     }
489    
490     }
491    
492     }
493     t->Fill();
494     }
495    
496    
497    
498    
499     //--------------------------------------------------------------------------------------------------
500     void HiInclusiveJetAnalyzer::fillL1Bits(const edm::Event &iEvent)
501     {
502     edm::Handle< L1GlobalTriggerReadoutRecord > L1GlobalTrigger;
503    
504     iEvent.getByLabel(L1gtReadout_, L1GlobalTrigger);
505     const TechnicalTriggerWord& technicalTriggerWordBeforeMask = L1GlobalTrigger->technicalTriggerWord();
506    
507     for (int i=0; i<64;i++)
508     {
509     jets_.l1TBit[i] = technicalTriggerWordBeforeMask.at(i);
510     }
511     jets_.nL1TBit = 64;
512    
513     int ntrigs = L1GlobalTrigger->decisionWord().size();
514     jets_.nL1ABit = ntrigs;
515    
516     for (int i=0; i != ntrigs; i++) {
517     bool accept = L1GlobalTrigger->decisionWord()[i];
518     //jets_.l1ABit[i] = (accept == true)? 1:0;
519     if(accept== true){
520     jets_.l1ABit[i] = 1;
521     }
522     else{
523     jets_.l1ABit[i] = 0;
524     }
525    
526     }
527     }
528    
529     //--------------------------------------------------------------------------------------------------
530     void HiInclusiveJetAnalyzer::fillHLTBits(const edm::Event &iEvent)
531     {
532     // Fill HLT trigger bits.
533     Handle<TriggerResults> triggerResultsHLT;
534     getProduct(hltResName_, triggerResultsHLT, iEvent);
535    
536     const TriggerResults *hltResults = triggerResultsHLT.product();
537     const TriggerNames & triggerNames = iEvent.triggerNames(*hltResults);
538    
539     jets_.nHLTBit = hltTrgNames_.size();
540    
541     for(size_t i=0;i<hltTrgNames_.size();i++){
542    
543     for(size_t j=0;j<triggerNames.size();++j) {
544    
545     if(triggerNames.triggerName(j) == hltTrgNames_[i]){
546    
547     //cout <<"hltTrgNames_(i) "<<hltTrgNames_[i]<<endl;
548     //cout <<"triggerName(j) "<<triggerNames.triggerName(j)<<endl;
549     //cout<<" result "<<triggerResultsHLT->accept(j)<<endl;
550     jets_.hltBit[i] = triggerResultsHLT->accept(j);
551     }
552    
553     }
554     }
555     }
556    
557     //--------------------------------------------------------------------------------------------------
558     template <typename TYPE>
559     inline void HiInclusiveJetAnalyzer::getProduct(const std::string name, edm::Handle<TYPE> &prod,
560     const edm::Event &event) const
561     {
562     // Try to access data collection from EDM file. We check if we really get just one
563     // product with the given name. If not we throw an exception.
564    
565     event.getByLabel(edm::InputTag(name),prod);
566     if (!prod.isValid())
567     throw edm::Exception(edm::errors::Configuration, "HiInclusiveJetAnalyzer::GetProduct()\n")
568     << "Collection with label '" << name << "' is not valid" << std::endl;
569     }
570    
571     //--------------------------------------------------------------------------------------------------
572     template <typename TYPE>
573     inline bool HiInclusiveJetAnalyzer::getProductSafe(const std::string name, edm::Handle<TYPE> &prod,
574     const edm::Event &event) const
575     {
576     // Try to safely access data collection from EDM file. We check if we really get just one
577     // product with the given name. If not, we return false.
578    
579     if (name.size()==0)
580     return false;
581    
582     try {
583     event.getByLabel(edm::InputTag(name),prod);
584     if (!prod.isValid())
585     return false;
586     } catch (...) {
587     return false;
588     }
589     return true;
590     }
591    
592    
593 mnguyen 1.8 int
594     HiInclusiveJetAnalyzer::getPFJetMuon(const pat::Jet& pfJet, const reco::PFCandidateCollection *pfCandidateColl)
595     {
596    
597     int pfMuonIndex = -1;
598     float ptMax = 0.;
599    
600    
601     for(unsigned icand=0;icand<pfCandidateColl->size(); icand++) {
602     const reco::PFCandidate pfCandidate = pfCandidateColl->at(icand);
603    
604     int id = pfCandidate.particleId();
605     if(abs(id) != 3) continue;
606 yilmaz 1.1
607 mnguyen 1.8 if(reco::deltaR(pfJet,pfCandidate)>0.5) continue;
608    
609     double pt = pfCandidate.pt();
610     if(pt>ptMax){
611     ptMax = pt;
612     pfMuonIndex = (int) icand;
613     }
614     }
615    
616     return pfMuonIndex;
617    
618     }
619    
620    
621     double
622     HiInclusiveJetAnalyzer::getPtRel(const reco::PFCandidate lep, const pat::Jet& jet )
623     {
624    
625     float lj_x = jet.p4().px();
626     float lj_y = jet.p4().py();
627     float lj_z = jet.p4().pz();
628    
629     // absolute values squared
630     float lj2 = lj_x*lj_x+lj_y*lj_y+lj_z*lj_z;
631     float lep2 = lep.px()*lep.px()+lep.py()*lep.py()+lep.pz()*lep.pz();
632    
633     // projection vec(mu) to lepjet axis
634     float lepXlj = lep.px()*lj_x+lep.py()*lj_y+lep.pz()*lj_z;
635    
636     // absolute value squared and normalized
637     float pLrel2 = lepXlj*lepXlj/lj2;
638    
639     // lep2 = pTrel2 + pLrel2
640     float pTrel2 = lep2-pLrel2;
641    
642     return (pTrel2 > 0) ? std::sqrt(pTrel2) : 0.0;
643     }
644 yilmaz 1.1
645     DEFINE_FWK_MODULE(HiInclusiveJetAnalyzer);