ROOT Tutorial II¶

  • TTree

TTrees¶

The TTree is the data structure ROOT provides to store large quantities of same types objects

  • Organised in branches, each one holding objects
  • Organised in independent events, e.g. collision events
  • Efficient disk space usage, optimised I/O runtime

Example¶

We’ll assume we’re recording the missing transverse energy, momentum (transverse momentum, theta, phi) for jets, where there could be between 0 and 5 jets in each event. This means we will define some C++ variables that will be used in the program.

image
In [1]:
Float_t met;    // Missing energy in the transverse direction.
Int_t njets;    // Necessary to keep track of the number of jets
Float_t pt[16]; // assume not write information for more than 16 jets. 
Float_t theta[16];
Float_t phi[16];

Write events to a TTree¶

In [2]:
TFile *f = new TFile("tree.root","recreate");    // Create a ROOT file, f.

Now we define the TTree object which will hold all of our variables and the data they represent.

  • A TTree object called tree.
  • The first argument is the name of the object as stored by ROOT.
  • The second argument is a short descriptor.
In [3]:
TTree *tree = new TTree("tree","A simple Tree with simple variables");

We now define the TBranch for the met variable.

  • The first argument is ROOT's internal name of the variable.
  • The second argument is the address of the actual variable we defined above
  • The third argument defines the type of the variable to be stored, and the F at the end signifies that this is a float.

The variable type must be 1 character:

  • I : a 32 bit signed integer (Int_t)
  • F : a 32 bit floating point (Float_t)
  • D : a 64 bit floating point (Double_t)
  • L : a 64 bit signed integer (Long64_t)
  • ...
In [4]:
tree->Branch("met", &met, "met/F");

Next we define the TBranches for each of the other variables, but the syntax is slightly different as these are acting as arrays with a varying number of entries for each event.

First we define njets where the syntax is the same as before, except we take care to identify this as an integer with the final I designation.

In [5]:
tree->Branch("njets", &njets, "njets/I");

We can now define the other variables, but we use a slightly different syntax for the third argument to identify the variable that will be used to count the number of entries per event.

In [6]:
tree->Branch("pt", &pt, "pt[njets]/F");
tree->Branch("theta", &theta, "theta[njets]/F");
tree->Branch("phi", &phi, "phi[njets]/F");
In [7]:
Int_t nevents = 100000;

for (Int_t i=0; i<nevents; i++) {

 // Generate random number between 10-60 (arbitrary)
    met = gRandom->Gaus(5,1) + 10;
    if(met < 0) continue;
    // Generate random number between 0-5, inclusive
    njets = gRandom->Integer(6);
    Double_t x, y, z;
    for (Int_t j=0; j<njets; j++) {
        pt[j] =  100 * gRandom->Exp(1);
        while(1){
            gRandom->Sphere(x,y,z,1);
            if(z>0) break;
        }
        TVector3 v(x,y,z);
        theta[j] = v.Theta() * 180./TMath::Pi();
        phi[j] = v.Phi() * 180./TMath::Pi();;
    }
    // After each event we need to *fill* the TTree
    if(njets > 0) tree->Fill();
}

After we've run over all the events, we "change directory" (cd) to the file object and write the tree to it. We can also print the tree, just as a visual identifier to ourselves that the program ran to completion.

In [8]:
f->cd();
tree->Write();
tree->Print();//print tree structure
f->Close();
******************************************************************************
*Tree    :tree      : A simple Tree with simple variables                    *
*Entries :    83122 : Total =         4674914 bytes  File  Size =    3519215 *
*        :          : Tree compression factor =   1.33                       *
******************************************************************************
*Br    0 :met       : met/F                                                  *
*Entries :    83122 : Total  Size=     333823 bytes  File Size  =     280745 *
*Baskets :       11 : Basket Size=      32000 bytes  Compression=   1.19     *
*............................................................................*
*Br    1 :njets     : njets/I                                                *
*Entries :    83122 : Total  Size=     333853 bytes  File Size  =      55164 *
*Baskets :       11 : Basket Size=      32000 bytes  Compression=   6.04     *
*............................................................................*
*Br    2 :pt        : pt[njets]/F                                            *
*Entries :    83122 : Total  Size=    1335605 bytes  File Size  =    1067327 *
*Baskets :       53 : Basket Size=      32000 bytes  Compression=   1.25     *
*............................................................................*
*Br    3 :theta     : theta[njets]/F                                         *
*Entries :    83122 : Total  Size=    1335776 bytes  File Size  =    1029845 *
*Baskets :       53 : Basket Size=      32000 bytes  Compression=   1.30     *
*............................................................................*
*Br    4 :phi       : phi[njets]/F                                           *
*Entries :    83122 : Total  Size=    1335662 bytes  File Size  =    1083658 *
*Baskets :       53 : Basket Size=      32000 bytes  Compression=   1.23     *
*............................................................................*

2. Put all codes into a file¶

  1. write2tree.cpp
#include "TTree.h"
#include "TFile.h"
#include "TRandom3.h"
#include "TVector3.h"

int main() {

    Float_t met; // Missing energy in the transverse direction.    
    Int_t njets; // Necessary to keep track of the number of jets
    Float_t pt[16];
    Float_t theta[16];
    Float_t phi[16];    

    TFile *f = new TFile("tree.root","recreate");
    TTree *tree = new TTree("tree","A simple Tree with simple variables");
    tree->Branch("met", &met, "met/F");
    tree->Branch("njets", &njets, "njets/I");
    tree->Branch("pt", &pt, "pt[njets]/F");
    tree->Branch("theta", &theta, "theta[njets]/F");
    tree->Branch("phi", &phi, "phi[njets]/F");
    Int_t nevents = 100000;

    for (Int_t i=0;i<nevents;i++) {

      // Generate random number between 10-60 (arbitrary)
      met = gRandom->Gaus(5,1) + 10;
      if(met < 0) continue;
      // Generate random number between 0-5, inclusive
      njets = gRandom->Integer(6);
      Double_t x, y, z;
      for (Int_t j=0; j<njets; j++) {
    pt[j] =  100 * gRandom->Exp(1);
        while(1){
      gRandom->Sphere(x,y,z,1);
      if(z>0) break;
        }
        TVector3 v(x,y,z);
        theta[j] = v.Theta() * 180./TMath::Pi();
        phi[j] = v.Phi() * 180./TMath::Pi();;
      }
      // After each event we need to *fill* the TTree
      if(njets > 0) tree->Fill();
    }

    f->cd();
    tree->Write();
    tree->Print();
    f->Close();
    return 0;
}

3. Compile¶

  • Copy the following commands to the terminal and execute them.

g++ write2tree.cpp $(root-config --glibs --cflags --libs) -o write2tree

./write2tree

In [9]:
!ls -lha *.root
-rw-r--r--  1 zhli  staff   3.4M Apr 12 09:06 out.root
-rw-r--r--  1 zhli  staff    16K Apr 12 08:11 output.root
-rw-r--r--  1 zhli  staff   3.4M Apr 12 09:12 tree.root

4. Data analysis of TTree¶

In [10]:
TFile *f = new TFile("tree.root"); //shell:       root -l tree.root
f->ls();                           //ROOT shell:  .ls
TFile**		tree.root	
 TFile*		tree.root	
  KEY: TTree	tree;1	A simple Tree with simple variables
In [11]:
TTree *tree = (TTree*)f->Get("tree"); //skip this line in ROOT shell
tree->Print();
******************************************************************************
*Tree    :tree      : A simple Tree with simple variables                    *
*Entries :    83122 : Total =         4674170 bytes  File  Size =    3519215 *
*        :          : Tree compression factor =   1.33                       *
******************************************************************************
*Br    0 :met       : met/F                                                  *
*Entries :    83122 : Total  Size=     333775 bytes  File Size  =     280745 *
*Baskets :       11 : Basket Size=      32000 bytes  Compression=   1.19     *
*............................................................................*
*Br    1 :njets     : njets/I                                                *
*Entries :    83122 : Total  Size=     333805 bytes  File Size  =      55164 *
*Baskets :       11 : Basket Size=      32000 bytes  Compression=   6.04     *
*............................................................................*
*Br    2 :pt        : pt[njets]/F                                            *
*Entries :    83122 : Total  Size=    1335389 bytes  File Size  =    1067327 *
*Baskets :       53 : Basket Size=      32000 bytes  Compression=   1.25     *
*............................................................................*
*Br    3 :theta     : theta[njets]/F                                         *
*Entries :    83122 : Total  Size=    1335560 bytes  File Size  =    1029845 *
*Baskets :       53 : Basket Size=      32000 bytes  Compression=   1.30     *
*............................................................................*
*Br    4 :phi       : phi[njets]/F                                           *
*Entries :    83122 : Total  Size=    1335446 bytes  File Size  =    1083658 *
*Baskets :       53 : Basket Size=      32000 bytes  Compression=   1.23     *
*............................................................................*
In [12]:
cout<<tree->GetEntriesFast()<<endl; //get number of entries
83122
In [13]:
tree->Show(0) //show first events
======> EVENT:0
 met             = 16.1056
 njets           = 1
 pt              = 72.3661
 theta           = 89.2434
 phi             = 121.201
In [14]:
tree->Show(100)
======> EVENT:100
 met             = 14.9587
 njets           = 2
 pt              = 133.223, 
                  64.0654
 theta           = 64.3745, 
                  59.7926
 phi             = -139.611, 
                  -44.2944
In [15]:
tree->Scan("met:njets:pt:theta:phi","abs(met-15)<2 && njets<3","",10,0);
***********************************************************************************
*    Row   * Instance *       met *     njets *        pt *     theta *       phi *
***********************************************************************************
*        0 *        0 * 16.105598 *         1 * 72.366081 * 89.243362 * 121.20128 *
*        1 *        0 * 14.014933 *         1 * 74.332588 * 51.230472 * 10.211810 *
*        3 *        0 * 16.216753 *         1 *  83.96875 * 83.571914 * -164.5037 *
*        4 *        0 * 14.706235 *         1 * 81.911666 * 20.178411 * -157.4595 *
*        6 *        0 * 14.212073 *         2 * 32.323730 * 64.427314 * -116.5691 *
*        6 *        1 * 14.212073 *         2 * 141.85296 * 54.492473 * -65.29646 *
*        9 *        0 * 15.920655 *         2 * 159.61865 * 32.048381 * 146.35559 *
*        9 *        1 * 15.920655 *         2 * 30.709390 * 41.741668 * -9.437234 *
***********************************************************************************
==> 8 selected entries
In [16]:
%jsroot on
In [17]:
//for juypter-notebook
TCanvas *c1 = new TCanvas;
In [18]:
tree->Draw("njets");//Draw a plot with default settings.
c1->Draw();
In [19]:
tree->Draw("njets >> hnjets(10,0,10)");//Saving the result of Draw to a histogram
hnjets->SetMinimum(10000);
c1->Draw();
In [20]:
tree->Draw("met>>hmet(100,0,30)");
c1->Draw();
In [21]:
tree->Draw("pt>>h1(200,0,200)");
c1->Draw();
In [22]:
tree->Draw("met:pt");
c1->Draw();
In [23]:
tree->Draw("met:pt>>h2(1000,0,200,100,10,20)","","colz");
c1->Draw();
In [24]:
tree->Draw("theta");
c1->Draw();
In [25]:
tree->Draw("theta:phi","","colz");// fill to a histogram
c1->Draw();

Draw() with selection¶

In [26]:
tree->Draw("met:pt>>(1000,0,200,100,10,20)","abs(met-15)<2 && njets<3","colz");//Draw with condition
c1->Draw();

Read events from a TTree¶

Without the 'recreate' argument, ROOT will assume this file exists to be read in.We will now Get the tree from the file and assign it to a new local variable.

In [27]:
TFile *f = new TFile("tree.root");
TTree *tree = (TTree*)f->Get("tree");

Just as we did in the write2tree_file.cpp example, we will define some local variables. These variables will actually get “filled” by the ROOT file when we loop over the events.

In [28]:
Float_t met;
Int_t njets;
Float_t pt[16];
Float_t theta[16];
Float_t phi[16];

We’ll now assign these local variables to specific TBranches in tree. Note that we’ll be using the address of each local variable when we precede the variable name with an ampersand &.

In [29]:
// Assign these variables to specific branch addresses
tree->SetBranchAddress("met", &met);
tree->SetBranchAddress("njets", &njets);
tree->SetBranchAddress("pt", &pt);
tree->SetBranchAddress("theta", &theta);
tree->SetBranchAddress("phi", &phi);

We’re ready now to loop over events! Each time we call tree->GetEntry(i), it pulls the ith values out of input_tree and “fills” the local variables with those values.

In [30]:
Long64_t nentries = tree->GetEntriesFast();
for(Long64_t jentry=0; jentry<nentries; jentry++) {
    // Get the values for the i`th event and fill 
    // all our local variables that were assigned to TBranches
       tree->GetEntry(jentry);

       // Print the number of jets in this event
       if(jentry<10) { 
           cout << njets << endl;;

           // Print out the momentum for each jet in this event
           for (Int_t j=0;j<njets;j++) {
               cout<< pt[j] <<", " << theta[j] <<", " << phi[j] <<endl;;
           }
       }
}
1
72.3661, 89.2434, 121.201
1
74.3326, 51.2305, 10.2118
3
69.593, 59.8529, -118.006
276.593, 53.1034, 153.543
222.736, 84.2762, 179.629
1
83.9688, 83.5719, -164.504
1
81.9117, 20.1784, -157.46
4
100.515, 56.4845, -158.523
235.163, 13.8857, 75.4596
99.309, 24.8166, -47.3553
29.56, 41.4123, 166.01
2
32.3237, 64.4273, -116.569
141.853, 54.4925, -65.2965
4
207.514, 84.9573, -119.977
195.033, 84.7252, -177.587
98.1751, 25.9072, -6.9883
21.7239, 13.3978, 4.79973
3
135.2, 70.839, 158.167
358.29, 77.4885, 62.8894
43.1928, 54.6431, 93.6167
2
159.619, 32.0484, 146.356
30.7094, 41.7417, -9.43723

Write events to a new ROOT file¶

You can save the following code as a script to run it in a root shell, or modify it into a source code(main function) then compile and execute .

  • script: root -l newROOT.cpp; root[0] .x newROOT.cpp
  • source code: make xxx
In [31]:
void newROOT()
{
    //input 
    TFile *f = new TFile("tree.root");
    TTree *tree = (TTree*)f->Get("tree");

    Float_t met;
    Int_t njets;
    Float_t pt[16];
    Float_t theta[16];
    Float_t phi[16];

    tree->SetBranchAddress("met", &met);
    tree->SetBranchAddress("njets", &njets);
    tree->SetBranchAddress("pt", &pt);
    tree->SetBranchAddress("theta", &theta);
    tree->SetBranchAddress("phi", &phi);

    //output 
    Float_t mass;
    Int_t hits;
    Float_t p[16];
    Float_t px[16];
    Float_t py[16];
    Float_t pz[16];

    TFile *fout = new TFile("out.root","recreate");
    TTree *tout = new TTree("tree","tree");

    tout->Branch("mass", &mass, "mass/F");
    tout->Branch("hits", &hits, "hits/I");
    tout->Branch("p", &p, "p[hits]/F");
    tout->Branch("px", &px, "px[hits]/F");
    tout->Branch("py", &py, "py[hits]/F");
    tout->Branch("pz", &pz, "pz[hits]/F");

    double rad = TMath::Pi()/180.;

    //events loop
    Long64_t nentries = tree->GetEntriesFast();
    for(Long64_t jentry=0; jentry<nentries; jentry++) {
        tree->GetEntry(jentry);
    
        mass = met;
        for(int i=0; i<njets; i++) {
            Double_t theta_r = theta[i]*rad;
            Double_t phi_r = phi[i]*rad;
            p[i] = pt[i];
            px[i] = pt[i]*sin(theta_r)*cos(phi_r);
            py[i] = pt[i]*sin(theta_r)*sin(phi_r);
            pz[i] = pt[i]*cos(theta_r);
        }
        hits = njets;
        if(hits>2) tout->Fill(); //fill to the tree
    
    if(jentry%10000 == 0) cout << ".";
        
    }
    //write to the ROOT file
    tout->Write();
    fout->Close();
    cout << endl;
}
In [32]:
newROOT();
.........
In [33]:
TFile *f = new TFile("out.root");
TTree *tree = (TTree*)f->Get("tree");
In [34]:
tree->Print();
******************************************************************************
*Tree    :tree      : tree                                                   *
*Entries :    49795 : Total =         4400943 bytes  File  Size =    3573672 *
*        :          : Tree compression factor =   1.23                       *
******************************************************************************
*Br    0 :mass      : mass/F                                                 *
*Entries :    49795 : Total  Size=     200158 bytes  File Size  =     168182 *
*Baskets :        7 : Basket Size=      32000 bytes  Compression=   1.19     *
*............................................................................*
*Br    1 :hits      : hits/I                                                 *
*Entries :    49795 : Total  Size=     200158 bytes  File Size  =      26211 *
*Baskets :        7 : Basket Size=      32000 bytes  Compression=   7.62     *
*............................................................................*
*Br    2 :p         : p[hits]/F                                              *
*Entries :    49795 : Total  Size=    1000108 bytes  File Size  =     829033 *
*Baskets :       38 : Basket Size=      32000 bytes  Compression=   1.20     *
*............................................................................*
*Br    3 :px        : px[hits]/F                                             *
*Entries :    49795 : Total  Size=    1000150 bytes  File Size  =     855750 *
*Baskets :       38 : Basket Size=      32000 bytes  Compression=   1.17     *
*............................................................................*
*Br    4 :py        : py[hits]/F                                             *
*Entries :    49795 : Total  Size=    1000150 bytes  File Size  =     855823 *
*Baskets :       38 : Basket Size=      32000 bytes  Compression=   1.17     *
*............................................................................*
*Br    5 :pz        : pz[hits]/F                                             *
*Entries :    49795 : Total  Size=    1000150 bytes  File Size  =     836699 *
*Baskets :       38 : Basket Size=      32000 bytes  Compression=   1.19     *
*............................................................................*
In [35]:
tree->Draw("hits>>(10,0,10)");
c1->Draw();
In [36]:
tree->Draw("px:py>>(500,-1000,1000,500,-1000,1000)","","colz");
c1->SetLogz();
c1->Draw();
In [ ]: