ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MPIAnalyzer/src/MPIntuple.cc
Revision: 1.16
Committed: Thu Oct 14 22:44:24 2010 UTC (14 years, 6 months ago) by naodell
Content type: text/plain
Branch: MAIN
Changes since 1.15: +120 -72 lines
Log Message:
Added various histograms

File Contents

# User Rev Content
1 naodell 1.1 // Original Author: "Nathaniel Odell"
2     // Secondary Author: Steven Won
3 andrey 1.2 // With contributions from: Andrey Pozdnyakov
4 naodell 1.1 // Created: Thurs April 22 2010
5 naodell 1.16 // $Id: MPIntuple.cc,v 1.15 2010/10/12 19:03:10 naodell Exp $
6 naodell 1.1
7     // system include files
8     #include <memory>
9     #include <string>
10    
11     // user include files
12     #include "FWCore/Framework/interface/Frameworkfwd.h"
13     #include "FWCore/Framework/interface/EDAnalyzer.h"
14     #include "FWCore/Framework/interface/ESHandle.h"
15     #include "FWCore/Framework/interface/EventSetup.h"
16     #include "FWCore/Framework/interface/Event.h"
17     #include "FWCore/Framework/interface/MakerMacros.h"
18     #include "FWCore/ParameterSet/interface/ParameterSet.h"
19     #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
20     #include "Geometry/Records/interface/CaloGeometryRecord.h"
21 naodell 1.10 #include "DataFormats/Math/interface/deltaPhi.h"
22 naodell 1.1
23 naodell 1.15 // Libraries for objects
24 naodell 1.1 #include "DataFormats/JetReco/interface/CaloJet.h"
25     #include "DataFormats/JetReco/interface/CaloJetCollection.h"
26     #include "DataFormats/JetReco/interface/PFJet.h"
27     #include "DataFormats/JetReco/interface/PFJetCollection.h"
28     #include "DataFormats/JetReco/interface/GenJet.h"
29     #include "DataFormats/JetReco/interface/GenJetCollection.h"
30     #include "DataFormats/JetReco/interface/Jet.h"
31     #include "DataFormats/VertexReco/interface/Vertex.h"
32     #include "DataFormats/VertexReco/interface/VertexFwd.h"
33 naodell 1.7 #include "DataFormats/BTauReco/interface/JetTag.h"
34 naodell 1.15 #include "DataFormats/TrackReco/interface/Track.h"
35 naodell 1.16 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
36 devildog 1.11 //#include "RecoVertex/PrimaryVertexProducer/interface/VertexHigherPtSquared.h"
37 naodell 1.7
38 naodell 1.10 //GenParticles
39 naodell 1.13 //#include "DataFormats/HepMCCandidate/interface/GenParticle.h"
40 naodell 1.1
41 naodell 1.4 // Need to get correctors
42     #include "JetMETCorrections/Objects/interface/JetCorrector.h"
43    
44     // ntuple storage classes
45     #include "TCJet.h"
46     #include "TCPrimaryVtx.h"
47    
48 andrey 1.2 // Need for HLT trigger info:
49     #include "FWCore/Common/interface/TriggerNames.h"
50     #include "DataFormats/Common/interface/TriggerResults.h"
51     #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
52    
53 naodell 1.4 //Root stuff
54 naodell 1.1 #include "TROOT.h"
55 naodell 1.15 #include "TH1.h"
56 naodell 1.16 #include "TH2.h"
57 naodell 1.15 #include "TProfile.h"
58 naodell 1.1 #include "TFile.h"
59     #include "TTree.h"
60     #include "TString.h"
61     #include "TObject.h"
62     #include "TObjArray.h"
63     #include "TClonesArray.h"
64     #include "TRefArray.h"
65     #include "TLorentzVector.h"
66     #include "TVector3.h"
67    
68     using namespace edm;
69     using namespace std;
70     using namespace reco;
71    
72     //
73     // class declaration
74     //
75    
76     class MPIntuple : public edm::EDAnalyzer {
77     public:
78     explicit MPIntuple(const edm::ParameterSet&);
79     ~MPIntuple();
80    
81    
82     private:
83     virtual void beginJob() ;
84     virtual void analyze(const edm::Event&, const edm::EventSetup&);
85     virtual void endJob() ;
86    
87 andrey 1.2 bool triggerDecision(edm::Handle<edm::TriggerResults> &hltR, int iTrigger);
88 devildog 1.11 double sumPtSquared(const Vertex & v);
89 andrey 1.2
90 naodell 1.1 // ----------member data ---------------------------
91    
92     edm::InputTag PFJetHandle_;
93     edm::InputTag GenJetHandle_;
94     edm::InputTag PrimaryVtxHandle_;
95 naodell 1.4 edm::InputTag triggerResultsTag_;
96 naodell 1.1 // Counting variables //
97    
98 naodell 1.8 int eventNumber, runNumber, lumiSection, bunchCross;
99 naodell 1.1
100 naodell 1.6 TTree* sTree;
101 naodell 1.1 TFile* ntupleFile;
102    
103 naodell 1.4 TClonesArray* jet_ak5PF;
104 naodell 1.1 TClonesArray* jetP4_ak5Gen;
105 naodell 1.5 TClonesArray* primaryVtx;
106 naodell 1.1
107     bool doGenJets_;
108     bool doPFJets_;
109 naodell 1.7 bool triggerHLT_;
110 naodell 1.8 bool isRealData;
111 naodell 1.1 string rootfilename;
112    
113 andrey 1.2 //Triggers
114 naodell 1.4 string hlTriggerResults_, hltName_, triggerName_;
115     TriggerNames triggerNames;
116 andrey 1.2 HLTConfigProvider hltConfig_;
117 naodell 1.4 vector<string> hlNames;
118 naodell 1.5 unsigned int triggerStatus;
119 andrey 1.2
120 naodell 1.15 //Histograms
121 naodell 1.16 TH1F * h1_nAssociatedTracks;
122 naodell 1.15 TH1F * h1_trackZ;
123 naodell 1.16 TH1F * h1_trackDxy;
124     TH1F * h1_allTrackDeltaZ_PV;
125     TH1F * h1_allTrackDeltaZ_LJ;
126     TH2F * h2_nAssTracksVsJetPt;
127 naodell 1.15 TProfile * p1_nVtcs;
128    
129    
130 andrey 1.2
131 naodell 1.1 };
132    
133     MPIntuple::MPIntuple(const edm::ParameterSet& iConfig):
134    
135     PFJetHandle_(iConfig.getUntrackedParameter<edm::InputTag>("PFJetTag")),
136     GenJetHandle_(iConfig.getUntrackedParameter<edm::InputTag>("GenJetTag")),
137     PrimaryVtxHandle_(iConfig.getUntrackedParameter<edm::InputTag>("PrimaryVtxTag")),
138     doGenJets_(iConfig.getUntrackedParameter<bool>("doGenJets")),
139     doPFJets_(iConfig.getUntrackedParameter<bool>("doPFJets")),
140 naodell 1.7 triggerHLT_(iConfig.getUntrackedParameter<bool>("triggerHLT")),
141 andrey 1.2 rootfilename(iConfig.getUntrackedParameter<string>("rootfilename")),
142     hlTriggerResults_(iConfig.getUntrackedParameter<string>("HLTriggerResults","TriggerResults")),
143     hltName_(iConfig.getUntrackedParameter<string>("hltName","HLT"))
144 naodell 1.1 {
145     }
146    
147     MPIntuple::~MPIntuple()
148     {
149    
150     }
151    
152     //
153     // member functions
154     //
155    
156     // ------------ method called to for each event ------------
157     void MPIntuple::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
158     {
159    
160 naodell 1.14
161 naodell 1.6 eventNumber = iEvent.id().event();
162     runNumber = iEvent.id().run();
163 naodell 1.1 lumiSection = (unsigned int)iEvent.getLuminosityBlock().luminosityBlock();
164 naodell 1.8 bunchCross = (unsigned int)iEvent.bunchCrossing();
165     isRealData = iEvent.isRealData();
166    
167 naodell 1.16 edm::Handle<reco::BeamSpot> beamSpotHandle;
168     iEvent.getByLabel("offlineBeamSpot", beamSpotHandle);
169     reco::BeamSpot vertexBeamSpot = *beamSpotHandle;
170    
171     int pfCount = 0;
172     int genCount = 0;
173     int vtxCount = 0;
174     float ljVertexZ = -999;
175     float primaryVertexZ = -999;
176    
177     //////////////////////////
178     //Get vertex information//
179     //////////////////////////
180    
181     Handle<reco::VertexCollection> primaryVtcs;
182     iEvent.getByLabel(PrimaryVtxHandle_, primaryVtcs);
183    
184     for(VertexCollection::const_iterator vtx_iter = primaryVtcs->begin(); vtx_iter!= primaryVtcs->end(); ++vtx_iter){
185    
186     reco::Vertex myVtx = reco::Vertex(*vtx_iter);
187     if(!myVtx.isValid() || myVtx.isFake()) continue;
188     TCPrimaryVtx* vtxCon = new ((*primaryVtx)[vtxCount]) TCPrimaryVtx;
189    
190     vtxCon->SetPosition(myVtx.x(), myVtx.y(), myVtx.z());
191     vtxCon->SetNDof(myVtx.ndof());
192     vtxCon->SetChi2(myVtx.chi2());
193     vtxCon->SetNTrks(myVtx.tracksSize());
194     vtxCon->SetSumPt2Trks(sumPtSquared(myVtx));
195    
196     if(vtxCount == 0) primaryVertexZ = myVtx.z();
197    
198     ++vtxCount;
199    
200     }
201    
202     p1_nVtcs->Fill(runNumber, vtxCount);
203    
204     ///////////////////////
205     //get jet information//
206     ///////////////////////
207 naodell 1.1
208     if(doPFJets_){
209 naodell 1.4
210 naodell 1.7 edm::Handle<reco::JetTagCollection> bTagHandle1;
211     iEvent.getByLabel("trackCountingHighEffBJetTags", bTagHandle1);
212     const reco::JetTagCollection & bTags1 = *(bTagHandle1.product());
213     reco::JetTagCollection::const_iterator jet_it_1;
214    
215     edm::Handle<reco::JetTagCollection> bTagHandle2;
216     iEvent.getByLabel("trackCountingHighPurBJetTags", bTagHandle2);
217     const reco::JetTagCollection & bTags2 = *(bTagHandle2.product());
218     reco::JetTagCollection::const_iterator jet_it_2;
219    
220 naodell 1.14
221 naodell 1.4 const JetCorrector* correctorL2 = JetCorrector::getJetCorrector ("ak5PFL2Relative",iSetup);
222     const JetCorrector* correctorL3 = JetCorrector::getJetCorrector ("ak5PFL3Absolute",iSetup);
223    
224 naodell 1.1 Handle<reco::PFJetCollection> PFJets;
225     iEvent.getByLabel(PFJetHandle_, PFJets);
226 naodell 1.16
227 naodell 1.1 for(PFJetCollection::const_iterator jet_iter = PFJets->begin(); jet_iter!= PFJets->end(); ++jet_iter){
228    
229     reco::PFJet myJet = reco::PFJet(*jet_iter);
230     if(myJet.et() > 5){
231 naodell 1.7
232     for(jet_it_1 = bTags1.begin(); jet_it_1 != bTags1.end(); jet_it_1++) {
233     if( (fabs(jet_it_1->first->eta() - myJet.eta()) < .005) && (deltaPhi(jet_it_1->first->phi(),myJet.phi()) < .005)) break;
234     }
235    
236     for(jet_it_2 = bTags2.begin(); jet_it_2 != bTags2.end(); jet_it_2++) {
237     if( (fabs(jet_it_2->first->eta() - myJet.eta()) < .005) && (deltaPhi(jet_it_2->first->phi(),myJet.phi()) < .005)) break;
238     }
239 naodell 1.6
240 naodell 1.16 TCJet* jetCon = new ((*jet_ak5PF)[pfCount]) TCJet;
241 naodell 1.6
242     jetCon->SetP4(myJet.px(), myJet.py(), myJet.pz(), myJet.energy());
243 naodell 1.15 jetCon->SetVtx(-999.0, -999.0, -999.0);
244 naodell 1.4
245 naodell 1.6 jetCon->SetChHadFrac(myJet.chargedHadronEnergyFraction());
246     jetCon->SetNeuHadFrac(myJet.neutralHadronEnergyFraction());
247     jetCon->SetChEmFrac(myJet.chargedEmEnergyFraction());
248     jetCon->SetNeuEmFrac(myJet.neutralEmEnergyFraction());
249 naodell 1.4
250 naodell 1.6 jetCon->SetNumConstit(myJet.chargedMultiplicity() + myJet.neutralMultiplicity());
251     jetCon->SetNumChPart(myJet.chargedMultiplicity());
252 naodell 1.4
253 naodell 1.6 jetCon->SetNumChPart(myJet.chargedMultiplicity());
254 naodell 1.4
255 naodell 1.7 if(jet_it_2 != bTags2.end()) jetCon->SetBDiscrTrkCountHiPure(jet_it_2->second);
256     if(jet_it_1 != bTags1.end()) jetCon->SetBDiscrTrkCountHiEff(jet_it_1->second);
257    
258     float scale2 = correctorL2->correction(myJet);
259     myJet.scaleEnergy(scale2);
260     float scale3 = correctorL3->correction(myJet);
261     myJet.scaleEnergy(scale3);
262    
263     //more corrections?
264    
265 naodell 1.6 jetCon->SetJetCorr(2, scale2);
266     jetCon->SetJetCorr(3, scale3);
267 naodell 1.7
268 naodell 1.16 /////////////////////////
269     //get associated tracks//
270     /////////////////////////
271 naodell 1.15
272     const reco::TrackRefVector &tracks = myJet.getTrackRefs();
273    
274 naodell 1.16 float sumTrackZ, sumTrackIP;
275     float meanTrackZ, meanTrackIP;
276     int nAssociatedTracks = 0;
277    
278     sumTrackZ = sumTrackIP = 0;
279     meanTrackZ = meanTrackIP = -999;
280 naodell 1.15
281 naodell 1.16 if(myJet.pt() > 20){
282 naodell 1.15
283     for(TrackRefVector::const_iterator iTrack = tracks.begin(); iTrack != tracks.end(); ++iTrack){
284    
285     const reco::Track& myTrack = **iTrack;
286    
287     sumTrackZ += myTrack.vz();
288 naodell 1.16 sumTrackIP += myTrack.dxy(vertexBeamSpot.position());
289 naodell 1.15 ++nAssociatedTracks;
290 naodell 1.16
291 naodell 1.15 }
292 naodell 1.16
293 naodell 1.15 meanTrackZ = sumTrackZ/nAssociatedTracks;
294 naodell 1.16 h1_nAssociatedTracks->Fill(nAssociatedTracks);
295     h2_nAssTracksVsJetPt->Fill(nAssociatedTracks, myJet.pt());
296    
297 naodell 1.15 for(TrackRefVector::const_iterator iTrack = tracks.begin(); iTrack != tracks.end(); ++iTrack){
298 naodell 1.16
299 naodell 1.15 const reco::Track& myTrack = **iTrack;
300 naodell 1.16
301 naodell 1.15 h1_trackZ->Fill(meanTrackZ - myTrack.vz());
302 naodell 1.16 h1_trackDxy->Fill(myTrack.dxy(vertexBeamSpot.position()));
303    
304     }
305    
306     if(pfCount == 0){
307    
308     ljVertexZ = meanTrackZ;
309    
310     }else{
311    
312     h1_allTrackDeltaZ_LJ->Fill(ljVertexZ - meanTrackZ);
313    
314 naodell 1.15 }
315 naodell 1.16
316     h1_allTrackDeltaZ_PV->Fill(primaryVertexZ - meanTrackZ);
317    
318 naodell 1.15 }
319 naodell 1.16
320     jetCon->SetVtx(0, 0, meanTrackZ);
321    
322     ++pfCount;
323 naodell 1.4 }
324     }
325 naodell 1.1 }
326    
327 naodell 1.10 if(!isRealData){
328 naodell 1.1
329     Handle<reco::GenJetCollection> GenJets;
330     iEvent.getByLabel(GenJetHandle_, GenJets);
331 naodell 1.16
332 naodell 1.1 for(GenJetCollection::const_iterator jet_iter = GenJets->begin(); jet_iter!= GenJets->end(); ++jet_iter){
333    
334     reco::GenJet myJet = reco::GenJet(*jet_iter);
335    
336 naodell 1.4 if(myJet.pt() > 5){
337 naodell 1.1
338 naodell 1.16 new ((*jetP4_ak5Gen)[genCount]) TLorentzVector(myJet.px(), myJet.py(), myJet.pz(), myJet.energy());
339 naodell 1.1
340 naodell 1.16 ++genCount;
341 naodell 1.1
342     }
343 naodell 1.13 }
344 naodell 1.1 }
345 naodell 1.10
346 naodell 1.16 ///////////////////////////
347     //get trigger information//
348     ///////////////////////////
349 andrey 1.2
350 naodell 1.7 if(triggerHLT_){
351 andrey 1.2
352 naodell 1.7 edm::Handle<TriggerResults> hltR;
353     triggerResultsTag_ = InputTag(hlTriggerResults_,"",hltName_);
354     iEvent.getByLabel(triggerResultsTag_,hltR);
355    
356     const TriggerNames & triggerNames = iEvent.triggerNames(*hltR);
357     hlNames=triggerNames.triggerNames();
358    
359 naodell 1.9 string MPI_TriggerNames[] = {"HLT_PixelTracks_Multiplicity70", "HLT_MinBiasBSC_NoBPTX", "HLT_PixelTracks_Multiplicity40","HLT_L1Tech_HCAL_HF", "HLT_IsoTrackHB_8E29", "HLT_IsoTrackHE_8E29", "HLT_L1Tech_RPC_TTU_RBst1_collisions", "HLT_L1_BscMinBiasOR_BptxPlusORMinus", "HLT_L1Tech_BSC_halo_forPhysicsBackground", "HLT_L1Tech_BSC_HighMultiplicity", "HLT_MinBiasPixel_DoubleIsoTrack5", "HLT_MinBiasPixel_DoubleTrack", "HLT_MinBiasPixel_SingleTrack", "HLT_ZeroBiasPixel_SingleTrack", "HLT_MinBiasBSC", "HLT_StoppedHSCP_8E29", "HLT_Jet15U_HcalNoiseFiltered", "HLT_QuadJet15U", "HLT_DiJetAve30U_8E29", "HLT_DiJetAve15U_8E29", "HLT_FwdJet20U", "HLT_Jet50U", "HLT_Jet30U", "HLT_Jet15U", "HLT_BTagMu_Jet10U", "HLT_DoubleJet15U_ForwardBackward", "HLT_BTagIP_Jet50U", "HLT_DoubleLooseIsoTau15", "HLT_SingleLooseIsoTau20", "HLT_HT100U", "HLT_MET100", "HLT_MET45"};
360 naodell 1.7
361     bool triggerPassed = false;
362     triggerStatus = 0x0;
363 naodell 1.4
364 naodell 1.7 for (uint i=0; i<hlNames.size(); ++i) {
365    
366     triggerPassed = triggerDecision(hltR, i);
367 naodell 1.4
368 naodell 1.7 if(triggerPassed){
369 andrey 1.3
370 naodell 1.7 for (uint j = 0; j != 32; ++j){
371    
372     if (hlNames[i] == MPI_TriggerNames[j])
373     {
374 naodell 1.16 //cout<<"trigger name: "<<hlNames[i]<<"list: "<<dec<<j+1<<endl;
375 naodell 1.7 triggerStatus |= 0x01 << j;
376    
377     }
378     }
379 naodell 1.4 }
380 naodell 1.7 }
381     }
382 naodell 1.16 //cout<< "total status: "<<hex<<triggerStatus<<endl;
383 naodell 1.14
384    
385    
386 naodell 1.16 if((pfCount > 3 || genCount > 3) && vtxCount > 0) sTree -> Fill();
387 naodell 1.15
388 naodell 1.4 jet_ak5PF->Clear();
389 naodell 1.5 primaryVtx->Clear();
390 naodell 1.4 jetP4_ak5Gen->Clear();
391 andrey 1.2
392 naodell 1.1 }
393    
394    
395     // ------------ method called once each job just before starting event loop ------------
396     void MPIntuple::beginJob()
397     {
398    
399     ntupleFile = new TFile(rootfilename.c_str(), "RECREATE");
400 naodell 1.5 sTree = new TTree("mpiTree", "Tree for Jets");
401 naodell 1.15
402 naodell 1.16 h1_trackZ = new TH1F("h1_trackZ", "#Delta z for associated tracks", 60, -0.3, 0.3);
403     h1_nAssociatedTracks = new TH1F("h1_nAssociatedTracks", "Number of tracks associated to jet", 20, -0.5, 19.5);
404     h1_trackDxy = new TH1F("h1_trackDxy", "dxy for all associated tracks", 50, -0.4, 0.4);
405     h1_allTrackDeltaZ_PV = new TH1F("h1_allTrackDeltaZ_PV", "#Delta z between jet and primary vertex", 60, -0.3, 0.3);
406     h1_allTrackDeltaZ_LJ = new TH1F("h1_allTrackDeltaZ_LJ", "#Delta z between leading jet and other jet vertices", 60, -0.3, 0.3);
407     h2_nAssTracksVsJetPt = new TH2F("h2_nAssTracksVsJetPt", "nTracks vs Jet p_{T}", 20, -0.5, 19.5, 50, 0, 100);
408 naodell 1.15 p1_nVtcs = new TProfile("p1_nVtcs", "<N> per run", 5200.0, 132300.0, 137500.0, 0.0, 6.0);
409    
410 naodell 1.4 jet_ak5PF = new TClonesArray("TCJet");
411 naodell 1.1 jetP4_ak5Gen = new TClonesArray("TLorentzVector");
412 naodell 1.5 primaryVtx = new TClonesArray("TCPrimaryVtx");
413 naodell 1.1
414 naodell 1.15
415 naodell 1.4 sTree->Branch("jet_ak5PF",&jet_ak5PF, 6400, 0);
416 naodell 1.1 sTree->Branch("jetP4_ak5Gen",&jetP4_ak5Gen, 6400, 0);
417 naodell 1.5 sTree->Branch("primaryVtx",&primaryVtx, 6400, 0);
418 naodell 1.1
419     sTree->Branch("eventNumber",&eventNumber, "eventNumber/I");
420     sTree->Branch("runNumber",&runNumber, "runNumber/I");
421     sTree->Branch("lumiSection",&lumiSection, "lumiSection/I");
422 naodell 1.6 sTree->Branch("triggerStatus",&triggerStatus, "triggerStatus/i");
423 naodell 1.8 sTree->Branch("isRealData",&isRealData, "isRealData/i");
424     sTree->Branch("bunchCross",&bunchCross, "bunchCross/i");
425 naodell 1.1
426     }
427    
428     // ------------ method called once each job just after ending the event loop ------------
429     void MPIntuple::endJob()
430     {
431 naodell 1.15
432     h1_trackZ->Write();
433 naodell 1.16 h1_nAssociatedTracks->Write();
434     h1_trackDxy->Write();
435     h1_allTrackDeltaZ_PV->Write();
436     h1_allTrackDeltaZ_LJ->Write();
437     h2_nAssTracksVsJetPt->Write();
438 naodell 1.15 p1_nVtcs->Write();
439    
440 naodell 1.1 ntupleFile->Write();
441     ntupleFile->Close();
442    
443     }
444    
445    
446 andrey 1.2 bool MPIntuple::triggerDecision(edm::Handle<edm::TriggerResults> &hltR, int iTrigger)
447     {
448     bool triggerPassed = false;
449     if(hltR->wasrun(iTrigger) &&
450     hltR->accept(iTrigger) &&
451     !hltR->error(iTrigger) ){
452     triggerPassed = true;
453     }
454     return triggerPassed;
455 naodell 1.16 }
456    
457     double MPIntuple::sumPtSquared(const Vertex & v)
458     {
459     double sum = 0.;
460     double pT;
461     for (Vertex::trackRef_iterator it = v.tracks_begin(); it != v.tracks_end(); it++) {
462     pT = (**it).pt();
463     double epT=(**it).ptError(); pT=pT>epT ? pT-epT : 0;
464    
465     sum += pT*pT;
466 andrey 1.2 }
467 naodell 1.16 return sum;
468     }
469 andrey 1.2
470 naodell 1.1 //define this as a plug-in
471     DEFINE_FWK_MODULE(MPIntuple);