ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MPIAnalyzer/src/MPIntuple.cc
Revision: 1.6
Committed: Sun May 23 21:56:11 2010 UTC (14 years, 11 months ago) by naodell
Content type: text/plain
Branch: MAIN
Changes since 1.5: +29 -32 lines
Log Message:
Corrections to triggerStatus variable.

File Contents

# User Rev Content
1 naodell 1.1 /*
2     Description: <one line class summary>
3    
4     Implementation:
5     n-tuple creator for jet studies
6     */
7     //
8     // Original Author: "Nathaniel Odell"
9     // Secondary Author: Steven Won
10 andrey 1.2 // With contributions from: Andrey Pozdnyakov
11 naodell 1.1 // Created: Thurs April 22 2010
12 naodell 1.6 // $Id: MPIntuple.cc,v 1.5 2010/05/23 02:33:41 naodell Exp $
13 naodell 1.1 //
14     //
15    
16    
17     // system include files
18     #include <memory>
19     #include <string>
20    
21     // user include files
22     #include "FWCore/Framework/interface/Frameworkfwd.h"
23     #include "FWCore/Framework/interface/EDAnalyzer.h"
24     #include "FWCore/Framework/interface/ESHandle.h"
25     #include "FWCore/Framework/interface/EventSetup.h"
26     #include "FWCore/Framework/interface/Event.h"
27     #include "FWCore/Framework/interface/MakerMacros.h"
28    
29     #include "FWCore/ParameterSet/interface/ParameterSet.h"
30    
31     #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
32     #include "Geometry/Records/interface/CaloGeometryRecord.h"
33    
34 naodell 1.4 // Jet and vertex functions
35 naodell 1.1 #include "DataFormats/JetReco/interface/CaloJet.h"
36     #include "DataFormats/JetReco/interface/CaloJetCollection.h"
37     #include "DataFormats/JetReco/interface/PFJet.h"
38     #include "DataFormats/JetReco/interface/PFJetCollection.h"
39     #include "DataFormats/JetReco/interface/GenJet.h"
40     #include "DataFormats/JetReco/interface/GenJetCollection.h"
41     #include "DataFormats/JetReco/interface/Jet.h"
42     #include "DataFormats/VertexReco/interface/Vertex.h"
43     #include "DataFormats/VertexReco/interface/VertexFwd.h"
44    
45 naodell 1.4 // Need to get correctors
46     #include "JetMETCorrections/Objects/interface/JetCorrector.h"
47    
48     // ntuple storage classes
49     #include "TCJet.h"
50     #include "TCPrimaryVtx.h"
51    
52 andrey 1.2 // Need for HLT trigger info:
53     #include "FWCore/Common/interface/TriggerNames.h"
54     #include "DataFormats/Common/interface/TriggerResults.h"
55     #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
56    
57 naodell 1.4 //Root stuff
58 naodell 1.1 #include "TROOT.h"
59     #include "TFile.h"
60     #include "TTree.h"
61     #include "TString.h"
62     #include "TObject.h"
63     #include "TObjArray.h"
64     #include "TClonesArray.h"
65     #include "TRefArray.h"
66     #include "TLorentzVector.h"
67     #include "TVector3.h"
68    
69     using namespace edm;
70     using namespace std;
71     using namespace reco;
72    
73     //
74     // class declaration
75     //
76    
77     class MPIntuple : public edm::EDAnalyzer {
78     public:
79     explicit MPIntuple(const edm::ParameterSet&);
80     ~MPIntuple();
81    
82    
83     private:
84     virtual void beginJob() ;
85     virtual void analyze(const edm::Event&, const edm::EventSetup&);
86     virtual void endJob() ;
87    
88 andrey 1.2 bool triggerDecision(edm::Handle<edm::TriggerResults> &hltR, int iTrigger);
89    
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
97    
98     // Counting variables //
99    
100     int eventNumber, runNumber, lumiSection;
101    
102 naodell 1.6 TTree* sTree;
103 naodell 1.1 TFile* ntupleFile;
104    
105 naodell 1.4 TClonesArray* jet_ak5PF;
106 naodell 1.1 TClonesArray* jetP4_ak5Gen;
107 naodell 1.5 TClonesArray* primaryVtx;
108 naodell 1.1
109     bool doGenJets_;
110     bool doPFJets_;
111 naodell 1.6 bool isMC;
112 naodell 1.1
113     string rootfilename;
114    
115 andrey 1.2 //Triggers
116 naodell 1.4 string hlTriggerResults_, hltName_, triggerName_;
117     TriggerNames triggerNames;
118 andrey 1.2 HLTConfigProvider hltConfig_;
119 naodell 1.4 vector<string> hlNames;
120 naodell 1.5 unsigned int triggerStatus;
121 andrey 1.2
122    
123 naodell 1.1 };
124    
125     MPIntuple::MPIntuple(const edm::ParameterSet& iConfig):
126    
127     PFJetHandle_(iConfig.getUntrackedParameter<edm::InputTag>("PFJetTag")),
128     GenJetHandle_(iConfig.getUntrackedParameter<edm::InputTag>("GenJetTag")),
129     PrimaryVtxHandle_(iConfig.getUntrackedParameter<edm::InputTag>("PrimaryVtxTag")),
130     doGenJets_(iConfig.getUntrackedParameter<bool>("doGenJets")),
131     doPFJets_(iConfig.getUntrackedParameter<bool>("doPFJets")),
132 andrey 1.2 rootfilename(iConfig.getUntrackedParameter<string>("rootfilename")),
133     hlTriggerResults_(iConfig.getUntrackedParameter<string>("HLTriggerResults","TriggerResults")),
134     hltName_(iConfig.getUntrackedParameter<string>("hltName","HLT"))
135 naodell 1.1 {
136 andrey 1.2 //edm::TriggerNames
137     // triggerNames(iConfig);
138 naodell 1.1
139     }
140    
141    
142     MPIntuple::~MPIntuple()
143     {
144    
145     }
146    
147     //
148     // member functions
149     //
150    
151     // ------------ method called to for each event ------------
152     void MPIntuple::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
153     {
154    
155 naodell 1.6 eventNumber = iEvent.id().event();
156     runNumber = iEvent.id().run();
157 naodell 1.1 lumiSection = (unsigned int)iEvent.getLuminosityBlock().luminosityBlock();
158 naodell 1.6 // isMC = iEvent.???();
159 naodell 1.1
160 naodell 1.6 int PFcount = 0;
161 naodell 1.1 int Gencount = 0;
162     int Vtxcount = 0;
163    
164     if(doPFJets_){
165 naodell 1.4
166     const JetCorrector* correctorL2 = JetCorrector::getJetCorrector ("ak5PFL2Relative",iSetup);
167     const JetCorrector* correctorL3 = JetCorrector::getJetCorrector ("ak5PFL3Absolute",iSetup);
168    
169 naodell 1.1 Handle<reco::PFJetCollection> PFJets;
170     iEvent.getByLabel(PFJetHandle_, PFJets);
171    
172     PFcount = 0;
173    
174     for(PFJetCollection::const_iterator jet_iter = PFJets->begin(); jet_iter!= PFJets->end(); ++jet_iter){
175    
176     reco::PFJet myJet = reco::PFJet(*jet_iter);
177 naodell 1.4
178     float scale2 = correctorL2->correction(myJet);
179     float scale3 = correctorL3->correction(myJet);
180    
181 naodell 1.1 if(myJet.et() > 5){
182 naodell 1.6
183     TCJet* jetCon = new ((*jet_ak5PF)[PFcount]) TCJet;
184    
185     jetCon->SetP4(myJet.px(), myJet.py(), myJet.pz(), myJet.energy());
186     jetCon->SetVtx(myJet.vx(), myJet.vy(), myJet.vz());
187 naodell 1.4
188 naodell 1.6 jetCon->SetChHadFrac(myJet.chargedHadronEnergyFraction());
189     jetCon->SetNeuHadFrac(myJet.neutralHadronEnergyFraction());
190     jetCon->SetChEmFrac(myJet.chargedEmEnergyFraction());
191     jetCon->SetNeuEmFrac(myJet.neutralEmEnergyFraction());
192 naodell 1.4
193 naodell 1.6 jetCon->SetNumConstit(myJet.chargedMultiplicity() + myJet.neutralMultiplicity());
194     jetCon->SetNumChPart(myJet.chargedMultiplicity());
195 naodell 1.4
196 naodell 1.6 jetCon->SetNumChPart(myJet.chargedMultiplicity());
197 naodell 1.4
198 naodell 1.6 jetCon->SetJetCorr(2, scale2);
199     jetCon->SetJetCorr(3, scale3);
200 naodell 1.4
201 naodell 1.6
202 naodell 1.1 ++PFcount;
203 naodell 1.4 }
204     }
205 naodell 1.1 }
206    
207     if(doGenJets_){
208    
209     Handle<reco::GenJetCollection> GenJets;
210     iEvent.getByLabel(GenJetHandle_, GenJets);
211    
212    
213     for(GenJetCollection::const_iterator jet_iter = GenJets->begin(); jet_iter!= GenJets->end(); ++jet_iter){
214    
215     reco::GenJet myJet = reco::GenJet(*jet_iter);
216    
217 naodell 1.4 if(myJet.pt() > 5){
218 naodell 1.1
219     new ((*jetP4_ak5Gen)[Gencount]) TLorentzVector(myJet.px(), myJet.py(), myJet.pz(), myJet.energy());
220    
221     ++Gencount;
222    
223     }
224     }
225     }
226 naodell 1.4
227    
228 naodell 1.5 Handle<reco::VertexCollection> primaryVtcs;
229     iEvent.getByLabel(PrimaryVtxHandle_, primaryVtcs);
230 naodell 1.1
231 naodell 1.5 for(VertexCollection::const_iterator vtx_iter = primaryVtcs->begin(); vtx_iter!= primaryVtcs->end(); ++vtx_iter){
232 naodell 1.1
233     reco::Vertex myVtx = reco::Vertex(*vtx_iter);
234    
235 naodell 1.6 TCPrimaryVtx* vtxCon = new ((*primaryVtx)[Vtxcount]) TCPrimaryVtx;
236 naodell 1.4
237 naodell 1.6 vtxCon->SetPosition(myVtx.x(), myVtx.y(), myVtx.z());
238     vtxCon->SetNDof(myVtx.ndof());
239     vtxCon->SetChi2(myVtx.chi2());
240    
241 naodell 1.1 ++Vtxcount;
242    
243     }
244    
245 andrey 1.2
246    
247     //---------- Filling HLT trigger bits! ------------
248    
249     edm::Handle<TriggerResults> hltR;
250     triggerResultsTag_ = InputTag(hlTriggerResults_,"",hltName_);
251     iEvent.getByLabel(triggerResultsTag_,hltR);
252 naodell 1.4
253     const TriggerNames & triggerNames = iEvent.triggerNames(*hltR);
254 andrey 1.2 hlNames=triggerNames.triggerNames();
255 naodell 1.4
256     string MPI_TriggerNames[32] = {"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"};
257    
258 andrey 1.2 bool triggerPassed = false;
259 naodell 1.5 triggerStatus = 0x0;
260 andrey 1.2
261 naodell 1.4 for (uint iT=0; iT<hlNames.size(); ++iT) {
262 andrey 1.3
263 naodell 1.4 triggerPassed = triggerDecision(hltR, iT);
264    
265     if(triggerPassed){
266    
267     for (int j = 0; j != 32; ++j){
268 andrey 1.3
269 naodell 1.4 if (hlNames[iT] == MPI_TriggerNames[j])
270     {
271     cout<<"trigger name: "<<hlNames[iT]<<" status: "<<triggerPassed<<endl;
272     triggerStatus |= 0x01 << j;
273     }
274     }
275 andrey 1.2 }
276 naodell 1.4 }
277    
278 andrey 1.2
279 naodell 1.4 if(PFcount > 0 || Gencount > 0 || Vtxcount > 0); sTree -> Fill();
280    
281     jet_ak5PF->Clear();
282 naodell 1.5 primaryVtx->Clear();
283 naodell 1.4 jetP4_ak5Gen->Clear();
284 andrey 1.2
285 naodell 1.1 }
286    
287    
288     // ------------ method called once each job just before starting event loop ------------
289     void MPIntuple::beginJob()
290     {
291    
292     ntupleFile = new TFile(rootfilename.c_str(), "RECREATE");
293 naodell 1.5 sTree = new TTree("mpiTree", "Tree for Jets");
294 naodell 1.1
295 naodell 1.4 jet_ak5PF = new TClonesArray("TCJet");
296 naodell 1.1 jetP4_ak5Gen = new TClonesArray("TLorentzVector");
297 naodell 1.5 primaryVtx = new TClonesArray("TCPrimaryVtx");
298 naodell 1.1
299 naodell 1.4 sTree->Branch("jet_ak5PF",&jet_ak5PF, 6400, 0);
300 naodell 1.1 sTree->Branch("jetP4_ak5Gen",&jetP4_ak5Gen, 6400, 0);
301 naodell 1.5 sTree->Branch("primaryVtx",&primaryVtx, 6400, 0);
302 naodell 1.1
303     sTree->Branch("eventNumber",&eventNumber, "eventNumber/I");
304     sTree->Branch("runNumber",&runNumber, "runNumber/I");
305     sTree->Branch("lumiSection",&lumiSection, "lumiSection/I");
306 naodell 1.6 sTree->Branch("triggerStatus",&triggerStatus, "triggerStatus/i");
307     sTree->Branch("isMC",&isMC,"isMC/O");
308 naodell 1.5
309 naodell 1.1
310     }
311    
312     // ------------ method called once each job just after ending the event loop ------------
313     void MPIntuple::endJob()
314     {
315     ntupleFile->Write();
316     ntupleFile->Close();
317    
318     }
319    
320    
321 andrey 1.2 bool MPIntuple::triggerDecision(edm::Handle<edm::TriggerResults> &hltR, int iTrigger)
322     {
323     bool triggerPassed = false;
324     if(hltR->wasrun(iTrigger) &&
325     hltR->accept(iTrigger) &&
326     !hltR->error(iTrigger) ){
327     triggerPassed = true;
328     }
329     return triggerPassed;
330     }
331    
332    
333 naodell 1.1 //define this as a plug-in
334     DEFINE_FWK_MODULE(MPIntuple);