ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/CmsHi/JetAnalysis/src/HiInclusiveJetAnalyzer.cc
Revision: 1.19
Committed: Fri May 11 19:07:29 2012 UTC (12 years, 11 months ago) by yilmaz
Content type: text/plain
Branch: MAIN
CVS Tags: HiForest_V02_05
Changes since 1.18: +28 -1 lines
Log Message:
cfg and setup update

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