Logo ROOT   6.10/09
Reference Guide
hsimpleReader.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_tree
3 /// \notebook
4 /// TTreeReader simplest example.
5 ///
6 /// Read data from hsimple.root (written by hsimple.C)
7 ///
8 /// \macro_code
9 ///
10 /// \author Anders Eie, 2013
11 
12 #include "TFile.h"
13 #include "TH1F.h"
14 #include "TTreeReader.h"
15 #include "TTreeReaderValue.h"
16 
17 void hsimpleReader() {
18  // Create a histogram for the values we read.
19  auto myHist = new TH1F("h1","ntuple",100,-4,4);
20 
21  // Open the file containing the tree.
22  auto myFile = TFile::Open("hsimple.root");
23  if (!myFile || myFile->IsZombie()) {
24  return;
25  }
26  // Create a TTreeReader for the tree, for instance by passing the
27  // TTree's name and the TDirectory / TFile it is in.
28  TTreeReader myReader("ntuple", myFile);
29 
30  // The branch "px" contains floats; access them as myPx.
31  TTreeReaderValue<Float_t> myPx(myReader, "px");
32  // The branch "py" contains floats, too; access those as myPy.
33  TTreeReaderValue<Float_t> myPy(myReader, "py");
34 
35  // Loop over all entries of the TTree or TChain.
36  while (myReader.Next()) {
37  // Just access the data as if myPx and myPy were iterators (note the '*'
38  // in front of them):
39  myHist->Fill(*myPx + *myPy);
40  }
41 
42  myHist->Draw();
43 }
TTreeReader is a simple, robust and fast interface to read values from a TTree, TChain or TNtuple...
Definition: TTreeReader.h:43
THist< 1, float, THistStatContent, THistStatUncertainty > TH1F
Definition: THist.hxx:311
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=1, Int_t netopt=0)
Create / open a file.
Definition: TFile.cxx:3909