ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/CmsHi/JetAnalysis/src/HiInclusiveJetAnalyzer.cc
Revision: 1.11
Committed: Thu May 3 20:29:42 2012 UTC (13 years ago) by mnguyen
Content type: text/plain
Branch: MAIN
Changes since 1.10: +161 -36 lines
Log Message:
fixing a bug in btagging extras and adding the option to save the gen particles from b

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 mnguyen 1.10 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
28 yilmaz 1.1 #include "DataFormats/JetReco/interface/GenJetCollection.h"
29     #include "DataFormats/VertexReco/interface/Vertex.h"
30     #include "DataFormats/Math/interface/deltaPhi.h"
31     #include "DataFormats/Math/interface/deltaR.h"
32    
33     #include "DataFormats/Common/interface/TriggerResults.h"
34     #include "FWCore/Common/interface/TriggerNames.h"
35     #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutSetup.h"
36     #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
37     #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerObjectMapRecord.h"
38     #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerObjectMap.h"
39     #include "L1Trigger/GlobalTrigger/interface/L1GlobalTrigger.h"
40    
41     using namespace std;
42     using namespace edm;
43     using namespace reco;
44    
45     HiInclusiveJetAnalyzer::HiInclusiveJetAnalyzer(const edm::ParameterSet& iConfig) {
46    
47    
48     jetTag_ = iConfig.getParameter<InputTag>("jetTag");
49 mnguyen 1.8 vtxTag_ = iConfig.getUntrackedParameter<edm::InputTag>("vtxTag",edm::InputTag("hiSelectedVertex"));
50 yilmaz 1.1
51     isMC_ = iConfig.getUntrackedParameter<bool>("isMC",false);
52 mnguyen 1.8
53 yilmaz 1.1 if(isMC_){
54     genjetTag_ = iConfig.getParameter<InputTag>("genjetTag");
55     eventInfoTag_ = iConfig.getParameter<InputTag>("eventInfoTag");
56     }
57     verbose_ = iConfig.getUntrackedParameter<bool>("verbose",false);
58    
59     useCentrality_ = iConfig.getUntrackedParameter<bool>("useCentrality",false);
60 yjlee 1.4 useVtx_ = iConfig.getUntrackedParameter<bool>("useVtx",false);
61 yilmaz 1.1 useJEC_ = iConfig.getUntrackedParameter<bool>("useJEC",true);
62 yilmaz 1.2 usePat_ = iConfig.getUntrackedParameter<bool>("usePAT",true);
63 yilmaz 1.1
64 mnguyen 1.8 doLifeTimeTagging_ = iConfig.getUntrackedParameter<bool>("doLifeTimeTagging",false);
65 mnguyen 1.10 doLifeTimeTaggingExtras_ = iConfig.getUntrackedParameter<bool>("doLifeTimeTaggingExtras",true);
66 mnguyen 1.11 saveBfragments_ = iConfig.getUntrackedParameter<bool>("saveBfragments",false);
67    
68 mnguyen 1.8 pfCandidateLabel_ = iConfig.getUntrackedParameter<edm::InputTag>("pfCandidateLabel",edm::InputTag("particleFlowTmp"));
69    
70 yilmaz 1.1 if(!isMC_){
71     L1gtReadout_ = iConfig.getParameter<edm::InputTag>("L1gtReadout");
72     hltResName_ = iConfig.getUntrackedParameter<string>("hltTrgResults","TriggerResults::HLT");
73    
74    
75     if (iConfig.exists("hltTrgNames"))
76     hltTrgNames_ = iConfig.getUntrackedParameter<vector<string> >("hltTrgNames");
77    
78     if (iConfig.exists("hltProcNames"))
79     hltProcNames_ = iConfig.getUntrackedParameter<vector<string> >("hltProcNames");
80     else {
81     hltProcNames_.push_back("FU");
82     hltProcNames_.push_back("HLT");
83     }
84     }
85    
86     cout<<" jet collection : "<<jetTag_<<endl;
87 yilmaz 1.9 if(isMC_){
88     cout<<" genjet collection : "<<genjetTag_<<endl;
89     genPtMin_ = iConfig.getUntrackedParameter<double>("genPtMin",0);
90     }
91 yilmaz 1.1
92    
93     }
94    
95    
96    
97     HiInclusiveJetAnalyzer::~HiInclusiveJetAnalyzer() { }
98    
99    
100    
101     void
102     HiInclusiveJetAnalyzer::beginRun(const edm::Run& run,
103     const edm::EventSetup & es) {}
104    
105     void
106     HiInclusiveJetAnalyzer::beginJob() {
107    
108     centrality_ = 0;
109    
110     //string jetTagName = jetTag_.label()+"_tree";
111     string jetTagTitle = jetTag_.label()+" Jet Analysis Tree";
112     t = fs1->make<TTree>("t",jetTagTitle.c_str());
113    
114 yjlee 1.6 // TTree* t= new TTree("t","Jet Response Analyzer");
115 yjlee 1.4 //t->Branch("run",&jets_.run,"run/I");
116 yilmaz 1.1 t->Branch("evt",&jets_.evt,"evt/I");
117 yjlee 1.4 //t->Branch("lumi",&jets_.lumi,"lumi/I");
118 yilmaz 1.1 t->Branch("b",&jets_.b,"b/F");
119 yjlee 1.4 if (useVtx_) {
120     t->Branch("vx",&jets_.vx,"vx/F");
121     t->Branch("vy",&jets_.vy,"vy/F");
122     t->Branch("vz",&jets_.vz,"vz/F");
123     }
124     if (useCentrality_) {
125     t->Branch("hf",&jets_.hf,"hf/F");
126     t->Branch("bin",&jets_.bin,"bin/I");
127     }
128 yjlee 1.7
129     t->Branch("nref",&jets_.nref,"nref/I");
130 yilmaz 1.1 t->Branch("rawpt",jets_.rawpt,"rawpt[nref]/F");
131     t->Branch("jtpt",jets_.jtpt,"jtpt[nref]/F");
132     t->Branch("jteta",jets_.jteta,"jteta[nref]/F");
133     t->Branch("jty",jets_.jty,"jty[nref]/F");
134     t->Branch("jtphi",jets_.jtphi,"jtphi[nref]/F");
135 yilmaz 1.2 t->Branch("jtpu",jets_.jtpu,"jtpu[nref]/F");
136 yilmaz 1.1
137 mnguyen 1.8 // b-jet discriminators
138     if (doLifeTimeTagging_) {
139    
140     t->Branch("discr_csvMva",jets_.discr_csvMva,"discr_csvMva[nref]/F");
141     t->Branch("discr_csvSimple",jets_.discr_csvSimple,"discr_csvSimple[nref]/F");
142     t->Branch("discr_muByIp3",jets_.discr_muByIp3,"discr_muByIp3[nref]/F");
143     t->Branch("discr_muByPt",jets_.discr_muByPt,"discr_muByPt[nref]/F");
144     t->Branch("discr_prob",jets_.discr_prob,"discr_prob[nref]/F");
145     t->Branch("discr_probb",jets_.discr_probb,"discr_probb[nref]/F");
146     t->Branch("discr_tcHighEff",jets_.discr_tcHighEff,"discr_tcHighEff[nref]/F");
147     t->Branch("discr_tcHighPur",jets_.discr_tcHighPur,"discr_tcHighPur[nref]/F");
148    
149 mnguyen 1.10 t->Branch("nsvtx", jets_.nsvtx, "nsvtx[nref]/I");
150     t->Branch("svtxntrk", jets_.svtxntrk, "svtxntrk[nref]/I");
151 mnguyen 1.8 t->Branch("svtxdl", jets_.svtxdl, "svtxdl[nref]/F");
152     t->Branch("svtxdls", jets_.svtxdls, "svtxdls[nref]/F");
153     t->Branch("svtxm", jets_.svtxm, "svtxm[nref]/F");
154     t->Branch("svtxpt", jets_.svtxpt, "svtxpt[nref]/F");
155    
156 mnguyen 1.10 t->Branch("nIPtrk",jets_.nIPtrk,"nIPtrk[nref]/I");
157     t->Branch("nselIPtrk",jets_.nselIPtrk,"nselIPtrk[nref]/I");
158 mnguyen 1.11
159 mnguyen 1.10 if (doLifeTimeTaggingExtras_) {
160 mnguyen 1.11 t->Branch("nIP",&jets_.nIP,"nIP/I");
161     t->Branch("ipJetIndex",jets_.ipJetIndex,"ipJetIndex[nIP]/I");
162     t->Branch("ipPt",jets_.ipPt,"ipPt[nIP]/F");
163     t->Branch("ipProb0",jets_.ipProb0,"ipProb0[nIP]/F");
164     t->Branch("ipProb1",jets_.ipProb1,"ipProb1[nIP]/F");
165     t->Branch("ip2d",jets_.ip2d,"ip2d[nIP]/F");
166     t->Branch("ip2dSig",jets_.ip2dSig,"ip2dSig[nIP]/F");
167     t->Branch("ip3d",jets_.ip3d,"ip3d[nIP]/F");
168     t->Branch("ip3dSig",jets_.ip3dSig,"ip3dSig[nIP]/F");
169     t->Branch("ipDist2Jet",jets_.ipDist2Jet,"ipDist2Jet[nIP]/F");
170     t->Branch("ipDist2JetSig",jets_.ipDist2JetSig,"ipDist2JetSig[nIP]/F");
171     t->Branch("ipClosest2Jet",jets_.ipClosest2Jet,"ipClosest2Jet[nIP]/F");
172    
173 mnguyen 1.10 }
174 mnguyen 1.8
175     t->Branch("mue", jets_.mue, "mue[nref]/F");
176     t->Branch("mupt", jets_.mupt, "mupt[nref]/F");
177     t->Branch("mueta", jets_.mueta, "mueta[nref]/F");
178     t->Branch("muphi", jets_.muphi, "muphi[nref]/F");
179     t->Branch("mudr", jets_.mudr, "mudr[nref]/F");
180     t->Branch("muptrel", jets_.muptrel, "muptrel[nref]/F");
181     t->Branch("muchg", jets_.muchg, "muchg[nref]/I");
182     }
183    
184 yilmaz 1.1 if(isMC_){
185 mnguyen 1.10 t->Branch("beamId1",&jets_.beamId1,"beamId1/I");
186     t->Branch("beamId2",&jets_.beamId2,"beamId2/I");
187    
188 yilmaz 1.1 t->Branch("pthat",&jets_.pthat,"pthat/F");
189    
190     // Only matched gen jets
191     t->Branch("refpt",jets_.refpt,"refpt[nref]/F");
192     t->Branch("refeta",jets_.refeta,"refeta[nref]/F");
193     t->Branch("refy",jets_.refy,"refy[nref]/F");
194     t->Branch("refphi",jets_.refphi,"refphi[nref]/F");
195     t->Branch("refdphijt",jets_.refdphijt,"refdphijt[nref]/F");
196     t->Branch("refdrjt",jets_.refdrjt,"refdrjt[nref]/F");
197     // matched parton
198     t->Branch("refparton_pt",jets_.refparton_pt,"refparton_pt[nref]/F");
199 mnguyen 1.8 t->Branch("refparton_flavor",jets_.refparton_flavor,"refparton_flavor[nref]/I");
200     t->Branch("refparton_flavorForB",jets_.refparton_flavorForB,"refparton_flavorForB[nref]/I");
201 yilmaz 1.1
202     // For all gen jets, matched or unmatched
203     t->Branch("ngen",&jets_.ngen,"ngen/I");
204     t->Branch("genmatchindex",jets_.genmatchindex,"genmatchindex[ngen]/I");
205     t->Branch("genpt",jets_.genpt,"genpt[ngen]/F");
206     t->Branch("geneta",jets_.geneta,"geneta[ngen]/F");
207     t->Branch("geny",jets_.geny,"geny[ngen]/F");
208     t->Branch("genphi",jets_.genphi,"genphi[ngen]/F");
209     t->Branch("gendphijt",jets_.gendphijt,"gendphijt[ngen]/F");
210     t->Branch("gendrjt",jets_.gendrjt,"gendrjt[ngen]/F");
211 mnguyen 1.11
212     if(saveBfragments_ ) {
213     t->Branch("bMult",&jets_.bMult,"bMult/I");
214     t->Branch("bJetIndex",jets_.bJetIndex,"bJetIndex[bMult]/I");
215     t->Branch("bStatus",jets_.bStatus,"bStatus[bMult]/I");
216     t->Branch("bVx",jets_.bVx,"bVx[bMult]/F");
217     t->Branch("bVy",jets_.bVy,"bVy[bMult]/F");
218     t->Branch("bVz",jets_.bVz,"bVz[bMult]/F");
219     t->Branch("bPt",jets_.bPt,"bPt[bMult]/F");
220     t->Branch("bEta",jets_.bEta,"bEta[bMult]/F");
221     t->Branch("bPhi",jets_.bPhi,"bPhi[bMult]/F");
222     t->Branch("bPdg",jets_.bPdg,"bPdg[bMult]/I");
223     t->Branch("bChg",jets_.bChg,"bChg[bMult]/I");
224     }
225 yilmaz 1.1 }
226 yjlee 1.7 /*
227 yilmaz 1.1 if(!isMC_){
228     t->Branch("nL1TBit",&jets_.nL1TBit,"nL1TBit/I");
229     t->Branch("l1TBit",jets_.l1TBit,"l1TBit[nL1TBit]/O");
230    
231     t->Branch("nL1ABit",&jets_.nL1ABit,"nL1ABit/I");
232     t->Branch("l1ABit",jets_.l1ABit,"l1ABit[nL1ABit]/O");
233    
234     t->Branch("nHLTBit",&jets_.nHLTBit,"nHLTBit/I");
235     t->Branch("hltBit",jets_.hltBit,"hltBit[nHLTBit]/O");
236    
237     }
238 yjlee 1.7 */
239 yilmaz 1.1 TH1D::SetDefaultSumw2();
240    
241    
242     }
243    
244    
245     void
246     HiInclusiveJetAnalyzer::analyze(const Event& iEvent,
247     const EventSetup& iSetup) {
248    
249     int event = iEvent.id().event();
250     int run = iEvent.id().run();
251     int lumi = iEvent.id().luminosityBlock();
252    
253     jets_.run = run;
254     jets_.evt = event;
255     jets_.lumi = lumi;
256    
257     LogDebug("HiInclusiveJetAnalyzer")<<"START event: "<<event<<" in run "<<run<<endl;
258    
259     int bin = -1;
260     double hf = 0.;
261     double b = 999.;
262    
263     if(useCentrality_){
264     if(!centrality_) centrality_ = new CentralityProvider(iSetup);
265     centrality_->newEvent(iEvent,iSetup); // make sure you do this first in every event
266     //double c = centrality_->centralityValue();
267     const reco::Centrality *cent = centrality_->raw();
268    
269     hf = cent->EtHFhitSum();
270    
271     bin = centrality_->getBin();
272     b = centrality_->bMean();
273     }
274    
275    
276    
277    
278    
279     // loop the events
280    
281     jets_.bin = bin;
282     jets_.hf = hf;
283    
284    
285     if (useVtx_) {
286     edm::Handle<vector<reco::Vertex> >vertex;
287     iEvent.getByLabel(vtxTag_, vertex);
288    
289     if(vertex->size()>0) {
290     jets_.vx=vertex->begin()->x();
291     jets_.vy=vertex->begin()->y();
292     jets_.vz=vertex->begin()->z();
293     }
294     }
295    
296 yilmaz 1.2 edm::Handle<pat::JetCollection> patjets;
297     if(usePat_)iEvent.getByLabel(jetTag_, patjets);
298    
299     edm::Handle<reco::JetView> jets;
300 yilmaz 1.1 iEvent.getByLabel(jetTag_, jets);
301    
302 mnguyen 1.8 edm::Handle<reco::PFCandidateCollection> pfCandidates;
303     if(doLifeTimeTagging_){
304     iEvent.getByLabel(pfCandidateLabel_,pfCandidates);
305     }
306 yilmaz 1.1 // FILL JRA TREE
307    
308     jets_.b = b;
309     jets_.nref = 0;
310    
311     if(!isMC_){
312     fillL1Bits(iEvent);
313     fillHLTBits(iEvent);
314     }
315    
316     for(unsigned int j = 0 ; j < jets->size(); ++j){
317 yilmaz 1.2 const reco::Jet& jet = (*jets)[j];
318 yilmaz 1.1
319 mnguyen 1.11
320 yilmaz 1.2 if (useJEC_ && usePat_){
321     jets_.rawpt[jets_.nref]=(*patjets)[j].correctedJet("Uncorrected").pt();
322     }
323 mnguyen 1.8
324     if(doLifeTimeTagging_){
325    
326     if(jetTag_.label()=="icPu5patJets"){
327     jets_.discr_csvMva[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5CombinedSecondaryVertexMVABJetTags");
328     jets_.discr_csvSimple[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5CombinedSecondaryVertexBJetTags");
329     jets_.discr_muByIp3[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5SoftMuonByIP3dBJetTags");
330     jets_.discr_muByPt[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5SoftMuonByPtBJetTags");
331     jets_.discr_prob[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5JetBProbabilityBJetTags");
332     jets_.discr_probb[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5JetProbabilityBJetTags");
333     jets_.discr_tcHighEff[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5TrackCountingHighEffBJetTags");
334     jets_.discr_tcHighPur[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5TrackCountingHighPurBJetTags");
335     }
336     else if(jetTag_.label()=="akPu3PFpatJets"){
337     jets_.discr_csvMva[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFCombinedSecondaryVertexMVABJetTags");
338     jets_.discr_csvSimple[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFCombinedSecondaryVertexBJetTags");
339     jets_.discr_muByIp3[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFSoftMuonByIP3dBJetTags");
340     jets_.discr_muByPt[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFSoftMuonByPtBJetTags");
341     jets_.discr_prob[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFJetBProbabilityBJetTags");
342     jets_.discr_probb[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFJetProbabilityBJetTags");
343     jets_.discr_tcHighEff[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFTrackCountingHighEffBJetTags");
344     jets_.discr_tcHighPur[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFTrackCountingHighPurBJetTags");
345     }
346     else{
347 mnguyen 1.10 cout<<" b-tagging variables not filled for this collection, turn of doLifeTimeTagging "<<endl;
348 mnguyen 1.8 }
349 mnguyen 1.10
350     const reco::SecondaryVertexTagInfo& tagInfoSV=*(*patjets)[j].tagInfoSecondaryVertex();
351    
352     jets_.nsvtx[jets_.nref] = tagInfoSV.nVertices();
353 mnguyen 1.8
354     if (tagInfoSV.nVertices()>0) {
355     jets_.svtxntrk[jets_.nref] = tagInfoSV.nVertexTracks(0);
356 mnguyen 1.10 // this is the 3d flight distance, for 2-D use (0,true)
357     Measurement1D m1D = tagInfoSV.flightDistance(0);
358 mnguyen 1.8 jets_.svtxdl[jets_.nref] = m1D.value();
359     jets_.svtxdls[jets_.nref] = m1D.significance();
360    
361 mnguyen 1.10 const reco::Vertex& svtx = tagInfoSV.secondaryVertex(0);
362 mnguyen 1.11 //cout<<" SV: vx: "<<svtx.x()<<" vy "<<svtx.y()<<" vz "<<svtx.z()<<endl;
363     //cout<<" PV: vx: "<<jet.vx()<<" vy "<<jet.vy()<<" vz "<<jet.vz()<<endl;
364 mnguyen 1.8 jets_.svtxm[jets_.nref] = svtx.p4().mass();
365 mnguyen 1.10 jets_.svtxpt[jets_.nref] = svtx.p4().pt();
366 mnguyen 1.11 //cout<<" chi2 "<<svtx.chi2()<<" ndof "<<svtx.ndof()<<endl;
367 mnguyen 1.8 }
368 mnguyen 1.10
369 mnguyen 1.8 const reco::TrackIPTagInfo& tagInfoIP=*(*patjets)[j].tagInfoTrackIP();
370 mnguyen 1.10
371     jets_.nIPtrk[jets_.nref] = tagInfoIP.tracks().size();
372     jets_.nselIPtrk[jets_.nref] = tagInfoIP.selectedTracks().size();
373 mnguyen 1.8
374 mnguyen 1.10 if (doLifeTimeTaggingExtras_) {
375 mnguyen 1.8
376 mnguyen 1.10 TrackRefVector selTracks=tagInfoIP.selectedTracks();
377    
378     GlobalPoint pv(tagInfoIP.primaryVertex()->position().x(),tagInfoIP.primaryVertex()->position().y(),tagInfoIP.primaryVertex()->position().z());
379    
380     for(int it=0;it<jets_.nselIPtrk[jets_.nref] ;it++)
381     {
382 mnguyen 1.11 jets_.ipJetIndex[jets_.nIP + it]= jets_.nref;
383 mnguyen 1.10 TrackIPTagInfo::TrackIPData data = tagInfoIP.impactParameterData()[it];
384 mnguyen 1.11 jets_.ipPt[jets_.nIP + it] = selTracks[it]->pt();
385     jets_.ipProb0[jets_.nIP + it] = tagInfoIP.probabilities(0)[it];
386     jets_.ipProb1[jets_.nIP + it] = tagInfoIP.probabilities(1)[it];
387     jets_.ip2d[jets_.nIP + it] = data.ip2d.value();
388     jets_.ip2dSig[jets_.nIP + it] = data.ip2d.significance();
389     jets_.ip3d[jets_.nIP + it] = data.ip3d.value();
390     jets_.ip3dSig[jets_.nIP + it] = data.ip3d.significance();
391     jets_.ipDist2Jet[jets_.nIP + it] = data.distanceToJetAxis.value();
392     jets_.ipDist2JetSig[jets_.nIP + it] = data.distanceToJetAxis.significance();
393     jets_.ipClosest2Jet[jets_.nIP + it] = (data.closestToJetAxis - pv).mag(); //decay length
394 mnguyen 1.10 }
395 mnguyen 1.11
396     jets_.nIP += jets_.nselIPtrk[jets_.nref];
397    
398 mnguyen 1.10 }
399    
400 mnguyen 1.8 const reco::PFCandidateCollection *pfCandidateColl = &(*pfCandidates);
401     int pfMuonIndex = getPFJetMuon(jet, pfCandidateColl);
402     if(pfMuonIndex >=0){
403     const reco::PFCandidate muon = pfCandidateColl->at(pfMuonIndex);
404     jets_.mupt[jets_.nref] = muon.pt();
405     jets_.mueta[jets_.nref] = muon.eta();
406     jets_.muphi[jets_.nref] = muon.phi();
407     jets_.mue[jets_.nref] = muon.energy();
408     jets_.mudr[jets_.nref] = reco::deltaR(jet,muon);
409     jets_.muptrel[jets_.nref] = getPtRel(muon, jet);
410     jets_.muchg[jets_.nref] = muon.charge();
411     }
412     }
413    
414 yilmaz 1.1 jets_.jtpt[jets_.nref] = jet.pt();
415     jets_.jteta[jets_.nref] = jet.eta();
416     jets_.jtphi[jets_.nref] = jet.phi();
417     jets_.jty[jets_.nref] = jet.eta();
418 yilmaz 1.2 jets_.jtpu[jets_.nref] = jet.pileup();
419 yilmaz 1.1
420 yilmaz 1.2 if(isMC_ && usePat_){
421 mnguyen 1.11
422     const reco::GenJet * genjet = (*patjets)[j].genJet();
423 mnguyen 1.8
424 yilmaz 1.2 if(genjet){
425     jets_.refpt[jets_.nref] = genjet->pt();
426     jets_.refeta[jets_.nref] = genjet->eta();
427     jets_.refphi[jets_.nref] = genjet->phi();
428     jets_.refy[jets_.nref] = genjet->eta();
429     jets_.refdphijt[jets_.nref] = reco::deltaPhi(jet.phi(), genjet->phi());
430     jets_.refdrjt[jets_.nref] = reco::deltaR(jet.eta(),jet.phi(),genjet->eta(),genjet->phi());
431     }else{
432     jets_.refpt[jets_.nref] = -999.;
433     jets_.refeta[jets_.nref] = -999.;
434     jets_.refphi[jets_.nref] = -999.;
435     jets_.refy[jets_.nref] = -999.;
436     jets_.refdphijt[jets_.nref] = -999.;
437     jets_.refdrjt[jets_.nref] = -999.;
438     }
439 yilmaz 1.1
440 mnguyen 1.11 jets_.refparton_flavorForB[jets_.nref] = (*patjets)[j].partonFlavour();
441    
442     // matched partons
443     const reco::GenParticle & parton = *(*patjets)[j].genParton();
444    
445     if((*patjets)[j].genParton()){
446     jets_.refparton_pt[jets_.nref] = parton.pt();
447     jets_.refparton_flavor[jets_.nref] = parton.pdgId();
448    
449     if(saveBfragments_ && abs(jets_.refparton_flavorForB[jets_.nref])==5){
450    
451     usedStringPts.clear();
452    
453     // uncomment this if you want to know the ugly truth about parton matching -matt
454     //if(jet.pt() > 50 &&abs(parton.pdgId())!=5 && parton.pdgId()!=21)
455     // cout<<" Identified as a b, but doesn't match b or gluon, id = "<<parton.pdgId()<<endl;
456    
457     jets_.bJetIndex[jets_.bMult] = jets_.nref;
458     jets_.bStatus[jets_.bMult] = parton.status();
459     jets_.bVx[jets_.bMult] = parton.vx();
460     jets_.bVy[jets_.bMult] = parton.vy();
461     jets_.bVz[jets_.bMult] = parton.vz();
462     jets_.bPt[jets_.bMult] = parton.pt();
463     jets_.bEta[jets_.bMult] = parton.eta();
464     jets_.bPhi[jets_.bMult] = parton.phi();
465     jets_.bPdg[jets_.bMult] = parton.pdgId();
466     jets_.bChg[jets_.bMult] = parton.charge();
467     jets_.bMult++;
468     saveDaughters(parton);
469     }
470    
471    
472 mnguyen 1.8 } else {
473 mnguyen 1.11 jets_.refparton_pt[jets_.nref] = -999;
474     jets_.refparton_flavor[jets_.nref] = -999;
475 mnguyen 1.8 }
476 mnguyen 1.11
477    
478 yilmaz 1.2 }
479 yilmaz 1.1
480     jets_.nref++;
481    
482    
483     }
484    
485    
486     if(isMC_){
487    
488 mnguyen 1.10 edm::Handle<HepMCProduct> hepMCProduct;
489     iEvent.getByLabel(eventInfoTag_,hepMCProduct);
490     const HepMC::GenEvent* MCEvt = hepMCProduct->GetEvent();
491     std::pair<HepMC::GenParticle*,HepMC::GenParticle*> beamParticles = MCEvt->beam_particles();
492     jets_.beamId1 = beamParticles.first->pdg_id();
493     jets_.beamId2 = beamParticles.second->pdg_id();
494    
495 yilmaz 1.1 edm::Handle<GenEventInfoProduct> hEventInfo;
496     iEvent.getByLabel(eventInfoTag_,hEventInfo);
497     //jets_.pthat = hEventInfo->binningValues()[0];
498    
499     // binning values and qscale appear to be equivalent, but binning values not always present
500     jets_.pthat = hEventInfo->qScale();
501    
502     edm::Handle<vector<reco::GenJet> >genjets;
503     iEvent.getByLabel(genjetTag_, genjets);
504    
505     jets_.ngen = 0;
506     for(unsigned int igen = 0 ; igen < genjets->size(); ++igen){
507     const reco::GenJet & genjet = (*genjets)[igen];
508    
509     float genjet_pt = genjet.pt();
510    
511     // threshold to reduce size of output in minbias PbPb
512 yilmaz 1.9 if(genjet_pt>genPtMin_){
513 yilmaz 1.1
514     jets_.genpt[jets_.ngen] = genjet_pt;
515     jets_.geneta[jets_.ngen] = genjet.eta();
516     jets_.genphi[jets_.ngen] = genjet.phi();
517     jets_.geny[jets_.ngen] = genjet.eta();
518    
519    
520     // find matching patJet if there is one
521    
522     jets_.gendrjt[jets_.ngen] = -1.0;
523     jets_.genmatchindex[jets_.ngen] = -1;
524    
525     for(unsigned int ijet = 0 ; ijet < jets->size(); ++ijet){
526     const pat::Jet& jet = (*jets)[ijet];
527 mnguyen 1.8 const reco::GenJet* matchedGenJet = (*patjets)[ijet].genJet();
528     if(matchedGenJet){
529     // poor man's matching, someone fix please
530     if(fabs(genjet.pt()-matchedGenJet->pt())<0.0001 &&
531     fabs(genjet.eta()-matchedGenJet->eta())<0.0001 &&
532     (fabs(genjet.phi()-matchedGenJet->phi())<0.0001 || fabs(genjet.phi()-matchedGenJet->phi() - 2.0*TMath::Pi()) < 0.0001 )){
533 yilmaz 1.1
534     jets_.genmatchindex[jets_.ngen] = (int)ijet;
535     jets_.gendphijt[jets_.ngen] = reco::deltaPhi(jet.phi(),genjet.phi());
536     jets_.gendrjt[jets_.ngen] = reco::deltaR(jet,genjet);
537    
538     break;
539     }
540     }
541     }
542     jets_.ngen++;
543     }
544    
545     }
546    
547     }
548     t->Fill();
549 mnguyen 1.10 memset(&jets_,0,sizeof jets_);
550 yilmaz 1.1 }
551    
552    
553    
554    
555     //--------------------------------------------------------------------------------------------------
556     void HiInclusiveJetAnalyzer::fillL1Bits(const edm::Event &iEvent)
557     {
558     edm::Handle< L1GlobalTriggerReadoutRecord > L1GlobalTrigger;
559    
560     iEvent.getByLabel(L1gtReadout_, L1GlobalTrigger);
561     const TechnicalTriggerWord& technicalTriggerWordBeforeMask = L1GlobalTrigger->technicalTriggerWord();
562    
563     for (int i=0; i<64;i++)
564     {
565     jets_.l1TBit[i] = technicalTriggerWordBeforeMask.at(i);
566     }
567     jets_.nL1TBit = 64;
568    
569     int ntrigs = L1GlobalTrigger->decisionWord().size();
570     jets_.nL1ABit = ntrigs;
571    
572     for (int i=0; i != ntrigs; i++) {
573     bool accept = L1GlobalTrigger->decisionWord()[i];
574     //jets_.l1ABit[i] = (accept == true)? 1:0;
575     if(accept== true){
576     jets_.l1ABit[i] = 1;
577     }
578     else{
579     jets_.l1ABit[i] = 0;
580     }
581    
582     }
583     }
584    
585     //--------------------------------------------------------------------------------------------------
586     void HiInclusiveJetAnalyzer::fillHLTBits(const edm::Event &iEvent)
587     {
588     // Fill HLT trigger bits.
589     Handle<TriggerResults> triggerResultsHLT;
590     getProduct(hltResName_, triggerResultsHLT, iEvent);
591    
592     const TriggerResults *hltResults = triggerResultsHLT.product();
593     const TriggerNames & triggerNames = iEvent.triggerNames(*hltResults);
594    
595     jets_.nHLTBit = hltTrgNames_.size();
596    
597     for(size_t i=0;i<hltTrgNames_.size();i++){
598    
599     for(size_t j=0;j<triggerNames.size();++j) {
600    
601     if(triggerNames.triggerName(j) == hltTrgNames_[i]){
602    
603     //cout <<"hltTrgNames_(i) "<<hltTrgNames_[i]<<endl;
604     //cout <<"triggerName(j) "<<triggerNames.triggerName(j)<<endl;
605     //cout<<" result "<<triggerResultsHLT->accept(j)<<endl;
606     jets_.hltBit[i] = triggerResultsHLT->accept(j);
607     }
608    
609     }
610     }
611     }
612    
613     //--------------------------------------------------------------------------------------------------
614     template <typename TYPE>
615     inline void HiInclusiveJetAnalyzer::getProduct(const std::string name, edm::Handle<TYPE> &prod,
616     const edm::Event &event) const
617     {
618     // Try to access data collection from EDM file. We check if we really get just one
619     // product with the given name. If not we throw an exception.
620    
621     event.getByLabel(edm::InputTag(name),prod);
622     if (!prod.isValid())
623     throw edm::Exception(edm::errors::Configuration, "HiInclusiveJetAnalyzer::GetProduct()\n")
624     << "Collection with label '" << name << "' is not valid" << std::endl;
625     }
626    
627     //--------------------------------------------------------------------------------------------------
628     template <typename TYPE>
629     inline bool HiInclusiveJetAnalyzer::getProductSafe(const std::string name, edm::Handle<TYPE> &prod,
630     const edm::Event &event) const
631     {
632     // Try to safely access data collection from EDM file. We check if we really get just one
633     // product with the given name. If not, we return false.
634    
635     if (name.size()==0)
636     return false;
637    
638     try {
639     event.getByLabel(edm::InputTag(name),prod);
640     if (!prod.isValid())
641     return false;
642     } catch (...) {
643     return false;
644     }
645     return true;
646     }
647    
648    
649 mnguyen 1.8 int
650     HiInclusiveJetAnalyzer::getPFJetMuon(const pat::Jet& pfJet, const reco::PFCandidateCollection *pfCandidateColl)
651     {
652    
653     int pfMuonIndex = -1;
654     float ptMax = 0.;
655    
656    
657     for(unsigned icand=0;icand<pfCandidateColl->size(); icand++) {
658     const reco::PFCandidate pfCandidate = pfCandidateColl->at(icand);
659    
660     int id = pfCandidate.particleId();
661     if(abs(id) != 3) continue;
662 yilmaz 1.1
663 mnguyen 1.8 if(reco::deltaR(pfJet,pfCandidate)>0.5) continue;
664    
665     double pt = pfCandidate.pt();
666     if(pt>ptMax){
667     ptMax = pt;
668     pfMuonIndex = (int) icand;
669     }
670     }
671    
672     return pfMuonIndex;
673    
674     }
675    
676    
677     double
678     HiInclusiveJetAnalyzer::getPtRel(const reco::PFCandidate lep, const pat::Jet& jet )
679     {
680    
681     float lj_x = jet.p4().px();
682     float lj_y = jet.p4().py();
683     float lj_z = jet.p4().pz();
684    
685     // absolute values squared
686     float lj2 = lj_x*lj_x+lj_y*lj_y+lj_z*lj_z;
687     float lep2 = lep.px()*lep.px()+lep.py()*lep.py()+lep.pz()*lep.pz();
688    
689     // projection vec(mu) to lepjet axis
690     float lepXlj = lep.px()*lj_x+lep.py()*lj_y+lep.pz()*lj_z;
691    
692     // absolute value squared and normalized
693     float pLrel2 = lepXlj*lepXlj/lj2;
694    
695     // lep2 = pTrel2 + pLrel2
696     float pTrel2 = lep2-pLrel2;
697    
698     return (pTrel2 > 0) ? std::sqrt(pTrel2) : 0.0;
699     }
700 mnguyen 1.11
701     // Recursive function, but this version gets called only the first time
702     void
703     HiInclusiveJetAnalyzer::saveDaughters(const reco::GenParticle &gen){
704    
705     for(unsigned i=0;i<gen.numberOfDaughters();i++){
706     const reco::Candidate & daughter = *gen.daughter(i);
707     double daughterPt = daughter.pt();
708     if(daughterPt<1.) continue;
709     double daughterEta = daughter.eta();
710     if(fabs(daughterEta)>3.) continue;
711     int daughterPdgId = daughter.pdgId();
712     int daughterStatus = daughter.status();
713     // Special case when b->b+string, both b and string contain all daughters, so only take the string
714     if(gen.pdgId()==daughterPdgId && gen.status()==3 && daughterStatus==2) continue;
715    
716     // cheesy way of finding strings which were already used
717     if(daughter.pdgId()==92){
718     for(unsigned ist=0;ist<usedStringPts.size();ist++){
719     if(fabs(daughter.pt() - usedStringPts[ist]) < 0.0001) return;
720     }
721     usedStringPts.push_back(daughter.pt());
722     }
723     jets_.bJetIndex[jets_.bMult] = jets_.nref;
724     jets_.bStatus[jets_.bMult] = daughterStatus;
725     jets_.bVx[jets_.bMult] = daughter.vx();
726     jets_.bVy[jets_.bMult] = daughter.vy();
727     jets_.bVz[jets_.bMult] = daughter.vz();
728     jets_.bPt[jets_.bMult] = daughterPt;
729     jets_.bEta[jets_.bMult] = daughterEta;
730     jets_.bPhi[jets_.bMult] = daughter.phi();
731     jets_.bPdg[jets_.bMult] = daughterPdgId;
732     jets_.bChg[jets_.bMult] = daughter.charge();
733     jets_.bMult++;
734     saveDaughters(daughter);
735     }
736     }
737    
738     // This version called for all subsequent calls
739     void
740     HiInclusiveJetAnalyzer::saveDaughters(const reco::Candidate &gen){
741    
742     for(unsigned i=0;i<gen.numberOfDaughters();i++){
743     const reco::Candidate & daughter = *gen.daughter(i);
744     double daughterPt = daughter.pt();
745     if(daughterPt<1.) continue;
746     double daughterEta = daughter.eta();
747     if(fabs(daughterEta)>3.) continue;
748     int daughterPdgId = daughter.pdgId();
749     int daughterStatus = daughter.status();
750     // Special case when b->b+string, both b and string contain all daughters, so only take the string
751     if(gen.pdgId()==daughterPdgId && gen.status()==3 && daughterStatus==2) continue;
752    
753     // cheesy way of finding strings which were already used
754     if(daughter.pdgId()==92){
755     for(unsigned ist=0;ist<usedStringPts.size();ist++){
756     if(fabs(daughter.pt() - usedStringPts[ist]) < 0.0001) return;
757     }
758     usedStringPts.push_back(daughter.pt());
759     }
760    
761     jets_.bJetIndex[jets_.bMult] = jets_.nref;
762     jets_.bStatus[jets_.bMult] = daughterStatus;
763     jets_.bVx[jets_.bMult] = daughter.vx();
764     jets_.bVy[jets_.bMult] = daughter.vy();
765     jets_.bVz[jets_.bMult] = daughter.vz();
766     jets_.bPt[jets_.bMult] = daughterPt;
767     jets_.bEta[jets_.bMult] = daughterEta;
768     jets_.bPhi[jets_.bMult] = daughter.phi();
769     jets_.bPdg[jets_.bMult] = daughterPdgId;
770     jets_.bChg[jets_.bMult] = daughter.charge();
771     jets_.bMult++;
772     saveDaughters(daughter);
773     }
774     }
775 yilmaz 1.1
776     DEFINE_FWK_MODULE(HiInclusiveJetAnalyzer);