Fork me on GitHub

source: git/examples/Example7.C@ b4786d3

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

added snowmass examples

  • Property mode set to 100644
File size: 2.8 KB
Line 
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");
32
33 // Book histograms
34 TH1 *histJetPT = new TH1F("jet_pt", "jet P_{T}", 100, 0.0, 1000.0);
35 TH1 *histElectronPT = new TH1F("ele_pt", "electron P_{T}", 100, 0.0, 1000.0);
36
37 // Loop over all events
38 for(Int_t entry = 0; entry < numberOfEntries; ++entry)
39 {
40 // Load selected branches with data from specified event
41 treeReader->ReadEntry(entry);
42
43 Double_t w =1.0;
44 // read MC event weight
45 if(branchWeight->GetEntries() > 0)
46 {
47 Weight *weight = (Weight*) branchWeight -> At(0);
48 w = weight->Weight;
49 }
50
51 // If event contains at least 1 jet
52 if(branchJet->GetEntries() > 0)
53 {
54 // Take first jet
55 Jet *jet = (Jet*) branchJet->At(0);
56
57 // 0 - Loose , 1 - Medium, 2 - Tight
58 Int_t wp = 1;
59
60 Bool_t BtagOk = ( jet->BTag & (1 << wp) );
61 Double_t pt = jet->PT;
62 Double_t eta = TMath::Abs(jet->Eta);
63
64 // require jet to be within acceptance and btagged with Medium WP
65 if (BtagOk && pt > 30. && eta < 5.)
66 histJetPT->Fill(jet->PT, w);
67 }
68
69
70 // If event contains at least 1 jet
71 if(branchElectron->GetEntries() > 0)
72 {
73 // Take first jet
74 Electron *electron = (Electron*) branchElectron->At(0);
75
76 Double_t pt = electron->PT;
77 Double_t eta = TMath::Abs(electron->Eta);
78
79 // 0.1 - Loose , 0.2 - Medium, 0.3 - Tight
80 Double_t IsoCut = 0.2;
81 Bool_t IsoOk = electron->IsolationVar < IsoCut;
82
83 // Plot jet transverse momentum
84 if (IsoOk && pt > 10. && eta < 5.)
85 histElectronPT->Fill(electron->PT, w);
86 }
87 }
88
89
90 //
91 TCanvas *cnv = new TCanvas("cnv", "cnv", 50, 50, 800, 500);
92 cnv->Divide(2, 1);
93 cnv->cd(1);
94 gStyle->SetOptStat(0);
95
96 histJetPT->Draw();
97
98 cnv->cd(2);
99 gStyle->SetOptStat(0);
100 // Show resulting histograms
101 histElectronPT->Draw();
102
103 cnv->Print("example7.png", "png");
104
105}
106
Note: See TracBrowser for help on using the repository browser.