ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/CmsHi/JetAnalysis/src/HiInclusiveJetAnalyzer.cc
Revision: 1.27
Committed: Tue Jan 15 14:04:25 2013 UTC (12 years, 3 months ago) by yilmaz
Content type: text/plain
Branch: MAIN
CVS Tags: HiForest_V02_53, HiForest_V02_52, HiForest_V02_51, HiForest_V02_50, HiForest_V02_49, HiForest_V02_48, HiForest_V02_47, HiForest_V02_46, HiForest_V02_45, HiForest_V02_44, HiForest_V02_43, HiForest_V02_42, HiForest_V02_41, HiForest_V02_40, HiForest_V02_39, HiForest_V02_38, HiForest_V02_37, HiForest_V02_36, HiForest_V02_35, HiForest_V02_34
Changes since 1.26: +6 -3 lines
Log Message:
up

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 yilmaz 1.25 #include "Geometry/Records/interface/CaloGeometryRecord.h"
23 yilmaz 1.1
24     #include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
25     #include "DataFormats/JetReco/interface/CaloJetCollection.h"
26     #include "SimDataFormats/HiGenData/interface/GenHIEvent.h"
27     #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
28 mnguyen 1.10 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
29 yilmaz 1.1 #include "DataFormats/JetReco/interface/GenJetCollection.h"
30     #include "DataFormats/VertexReco/interface/Vertex.h"
31     #include "DataFormats/Math/interface/deltaPhi.h"
32     #include "DataFormats/Math/interface/deltaR.h"
33    
34     #include "DataFormats/Common/interface/TriggerResults.h"
35     #include "FWCore/Common/interface/TriggerNames.h"
36     #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutSetup.h"
37     #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
38     #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerObjectMapRecord.h"
39     #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerObjectMap.h"
40     #include "L1Trigger/GlobalTrigger/interface/L1GlobalTrigger.h"
41    
42 yilmaz 1.25 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
43     #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
44    
45 yilmaz 1.1 using namespace std;
46     using namespace edm;
47     using namespace reco;
48    
49 yilmaz 1.25 HiInclusiveJetAnalyzer::HiInclusiveJetAnalyzer(const edm::ParameterSet& iConfig) :
50     geo(0)
51     {
52 yilmaz 1.1
53     jetTag_ = iConfig.getParameter<InputTag>("jetTag");
54 yilmaz 1.18 matchTag_ = iConfig.getUntrackedParameter<InputTag>("matchTag",jetTag_);
55    
56 mnguyen 1.8 vtxTag_ = iConfig.getUntrackedParameter<edm::InputTag>("vtxTag",edm::InputTag("hiSelectedVertex"));
57 yilmaz 1.17 trackTag_ = iConfig.getParameter<InputTag>("trackTag");
58 yilmaz 1.22 useQuality_ = iConfig.getUntrackedParameter<bool>("useQuality",1);
59     trackQuality_ = iConfig.getUntrackedParameter<string>("trackQuality","highPurity");
60 yilmaz 1.1
61     isMC_ = iConfig.getUntrackedParameter<bool>("isMC",false);
62 yilmaz 1.21 fillGenJets_ = iConfig.getUntrackedParameter<bool>("fillGenJets",false);
63    
64 yilmaz 1.14 doTrigger_ = iConfig.getUntrackedParameter<bool>("doTrigger",false);
65 yilmaz 1.17
66     rParam = iConfig.getParameter<double>("rParam");
67 yilmaz 1.19 hardPtMin_ = iConfig.getUntrackedParameter<double>("hardPtMin",4);
68    
69 yilmaz 1.1 if(isMC_){
70     genjetTag_ = iConfig.getParameter<InputTag>("genjetTag");
71     eventInfoTag_ = iConfig.getParameter<InputTag>("eventInfoTag");
72     }
73     verbose_ = iConfig.getUntrackedParameter<bool>("verbose",false);
74    
75     useCentrality_ = iConfig.getUntrackedParameter<bool>("useCentrality",false);
76 yjlee 1.4 useVtx_ = iConfig.getUntrackedParameter<bool>("useVtx",false);
77 yilmaz 1.1 useJEC_ = iConfig.getUntrackedParameter<bool>("useJEC",true);
78 yilmaz 1.2 usePat_ = iConfig.getUntrackedParameter<bool>("usePAT",true);
79 yilmaz 1.1
80 mnguyen 1.8 doLifeTimeTagging_ = iConfig.getUntrackedParameter<bool>("doLifeTimeTagging",false);
81 mnguyen 1.10 doLifeTimeTaggingExtras_ = iConfig.getUntrackedParameter<bool>("doLifeTimeTaggingExtras",true);
82 mnguyen 1.11 saveBfragments_ = iConfig.getUntrackedParameter<bool>("saveBfragments",false);
83 yilmaz 1.27 skipCorrections_ = iConfig.getUntrackedParameter<bool>("skipCorrections",false);
84 mnguyen 1.11
85 mnguyen 1.8 pfCandidateLabel_ = iConfig.getUntrackedParameter<edm::InputTag>("pfCandidateLabel",edm::InputTag("particleFlowTmp"));
86    
87 yilmaz 1.25 EBSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("EBRecHitSrc",edm::InputTag("ecalRecHit","EcalRecHitsEB"));
88     EESrc_ = iConfig.getUntrackedParameter<edm::InputTag>("EERecHitSrc",edm::InputTag("ecalRecHit","EcalRecHitsEE"));
89     HcalRecHitHFSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("hcalHFRecHitSrc",edm::InputTag("hfreco"));
90     HcalRecHitHBHESrc_ = iConfig.getUntrackedParameter<edm::InputTag>("hcalHBHERecHitSrc",edm::InputTag("hbhereco"));
91    
92 yilmaz 1.26 genParticleSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("genParticles",edm::InputTag("hiGenParticles"));
93    
94 yilmaz 1.14 if(doTrigger_){
95 yilmaz 1.1 L1gtReadout_ = iConfig.getParameter<edm::InputTag>("L1gtReadout");
96     hltResName_ = iConfig.getUntrackedParameter<string>("hltTrgResults","TriggerResults::HLT");
97    
98    
99     if (iConfig.exists("hltTrgNames"))
100     hltTrgNames_ = iConfig.getUntrackedParameter<vector<string> >("hltTrgNames");
101    
102     if (iConfig.exists("hltProcNames"))
103     hltProcNames_ = iConfig.getUntrackedParameter<vector<string> >("hltProcNames");
104     else {
105     hltProcNames_.push_back("FU");
106     hltProcNames_.push_back("HLT");
107     }
108     }
109    
110     cout<<" jet collection : "<<jetTag_<<endl;
111 yilmaz 1.12 doSubEvent_ = 0;
112    
113 yilmaz 1.9 if(isMC_){
114     cout<<" genjet collection : "<<genjetTag_<<endl;
115     genPtMin_ = iConfig.getUntrackedParameter<double>("genPtMin",0);
116 yilmaz 1.12 doSubEvent_ = iConfig.getUntrackedParameter<bool>("doSubEvent",1);
117 yilmaz 1.9 }
118 yilmaz 1.1
119    
120     }
121    
122    
123    
124     HiInclusiveJetAnalyzer::~HiInclusiveJetAnalyzer() { }
125    
126    
127    
128     void
129     HiInclusiveJetAnalyzer::beginRun(const edm::Run& run,
130     const edm::EventSetup & es) {}
131    
132     void
133     HiInclusiveJetAnalyzer::beginJob() {
134    
135     centrality_ = 0;
136    
137     //string jetTagName = jetTag_.label()+"_tree";
138     string jetTagTitle = jetTag_.label()+" Jet Analysis Tree";
139     t = fs1->make<TTree>("t",jetTagTitle.c_str());
140    
141 yjlee 1.6 // TTree* t= new TTree("t","Jet Response Analyzer");
142 yjlee 1.4 //t->Branch("run",&jets_.run,"run/I");
143 yilmaz 1.1 t->Branch("evt",&jets_.evt,"evt/I");
144 yjlee 1.4 //t->Branch("lumi",&jets_.lumi,"lumi/I");
145 yilmaz 1.1 t->Branch("b",&jets_.b,"b/F");
146 yjlee 1.4 if (useVtx_) {
147     t->Branch("vx",&jets_.vx,"vx/F");
148     t->Branch("vy",&jets_.vy,"vy/F");
149     t->Branch("vz",&jets_.vz,"vz/F");
150     }
151     if (useCentrality_) {
152     t->Branch("hf",&jets_.hf,"hf/F");
153     t->Branch("bin",&jets_.bin,"bin/I");
154     }
155 yjlee 1.7
156     t->Branch("nref",&jets_.nref,"nref/I");
157 yilmaz 1.1 t->Branch("rawpt",jets_.rawpt,"rawpt[nref]/F");
158 yilmaz 1.27 if(!skipCorrections_) t->Branch("jtpt",jets_.jtpt,"jtpt[nref]/F");
159 yilmaz 1.1 t->Branch("jteta",jets_.jteta,"jteta[nref]/F");
160     t->Branch("jty",jets_.jty,"jty[nref]/F");
161     t->Branch("jtphi",jets_.jtphi,"jtphi[nref]/F");
162 yilmaz 1.2 t->Branch("jtpu",jets_.jtpu,"jtpu[nref]/F");
163 yilmaz 1.23 t->Branch("jtm",jets_.jtm,"jtm[nref]/F");
164 yilmaz 1.1
165 yilmaz 1.17 // jet ID information, jet composition
166     t->Branch("discr_fr01", jets_.discr_fr01,"discr_fr01[nref]/F");
167    
168     t->Branch("trackMax", jets_.trackMax,"trackMax[nref]/F");
169     t->Branch("trackSum", jets_.trackSum,"trackSum[nref]/F");
170     t->Branch("trackN", jets_.trackN,"trackN[nref]/I");
171 yilmaz 1.19 t->Branch("trackHardSum", jets_.trackHardSum,"trackHardSum[nref]/F");
172     t->Branch("trackHardN", jets_.trackHardN,"trackHardN[nref]/I");
173 yilmaz 1.17
174     t->Branch("chargedMax", jets_.chargedMax,"chargedMax[nref]/F");
175     t->Branch("chargedSum", jets_.chargedSum,"chargedSum[nref]/F");
176     t->Branch("chargedN", jets_.chargedN,"chargedN[nref]/I");
177 yilmaz 1.19 t->Branch("chargedHardSum", jets_.chargedHardSum,"chargedHardSum[nref]/F");
178     t->Branch("chargedHardN", jets_.chargedHardN,"chargedHardN[nref]/I");
179 yilmaz 1.17
180     t->Branch("photonMax", jets_.photonMax,"photonMax[nref]/F");
181     t->Branch("photonSum", jets_.photonSum,"photonSum[nref]/F");
182     t->Branch("photonN", jets_.photonN,"photonN[nref]/I");
183 yilmaz 1.19 t->Branch("photonHardSum", jets_.photonHardSum,"photonHardSum[nref]/F");
184     t->Branch("photonHardN", jets_.photonHardN,"photonHardN[nref]/I");
185 yilmaz 1.17
186     t->Branch("neutralMax", jets_.neutralMax,"neutralMax[nref]/F");
187     t->Branch("neutralSum", jets_.neutralSum,"neutralSum[nref]/F");
188     t->Branch("neutralN", jets_.neutralN,"neutralN[nref]/I");
189    
190 yilmaz 1.25 t->Branch("hcalSum", jets_.hcalSum,"hcalSum[nref]/F");
191     t->Branch("ecalSum", jets_.ecalSum,"ecalSum[nref]/F");
192    
193 yilmaz 1.17 t->Branch("eMax", jets_.eMax,"eMax[nref]/F");
194     t->Branch("eSum", jets_.eSum,"eSum[nref]/F");
195     t->Branch("eN", jets_.eN,"eN[nref]/I");
196    
197     t->Branch("muMax", jets_.muMax,"muMax[nref]/F");
198     t->Branch("muSum", jets_.muSum,"muSum[nref]/F");
199     t->Branch("muN", jets_.muN,"muN[nref]/I");
200    
201 yilmaz 1.27 if(!skipCorrections_) t->Branch("matchedPt", jets_.matchedPt,"matchedPt[nref]/F");
202     t->Branch("matchedRawPt", jets_.matchedRawPt,"matchedRawPt[nref]/F");
203 yilmaz 1.18 t->Branch("matchedR", jets_.matchedR,"matchedR[nref]/F");
204 yilmaz 1.17
205 mnguyen 1.8 // b-jet discriminators
206     if (doLifeTimeTagging_) {
207    
208     t->Branch("discr_csvMva",jets_.discr_csvMva,"discr_csvMva[nref]/F");
209     t->Branch("discr_csvSimple",jets_.discr_csvSimple,"discr_csvSimple[nref]/F");
210     t->Branch("discr_muByIp3",jets_.discr_muByIp3,"discr_muByIp3[nref]/F");
211     t->Branch("discr_muByPt",jets_.discr_muByPt,"discr_muByPt[nref]/F");
212     t->Branch("discr_prob",jets_.discr_prob,"discr_prob[nref]/F");
213     t->Branch("discr_probb",jets_.discr_probb,"discr_probb[nref]/F");
214     t->Branch("discr_tcHighEff",jets_.discr_tcHighEff,"discr_tcHighEff[nref]/F");
215     t->Branch("discr_tcHighPur",jets_.discr_tcHighPur,"discr_tcHighPur[nref]/F");
216    
217 mnguyen 1.10 t->Branch("nsvtx", jets_.nsvtx, "nsvtx[nref]/I");
218     t->Branch("svtxntrk", jets_.svtxntrk, "svtxntrk[nref]/I");
219 mnguyen 1.8 t->Branch("svtxdl", jets_.svtxdl, "svtxdl[nref]/F");
220     t->Branch("svtxdls", jets_.svtxdls, "svtxdls[nref]/F");
221     t->Branch("svtxm", jets_.svtxm, "svtxm[nref]/F");
222     t->Branch("svtxpt", jets_.svtxpt, "svtxpt[nref]/F");
223    
224 mnguyen 1.10 t->Branch("nIPtrk",jets_.nIPtrk,"nIPtrk[nref]/I");
225     t->Branch("nselIPtrk",jets_.nselIPtrk,"nselIPtrk[nref]/I");
226 mnguyen 1.11
227 mnguyen 1.10 if (doLifeTimeTaggingExtras_) {
228 mnguyen 1.11 t->Branch("nIP",&jets_.nIP,"nIP/I");
229     t->Branch("ipJetIndex",jets_.ipJetIndex,"ipJetIndex[nIP]/I");
230     t->Branch("ipPt",jets_.ipPt,"ipPt[nIP]/F");
231     t->Branch("ipProb0",jets_.ipProb0,"ipProb0[nIP]/F");
232     t->Branch("ipProb1",jets_.ipProb1,"ipProb1[nIP]/F");
233     t->Branch("ip2d",jets_.ip2d,"ip2d[nIP]/F");
234     t->Branch("ip2dSig",jets_.ip2dSig,"ip2dSig[nIP]/F");
235     t->Branch("ip3d",jets_.ip3d,"ip3d[nIP]/F");
236     t->Branch("ip3dSig",jets_.ip3dSig,"ip3dSig[nIP]/F");
237     t->Branch("ipDist2Jet",jets_.ipDist2Jet,"ipDist2Jet[nIP]/F");
238     t->Branch("ipDist2JetSig",jets_.ipDist2JetSig,"ipDist2JetSig[nIP]/F");
239     t->Branch("ipClosest2Jet",jets_.ipClosest2Jet,"ipClosest2Jet[nIP]/F");
240    
241 mnguyen 1.10 }
242 mnguyen 1.8
243     t->Branch("mue", jets_.mue, "mue[nref]/F");
244     t->Branch("mupt", jets_.mupt, "mupt[nref]/F");
245     t->Branch("mueta", jets_.mueta, "mueta[nref]/F");
246     t->Branch("muphi", jets_.muphi, "muphi[nref]/F");
247     t->Branch("mudr", jets_.mudr, "mudr[nref]/F");
248     t->Branch("muptrel", jets_.muptrel, "muptrel[nref]/F");
249     t->Branch("muchg", jets_.muchg, "muchg[nref]/I");
250     }
251 ylai 1.16
252 mnguyen 1.8
253 yilmaz 1.1 if(isMC_){
254 mnguyen 1.10 t->Branch("beamId1",&jets_.beamId1,"beamId1/I");
255     t->Branch("beamId2",&jets_.beamId2,"beamId2/I");
256    
257 yilmaz 1.1 t->Branch("pthat",&jets_.pthat,"pthat/F");
258    
259     // Only matched gen jets
260     t->Branch("refpt",jets_.refpt,"refpt[nref]/F");
261     t->Branch("refeta",jets_.refeta,"refeta[nref]/F");
262     t->Branch("refy",jets_.refy,"refy[nref]/F");
263     t->Branch("refphi",jets_.refphi,"refphi[nref]/F");
264     t->Branch("refdphijt",jets_.refdphijt,"refdphijt[nref]/F");
265     t->Branch("refdrjt",jets_.refdrjt,"refdrjt[nref]/F");
266     // matched parton
267     t->Branch("refparton_pt",jets_.refparton_pt,"refparton_pt[nref]/F");
268 mnguyen 1.8 t->Branch("refparton_flavor",jets_.refparton_flavor,"refparton_flavor[nref]/I");
269     t->Branch("refparton_flavorForB",jets_.refparton_flavorForB,"refparton_flavorForB[nref]/I");
270 yilmaz 1.1
271 yilmaz 1.25 t->Branch("genChargedSum", jets_.genChargedSum,"genChargedSum[nref]/F");
272     t->Branch("genHardSum", jets_.genHardSum,"genHardSum[nref]/F");
273     t->Branch("signalChargedSum", jets_.signalChargedSum,"signalChargedSum[nref]/F");
274     t->Branch("signalHardSum", jets_.signalHardSum,"signalHardSum[nref]/F");
275    
276     if(doSubEvent_){
277     t->Branch("subid",jets_.subid,"subid[nref]/I");
278     }
279    
280 yilmaz 1.21 if(fillGenJets_){
281     // For all gen jets, matched or unmatched
282     t->Branch("ngen",&jets_.ngen,"ngen/I");
283     t->Branch("genmatchindex",jets_.genmatchindex,"genmatchindex[ngen]/I");
284     t->Branch("genpt",jets_.genpt,"genpt[ngen]/F");
285     t->Branch("geneta",jets_.geneta,"geneta[ngen]/F");
286     t->Branch("geny",jets_.geny,"geny[ngen]/F");
287     t->Branch("genphi",jets_.genphi,"genphi[ngen]/F");
288     t->Branch("gendphijt",jets_.gendphijt,"gendphijt[ngen]/F");
289     t->Branch("gendrjt",jets_.gendrjt,"gendrjt[ngen]/F");
290    
291     if(doSubEvent_){
292     t->Branch("gensubid",jets_.gensubid,"gensubid[ngen]/I");
293     }
294 yilmaz 1.12 }
295 yilmaz 1.21
296 mnguyen 1.11 if(saveBfragments_ ) {
297     t->Branch("bMult",&jets_.bMult,"bMult/I");
298     t->Branch("bJetIndex",jets_.bJetIndex,"bJetIndex[bMult]/I");
299     t->Branch("bStatus",jets_.bStatus,"bStatus[bMult]/I");
300     t->Branch("bVx",jets_.bVx,"bVx[bMult]/F");
301     t->Branch("bVy",jets_.bVy,"bVy[bMult]/F");
302     t->Branch("bVz",jets_.bVz,"bVz[bMult]/F");
303     t->Branch("bPt",jets_.bPt,"bPt[bMult]/F");
304     t->Branch("bEta",jets_.bEta,"bEta[bMult]/F");
305     t->Branch("bPhi",jets_.bPhi,"bPhi[bMult]/F");
306     t->Branch("bPdg",jets_.bPdg,"bPdg[bMult]/I");
307     t->Branch("bChg",jets_.bChg,"bChg[bMult]/I");
308     }
309 yilmaz 1.13
310 yilmaz 1.1 }
311 yjlee 1.7 /*
312 yilmaz 1.1 if(!isMC_){
313     t->Branch("nL1TBit",&jets_.nL1TBit,"nL1TBit/I");
314     t->Branch("l1TBit",jets_.l1TBit,"l1TBit[nL1TBit]/O");
315    
316     t->Branch("nL1ABit",&jets_.nL1ABit,"nL1ABit/I");
317     t->Branch("l1ABit",jets_.l1ABit,"l1ABit[nL1ABit]/O");
318    
319     t->Branch("nHLTBit",&jets_.nHLTBit,"nHLTBit/I");
320     t->Branch("hltBit",jets_.hltBit,"hltBit[nHLTBit]/O");
321    
322     }
323 yjlee 1.7 */
324 yilmaz 1.1 TH1D::SetDefaultSumw2();
325    
326    
327     }
328    
329    
330     void
331     HiInclusiveJetAnalyzer::analyze(const Event& iEvent,
332     const EventSetup& iSetup) {
333    
334     int event = iEvent.id().event();
335     int run = iEvent.id().run();
336     int lumi = iEvent.id().luminosityBlock();
337    
338     jets_.run = run;
339     jets_.evt = event;
340     jets_.lumi = lumi;
341    
342     LogDebug("HiInclusiveJetAnalyzer")<<"START event: "<<event<<" in run "<<run<<endl;
343    
344     int bin = -1;
345     double hf = 0.;
346     double b = 999.;
347    
348 yilmaz 1.25 if(!geo){
349     edm::ESHandle<CaloGeometry> pGeo;
350     iSetup.get<CaloGeometryRecord>().get(pGeo);
351     geo = pGeo.product();
352     }
353 yilmaz 1.1 if(useCentrality_){
354     if(!centrality_) centrality_ = new CentralityProvider(iSetup);
355     centrality_->newEvent(iEvent,iSetup); // make sure you do this first in every event
356     //double c = centrality_->centralityValue();
357     const reco::Centrality *cent = centrality_->raw();
358    
359     hf = cent->EtHFhitSum();
360    
361     bin = centrality_->getBin();
362     b = centrality_->bMean();
363     }
364    
365     // loop the events
366    
367     jets_.bin = bin;
368     jets_.hf = hf;
369    
370 yilmaz 1.25 reco::Vertex::Point vtx(0,0,0);
371     if (useVtx_) {
372     edm::Handle<vector<reco::Vertex> >vertex;
373     iEvent.getByLabel(vtxTag_, vertex);
374    
375     if(vertex->size()>0) {
376     jets_.vx=vertex->begin()->x();
377     jets_.vy=vertex->begin()->y();
378     jets_.vz=vertex->begin()->z();
379     vtx = vertex->begin()->position();
380     }
381     }
382 yilmaz 1.1
383 yilmaz 1.2 edm::Handle<pat::JetCollection> patjets;
384     if(usePat_)iEvent.getByLabel(jetTag_, patjets);
385 yilmaz 1.18
386     edm::Handle<pat::JetCollection> matchedjets;
387     iEvent.getByLabel(matchTag_, matchedjets);
388 yilmaz 1.2
389     edm::Handle<reco::JetView> jets;
390 yilmaz 1.1 iEvent.getByLabel(jetTag_, jets);
391    
392 mnguyen 1.8 edm::Handle<reco::PFCandidateCollection> pfCandidates;
393 yilmaz 1.17 iEvent.getByLabel(pfCandidateLabel_,pfCandidates);
394    
395     edm::Handle<reco::TrackCollection> tracks;
396     iEvent.getByLabel(trackTag_,tracks);
397    
398 yilmaz 1.25 edm::Handle<edm::SortedCollection<EcalRecHit,edm::StrictWeakOrdering<EcalRecHit> > > ebHits;
399     edm::Handle<edm::SortedCollection<EcalRecHit,edm::StrictWeakOrdering<EcalRecHit> > > eeHits;
400    
401     edm::Handle<HFRecHitCollection> hfHits;
402     edm::Handle<HBHERecHitCollection> hbheHits;
403    
404     iEvent.getByLabel(HcalRecHitHBHESrc_,hbheHits);
405     iEvent.getByLabel(HcalRecHitHFSrc_,hfHits);
406     iEvent.getByLabel(EBSrc_,ebHits);
407     iEvent.getByLabel(EESrc_,eeHits);
408    
409 yilmaz 1.26 edm::Handle<reco::GenParticleCollection> genparts;
410     iEvent.getByLabel(genParticleSrc_,genparts);
411    
412 yilmaz 1.25
413 yilmaz 1.1 // FILL JRA TREE
414     jets_.b = b;
415     jets_.nref = 0;
416 yilmaz 1.24
417 yilmaz 1.1
418 yilmaz 1.15 if(doTrigger_){
419 yilmaz 1.1 fillL1Bits(iEvent);
420     fillHLTBits(iEvent);
421     }
422    
423 yilmaz 1.24 for(unsigned int j = 0; j < jets->size(); ++j){
424     const reco::Jet& jet = (*jets)[j];
425 yilmaz 1.2 if (useJEC_ && usePat_){
426     jets_.rawpt[jets_.nref]=(*patjets)[j].correctedJet("Uncorrected").pt();
427     }
428 mnguyen 1.8
429     if(doLifeTimeTagging_){
430    
431     if(jetTag_.label()=="icPu5patJets"){
432     jets_.discr_csvMva[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5CombinedSecondaryVertexMVABJetTags");
433     jets_.discr_csvSimple[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5CombinedSecondaryVertexBJetTags");
434     jets_.discr_muByIp3[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5SoftMuonByIP3dBJetTags");
435     jets_.discr_muByPt[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5SoftMuonByPtBJetTags");
436     jets_.discr_prob[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5JetBProbabilityBJetTags");
437     jets_.discr_probb[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5JetProbabilityBJetTags");
438     jets_.discr_tcHighEff[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5TrackCountingHighEffBJetTags");
439     jets_.discr_tcHighPur[jets_.nref] = (*patjets)[j].bDiscriminator("icPu5TrackCountingHighPurBJetTags");
440     }
441     else if(jetTag_.label()=="akPu3PFpatJets"){
442     jets_.discr_csvMva[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFCombinedSecondaryVertexMVABJetTags");
443     jets_.discr_csvSimple[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFCombinedSecondaryVertexBJetTags");
444     jets_.discr_muByIp3[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFSoftMuonByIP3dBJetTags");
445     jets_.discr_muByPt[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFSoftMuonByPtBJetTags");
446     jets_.discr_prob[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFJetBProbabilityBJetTags");
447     jets_.discr_probb[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFJetProbabilityBJetTags");
448     jets_.discr_tcHighEff[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFTrackCountingHighEffBJetTags");
449     jets_.discr_tcHighPur[jets_.nref] = (*patjets)[j].bDiscriminator("akPu3PFTrackCountingHighPurBJetTags");
450     }
451     else{
452 mnguyen 1.10 cout<<" b-tagging variables not filled for this collection, turn of doLifeTimeTagging "<<endl;
453 mnguyen 1.8 }
454 mnguyen 1.10
455     const reco::SecondaryVertexTagInfo& tagInfoSV=*(*patjets)[j].tagInfoSecondaryVertex();
456    
457     jets_.nsvtx[jets_.nref] = tagInfoSV.nVertices();
458 mnguyen 1.8
459     if (tagInfoSV.nVertices()>0) {
460     jets_.svtxntrk[jets_.nref] = tagInfoSV.nVertexTracks(0);
461 mnguyen 1.10 // this is the 3d flight distance, for 2-D use (0,true)
462     Measurement1D m1D = tagInfoSV.flightDistance(0);
463 mnguyen 1.8 jets_.svtxdl[jets_.nref] = m1D.value();
464     jets_.svtxdls[jets_.nref] = m1D.significance();
465    
466 mnguyen 1.10 const reco::Vertex& svtx = tagInfoSV.secondaryVertex(0);
467 mnguyen 1.11 //cout<<" SV: vx: "<<svtx.x()<<" vy "<<svtx.y()<<" vz "<<svtx.z()<<endl;
468     //cout<<" PV: vx: "<<jet.vx()<<" vy "<<jet.vy()<<" vz "<<jet.vz()<<endl;
469 mnguyen 1.8 jets_.svtxm[jets_.nref] = svtx.p4().mass();
470 mnguyen 1.10 jets_.svtxpt[jets_.nref] = svtx.p4().pt();
471 mnguyen 1.11 //cout<<" chi2 "<<svtx.chi2()<<" ndof "<<svtx.ndof()<<endl;
472 mnguyen 1.8 }
473 mnguyen 1.10
474 mnguyen 1.8 const reco::TrackIPTagInfo& tagInfoIP=*(*patjets)[j].tagInfoTrackIP();
475 mnguyen 1.10
476     jets_.nIPtrk[jets_.nref] = tagInfoIP.tracks().size();
477     jets_.nselIPtrk[jets_.nref] = tagInfoIP.selectedTracks().size();
478 mnguyen 1.8
479 mnguyen 1.10 if (doLifeTimeTaggingExtras_) {
480 mnguyen 1.8
481 mnguyen 1.10 TrackRefVector selTracks=tagInfoIP.selectedTracks();
482    
483     GlobalPoint pv(tagInfoIP.primaryVertex()->position().x(),tagInfoIP.primaryVertex()->position().y(),tagInfoIP.primaryVertex()->position().z());
484    
485     for(int it=0;it<jets_.nselIPtrk[jets_.nref] ;it++)
486     {
487 mnguyen 1.11 jets_.ipJetIndex[jets_.nIP + it]= jets_.nref;
488 mnguyen 1.10 TrackIPTagInfo::TrackIPData data = tagInfoIP.impactParameterData()[it];
489 mnguyen 1.11 jets_.ipPt[jets_.nIP + it] = selTracks[it]->pt();
490     jets_.ipProb0[jets_.nIP + it] = tagInfoIP.probabilities(0)[it];
491     jets_.ipProb1[jets_.nIP + it] = tagInfoIP.probabilities(1)[it];
492     jets_.ip2d[jets_.nIP + it] = data.ip2d.value();
493     jets_.ip2dSig[jets_.nIP + it] = data.ip2d.significance();
494     jets_.ip3d[jets_.nIP + it] = data.ip3d.value();
495     jets_.ip3dSig[jets_.nIP + it] = data.ip3d.significance();
496     jets_.ipDist2Jet[jets_.nIP + it] = data.distanceToJetAxis.value();
497     jets_.ipDist2JetSig[jets_.nIP + it] = data.distanceToJetAxis.significance();
498     jets_.ipClosest2Jet[jets_.nIP + it] = (data.closestToJetAxis - pv).mag(); //decay length
499 mnguyen 1.10 }
500 mnguyen 1.11
501     jets_.nIP += jets_.nselIPtrk[jets_.nref];
502    
503 mnguyen 1.10 }
504    
505 mnguyen 1.8 const reco::PFCandidateCollection *pfCandidateColl = &(*pfCandidates);
506     int pfMuonIndex = getPFJetMuon(jet, pfCandidateColl);
507     if(pfMuonIndex >=0){
508     const reco::PFCandidate muon = pfCandidateColl->at(pfMuonIndex);
509     jets_.mupt[jets_.nref] = muon.pt();
510     jets_.mueta[jets_.nref] = muon.eta();
511     jets_.muphi[jets_.nref] = muon.phi();
512     jets_.mue[jets_.nref] = muon.energy();
513     jets_.mudr[jets_.nref] = reco::deltaR(jet,muon);
514     jets_.muptrel[jets_.nref] = getPtRel(muon, jet);
515     jets_.muchg[jets_.nref] = muon.charge();
516 yilmaz 1.17 }else{
517 ylai 1.16 jets_.mupt[jets_.nref] = 0.0;
518     jets_.mueta[jets_.nref] = 0.0;
519     jets_.muphi[jets_.nref] = 0.0;
520     jets_.mue[jets_.nref] = 0.0;
521     jets_.mudr[jets_.nref] = 9.9;
522     jets_.muptrel[jets_.nref] = 0.0;
523     jets_.muchg[jets_.nref] = 0;
524     }
525    
526 mnguyen 1.8 }
527    
528 yilmaz 1.17 // Jet ID variables
529    
530     jets_.muMax[jets_.nref] = 0;
531     jets_.muSum[jets_.nref] = 0;
532     jets_.muN[jets_.nref] = 0;
533    
534     jets_.eMax[jets_.nref] = 0;
535     jets_.eSum[jets_.nref] = 0;
536     jets_.eN[jets_.nref] = 0;
537    
538     jets_.neutralMax[jets_.nref] = 0;
539     jets_.neutralSum[jets_.nref] = 0;
540     jets_.neutralN[jets_.nref] = 0;
541    
542     jets_.photonMax[jets_.nref] = 0;
543     jets_.photonSum[jets_.nref] = 0;
544     jets_.photonN[jets_.nref] = 0;
545 yilmaz 1.19 jets_.photonHardSum[jets_.nref] = 0;
546     jets_.photonHardN[jets_.nref] = 0;
547 yilmaz 1.17
548     jets_.chargedMax[jets_.nref] = 0;
549     jets_.chargedSum[jets_.nref] = 0;
550     jets_.chargedN[jets_.nref] = 0;
551 yilmaz 1.19 jets_.chargedHardSum[jets_.nref] = 0;
552     jets_.chargedHardN[jets_.nref] = 0;
553 yilmaz 1.17
554     jets_.trackMax[jets_.nref] = 0;
555     jets_.trackSum[jets_.nref] = 0;
556     jets_.trackN[jets_.nref] = 0;
557 yilmaz 1.19 jets_.trackHardSum[jets_.nref] = 0;
558     jets_.trackHardN[jets_.nref] = 0;
559 yilmaz 1.17
560 yilmaz 1.25 jets_.hcalSum[jets_.nref] = 0;
561     jets_.ecalSum[jets_.nref] = 0;
562    
563     jets_.genChargedSum[jets_.nref] = 0;
564     jets_.genHardSum[jets_.nref] = 0;
565    
566     jets_.signalChargedSum[jets_.nref] = 0;
567     jets_.signalHardSum[jets_.nref] = 0;
568    
569     jets_.subid[jets_.nref] = -1;
570    
571 yilmaz 1.17 for(unsigned int icand = 0; icand < tracks->size(); ++icand){
572     const reco::Track& track = (*tracks)[icand];
573 yilmaz 1.22 if(useQuality_ ){
574     bool goodtrack = track.quality(reco::TrackBase::qualityByName(trackQuality_));
575     if(!goodtrack) continue;
576     }
577 yilmaz 1.18
578 yilmaz 1.17 double dr = deltaR(jet,track);
579     if(dr < rParam){
580     double ptcand = track.pt();
581     jets_.trackSum[jets_.nref] += ptcand;
582     jets_.trackN[jets_.nref] += 1;
583 yilmaz 1.19
584     if(ptcand > hardPtMin_){
585     jets_.trackHardSum[jets_.nref] += ptcand;
586     jets_.trackHardN[jets_.nref] += 1;
587    
588     }
589 yilmaz 1.17 if(ptcand > jets_.trackMax[jets_.nref]) jets_.trackMax[jets_.nref] = ptcand;
590    
591     }
592     }
593    
594     for(unsigned int icand = 0; icand < pfCandidates->size(); ++icand){
595     const reco::PFCandidate& track = (*pfCandidates)[icand];
596     double dr = deltaR(jet,track);
597     if(dr < rParam){
598     double ptcand = track.pt();
599     int pfid = track.particleId();
600    
601     switch(pfid){
602    
603     case 1:
604     jets_.chargedSum[jets_.nref] += ptcand;
605     jets_.chargedN[jets_.nref] += 1;
606 yilmaz 1.19 if(ptcand > hardPtMin_){
607     jets_.chargedHardSum[jets_.nref] += ptcand;
608     jets_.chargedHardN[jets_.nref] += 1;
609     }
610 yilmaz 1.17 if(ptcand > jets_.chargedMax[jets_.nref]) jets_.chargedMax[jets_.nref] = ptcand;
611     break;
612    
613     case 2:
614     jets_.eSum[jets_.nref] += ptcand;
615     jets_.eN[jets_.nref] += 1;
616     if(ptcand > jets_.eMax[jets_.nref]) jets_.eMax[jets_.nref] = ptcand;
617     break;
618    
619     case 3:
620     jets_.muSum[jets_.nref] += ptcand;
621     jets_.muN[jets_.nref] += 1;
622     if(ptcand > jets_.muMax[jets_.nref]) jets_.muMax[jets_.nref] = ptcand;
623     break;
624    
625     case 4:
626     jets_.photonSum[jets_.nref] += ptcand;
627     jets_.photonN[jets_.nref] += 1;
628 yilmaz 1.19 if(ptcand > hardPtMin_){
629     jets_.photonHardSum[jets_.nref] += ptcand;
630     jets_.photonHardN[jets_.nref] += 1;
631     }
632 yilmaz 1.17 if(ptcand > jets_.photonMax[jets_.nref]) jets_.photonMax[jets_.nref] = ptcand;
633     break;
634    
635     case 5:
636     jets_.neutralSum[jets_.nref] += ptcand;
637     jets_.neutralN[jets_.nref] += 1;
638     if(ptcand > jets_.neutralMax[jets_.nref]) jets_.neutralMax[jets_.nref] = ptcand;
639     break;
640    
641 mnguyen 1.20 default:
642     break;
643    
644 yilmaz 1.17 }
645     }
646     }
647    
648 yilmaz 1.25 // Calorimeter fractions
649    
650     for(unsigned int i = 0; i < hbheHits->size(); ++i){
651     const HBHERecHit & hit= (*hbheHits)[i];
652     math::XYZPoint pos = getPosition(hit.id(),vtx);
653     double dr = deltaR(jet.eta(),jet.phi(),pos.eta(),pos.phi());
654     if(dr < rParam){
655     jets_.hcalSum[jets_.nref] += getEt(pos,hit.energy());
656     }
657     }
658    
659     for(unsigned int i = 0; i < hfHits->size(); ++i){
660     const HFRecHit & hit= (*hfHits)[i];
661     math::XYZPoint pos = getPosition(hit.id(),vtx);
662     double dr = deltaR(jet.eta(),jet.phi(),pos.eta(),pos.phi());
663     if(dr < rParam){
664     jets_.hcalSum[jets_.nref] += getEt(pos,hit.energy());
665     }
666     }
667    
668    
669     for(unsigned int i = 0; i < ebHits->size(); ++i){
670     const EcalRecHit & hit= (*ebHits)[i];
671     math::XYZPoint pos = getPosition(hit.id(),vtx);
672     double dr = deltaR(jet.eta(),jet.phi(),pos.eta(),pos.phi());
673     if(dr < rParam){
674     jets_.ecalSum[jets_.nref] += getEt(pos,hit.energy());
675     }
676     }
677    
678     for(unsigned int i = 0; i < eeHits->size(); ++i){
679     const EcalRecHit & hit= (*eeHits)[i];
680     math::XYZPoint pos = getPosition(hit.id(),vtx);
681     double dr = deltaR(jet.eta(),jet.phi(),pos.eta(),pos.phi());
682     if(dr < rParam){
683     jets_.ecalSum[jets_.nref] += getEt(pos,hit.energy());
684     }
685     }
686    
687     // Alternative reconstruction matching (PF for calo, calo for PF)
688    
689 yilmaz 1.18 double drMin = 100;
690 mnguyen 1.20 for(unsigned int imatch = 0 ; imatch < matchedjets->size(); ++imatch){
691 yilmaz 1.27 const pat::Jet& mjet = (*matchedjets)[imatch];
692 yilmaz 1.18
693     double dr = deltaR(jet,mjet);
694     if(dr < drMin){
695     jets_.matchedPt[jets_.nref] = mjet.pt();
696 yilmaz 1.27 jets_.matchedRawPt[jets_.nref] = mjet.correctedJet("Uncorrected").pt();
697 yilmaz 1.18 jets_.matchedR[jets_.nref] = dr;
698 yilmaz 1.22 drMin = dr;
699 yilmaz 1.18 }
700     }
701 mnguyen 1.20
702 yilmaz 1.17 // if(etrk.quality(reco::TrackBase::qualityByName(qualityString_))) pev_.trkQual[pev_.nTrk]=1;
703    
704 ylai 1.16 /////////////////////////////////////////////////////////////////
705     // Jet core pt^2 discriminant for fake jets
706     // Edited by Yue Shi Lai <ylai@mit.edu>
707    
708     // Initial value is 0
709     jets_.discr_fr01[jets_.nref] = 0;
710     // Start with no directional adaption, i.e. the fake rejection
711     // axis is the jet axis
712     float pseudorapidity_adapt = jets_.jteta[jets_.nref];
713     float azimuth_adapt = jets_.jtphi[jets_.nref];
714    
715     // Unadapted discriminant with adaption search
716     for (size_t iteration = 0; iteration < 2; iteration++) {
717     float pseudorapidity_adapt_new = pseudorapidity_adapt;
718     float azimuth_adapt_new = azimuth_adapt;
719     float max_weighted_perp = 0;
720     float perp_square_sum = 0;
721    
722     for (size_t index_pf_candidate = 0;
723     index_pf_candidate < pfCandidates->size();
724     index_pf_candidate++) {
725     const reco::PFCandidate &p =
726     (*pfCandidates)[index_pf_candidate];
727    
728     switch (p.particleId()) {
729 mnguyen 1.20 //case 1: // Charged hadron
730     //case 3: // Muon
731 ylai 1.16 case 4: // Photon
732     {
733     const float dpseudorapidity =
734     p.eta() - pseudorapidity_adapt;
735     const float dazimuth =
736     reco::deltaPhi(p.phi(), azimuth_adapt);
737     // The Gaussian scale factor is 0.5 / (0.1 * 0.1)
738     // = 50
739     const float angular_weight =
740     exp(-50.0F * (dpseudorapidity * dpseudorapidity +
741     dazimuth * dazimuth));
742     const float weighted_perp =
743     angular_weight * p.pt() * p.pt();
744     const float weighted_perp_square =
745     weighted_perp * p.pt();
746    
747     perp_square_sum += weighted_perp_square;
748     if (weighted_perp >= max_weighted_perp) {
749     pseudorapidity_adapt_new = p.eta();
750     azimuth_adapt_new = p.phi();
751     max_weighted_perp = weighted_perp;
752     }
753     }
754 mnguyen 1.20 default:
755     break;
756 ylai 1.16 }
757     }
758     // Update the fake rejection value
759     jets_.discr_fr01[jets_.nref] = std::max(
760     jets_.discr_fr01[jets_.nref], perp_square_sum);
761     // Update the directional adaption
762     pseudorapidity_adapt = pseudorapidity_adapt_new;
763     azimuth_adapt = azimuth_adapt_new;
764     }
765    
766 yilmaz 1.1 jets_.jtpt[jets_.nref] = jet.pt();
767     jets_.jteta[jets_.nref] = jet.eta();
768     jets_.jtphi[jets_.nref] = jet.phi();
769     jets_.jty[jets_.nref] = jet.eta();
770 yilmaz 1.2 jets_.jtpu[jets_.nref] = jet.pileup();
771 yilmaz 1.23 jets_.jtm[jets_.nref] = jet.mass();
772 yilmaz 1.1
773 yilmaz 1.2 if(isMC_ && usePat_){
774 mnguyen 1.11
775 yilmaz 1.26 for(UInt_t i = 0; i < genparts->size(); ++i){
776     const reco::GenParticle& p = (*genparts)[i];
777     if (p.status()!=1) continue;
778     if (p.charge()==0) continue;
779     double dr = deltaR(jet,p);
780     if(dr < rParam){
781     double ppt = p.pt();
782     jets_.genChargedSum[jets_.nref] += ppt;
783     if(ppt > hardPtMin_) jets_.genHardSum[jets_.nref] += ppt;
784     if(p.collisionId() == 0){
785     jets_.signalChargedSum[jets_.nref] += ppt;
786     if(ppt > hardPtMin_) jets_.signalHardSum[jets_.nref] += ppt;
787     }
788    
789     }
790     }
791    
792 mnguyen 1.11 const reco::GenJet * genjet = (*patjets)[j].genJet();
793 mnguyen 1.8
794 yilmaz 1.2 if(genjet){
795     jets_.refpt[jets_.nref] = genjet->pt();
796     jets_.refeta[jets_.nref] = genjet->eta();
797     jets_.refphi[jets_.nref] = genjet->phi();
798     jets_.refy[jets_.nref] = genjet->eta();
799     jets_.refdphijt[jets_.nref] = reco::deltaPhi(jet.phi(), genjet->phi());
800     jets_.refdrjt[jets_.nref] = reco::deltaR(jet.eta(),jet.phi(),genjet->eta(),genjet->phi());
801 yilmaz 1.25
802     if(doSubEvent_){
803     const GenParticle* gencon = genjet->getGenConstituent(0);
804     jets_.subid[jets_.nref] = gencon->collisionId();
805     }
806    
807 yilmaz 1.2 }else{
808     jets_.refpt[jets_.nref] = -999.;
809     jets_.refeta[jets_.nref] = -999.;
810     jets_.refphi[jets_.nref] = -999.;
811     jets_.refy[jets_.nref] = -999.;
812     jets_.refdphijt[jets_.nref] = -999.;
813     jets_.refdrjt[jets_.nref] = -999.;
814     }
815 yilmaz 1.1
816 mnguyen 1.11 jets_.refparton_flavorForB[jets_.nref] = (*patjets)[j].partonFlavour();
817    
818     // matched partons
819     const reco::GenParticle & parton = *(*patjets)[j].genParton();
820    
821     if((*patjets)[j].genParton()){
822     jets_.refparton_pt[jets_.nref] = parton.pt();
823     jets_.refparton_flavor[jets_.nref] = parton.pdgId();
824    
825     if(saveBfragments_ && abs(jets_.refparton_flavorForB[jets_.nref])==5){
826    
827     usedStringPts.clear();
828    
829     // uncomment this if you want to know the ugly truth about parton matching -matt
830     //if(jet.pt() > 50 &&abs(parton.pdgId())!=5 && parton.pdgId()!=21)
831     // cout<<" Identified as a b, but doesn't match b or gluon, id = "<<parton.pdgId()<<endl;
832    
833     jets_.bJetIndex[jets_.bMult] = jets_.nref;
834     jets_.bStatus[jets_.bMult] = parton.status();
835     jets_.bVx[jets_.bMult] = parton.vx();
836     jets_.bVy[jets_.bMult] = parton.vy();
837     jets_.bVz[jets_.bMult] = parton.vz();
838     jets_.bPt[jets_.bMult] = parton.pt();
839     jets_.bEta[jets_.bMult] = parton.eta();
840     jets_.bPhi[jets_.bMult] = parton.phi();
841     jets_.bPdg[jets_.bMult] = parton.pdgId();
842     jets_.bChg[jets_.bMult] = parton.charge();
843     jets_.bMult++;
844     saveDaughters(parton);
845     }
846    
847    
848 mnguyen 1.8 } else {
849 mnguyen 1.11 jets_.refparton_pt[jets_.nref] = -999;
850     jets_.refparton_flavor[jets_.nref] = -999;
851 mnguyen 1.8 }
852 mnguyen 1.11
853    
854 yilmaz 1.2 }
855 yilmaz 1.1
856     jets_.nref++;
857    
858    
859     }
860    
861    
862     if(isMC_){
863    
864 mnguyen 1.10 edm::Handle<HepMCProduct> hepMCProduct;
865     iEvent.getByLabel(eventInfoTag_,hepMCProduct);
866     const HepMC::GenEvent* MCEvt = hepMCProduct->GetEvent();
867 yilmaz 1.24
868     std::pair<HepMC::GenParticle*,HepMC::GenParticle*> beamParticles = MCEvt->beam_particles();
869     if(beamParticles.first != 0)jets_.beamId1 = beamParticles.first->pdg_id();
870     if(beamParticles.second != 0)jets_.beamId2 = beamParticles.second->pdg_id();
871 mnguyen 1.10
872 yilmaz 1.1 edm::Handle<GenEventInfoProduct> hEventInfo;
873     iEvent.getByLabel(eventInfoTag_,hEventInfo);
874     //jets_.pthat = hEventInfo->binningValues()[0];
875    
876     // binning values and qscale appear to be equivalent, but binning values not always present
877     jets_.pthat = hEventInfo->qScale();
878    
879     edm::Handle<vector<reco::GenJet> >genjets;
880     iEvent.getByLabel(genjetTag_, genjets);
881    
882     jets_.ngen = 0;
883     for(unsigned int igen = 0 ; igen < genjets->size(); ++igen){
884     const reco::GenJet & genjet = (*genjets)[igen];
885    
886     float genjet_pt = genjet.pt();
887    
888     // threshold to reduce size of output in minbias PbPb
889 yilmaz 1.9 if(genjet_pt>genPtMin_){
890 yilmaz 1.1
891     jets_.genpt[jets_.ngen] = genjet_pt;
892     jets_.geneta[jets_.ngen] = genjet.eta();
893     jets_.genphi[jets_.ngen] = genjet.phi();
894     jets_.geny[jets_.ngen] = genjet.eta();
895    
896 yilmaz 1.12 if(doSubEvent_){
897     const GenParticle* gencon = genjet.getGenConstituent(0);
898     jets_.gensubid[jets_.ngen] = gencon->collisionId();
899     }
900 yilmaz 1.1
901     // find matching patJet if there is one
902    
903     jets_.gendrjt[jets_.ngen] = -1.0;
904     jets_.genmatchindex[jets_.ngen] = -1;
905    
906     for(unsigned int ijet = 0 ; ijet < jets->size(); ++ijet){
907     const pat::Jet& jet = (*jets)[ijet];
908 mnguyen 1.8 const reco::GenJet* matchedGenJet = (*patjets)[ijet].genJet();
909     if(matchedGenJet){
910     // poor man's matching, someone fix please
911     if(fabs(genjet.pt()-matchedGenJet->pt())<0.0001 &&
912     fabs(genjet.eta()-matchedGenJet->eta())<0.0001 &&
913     (fabs(genjet.phi()-matchedGenJet->phi())<0.0001 || fabs(genjet.phi()-matchedGenJet->phi() - 2.0*TMath::Pi()) < 0.0001 )){
914 yilmaz 1.1
915     jets_.genmatchindex[jets_.ngen] = (int)ijet;
916     jets_.gendphijt[jets_.ngen] = reco::deltaPhi(jet.phi(),genjet.phi());
917     jets_.gendrjt[jets_.ngen] = reco::deltaR(jet,genjet);
918    
919     break;
920     }
921     }
922     }
923     jets_.ngen++;
924     }
925    
926     }
927    
928     }
929 yilmaz 1.24
930 yilmaz 1.1 t->Fill();
931 mnguyen 1.10 memset(&jets_,0,sizeof jets_);
932 yilmaz 1.24
933 yilmaz 1.1 }
934    
935    
936    
937    
938     //--------------------------------------------------------------------------------------------------
939     void HiInclusiveJetAnalyzer::fillL1Bits(const edm::Event &iEvent)
940     {
941     edm::Handle< L1GlobalTriggerReadoutRecord > L1GlobalTrigger;
942    
943     iEvent.getByLabel(L1gtReadout_, L1GlobalTrigger);
944     const TechnicalTriggerWord& technicalTriggerWordBeforeMask = L1GlobalTrigger->technicalTriggerWord();
945    
946     for (int i=0; i<64;i++)
947     {
948     jets_.l1TBit[i] = technicalTriggerWordBeforeMask.at(i);
949     }
950     jets_.nL1TBit = 64;
951    
952     int ntrigs = L1GlobalTrigger->decisionWord().size();
953     jets_.nL1ABit = ntrigs;
954    
955     for (int i=0; i != ntrigs; i++) {
956     bool accept = L1GlobalTrigger->decisionWord()[i];
957     //jets_.l1ABit[i] = (accept == true)? 1:0;
958     if(accept== true){
959     jets_.l1ABit[i] = 1;
960     }
961     else{
962     jets_.l1ABit[i] = 0;
963     }
964    
965     }
966     }
967    
968     //--------------------------------------------------------------------------------------------------
969     void HiInclusiveJetAnalyzer::fillHLTBits(const edm::Event &iEvent)
970     {
971     // Fill HLT trigger bits.
972     Handle<TriggerResults> triggerResultsHLT;
973     getProduct(hltResName_, triggerResultsHLT, iEvent);
974    
975     const TriggerResults *hltResults = triggerResultsHLT.product();
976     const TriggerNames & triggerNames = iEvent.triggerNames(*hltResults);
977    
978     jets_.nHLTBit = hltTrgNames_.size();
979    
980     for(size_t i=0;i<hltTrgNames_.size();i++){
981    
982     for(size_t j=0;j<triggerNames.size();++j) {
983    
984     if(triggerNames.triggerName(j) == hltTrgNames_[i]){
985    
986     //cout <<"hltTrgNames_(i) "<<hltTrgNames_[i]<<endl;
987     //cout <<"triggerName(j) "<<triggerNames.triggerName(j)<<endl;
988     //cout<<" result "<<triggerResultsHLT->accept(j)<<endl;
989     jets_.hltBit[i] = triggerResultsHLT->accept(j);
990     }
991    
992     }
993     }
994     }
995    
996     //--------------------------------------------------------------------------------------------------
997     template <typename TYPE>
998     inline void HiInclusiveJetAnalyzer::getProduct(const std::string name, edm::Handle<TYPE> &prod,
999     const edm::Event &event) const
1000     {
1001     // Try to access data collection from EDM file. We check if we really get just one
1002     // product with the given name. If not we throw an exception.
1003    
1004     event.getByLabel(edm::InputTag(name),prod);
1005     if (!prod.isValid())
1006     throw edm::Exception(edm::errors::Configuration, "HiInclusiveJetAnalyzer::GetProduct()\n")
1007     << "Collection with label '" << name << "' is not valid" << std::endl;
1008     }
1009    
1010     //--------------------------------------------------------------------------------------------------
1011     template <typename TYPE>
1012     inline bool HiInclusiveJetAnalyzer::getProductSafe(const std::string name, edm::Handle<TYPE> &prod,
1013     const edm::Event &event) const
1014     {
1015     // Try to safely access data collection from EDM file. We check if we really get just one
1016     // product with the given name. If not, we return false.
1017    
1018     if (name.size()==0)
1019     return false;
1020    
1021     try {
1022     event.getByLabel(edm::InputTag(name),prod);
1023     if (!prod.isValid())
1024     return false;
1025     } catch (...) {
1026     return false;
1027     }
1028     return true;
1029     }
1030    
1031    
1032 mnguyen 1.8 int
1033     HiInclusiveJetAnalyzer::getPFJetMuon(const pat::Jet& pfJet, const reco::PFCandidateCollection *pfCandidateColl)
1034     {
1035    
1036     int pfMuonIndex = -1;
1037     float ptMax = 0.;
1038    
1039    
1040     for(unsigned icand=0;icand<pfCandidateColl->size(); icand++) {
1041     const reco::PFCandidate pfCandidate = pfCandidateColl->at(icand);
1042    
1043     int id = pfCandidate.particleId();
1044     if(abs(id) != 3) continue;
1045 yilmaz 1.1
1046 mnguyen 1.8 if(reco::deltaR(pfJet,pfCandidate)>0.5) continue;
1047    
1048     double pt = pfCandidate.pt();
1049     if(pt>ptMax){
1050     ptMax = pt;
1051     pfMuonIndex = (int) icand;
1052     }
1053     }
1054    
1055     return pfMuonIndex;
1056    
1057     }
1058    
1059    
1060     double
1061     HiInclusiveJetAnalyzer::getPtRel(const reco::PFCandidate lep, const pat::Jet& jet )
1062     {
1063    
1064     float lj_x = jet.p4().px();
1065     float lj_y = jet.p4().py();
1066     float lj_z = jet.p4().pz();
1067    
1068     // absolute values squared
1069     float lj2 = lj_x*lj_x+lj_y*lj_y+lj_z*lj_z;
1070     float lep2 = lep.px()*lep.px()+lep.py()*lep.py()+lep.pz()*lep.pz();
1071    
1072     // projection vec(mu) to lepjet axis
1073     float lepXlj = lep.px()*lj_x+lep.py()*lj_y+lep.pz()*lj_z;
1074    
1075     // absolute value squared and normalized
1076     float pLrel2 = lepXlj*lepXlj/lj2;
1077    
1078     // lep2 = pTrel2 + pLrel2
1079     float pTrel2 = lep2-pLrel2;
1080    
1081     return (pTrel2 > 0) ? std::sqrt(pTrel2) : 0.0;
1082     }
1083 mnguyen 1.11
1084     // Recursive function, but this version gets called only the first time
1085     void
1086     HiInclusiveJetAnalyzer::saveDaughters(const reco::GenParticle &gen){
1087    
1088     for(unsigned i=0;i<gen.numberOfDaughters();i++){
1089     const reco::Candidate & daughter = *gen.daughter(i);
1090     double daughterPt = daughter.pt();
1091     if(daughterPt<1.) continue;
1092     double daughterEta = daughter.eta();
1093     if(fabs(daughterEta)>3.) continue;
1094     int daughterPdgId = daughter.pdgId();
1095     int daughterStatus = daughter.status();
1096     // Special case when b->b+string, both b and string contain all daughters, so only take the string
1097     if(gen.pdgId()==daughterPdgId && gen.status()==3 && daughterStatus==2) continue;
1098    
1099     // cheesy way of finding strings which were already used
1100     if(daughter.pdgId()==92){
1101     for(unsigned ist=0;ist<usedStringPts.size();ist++){
1102     if(fabs(daughter.pt() - usedStringPts[ist]) < 0.0001) return;
1103     }
1104     usedStringPts.push_back(daughter.pt());
1105     }
1106     jets_.bJetIndex[jets_.bMult] = jets_.nref;
1107     jets_.bStatus[jets_.bMult] = daughterStatus;
1108     jets_.bVx[jets_.bMult] = daughter.vx();
1109     jets_.bVy[jets_.bMult] = daughter.vy();
1110     jets_.bVz[jets_.bMult] = daughter.vz();
1111     jets_.bPt[jets_.bMult] = daughterPt;
1112     jets_.bEta[jets_.bMult] = daughterEta;
1113     jets_.bPhi[jets_.bMult] = daughter.phi();
1114     jets_.bPdg[jets_.bMult] = daughterPdgId;
1115     jets_.bChg[jets_.bMult] = daughter.charge();
1116     jets_.bMult++;
1117     saveDaughters(daughter);
1118     }
1119     }
1120    
1121     // This version called for all subsequent calls
1122     void
1123     HiInclusiveJetAnalyzer::saveDaughters(const reco::Candidate &gen){
1124    
1125     for(unsigned i=0;i<gen.numberOfDaughters();i++){
1126     const reco::Candidate & daughter = *gen.daughter(i);
1127     double daughterPt = daughter.pt();
1128     if(daughterPt<1.) continue;
1129     double daughterEta = daughter.eta();
1130     if(fabs(daughterEta)>3.) continue;
1131     int daughterPdgId = daughter.pdgId();
1132     int daughterStatus = daughter.status();
1133     // Special case when b->b+string, both b and string contain all daughters, so only take the string
1134     if(gen.pdgId()==daughterPdgId && gen.status()==3 && daughterStatus==2) continue;
1135    
1136     // cheesy way of finding strings which were already used
1137     if(daughter.pdgId()==92){
1138     for(unsigned ist=0;ist<usedStringPts.size();ist++){
1139     if(fabs(daughter.pt() - usedStringPts[ist]) < 0.0001) return;
1140     }
1141     usedStringPts.push_back(daughter.pt());
1142     }
1143    
1144     jets_.bJetIndex[jets_.bMult] = jets_.nref;
1145     jets_.bStatus[jets_.bMult] = daughterStatus;
1146     jets_.bVx[jets_.bMult] = daughter.vx();
1147     jets_.bVy[jets_.bMult] = daughter.vy();
1148     jets_.bVz[jets_.bMult] = daughter.vz();
1149     jets_.bPt[jets_.bMult] = daughterPt;
1150     jets_.bEta[jets_.bMult] = daughterEta;
1151     jets_.bPhi[jets_.bMult] = daughter.phi();
1152     jets_.bPdg[jets_.bMult] = daughterPdgId;
1153     jets_.bChg[jets_.bMult] = daughter.charge();
1154     jets_.bMult++;
1155     saveDaughters(daughter);
1156     }
1157     }
1158 yilmaz 1.25
1159     double HiInclusiveJetAnalyzer::getEt(math::XYZPoint pos, double energy){
1160     double et = energy*sin(pos.theta());
1161     return et;
1162     }
1163    
1164     math::XYZPoint HiInclusiveJetAnalyzer::getPosition(const DetId &id, reco::Vertex::Point vtx){
1165     const GlobalPoint& pos=geo->getPosition(id);
1166     math::XYZPoint posV(pos.x() - vtx.x(),pos.y() - vtx.y(),pos.z() - vtx.z());
1167     return posV;
1168     }
1169    
1170    
1171 yilmaz 1.1
1172     DEFINE_FWK_MODULE(HiInclusiveJetAnalyzer);