Fork me on GitHub

source: git/examples/Example7.C@ 4006893

Last change on this file since 4006893 was c6f9311, checked in by Michele Selvaggi <michele.selvaggi@…>, 4 years ago

fixed Example7.C

  • Property mode set to 100644
File size: 3.1 KB
RevLine 
[f8e61b2]1/*
2Simple macro showing how to access branches from the delphes output root file for Snowmass studies,
3loop over events, and plot simple quantities such as the leading electron and jet pt.
4
5root -l examples/Example7.C'("delphes_output.root")'
6*/
7
8#ifdef __CLING__
9R__LOAD_LIBRARY(libDelphes)
10#include "classes/DelphesClasses.h"
11#include "external/ExRootAnalysis/ExRootTreeReader.h"
12#endif
13
14//------------------------------------------------------------------------------
15
16void Example7(const char *inputFile)
17{
18 gSystem->Load("libDelphes");
19
20 // Create chain of root trees
21 TChain chain("Delphes");
22 chain.Add(inputFile);
23
24 // Create object of class ExRootTreeReader
25 ExRootTreeReader *treeReader = new ExRootTreeReader(&chain);
26 Long64_t numberOfEntries = treeReader->GetEntries();
27
28 // Get pointers to branches used in this analysis
29 TClonesArray *branchJet = treeReader->UseBranch("JetPUPPITight");
30 TClonesArray *branchElectron = treeReader->UseBranch("ElectronMedium");
31 TClonesArray *branchWeight = treeReader->UseBranch("Weight");
[b03c0d0]32 TClonesArray *branchEvent = treeReader->UseBranch("Event");
[f8e61b2]33
34 // Book histograms
35 TH1 *histJetPT = new TH1F("jet_pt", "jet P_{T}", 100, 0.0, 1000.0);
36 TH1 *histElectronPT = new TH1F("ele_pt", "electron P_{T}", 100, 0.0, 1000.0);
37
38 // Loop over all events
39 for(Int_t entry = 0; entry < numberOfEntries; ++entry)
40 {
41 // Load selected branches with data from specified event
42 treeReader->ReadEntry(entry);
[b03c0d0]43
44 // main MC event weight
45 HepMCEvent *event = (HepMCEvent*) branchEvent -> At(0);
46 Double_t w = event->Weight;
47
48 // read lhe event weights
[c6f9311]49 for(Int_t i = 0; i < branchWeight->GetEntriesFast(); ++i)
[f8e61b2]50 {
[c6f9311]51 Weight *weight = (Weight*) branchWeight -> At(i);
[b03c0d0]52 Double_t lhe_weight = weight->Weight;
53
54 //cout<<lhe_weight<<endl;
55 // do stuff ...
[f8e61b2]56 }
[b03c0d0]57
[f8e61b2]58 // If event contains at least 1 jet
59 if(branchJet->GetEntries() > 0)
60 {
61 // Take first jet
62 Jet *jet = (Jet*) branchJet->At(0);
63
64 // 0 - Loose , 1 - Medium, 2 - Tight
65 Int_t wp = 1;
66
67 Bool_t BtagOk = ( jet->BTag & (1 << wp) );
68 Double_t pt = jet->PT;
69 Double_t eta = TMath::Abs(jet->Eta);
70
71 // require jet to be within acceptance and btagged with Medium WP
72 if (BtagOk && pt > 30. && eta < 5.)
73 histJetPT->Fill(jet->PT, w);
74 }
75
76 // If event contains at least 1 jet
77 if(branchElectron->GetEntries() > 0)
78 {
79 // Take first jet
80 Electron *electron = (Electron*) branchElectron->At(0);
81
82 Double_t pt = electron->PT;
83 Double_t eta = TMath::Abs(electron->Eta);
84
85 // 0.1 - Loose , 0.2 - Medium, 0.3 - Tight
86 Double_t IsoCut = 0.2;
87 Bool_t IsoOk = electron->IsolationVar < IsoCut;
88
89 // Plot jet transverse momentum
90 if (IsoOk && pt > 10. && eta < 5.)
91 histElectronPT->Fill(electron->PT, w);
92 }
93 }
94
95
96 //
97 TCanvas *cnv = new TCanvas("cnv", "cnv", 50, 50, 800, 500);
98 cnv->Divide(2, 1);
99 cnv->cd(1);
100 gStyle->SetOptStat(0);
101
102 histJetPT->Draw();
103
104 cnv->cd(2);
105 gStyle->SetOptStat(0);
106 // Show resulting histograms
107 histElectronPT->Draw();
108
109 cnv->Print("example7.png", "png");
110
111}
112
Note: See TracBrowser for help on using the repository browser.