How to read a Tree ?

There are 2 ways of reading a Tree and to make selections. Both ways will be illustrated by using the file generated in the previous section "How to Create and Fill a Tree".

The simple way via TTree::Draw

Use this function if you want to quickly histogram a few variables with a very simple selection mechanism. The loop on the Tree is automatically done by ROOT.

void Draw(Text_t *varexp, Text_t *selection, 
          Option_t *option,Int_t nevents, Int_t firstevent);

where the parameters have the following meaning:

  • varexp is an expression with the forms e1, e1:e2 or e1:e2:e3 to generate 1-D, 2-D or 3-D distributions respectively and where e1,e2,e3 are expressions with combinations of the original Tree variables. Examples:
    varexp = x         : draw a 1-Dim distribution of column named x
           = sqrt(x)   : draw distribution of sqrt(x)
           = x*y/z
           = y:sqrt(x) : 2-Dim distribution of y versus sqrt(x)

     

    • selection is an expression with a combination of the Tree variables. In a selection all the C++ operators are authorized. Example:
            selection = "x<y && sqrt(z)>3.2"
    • option is the drawing option. See TH1::Draw for the list of all drawing options. If option contains the string "goff", no graphics is generated.
    • nevents is the number of events to process (default is all).
    • first is the first event to process (default is 0).

    By default the temporary histogram created is called htemp. If varexp contains >>hnew, the new histogram created is called hnew and it is kept in the current directory. Example:

        tree.Draw("sqrt(x)>>hsqrt","y>0");
    will draw sqrt(x) and save the histogram as "hsqrt" in the current directory.

    By default, the specified histogram is reset. To continue to append data to an existing histogram, use "+" in front of the histogram name.

     
        tree.Draw("sqrt(x)>>+hsqrt","y>0");
    will not reset hsqrt, but will continue filling.

    Look at the tutorial macro ntuple1 for a more detailed example.

    The kind of expressions you can specify in this function is restricted on purpose (see class TFormula). Our previous experience with the query processor of PAW indicates that we better restrict this kind of selection mechanism to trivial cases if we want our users to trust their results.

    Note that this method does not require access to the original class code, if you only want to perform the kind of operations described above.

    The general way

    Again, at the difference of PAW, we believe that the ROOT user must be responsible to organize himself the loop on all events (rows) of a Tree. In a loop, in general, one wants to generate many histograms, or correlate variables from one Tree with variables of another Tree, etc. We reproduce below a slightly simplified version of the macro eventb. This macro illustrates:

    • How to loop on a tree.
    • How to read a sub-branch only (the branch containing the number of tracks).
    • When the selection is positive, how to read all the branches.
    • How to histogram the result of a member function of the event object.
    {
       gROOT->Reset();
     
    //   Connect file generated in $ROOTSYS/test
       TFile f("Event.root");
     
    //   Create an histogram
       TH1F *hnseg = new TH1F("hnseg","Number of segments for selected tracks",
                               5000,0,5000);
     
    //   Start main loop on all events
       TClonesArray *tracks = 0;
       Event *event = new Event(); //object must be created before
                                   //setting the branch address
     
       TBranch *bntrack = T.GetBranch("fNtrack");
       TBranch *branch  = T.GetBranch("event");
       branch->SetAddress(&event);
       Int_t nevent = T.GetEntries();
       Int_t nselected = 0;
       for (Int_t i=0;i<nevent;i++) {
          bntrack->GetEvent(i);
          if (event->GetNtrack() > 587)continue;
          nselected++;
          T.GetEvent(i);                  //Read complete accepted event
                                          //in memory.
          hnseg->Fill(event->GetNseg());  //Fill histogram with number of
                                          //segments.
          tracks = event->GetTracks();    //Get pointer to the
                                          //TClonesArray object.
          tracks->Clear();                //Clear it.
       }
    //  Delete ClonesArray and histogram objects
       event->Finish();
     
    //  Draw the histogram
       hnseg->Draw();
    }

    This method requires access to the class implementation when one of the branches is an object.

    Note that in the above example, the object of class Event must be created via the new operator. It cannot be a static object because a new object will be recreated for each event.

    ROOT includes a function TTree::MakeCode to automatically generate the code for a skeleton analysis function.