ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/csander/HEPTutorial/MyAnalysis.C
Revision: 1.20
Committed: Thu May 31 20:00:45 2012 UTC (12 years, 11 months ago) by csander
Content type: text/plain
Branch: MAIN
Changes since 1.19: +69 -24 lines
Log Message:
more event lists

File Contents

# User Rev Content
1 csander 1.1 #define MyAnalysis_cxx
2     // The class definition in MyAnalysis.h has been generated automatically
3     // by the ROOT utility TTree::MakeSelector(). This class is derived
4     // from the ROOT class TSelector. For more information on the TSelector
5     // framework see $ROOTSYS/README/README.SELECTOR or the ROOT User Manual.
6    
7     // The following methods are defined in this file:
8     // Begin(): called every time a loop on the tree starts,
9     // a convenient place to create your histograms.
10     // SlaveBegin(): called after Begin(), when on PROOF called only on the
11     // slave servers.
12     // Process(): called for each event, in this function you decide what
13     // to read and fill your histograms.
14     // SlaveTerminate: called at the end of the loop on the tree, when on PROOF
15     // called only on the slave servers.
16     // Terminate(): called at the end of the loop on the tree,
17     // a convenient place to draw/fit your histograms.
18     //
19     // To use this file, try the following session on your Tree T:
20     //
21     // Root > T->Process("MyAnalysis.C")
22     // Root > T->Process("MyAnalysis.C","some options")
23     // Root > T->Process("MyAnalysis.C+")
24     //
25    
26     #include "MyAnalysis.h"
27 csander 1.3 #include <iostream>
28 csander 1.2 #include <TH1F.h>
29 csander 1.9 #include <TLatex.h>
30 csander 1.1
31 csander 1.3 using namespace std;
32    
33     void MyAnalysis::BuildEvent() {
34 csander 1.19 Jets.clear();
35     for (int i = 0; i < NJet; ++i) {
36     MyJet jet(Jet_Px[i], Jet_Py[i], Jet_Pz[i], Jet_E[i]);
37     jet.SetBTagDiscriminator(Jet_btag[i]);
38     jet.SetJetID(Jet_ID[i]);
39 csander 1.20 jet.SetMatchedGenJet(Jet_Px_gen[i], Jet_Py_gen[i], Jet_Pz_gen[i], Jet_E_gen[i]);
40 csander 1.19 jet.SetFlavor(Jet_Flavor_gen[i]);
41     Jets.push_back(jet);
42     }
43     Muons.clear();
44     for (int i = 0; i < NMuon; ++i) {
45     MyMuon muon(Muon_Px[i], Muon_Py[i], Muon_Pz[i], Muon_E[i]);
46     muon.SetIsolation(Muon_Iso[i]);
47     muon.SetCharge(Muon_Charge[i]);
48     muon.SetChargeGen(Muon_Charge_gen[i]);
49 csander 1.20 muon.SetMatchedGenMuon(Muon_Px_gen[i], Muon_Py_gen[i], Muon_Pz_gen[i], Muon_E_gen[i]);
50 csander 1.19 Muons.push_back(muon);
51     }
52     Electrons.clear();
53     for (int i = 0; i < NElectron; ++i) {
54 csander 1.20 MyElectron electron(Electron_Px[i], Electron_Py[i], Electron_Pz[i], Electron_E[i]);
55 csander 1.19 electron.SetIsolation(Electron_Iso[i]);
56     electron.SetCharge(Electron_Charge[i]);
57     electron.SetChargeGen(Electron_Charge_gen[i]);
58 csander 1.20 electron.SetMatchedGenElectron(Electron_Px_gen[i], Electron_Py_gen[i], Electron_Pz_gen[i], Electron_E_gen[i]);
59 csander 1.19 Electrons.push_back(electron);
60     }
61     Photons.clear();
62     for (int i = 0; i < NPhoton; ++i) {
63     MyPhoton photon(Photon_Px[i], Photon_Py[i], Photon_Pz[i], Photon_E[i]);
64     photon.SetIsolation(Photon_Iso[i]);
65 csander 1.20 photon.SetMatchedGenPhoton(Photon_Px_gen[i], Photon_Py_gen[i], Photon_Pz_gen[i], Photon_E_gen[i]);
66 csander 1.19 Photons.push_back(photon);
67     }
68     Taus.clear();
69     for (int i = 0; i < NTau; ++i) {
70     MyTau tau(Tau_Px[i], Tau_Py[i], Tau_Pz[i], Tau_E[i]);
71     tau.SetTauID(Jet_ID[i]);
72 csander 1.20 tau.SetMatchedGenJet(Tau_Px_gen[i], Tau_Py_gen[i], Tau_Pz_gen[i], Tau_E_gen[i]);
73 csander 1.19 Taus.push_back(tau);
74     }
75     met.SetXYZM(MET_px, MET_py, 0., 0.);
76     genmet.SetXYZM(MET_px_gen, MET_py_gen, 0., 0.);
77    
78     if (EventWeight < 1.e-5)
79     EventWeight = 1.;
80     EventWeight *= weight_factor;
81 csander 1.9
82 csander 1.3 }
83    
84 csander 1.2 void MyAnalysis::Begin(TTree * /*tree*/) {
85 csander 1.19 // The Begin() function is called at the start of the query.
86     // When running with PROOF Begin() is only called on the client.
87     // The tree argument is deprecated (on PROOF 0 is passed).
88 csander 1.1
89 csander 1.19 TString option = GetOption();
90 csander 1.1
91     }
92    
93 csander 1.2 void MyAnalysis::SlaveBegin(TTree * /*tree*/) {
94 csander 1.19 // The SlaveBegin() function is called after the Begin() function.
95     // When running with PROOF SlaveBegin() is called on each slave server.
96     // The tree argument is deprecated (on PROOF 0 is passed).
97    
98     TString option = GetOption();
99    
100     h_Mmumu = new TH1F("Mmumu", "Invariant di-muon mass", 100, 0, 1000);
101     h_Mmumu->SetXTitle("m_{#mu#mu}");
102     h_Mmumu->Sumw2();
103    
104     h_NJet = new TH1F("NJet", "Number of jets", 20, 0, 20);
105     h_NJet->SetXTitle("NJet");
106     h_NJet->Sumw2();
107    
108     h_NBJet = new TH1F("NBJet", "Number of b-jets", 5, 0, 5);
109     h_NBJet->SetXTitle("NBJet");
110     h_NBJet->Sumw2();
111    
112 csander 1.20 h_Jet1_Pt = new TH1F("Jet1Pt", "Pt of jet1", 50, 0, 1000);
113 csander 1.19 h_Jet1_Pt->SetXTitle("p_{T} [GeV]");
114     h_Jet1_Pt->Sumw2();
115    
116 csander 1.20 h_Jet2_Pt = new TH1F("Jet2Pt", "Pt of jet2", 50, 0, 1000);
117 csander 1.19 h_Jet2_Pt->SetXTitle("p_{T} [GeV]");
118     h_Jet2_Pt->Sumw2();
119    
120 csander 1.20 h_Jet3_Pt = new TH1F("Jet3Pt", "Pt of jet3", 50, 0, 1000);
121 csander 1.19 h_Jet3_Pt->SetXTitle("p_{T} [GeV]");
122     h_Jet3_Pt->Sumw2();
123    
124     h_Jet1_Eta = new TH1F("Jet1Eta", "#eta of jet1", 50, -5, 5);
125     h_Jet1_Eta->SetXTitle("#eta");
126     h_Jet1_Eta->Sumw2();
127    
128     h_Jet2_Eta = new TH1F("Jet2Eta", "#eta of jet2", 50, -5, 5);
129     h_Jet2_Eta->SetXTitle("#eta");
130     h_Jet2_Eta->Sumw2();
131    
132     h_Jet3_Eta = new TH1F("Jet3Eta", "#eta of jet3", 50, -5, 5);
133     h_Jet3_Eta->SetXTitle("#eta");
134     h_Jet3_Eta->Sumw2();
135    
136 csander 1.20 h_BJet1_Pt = new TH1F("BJet1Pt", "Pt of b-jet1", 50, 0, 1000);
137 csander 1.19 h_BJet1_Pt->SetXTitle("p_{T} [GeV]");
138     h_BJet1_Pt->Sumw2();
139    
140 csander 1.20 h_BJet2_Pt = new TH1F("BJet2Pt", "Pt of b-jet2", 50, 0, 1000);
141 csander 1.19 h_BJet2_Pt->SetXTitle("p_{T} [GeV]");
142     h_BJet2_Pt->Sumw2();
143    
144 csander 1.20 h_BJet1_Eta = new TH1F("BJet1Eta", "#eta of b-jet1", 50, -5, 5);
145 csander 1.19 h_BJet1_Eta->SetXTitle("#eta");
146     h_BJet1_Eta->Sumw2();
147    
148 csander 1.20 h_BJet2_Eta = new TH1F("BJet2Eta", "#eta of b-jet2", 50, -5, 5);
149 csander 1.19 h_BJet2_Eta->SetXTitle("#eta");
150     h_BJet2_Eta->Sumw2();
151    
152 csander 1.20 h_MET = new TH1F("MET", "MET", 50, 0, 2000);
153 csander 1.19 h_MET->SetXTitle("MET");
154     h_MET->Sumw2();
155    
156 csander 1.20 h_MHT = new TH1F("MHT", "MHT", 50, 0, 2000);
157 csander 1.19 h_MHT->SetXTitle("MHT");
158     h_MHT->Sumw2();
159    
160 csander 1.20 h_HT = new TH1F("HT", "HT", 50, 0, 5000);
161 csander 1.19 h_HT->SetXTitle("HT");
162     h_HT->Sumw2();
163 csander 1.20
164     h_DeltaPhiJet1Jet2 = new TH1F("DeltaPhiJet1Jet2", "DeltaPhiJet1Jet2", 100, -TMath::Pi(), TMath::Pi());
165     h_DeltaPhiJet1Jet2->SetXTitle("DeltaPhiJet1Jet2");
166     h_DeltaPhiJet1Jet2->Sumw2();
167    
168     h_Jet1OverHT = new TH1F("h_Jet1OverHT", "h_Jet1OverHT", 50, 0, 1.0);
169     h_Jet1OverHT->SetXTitle("h_Jet1OverHT");
170     h_Jet1OverHT->Sumw2();
171    
172     h_Jet2OverHT = new TH1F("h_Jet2OverHT", "h_Jet2OverHT", 50, 0, 1.0);
173     h_Jet2OverHT->SetXTitle("h_Jet2OverHT");
174     h_Jet2OverHT->Sumw2();
175    
176     h_Jet12OverHT = new TH1F("h_Jet12OverHT", "h_Jet12OverHT", 50, 0, 1.0);
177     h_Jet12OverHT->SetXTitle("h_Jet12OverHT");
178     h_Jet12OverHT->Sumw2();
179    
180     h_MHTOverHT = new TH1F("h_MHTOverHT", "h_MHTOverHT", 50, 0, 1.0);
181     h_MHTOverHT->SetXTitle("h_MHTOverHT");
182     h_MHTOverHT->Sumw2();
183    
184     h_Muon1_Pt = new TH1F("Muon1Pt", "Pt of muon1", 50, 0, 1000);
185     h_Muon1_Pt->SetXTitle("p_{T} [GeV]");
186     h_Muon1_Pt->Sumw2();
187    
188     h_Electron1_Pt = new TH1F("Electron1Pt", "Pt of electron1", 50, 0, 1000);
189     h_Electron1_Pt->SetXTitle("p_{T} [GeV]");
190     h_Electron1_Pt->Sumw2();
191    
192 csander 1.1 }
193    
194 csander 1.2 Bool_t MyAnalysis::Process(Long64_t entry) {
195 csander 1.19 // The Process() function is called for each entry in the tree (or possibly
196     // keyed object in the case of PROOF) to be processed. The entry argument
197     // specifies which entry in the currently loaded tree is to be processed.
198     // It can be passed to either MyAnalysis::GetEntry() or TBranch::GetEntry()
199     // to read either all or the required parts of the data. When processing
200     // keyed objects with PROOF, the object is already loaded and is available
201     // via the fObject pointer.
202     //
203     // This function should contain the "body" of the analysis. It can contain
204     // simple or elaborate selection criteria, run algorithms on the data
205     // of the event and typically fill histograms.
206     //
207     // The processing can be stopped by calling Abort().
208     //
209     // Use fStatus to set the return value of TTree::Process().
210     //
211     // The return value is currently not used.
212    
213     ++TotalEvents;
214    
215     GetEntry(entry);
216    
217     Event event(EventNumber, LumiNumber, RunNumber);
218    
219     if (TotalEvents % 10000 == 0)
220     cout << "Next event -----> " << TotalEvents << endl;
221    
222     if (EventList != 0) {
223     if (!EventList->IsInList(event))
224     return 0;
225     }
226    
227     BuildEvent();
228    
229     // cout << "Jets: " << endl;
230     // for (vector<MyJet>::iterator it = Jets.begin(); it != Jets.end(); ++it) {
231     // cout << "pt, eta, phi, btag: " << it->Pt() << ", " << it->Eta() << ", " << it->Phi() << ", " << it->IsBTagged()
232     // << endl;
233     // }
234     // cout << "Muons: " << endl;
235     // for (vector<MyMuon>::iterator it = Muons.begin(); it != Muons.end(); ++it) {
236     // cout << "pt, eta, phi, iso, charge: " << it->Pt() << ", " << it->Eta() << ", " << it->Phi() << ", "
237     // << it->GetIsolation() << ", " << it->GetCharge() << endl;
238     // }
239     // cout << "Electrons: " << endl;
240     // for (vector<MyElectron>::iterator it = Electrons.begin(); it != Electrons.end(); ++it) {
241     // cout << "pt, eta, phi, iso, charge: " << it->Pt() << ", " << it->Eta() << ", " << it->Phi() << ", "
242     // << it->GetIsolation() << ", " << it->GetCharge() << endl;
243     // }
244     // cout << "Photons: " << endl;
245     // for (vector<MyPhoton>::iterator it = Photons.begin(); it != Photons.end(); ++it) {
246     // cout << "pt, eta, phi, iso: " << it->Pt() << ", " << it->Eta() << ", " << it->Phi() << ", " << it->GetIsolation()
247     // << endl;
248     // }
249    
250     //////////////////////////////
251     int N_BJet = 0;
252     int N_Jet = 0;
253     double HT = 0;
254     TLorentzVector MHT(0, 0, 0, 0);
255    
256     for (vector<MyJet>::iterator it = Jets.begin(); it != Jets.end(); ++it) {
257     ++N_Jet;
258    
259     if (!it->GetJetID() && it->Pt() > 30)
260     return 0;
261    
262     if (it->Pt() > 30)
263     MHT += *it;
264    
265     if (it->Pt() > 50 && fabs(it->Eta()) < 2.5)
266     HT += it->Pt();
267    
268     if (N_Jet == 1) {
269     h_Jet1_Pt->Fill(it->Pt(), EventWeight);
270     h_Jet1_Eta->Fill(it->Eta(), EventWeight);
271     } else if (N_Jet == 2) {
272     h_Jet2_Pt->Fill(it->Pt(), EventWeight);
273     h_Jet2_Eta->Fill(it->Eta(), EventWeight);
274     } else if (N_Jet == 3) {
275     h_Jet3_Pt->Fill(it->Pt(), EventWeight);
276     h_Jet3_Eta->Fill(it->Eta(), EventWeight);
277     }
278     if (it->IsBTagged()) {
279     ++N_BJet;
280     if (N_BJet == 1) {
281     h_BJet1_Pt->Fill(it->Pt(), EventWeight);
282     h_BJet1_Eta->Fill(it->Eta(), EventWeight);
283     } else if (N_BJet == 2) {
284     h_BJet2_Pt->Fill(it->Pt(), EventWeight);
285     h_BJet2_Eta->Fill(it->Eta(), EventWeight);
286     }
287     }
288     }
289     h_NBJet->Fill(N_BJet, EventWeight);
290     h_NJet->Fill(N_Jet, EventWeight);
291    
292     h_MET->Fill(met.Pt(), EventWeight);
293     h_MHT->Fill(MHT.Pt(), EventWeight);
294     h_HT->Fill(HT, EventWeight);
295     //////////////////////////////
296    
297     //////////////////////////////
298 csander 1.20 if (Jets.size() > 1) {
299     h_DeltaPhiJet1Jet2->Fill(Jets.at(0).DeltaPhi(Jets.at(1)), EventWeight);
300     h_Jet2OverHT->Fill(Jets.at(1).Pt() / HT, EventWeight);
301     h_Jet12OverHT->Fill((Jets.at(0).Pt() + Jets.at(1).Pt()) / HT, EventWeight);
302     }
303     h_Jet1OverHT->Fill(Jets.at(0).Pt() / HT, EventWeight);
304     h_MHTOverHT->Fill(MHT.Pt() / HT, EventWeight);
305     //////////////////////////////
306    
307     //////////////////////////////
308     if (NMuon > 0) {
309     if (Muons.at(0).IsIsolated()) {
310     h_Muon1_Pt->Fill(Muons.at(0).Pt(), EventWeight);
311     }
312     }
313     if (NElectron > 0) {
314     if (Electrons.at(0).IsIsolated()) {
315     h_Electron1_Pt->Fill(Electrons.at(0).Pt(), EventWeight);
316     }
317     }
318     //////////////////////////////
319    
320     //////////////////////////////
321 csander 1.19 if (NMuon > 1) {
322 csander 1.20 if (Muons.at(0).GetCharge() * Muons.at(1).GetCharge() < 1 && Muons.at(0).IsIsolated() && Muons.at(1).IsIsolated()) {
323 csander 1.19 h_Mmumu->Fill((Muons.at(0) + Muons.at(1)).M(), EventWeight);
324     }
325     }
326     //////////////////////////////
327 csander 1.8
328 csander 1.19 return kTRUE;
329 csander 1.1 }
330    
331 csander 1.2 void MyAnalysis::SlaveTerminate() {
332 csander 1.19 // The SlaveTerminate() function is called after all entries or objects
333     // have been processed. When running with PROOF SlaveTerminate() is called
334     // on each slave server.
335 csander 1.1
336     }
337    
338 csander 1.2 void MyAnalysis::Terminate() {
339 csander 1.19 // The Terminate() function is the last function to be called during
340     // a query. It always runs on the client, it can be used to present
341     // the results graphically or save the results to file.
342 csander 1.1
343     }