Fork me on GitHub

source: git/examples/Example1.C@ a7c3c84

ImprovedOutputFile Timing dual_readout llp
Last change on this file since a7c3c84 was 9047a02, checked in by Michele Selvaggi <michele.selvaggi@…>, 6 years ago

added accessing to weight in Example1 (commented)

  • Property mode set to 100644
File size: 2.3 KB
Line 
1/*
2Simple macro showing how to access branches from the delphes output root file,
3loop over events, and plot simple quantities such as the jet pt and the di-electron invariant
4mass.
5
6root -l examples/Example1.C'("delphes_output.root")'
7*/
8
9#ifdef __CLING__
10R__LOAD_LIBRARY(libDelphes)
11#include "classes/DelphesClasses.h"
12#include "external/ExRootAnalysis/ExRootTreeReader.h"
13#endif
14
15//------------------------------------------------------------------------------
16
17void Example1(const char *inputFile)
18{
19 gSystem->Load("libDelphes");
20
21 // Create chain of root trees
22 TChain chain("Delphes");
23 chain.Add(inputFile);
24
25 // Create object of class ExRootTreeReader
26 ExRootTreeReader *treeReader = new ExRootTreeReader(&chain);
27 Long64_t numberOfEntries = treeReader->GetEntries();
28
29 // Get pointers to branches used in this analysis
30 TClonesArray *branchJet = treeReader->UseBranch("Jet");
31 TClonesArray *branchElectron = treeReader->UseBranch("Electron");
32 TClonesArray *branchEvent = treeReader->UseBranch("Event");
33
34 // Book histograms
35 TH1 *histJetPT = new TH1F("jet_pt", "jet P_{T}", 100, 0.0, 100.0);
36 TH1 *histMass = new TH1F("mass", "M_{inv}(e_{1}, e_{2})", 100, 40.0, 140.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);
43
44 //HepMCEvent *event = (HepMCEvent*) branchEvent -> At(0);
45 //LHEFEvent *event = (LHEFEvent*) branchEvent -> At(0);
46 //Float_t weight = event->Weight;
47
48 // If event contains at least 1 jet
49 if(branchJet->GetEntries() > 0)
50 {
51 // Take first jet
52 Jet *jet = (Jet*) branchJet->At(0);
53
54 // Plot jet transverse momentum
55 histJetPT->Fill(jet->PT, weight);
56
57 // Print jet transverse momentum
58 cout << "Jet pt: "<<jet->PT << endl;
59 }
60
61 Electron *elec1, *elec2;
62
63 // If event contains at least 2 electrons
64 if(branchElectron->GetEntries() > 1)
65 {
66 // Take first two electrons
67 elec1 = (Electron *) branchElectron->At(0);
68 elec2 = (Electron *) branchElectron->At(1);
69
70 // Plot their invariant mass
71 histMass->Fill(((elec1->P4()) + (elec2->P4())).M(), weight);
72 }
73 }
74
75 // Show resulting histograms
76 histJetPT->Draw();
77 histMass->Draw();
78}
79
Note: See TracBrowser for help on using the repository browser.