This example is a variant of hsimple.C but using a TTree instead of a TNtuple.
It shows:
- how to fill a Tree with a few simple variables.
- how to read this Tree
- how to browse and analyze the Tree via the TBrowser and TTreeViewer This example can be run in many different ways:
- Using the Cling interpreter
- Using the automatic compiler interface
.L tree1.C or .L tree1.C++
tree1()
One can also run the write and read parts in two separate sessions. For example following one of the sessions above, one can start the session:
void tree1w()
{
TFile f(
"tree1.root",
"recreate");
TTree t1(
"t1",
"a simple Tree with simple variables");
t1.Branch(
"px",&px,
"px/F");
t1.Branch(
"py",&py,
"py/F");
t1.Branch(
"pz",&pz,
"pz/F");
t1.Branch(
"random",&random,
"random/D");
t1.Branch(
"ev",&ev,
"ev/I");
for (
Int_t i=0;i<10000;i++) {
pz = px*px + py*py;
ev = i;
}
}
void tree1r()
{
t1->SetBranchAddress(
"px",&px);
t1->SetBranchAddress(
"py",&py);
t1->SetBranchAddress(
"pz",&pz);
t1->SetBranchAddress(
"random",&random);
t1->SetBranchAddress(
"ev",&ev);
TH1F *hpx =
new TH1F(
"hpx",
"px distribution",100,-3,3);
TH2F *hpxpy =
new TH2F(
"hpxpy",
"py vs px",30,-3,3,30,-3,3);
}
if (
gROOT->IsBatch())
return;
t1->ResetBranchAddresses();
}
void tree1() {
tree1w();
tree1r();
}
R__EXTERN TRandom * gRandom
Using a TBrowser one can browse all ROOT objects.
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
1-D histogram with a float per channel (see TH1 documentation)}
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
2-D histogram with a float per channel (see TH1 documentation)}
Int_t Fill(Double_t)
Invalid Fill method.
virtual void Rannor(Float_t &a, Float_t &b)
Return 2 numbers distributed following a gaussian with mean=0 and sigma=1.
virtual Double_t Rndm()
Machine independent random number generator.
A TTree represents a columnar dataset.
- Author
- Rene Brun
Definition in file tree1.C.