Fork me on GitHub

source: git/examples/Example1.C@ a190d94

ImprovedOutputFile Timing dual_readout llp
Last change on this file since a190d94 was fbb617b, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

improve formatting and remove trailing spaces

  • Property mode set to 100644
File size: 1.9 KB
RevLine 
[e7e90df]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
[fbb617b]4mass.
[e7e90df]5
6root -l examples/Example1.C'("delphes_output.root")'
7*/
8
9//------------------------------------------------------------------------------
10
11void Example1(const char *inputFile)
12{
13 gSystem->Load("libDelphes");
14
15 // Create chain of root trees
16 TChain chain("Delphes");
17 chain.Add(inputFile);
[fbb617b]18
[e7e90df]19 // Create object of class ExRootTreeReader
20 ExRootTreeReader *treeReader = new ExRootTreeReader(&chain);
21 Long64_t numberOfEntries = treeReader->GetEntries();
[fbb617b]22
[e7e90df]23 // Get pointers to branches used in this analysis
24 TClonesArray *branchJet = treeReader->UseBranch("Jet");
25 TClonesArray *branchElectron = treeReader->UseBranch("Electron");
[fbb617b]26
[e7e90df]27 // Book histograms
28 TH1 *histJetPT = new TH1F("jet_pt", "jet P_{T}", 100, 0.0, 100.0);
29 TH1 *histMass = new TH1F("mass", "M_{inv}(e_{1}, e_{2})", 100, 40.0, 140.0);
30
31 // Loop over all events
32 for(Int_t entry = 0; entry < numberOfEntries; ++entry)
33 {
34 // Load selected branches with data from specified event
35 treeReader->ReadEntry(entry);
[fbb617b]36
[e7e90df]37 // If event contains at least 1 jet
38 if(branchJet->GetEntries() > 0)
39 {
40 // Take first jet
41 Jet *jet = (Jet*) branchJet->At(0);
[fbb617b]42
[e7e90df]43 // Plot jet transverse momentum
44 histJetPT->Fill(jet->PT);
[fbb617b]45
[e7e90df]46 // Print jet transverse momentum
47 cout << "Jet pt: "<<jet->PT << endl;
48 }
49
50 Electron *elec1, *elec2;
51
52 // If event contains at least 2 electrons
53 if(branchElectron->GetEntries() > 1)
54 {
55 // Take first two electrons
56 elec1 = (Electron *) branchElectron->At(0);
57 elec2 = (Electron *) branchElectron->At(1);
58
59 // Plot their invariant mass
60 histMass->Fill(((elec1->P4()) + (elec2->P4())).M());
61 }
62 }
63
64 // Show resulting histograms
65 histJetPT->Draw();
66 histMass->Draw();
67}
68
Note: See TracBrowser for help on using the repository browser.