swh:1:snp:af87cd67498ef4fe47c76ed3e7caffe5b61facaf
Tip revision: 86108d5adde6979551fcbf3d95012f72653dd102 authored by Pere Mato on 03 March 2016, 09:36:03 UTC
Update ROOT version files to v6.06/02.
Update ROOT version files to v6.06/02.
Tip revision: 86108d5
hsimpleReader.C
/// \file
/// \ingroup tutorial_tree
///
/// TTreeReader simplest example.
///
/// Read data from hsimple.root (written by hsimple.C)
/// \macro_code
/// \author Anders Eie, 2013
#include "TFile.h"
#include "TH1F.h"
#include "TTreeReader.h"
#include "TTreeReaderValue.h"
void hsimpleReader() {
// Create a histogram for the values we read.
auto myHist = new TH1F("h1","ntuple",100,-4,4);
// Open the file containing the tree.
auto myFile = TFile::Open("hsimple.root");
if (!myFile || myFile->IsZombie()) {
return;
}
// Create a TTreeReader for the tree, for instance by passing the
// TTree's name and the TDirectory / TFile it is in.
TTreeReader myReader("ntuple", myFile);
// The branch "px" contains floats; access them as myPx.
TTreeReaderValue<Float_t> myPx(myReader, "px");
// The branch "py" contains floats, too; access those as myPy.
TTreeReaderValue<Float_t> myPy(myReader, "py");
// Loop over all entries of the TTree or TChain.
while (myReader.Next()) {
// Just access the data as if myPx and myPy were iterators (note the '*'
// in front of them):
myHist->Fill(*myPx + *myPy);
}
myHist->Draw();
}