ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/CmsHi/JetAnalysis/src/HiInclusiveJetAnalyzer.cc
Revision: 1.24
Committed: Fri May 25 23:22:28 2012 UTC (12 years, 11 months ago) by yilmaz
Content type: text/plain
Branch: MAIN
CVS Tags: HiForest_V02_23, HiForest_V02_22, HiForest_V02_21, HiForest_V02_20, HiForest_V02_19, HiForest_V02_18, HiForest_V02_17, HiForest_V02_16, HiForest_V02_15, HiForest_V02_14, HiForest_V02_13, HiForest_V02_12, HiForest_V02_11, HiForest_V02_10
Changes since 1.23: +9 -19 lines
Log Message:
Bug fix for Hydjet events

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