library: libTree
#include "TTree.h"

TTree


class description - source file - inheritance tree (.pdf)

class TTree : public TNamed, public TAttLine, public TAttFill, public TAttMarker

Inheritance Chart:
TObject
<-
TNamed
TAttLine
TAttFill
TAttMarker
<-
TTree
<-
TChain
TChainProof
THbookTree
TNtuple
TNtupleD
TTreeSQL

    protected:
void AddClone(TTree*) virtual TBranch* BranchImp(const char* branchname, const char* classname, TClass* ptrClass, void* addobj, Int_t bufsize, Int_t splitlevel) virtual TBranch* BranchImp(const char* branchname, TClass* ptrClass, void* addobj, Int_t bufsize, Int_t splitlevel) virtual TFile* ChangeFile(TFile* file) virtual Bool_t CheckBranchAddressType(TBranch* branch, TClass* ptrClass, EDataType datatype, Bool_t ptr) const char* GetNameByIndex(TString& varexp, Int_t* index, Int_t colindex) const virtual void KeepCircular() virtual void MakeIndex(TString& varexp, Int_t* index) public:
TTree() TTree(const char* name, const char* title, Int_t splitlevel = 99) virtual ~TTree() virtual TFriendElement* AddFriend(const char* treename, const char* filename = "") virtual TFriendElement* AddFriend(const char* treename, TFile* file) virtual TFriendElement* AddFriend(TTree* tree, const char* alias = "", Bool_t warn = kFALSE) virtual void AddTotBytes(Int_t tot) virtual void AddZipBytes(Int_t zip) virtual Long64_t AutoSave(Option_t* option = "") virtual Int_t Branch(TCollection* list, Int_t bufsize = 32000, Int_t splitlevel = 99, const char* name = "") virtual Int_t Branch(TList* list, Int_t bufsize = 32000, Int_t splitlevel = 99) virtual Int_t Branch(const char* folder, Int_t bufsize = 32000, Int_t splitlevel = 99) virtual TBranch* Branch(const char* name, void* address, const char* leaflist, Int_t bufsize = 32000) TBranch* Branch(const char* name, const char* classname, void*** addobj, Int_t bufsize = 32000, Int_t splitlevel = 99) TBranch* Branch(const char* name, void*** addobj, Int_t bufsize = 32000, Int_t splitlevel = 99) virtual TBranch* BranchOld(const char* name, const char* classname, void* addobj, Int_t bufsize = 32000, Int_t splitlevel = 1) virtual TBranch* BranchRef() virtual TBranch* Bronch(const char* name, const char* classname, void* addobj, Int_t bufsize = 32000, Int_t splitlevel = 99) virtual void Browse(TBrowser* b) virtual Int_t BuildIndex(const char* majorname, const char* minorname = "0") TStreamerInfo* BuildStreamerInfo(TClass* cl, void* pointer = 0) static TClass* Class() virtual TTree* CloneTree(Long64_t nentries = -1, Option_t* option = "") virtual void CopyAddresses(TTree* tree) virtual Long64_t CopyEntries(TTree* tree, Long64_t nentries = -1) virtual TTree* CopyTree(const char* selection, Option_t* option = "", Long64_t nentries = 1000000000, Long64_t firstentry = 0) virtual TBasket* CreateBasket(TBranch* branch) Int_t Debug() const virtual void Delete(Option_t* option = "") virtual void Draw(Option_t* opt) virtual Long64_t Draw(const char* varexp, const TCut& selection, Option_t* option = "", Long64_t nentries = 1000000000, Long64_t firstentry = 0) virtual Long64_t Draw(const char* varexp, const char* selection, Option_t* option = "", Long64_t nentries = 1000000000, Long64_t firstentry = 0) virtual void DropBaskets() virtual void DropBuffers(Int_t nbytes) virtual Int_t Fill() virtual TBranch* FindBranch(const char* name) virtual TLeaf* FindLeaf(const char* name) virtual Long64_t Fit(const char* funcname, const char* varexp, const char* selection = "", Option_t* option = "", Option_t* goption = "", Long64_t nentries = 1000000000, Long64_t firstentry = 0) virtual const char* GetAlias(const char* aliasName) const virtual TBranch* GetBranch(const char* name) virtual TBranchRef* GetBranchRef() const virtual Bool_t GetBranchStatus(const char* branchname) const static Int_t GetBranchStyle() virtual Long64_t GetChainEntryNumber(Long64_t entry) const virtual Long64_t GetChainOffset() const TFile* GetCurrentFile() const Long64_t GetDebugMax() const Long64_t GetDebugMin() const TDirectory* GetDirectory() const virtual Long64_t GetEntries() const virtual Long64_t GetEntriesFast() const virtual Long64_t GetEntriesFriend() const virtual Int_t GetEntry(Long64_t entry = 0, Int_t getall = 0) virtual Long64_t GetEntryNumber(Long64_t entry) const virtual Long64_t GetEntryNumberWithBestIndex(Int_t major, Int_t minor = 0) const virtual Long64_t GetEntryNumberWithIndex(Int_t major, Int_t minor = 0) const virtual Int_t GetEntryWithIndex(Int_t major, Int_t minor = 0) virtual Long64_t GetEstimate() const Int_t GetEvent(Long64_t entry = 0, Int_t getall = 0) TEventList* GetEventList() const virtual Int_t GetFileNumber() const virtual const char* GetFriendAlias(TTree*) const TH1* GetHistogram() virtual Int_t* GetIndex() virtual Double_t* GetIndexValues() virtual TIterator* GetIteratorOnAllLeaves(Bool_t dir = kIterForward) virtual TLeaf* GetLeaf(const char* name) virtual TSeqCollection* GetListOfAliases() const virtual TObjArray* GetListOfBranches() virtual TList* GetListOfClones() virtual TList* GetListOfFriends() const virtual TObjArray* GetListOfLeaves() Int_t GetMakeClass() const virtual Long64_t GetMaxEntryLoop() const virtual Double_t GetMaximum(const char* columname) static Long64_t GetMaxTreeSize() virtual Long64_t GetMaxVirtualSize() const virtual Double_t GetMinimum(const char* columname) virtual Int_t GetNbranches() TObject* GetNotify() const virtual Int_t GetPacketSize() const TVirtualTreePlayer* GetPlayer() virtual Long64_t GetReadEntry() const virtual Long64_t GetReadEvent() const virtual Int_t GetScanField() const TTreeFormula* GetSelect() virtual Long64_t GetSelectedRows() virtual Int_t GetTimerInterval() const virtual Long64_t GetTotBytes() const virtual TTree* GetTree() const virtual TVirtualIndex* GetTreeIndex() const virtual Int_t GetTreeNumber() const virtual Int_t GetUpdate() const virtual TList* GetUserInfo() virtual Double_t* GetV1() virtual Double_t* GetV2() virtual Double_t* GetV3() virtual Double_t* GetV4() TTreeFormula* GetVar1() TTreeFormula* GetVar2() TTreeFormula* GetVar3() TTreeFormula* GetVar4() virtual Double_t* GetW() virtual Double_t GetWeight() const virtual Long64_t GetZipBytes() const virtual void IncrementTotalBuffers(Int_t nbytes) virtual TClass* IsA() const virtual Bool_t IsFolder() const virtual Int_t LoadBaskets(Long64_t maxmemory = 2000000000) virtual Long64_t LoadTree(Long64_t entry) virtual Long64_t LoadTreeFriend(Long64_t entry, TTree* T) virtual Int_t MakeClass(const char* classname = "0", Option_t* option = "") virtual Int_t MakeCode(const char* filename = "0") virtual Int_t MakeProxy(const char* classname, const char* macrofilename = "0", const char* cutfilename = "0", const char* option = "0", Int_t maxUnrolling = 3) virtual Int_t MakeSelector(const char* selector = "0") Bool_t MemoryFull(Int_t nbytes) virtual Long64_t Merge(TCollection* list) static TTree* MergeTrees(TList* list) virtual Bool_t Notify() TPrincipal* Principal(const char* varexp = "", const char* selection = "", Option_t* option = "np", Long64_t nentries = 1000000000, Long64_t firstentry = 0) virtual void Print(Option_t* option = "") const virtual Long64_t Process(const char* filename, Option_t* option = "", Long64_t nentries = 1000000000, Long64_t firstentry = 0) virtual Long64_t Process(TSelector* selector, Option_t* option = "", Long64_t nentries = 1000000000, Long64_t firstentry = 0) virtual Long64_t Project(const char* hname, const char* varexp, const char* selection = "", Option_t* option = "", Long64_t nentries = 1000000000, Long64_t firstentry = 0) virtual TSQLResult* Query(const char* varexp = "", const char* selection = "", Option_t* option = "", Long64_t nentries = 1000000000, Long64_t firstentry = 0) virtual Long64_t ReadFile(const char* filename, const char* branchDescriptor = "") virtual void Refresh() virtual void RemoveFriend(TTree*) virtual void Reset(Option_t* option = "") virtual void ResetBranchAddresses() virtual Long64_t Scan(const char* varexp = "", const char* selection = "", Option_t* option = "", Long64_t nentries = 1000000000, Long64_t firstentry = 0) virtual Bool_t SetAlias(const char* aliasName, const char* aliasFormula) virtual void SetAutoSave(Long64_t autos = 10000000) virtual void SetBasketSize(const char* bname, Int_t buffsize = 16000) virtual void SetBranchAddress(const char* bname, void* add, TClass* realClass, EDataType datatype, Bool_t ptr) void SetBranchAddress(const char* bname, void** add) virtual void SetBranchStatus(const char* bname, Bool_t status = 1, UInt_t* found = 0) static void SetBranchStyle(Int_t style = 1) virtual void SetChainOffset(Int_t offset = 0) virtual void SetCircular(Long64_t maxEntries) virtual void SetDebug(Int_t level = 1, Long64_t min = 0, Long64_t max = 9999999) virtual void SetDirectory(TDirectory* dir) virtual Long64_t SetEntries(Long64_t n = -1) virtual void SetEstimate(Long64_t nentries = 10000) virtual void SetEventList(TEventList* list) virtual void SetFileNumber(Int_t number = 0) virtual void SetMakeClass(Int_t make) virtual void SetMaxEntryLoop(Long64_t maxev = 1000000000) static void SetMaxTreeSize(Long64_t maxsize = 1900000000) virtual void SetMaxVirtualSize(Long64_t size = 0) virtual void SetName(const char* name) virtual void SetNotify(TObject* obj) virtual void SetObject(const char* name, const char* title) virtual void SetScanField(Int_t n = 50) virtual void SetTimerInterval(Int_t msec = 333) virtual void SetTreeIndex(TVirtualIndex* index) virtual void SetUpdate(Int_t freq = 0) virtual void SetWeight(Double_t w = 1, Option_t* option = "") virtual void Show(Long64_t entry = -1, Int_t lenmax = 20) virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual void StartViewer() virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b) virtual Long64_t UnbinnedFit(const char* funcname, const char* varexp, const char* selection = "", Option_t* option = "", Long64_t nentries = 1000000000, Long64_t firstentry = 0) virtual void UseCurrentStyle()

Data Members


    protected:
Long64_t fEntries Number of entries Long64_t fTotBytes Total number of bytes in all branches before compression Long64_t fZipBytes Total number of bytes in all branches after compression Long64_t fSavedBytes Number of autosaved bytes Double_t fWeight Tree weight (see TTree::SetWeight) Int_t fTimerInterval Timer interval in milliseconds Int_t fScanField Number of runs before prompting in Scan Int_t fUpdate Update frequency for EntryLoop Long64_t fMaxEntries Maximum number of entries in case of circular buffers Long64_t fMaxEntryLoop Maximum number of entries to process Long64_t fMaxVirtualSize Maximum total size of buffers kept in memory Long64_t fAutoSave Autosave tree when fAutoSave bytes produced Long64_t fEstimate Number of entries to estimate histogram limits Long64_t fChainOffset ! Offset of 1st entry of this Tree in a TChain Long64_t fReadEntry ! Number of the entry being processed Long64_t fTotalBuffers ! Total number of bytes in branch buffers Int_t fPacketSize ! Number of entries in one packet for parallel root Int_t fNfill ! Local for EntryLoop Int_t fDebug ! Debug level Long64_t fDebugMin ! First entry number to debug Long64_t fDebugMax ! Last entry number to debug Int_t fMakeClass ! not zero when processing code generated by MakeClass Int_t fFileNumber ! current file number (if file extensions) TObject* fNotify ! Object to be notified when loading a Tree TDirectory* fDirectory ! Pointer to directory holding this tree TObjArray fBranches List of Branches TObjArray fLeaves Direct pointers to individual branch leaves TList* fAliases List of aliases for expressions based on the tree branches. TEventList* fEventList ! Pointer to event selection list (if one) TArrayD fIndexValues Sorted index values TArrayI fIndex Index of sorted values TVirtualIndex* fTreeIndex Pointer to the tree Index (if any) TList* fFriends pointer to list of friend elements TList* fUserInfo pointer to a list of user objects associated to this Tree TVirtualTreePlayer* fPlayer ! Pointer to current Tree player TList* fClones ! List of cloned trees which share our addresses TBranchRef* fBranchRef Branch supporting the TRefTable (if any) UInt_t fFriendLockStatus ! Record which method is locking the friend recursion static Int_t fgBranchStyle Old/New branch style static Long64_t fgMaxTreeSize Maximum size of a file containg a Tree public:
static const TTree::ELockStatusBits kFindBranch static const TTree::ELockStatusBits kFindLeaf static const TTree::ELockStatusBits kGetAlias static const TTree::ELockStatusBits kGetBranch static const TTree::ELockStatusBits kGetEntry static const TTree::ELockStatusBits kGetEntryWithIndex static const TTree::ELockStatusBits kGetFriendAlias static const TTree::ELockStatusBits kGetLeaf static const TTree::ELockStatusBits kLoadTree static const TTree::ELockStatusBits kPrint static const TTree::ELockStatusBits kRemoveFriend static const TTree::ELockStatusBits kSetBranchStatus static const enum TTree:: kForceRead

Class Description

                                                                      
 TTree                                                                
                                                                      
  a TTree object has a header with a name and a title.
  It consists of a list of independent branches (TBranch). Each branch
  has its own definition and list of buffers. Branch buffers may be
  automatically written to disk or kept in memory until the Tree attribute
  fMaxVirtualSize is reached.
  Variables of one branch are written to the same buffer.
  A branch buffer is automatically compressed if the file compression
  attribute is set (default).

  Branches may be written to different files (see TBranch::SetFile).

  The ROOT user can decide to make one single branch and serialize one
  object into one single I/O buffer or to make several branches.
  Making one single branch and one single buffer can be the right choice
  when one wants to process only a subset of all entries in the tree.
  (you know for example the list of entry numbers you want to process).
  Making several branches is particularly interesting in the data analysis
  phase, when one wants to histogram some attributes of an object (entry)
  without reading all the attributes.
/* */

  ==> TTree *tree = new TTree(name, title, maxvirtualsize)
     Creates a Tree with name and title. Maxvirtualsize is by default 64Mbytes,
     maxvirtualsize = 64000000(default) means: Keeps as many buffers in memory until
     the sum of all buffers is greater than 64 Megabyte. When this happens,
     memory buffers are written to disk and deleted until the size of all
     buffers is again below the threshold.
     maxvirtualsize = 0 means: keep only one buffer in memory.

     Various kinds of branches can be added to a tree:
       A - simple structures or list of variables. (may be for C or Fortran structures)
       B - any object (inheriting from TObject). (we expect this option be the most frequent)
       C - a ClonesArray. (a specialized object for collections of same class objects)

  ==> Case A
      ======
     TBranch *branch = tree->Branch(branchname,address, leaflist, bufsize)
       * address is the address of the first item of a structure
       * leaflist is the concatenation of all the variable names and types
         separated by a colon character :
         The variable name and the variable type are separated by a slash (/).
         The variable type may be 0,1 or 2 characters. If no type is given,
         the type of the variable is assumed to be the same as the previous
         variable. If the first variable does not have a type, it is assumed
         of type F by default. The list of currently supported types is given below:
            - C : a character string terminated by the 0 character
            - B : an 8 bit signed integer (Char_t)
            - b : an 8 bit unsigned integer (UChar_t)
            - S : a 16 bit signed integer (Short_t)
            - s : a 16 bit unsigned integer (UShort_t)
            - I : a 32 bit signed integer (Int_t)
            - i : a 32 bit unsigned integer (UInt_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)
            - l : a 64 bit unsigned integer (ULong64_t)
            - O : a boolean (Bool_t)

  ==> Case B
      ======
     TBranch *branch = tree->Branch(branchname,className,object, bufsize, splitlevel)
          object is the address of a pointer to an existing object (derived from TObject).
        if splitlevel=0, the object is serialized in the branch buffer.
        if splitlevel=1 (default), this branch will automatically be split
          into subbranches, with one subbranch for each data member or object
          of the object itself. In case the object member is a TClonesArray,
          the mechanism described in case C is applied to this array.
        if splitlevel=2 ,this branch will automatically be split
          into subbranches, with one subbranch for each data member or object
          of the object itself. In case the object member is a TClonesArray,
          it is processed as a TObject*, only one branch.

  ==> Case C
      ======
     TBranch *branch = tree->Branch(branchname,clonesarray, bufsize, splitlevel)
         clonesarray is the address of a pointer to a TClonesArray.
         The TClonesArray is a direct access list of objects of the same class.
         For example, if the TClonesArray is an array of TTrack objects,
         this function will create one subbranch for each data member of
         the object TTrack.


  ==> branch->SetAddress(Void *address)
      In case of dynamic structures changing with each entry for example, one must
      redefine the branch address before filling the branch again.
      This is done via the TBranch::SetAddress member function.

  ==> tree->Fill()
      loops on all defined branches and for each branch invokes the Fill function.

         See also the class TNtuple (a simple Tree with branches of floats)

       Adding a Branch to an Existing Tree
       ===================================
 You may want to add a branch to an existing tree. For example,
 if one variable in the tree was computed with a certain algorithm,
 you may want to try another algorithm and compare the results.
 One solution is to add a new branch, fill it, and save the tree.
 The code below adds a simple branch to an existing tree.
 Note the kOverwrite option in the Write method, it overwrites the
 existing tree. If it is not specified, two copies of the tree headers
 are saved.

 void tree3AddBranch(){
   TFile f("tree3.root","update");

   Float_t new_v;
   TTree *t3 = (TTree*)f->Get("t3");
   TBranch *newBranch = t3->Branch("new_v",&new_v,"new_v/F");

   //read the number of entries in the t3
   Long64_t nentries = t3->GetEntries();

   for (Long64_t i = 0; i < nentries; i++){
     new_v= gRandom->Gaus(0,1);
     newBranch->Fill();
   }
   // save only the new version of the tree
   t3->Write("",TObject::kOverwrite);
 }
 Adding a branch is often not possible because the tree is in a read-only
 file and you do not have permission to save the modified tree with the
 new branch. Even if you do have the permission, you risk loosing the
 original tree with an unsuccessful attempt to save  the modification.
 Since trees are usually large, adding a branch could extend it over the
 2GB  limit. In this case, the attempt to write the tree fails, and the
 original data is erased.
 In addition, adding a branch to a tree enlarges the tree and increases
 the amount of memory needed to read an entry, and therefore decreases
 the performance.
 For these reasons, ROOT offers the concept of friends for trees (and chains).
 We encourage you to use TTree::AddFriend rather than adding a branch manually.

/* */
  =============================================================================
______________________________________________________________________________
*-*-*-*-*-*-*A simple example with histograms and a tree*-*-*-*-*-*-*-*-*-*
*-*          ===========================================

  This program creates :
    - a one dimensional histogram
    - a two dimensional histogram
    - a profile histogram
    - a tree

  These objects are filled with some random numbers and saved on a file.

-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

 #include "TFile.h"
 #include "TH1.h"
 #include "TH2.h"
 #include "TProfile.h"
 #include "TRandom.h"
 #include "TTree.h"


 //______________________________________________________________________________
 main(int argc, char **argv)
 {
 // Create a new ROOT binary machine independent file.
 // Note that this file may contain any kind of ROOT objects, histograms,trees
 // pictures, graphics objects, detector geometries, tracks, events, etc..
 // This file is now becoming the current directory.
   TFile hfile("htree.root","RECREATE","Demo ROOT file with histograms & trees");

 // Create some histograms and a profile histogram
   TH1F *hpx   = new TH1F("hpx","This is the px distribution",100,-4,4);
   TH2F *hpxpy = new TH2F("hpxpy","py ps px",40,-4,4,40,-4,4);
   TProfile *hprof = new TProfile("hprof","Profile of pz versus px",100,-4,4,0,20);

 // Define some simple structures
   typedef struct {Float_t x,y,z;} POINT;
   typedef struct {
      Int_t ntrack,nseg,nvertex;
      UInt_t flag;
      Float_t temperature;
   } EVENTN;
   static POINT point;
   static EVENTN eventn;

 // Create a ROOT Tree
   TTree *tree = new TTree("T","An example of ROOT tree with a few branches");
   tree->Branch("point",&point,"x:y:z");
   tree->Branch("eventn",&eventn,"ntrack/I:nseg:nvertex:flag/i:temperature/F");
   tree->Branch("hpx","TH1F",&hpx,128000,0);

   Float_t px,py,pz;
   static Float_t p[3];

 //--------------------Here we start a loop on 1000 events
   for ( Int_t i=0; i<1000; i++) {
      gRandom->Rannor(px,py);
      pz = px*px + py*py;
      Float_t random = gRandom->::Rndm(1);

 //         Fill histograms
      hpx->Fill(px);
      hpxpy->Fill(px,py,1);
      hprof->Fill(px,pz,1);

 //         Fill structures
      p[0] = px;
      p[1] = py;
      p[2] = pz;
      point.x = 10*(random-1);;
      point.y = 5*random;
      point.z = 20*random;
      eventn.ntrack  = Int_t(100*random);
      eventn.nseg    = Int_t(2*eventn.ntrack);
      eventn.nvertex = 1;
      eventn.flag    = Int_t(random+0.5);
      eventn.temperature = 20+random;

 //        Fill the tree. For each event, save the 2 structures and 3 objects
 //      In this simple example, the objects hpx, hprof and hpxpy are slightly
 //      different from event to event. We expect a big compression factor!
      tree->Fill();
   }
  //--------------End of the loop

   tree->Print();

 // Save all objects in this file
   hfile.Write();

 // Close the file. Note that this is automatically done when you leave
 // the application.
   hfile.Close();

   return 0;
 }
                                                                      


TTree(): TNamed(),fFriendLockStatus(0)
*-*-*-*-*-*-*-*-*-*-*Default Tree constructor*-*-*-*-*-*-*-*-*-*-*-*-*-*
*-*                  ========================

TTree(const char *name,const char *title, Int_t splitlevel) :TNamed(name,title),fFriendLockStatus(0)
*-*-*-*-*-*-*-*-*-*Normal Tree constructor*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
*-*                ======================

   The Tree is created in the current directory
   Use the various functions Branch below to add branches to this Tree.

 If the first character of title is a "/", the function assumes a folder name.
 In this case, it creates automatically branches following the folder hierarchy.
 splitlevel may be used in this case to control the split level.

~TTree()
*-*-*-*-*-*-*-*-*-*-*Tree destructor*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
*-*                  =================

void AddClone(TTree *clone)
 Add a cloned tree to our list of tree to be notify whenever we changes our
 addresses and are being deleted.

TFriendElement* AddFriend(const char *treename, const char *filename)
 Add a TFriendElement to the list of friends.
 This function:
   -opens a file if filename is specified
   -reads a Tree with name treename from the file (current directory)
   -adds the Tree to the list of friends
 see other AddFriend functions

 A TFriendElement TF describes a TTree object TF in a file.
 When a TFriendElement TF is added to the the list of friends of an
 existing TTree T, any variable from TF can be referenced in a query
 to T.

   A tree keeps a list of friends. In the context of a tree (or a chain),
 friendship means unrestricted access to the friends data. In this way
 it is much like adding another branch to the tree without taking the risk
 of damaging it. To add a friend to the list, you can use the TTree::AddFriend
 method.  The tree in the diagram below has two friends (friend_tree1 and
 friend_tree2) and now has access to the variables a,b,c,i,j,k,l and m.

/* */

 The AddFriend method has two parameters, the first is the tree name and the
 second is the name of the ROOT file where the friend tree is saved.
 AddFriend automatically opens the friend file. If no file name is given,
 the tree called ft1 is assumed to be in the same file as the original tree.

 tree.AddFriend("ft1","friendfile1.root");
 If the friend tree has the same name as the original tree, you can give it
 an alias sin the context of the friendship:

 tree.AddFriend("tree1 = tree","friendfile1.root");
 Once the tree has friends, we can use TTree::Draw as if the friend's
 variables were in the original tree. To specify which tree to use in
 the Draw method, use the syntax:

 <treeName>.<branchname>.<varname>
 If the variablename is enough to uniquely identify the variable, you can
 leave out the tree and/or branch name.
 For example, these commands generate a 3-d scatter plot of variable "var"
 in the TTree tree versus variable v1 in TTree ft1 versus variable v2 in
 TTree ft2.

 tree.AddFriend("ft1","friendfile1.root");
 tree.AddFriend("ft2","friendfile2.root");
 tree.Draw("var:ft1.v1:ft2.v2");

/* */

 The picture illustrates the access of the tree and its friends with a
 Draw command.
 When AddFriend is called, the ROOT file is automatically opened and the
 friend tree (ft1) is read into memory. The new friend (ft1) is added to
 the list of friends of tree.
 The number of entries in the friend must be equal or greater to the number
 of entries of the original tree. If the friend tree has fewer entries a
 warning is given and the missing entries are not included in the histogram.
 To retrieve the list of friends from a tree use TTree::GetListOfFriends.
 When the tree is written to file (TTree::Write), the friends list is saved
 with it. And when the tree is retrieved, the trees on the friends list are
 also retrieved and the friendship restored.
 When a tree is deleted, the elements of the friend list are also deleted.
 It is possible to declare a friend tree that has the same internal
 structure (same branches and leaves) as the original tree, and compare the
 same values by specifying the tree.

  tree.Draw("var:ft1.var:ft2.var")

TFriendElement* AddFriend(const char *treename, TFile *file)
 Add a TFriendElement to the list of friends. The TFile is managed by
 the user (e.g. the user must delete the file).
 For complete description see AddFriend(const char *, const char *).
 This function:
   -reads a Tree with name treename from the file
   -adds the Tree to the list of friends

TFriendElement* AddFriend(TTree *tree, const char* alias, Bool_t warn)
 Add a TFriendElement to the list of friends. The TTree is managed by
 the user (e.g. the user must delete the file).
 For complete description see AddFriend(const char *, const char *).

Long64_t AutoSave(Option_t *option)
*-*-*-*-*-*-*-*-*-*-*AutoSave tree header every fAutoSave bytes*-*-*-*-*-*
*-*                  ==========================================

   When large Trees are produced, it is safe to activate the AutoSave
   procedure. Some branches may have buffers holding many entries.
   AutoSave is automatically called by TTree::Fill when the number of bytes
   generated since the previous AutoSave is greater than fAutoSave bytes.
   This function may also be invoked by the user, for example every
   N entries.
   Each AutoSave generates a new key on the file.
   Once the key with the tree header has been written, the previous cycle
   (if any) is deleted.

   Note that calling TTree::AutoSave too frequently (or similarly calling
   TTree::SetAutoSave with a small value) is an expensive operation.
   You should make tests for your own application to find a compromize
   between speed and the quantity of information you may loose in case of
   a job crash.

   In case your program crashes before closing the file holding this tree,
   the file will be automatically recovered when you will connect the file
   in UPDATE mode.
   The Tree will be recovered at the status corresponding to the last AutoSave.

   if option contains "SaveSelf", gDirectory->SaveSelf() is called.
   This allows another process to analyze the Tree while the Tree is being filled.

   By default the previous header is deleted after having written the new header.
   if option contains "Overwrite", the previous Tree header is deleted
   before written the new header. This option is slightly faster, but
   the default option is safer in case of a problem (disk quota exceeded)
   when writing the new header.

   The function returns the number of bytes written to the file.
   if the number of bytes is null, an error has occured while writing
   the header to the file.

   How to write a Tree in one process and view it from another process
   ===================================================================
   The following two scripts illustrate how to do this.
   The script treew.C is executed by process1, treer.C by process2

   ----- script treew.C
   void treew() {
     TFile f("test.root","recreate");
     TNtuple *ntuple = new TNtuple("ntuple","Demo","px:py:pz:random:i");
     Float_t px, py, pz;
     for ( Int_t i=0; i<10000000; i++) {
        gRandom->Rannor(px,py);
        pz = px*px + py*py;
        Float_t random = gRandom->Rndm(1);
        ntuple->Fill(px,py,pz,random,i);
        if (i%1000 == 1) ntuple->AutoSave("SaveSelf");
     }
   }

   ----- script treer.C
   void treer() {
      TFile f("test.root");
      TTree *ntuple = (TTree*)f.Get("ntuple");
      TCanvas c1;
      Int_t first = 0;
      while(1) {
         if (first == 0) ntuple->Draw("px>>hpx", "","",10000000,first);
         else            ntuple->Draw("px>>+hpx","","",10000000,first);
         first = (Int_t)ntuple->GetEntries();
         c1.Update();
         gSystem->Sleep(1000); //sleep 1 second
         ntuple->Refresh();
      }
   }

TBranch* BranchImp(const char *branchname, const char *classname, TClass *ptrClass, void *addobj, Int_t bufsize, Int_t splitlevel)
 Same as TTree::Branch with added check that the address passed in addobj
 corresponding to className.  See TTree::Branch for other details.

TBranch* BranchImp(const char *branchname, TClass *ptrClass, void *addobj, Int_t bufsize, Int_t splitlevel)
 Same as TTree::Branch but automatic detection of the class name
 See TTree::Branch for other details.

Int_t Branch(TList *li, Int_t bufsize, Int_t splitlevel)
   Deprecated function. Use next function instead.

Int_t Branch(TCollection *li, Int_t bufsize, Int_t splitlevel, const char *name)
   This function creates one branch for each element in the collection.
   Each entry in the collection becomes a top level branch if the
   corresponding class is not a collection. If it is a collection, the entry
   in the collection becomes in turn top level branches, etc.
   The splitlevel is decreased by 1 everytime a new collection is found.
   For example if list is a TObjArray*
     - if splitlevel = 1, one top level branch is created for each element
        of the TObjArray.
     - if splitlevel = 2, one top level branch is created for each array element.
       if, in turn, one of the array elements is a TCollection, one top level
       branch will be created for each element of this collection.

   In case a collection element is a TClonesArray, the special Tree constructor
   for TClonesArray is called.
   The collection itself cannot be a TClonesArray.

   The function returns the total number of branches created.

   If name is given, all branch names will be prefixed with name_.

 IMPORTANT NOTE1: This function should not be called with splitlevel < 1.

 IMPORTANT NOTE2: The branches created by this function will have names
 corresponding to the collection or object names. It is important
 to give names to collections to avoid misleading branch names or
 identical branch names. By default collections have a name equal to
 the corresponding class name, eg the default name for a TList is "TList".

 Example--------------------------------------------------------------:

Int_t Branch(const char *foldername, Int_t bufsize, Int_t splitlevel)
   This function creates one branch for each element in the folder.
   The function returns the total number of branches created.

TBranch* Branch(const char *name, void *address, const char *leaflist,Int_t bufsize)
*-*-*-*-*-*-*-*-*-*-*Create a new TTree Branch*-*-*-*-*-*-*-*-*-*-*-*-*
*-*                  =========================

     This Branch constructor is provided to support non-objects in
     a Tree. The variables described in leaflist may be simple variables
     or structures.
    See the two following constructors for writing objects in a Tree.

    By default the branch buffers are stored in the same file as the Tree.
    use TBranch::SetFile to specify a different file

       * address is the address of the first item of a structure
         or the address of a pointer to an object (see example).
       * leaflist is the concatenation of all the variable names and types
         separated by a colon character :
         The variable name and the variable type are separated by a slash (/).
         The variable type may be 0,1 or 2 characters. If no type is given,
         the type of the variable is assumed to be the same as the previous
         variable. If the first variable does not have a type, it is assumed
         of type F by default. The list of currently supported types is given below:
            - C : a character string terminated by the 0 character
            - B : an 8 bit signed integer (Char_t)
            - b : an 8 bit unsigned integer (UChar_t)
            - S : a 16 bit signed integer (Short_t)
            - s : a 16 bit unsigned integer (UShort_t)
            - I : a 32 bit signed integer (Int_t)
            - i : a 32 bit unsigned integer (UInt_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)
            - l : a 64 bit unsigned integer (ULong64_t)
            - O : a boolean (Bool_t)

         By default, a variable will be copied to the buffer with the number of
         bytes specified in the type descriptor character. However, if the type
         consists of 2 characters, the second character is an integer that
         specifies the number of bytes to be used when copying the variable
         to the output buffer. Example:
             X         ; variable X, type Float_t
             Y/I       : variable Y, type Int_t
             Y/I2      ; variable Y, type Int_t converted to a 16 bits integer


TBranch* Branch(const char *name, const char *classname, void *addobj, Int_t bufsize, Int_t splitlevel)
 create a new branch with the object of class classname at address addobj.

 WARNING:
 Starting with Root version 3.01, the Branch function uses the new style
 branches (TBranchElement). To get the old behaviour, you can:
   - call BranchOld or
   - call TTree::SetBranchStyle(0)

 Note that with the new style, classname does not need to derive from TObject.
 It must derived from TObject if the branch style has been set to 0 (old)

 Use splitlevel < 0 instead of splitlevel=0 when the class
 has a custom Streamer

 Note: if the split level is set to the default (99),  TTree::Branch will
 not issue a warning if the class can not be split.

TBranch* BranchOld(const char *name, const char *classname, void *addobj, Int_t bufsize, Int_t splitlevel)
*-*-*-*-*-*-*-*-*-*-*Create a new TTree BranchObject*-*-*-*-*-*-*-*-*-*-*-*
*-*                  ===============================

    Build a TBranchObject for an object of class classname.
    addobj is the address of a pointer to an object of class classname.
    IMPORTANT: classname must derive from TObject.
    The class dictionary must be available (ClassDef in class header).

    This option requires access to the library where the corresponding class
    is defined. Accessing one single data member in the object implies
    reading the full object.
    See the next Branch constructor for a more efficient storage
    in case the entry consists of arrays of identical objects.

    By default the branch buffers are stored in the same file as the Tree.
    use TBranch::SetFile to specify a different file

      IMPORTANT NOTE about branch names
    In case two or more master branches contain subbranches with
    identical names, one must add a "." (dot) character at the end
    of the master branch name. This will force the name of the subbranch
    to be master.subbranch instead of simply subbranch.
    This situation happens when the top level object (say event)
    has two or more members referencing the same class.
    For example, if a Tree has two branches B1 and B2 corresponding
    to objects of the same class MyClass, one can do:
       tree.Branch("B1.","MyClass",&b1,8000,1);
       tree.Branch("B2.","MyClass",&b2,8000,1);
    if MyClass has 3 members a,b,c, the two instructions above will generate
    subbranches called B1.a, B1.b ,B1.c, B2.a, B2.b, B2.c

TBranch* BranchRef()
 Build the optional branch supporting the TRefTable
 This branch will keep all the information to find the branches
 containing referenced objects

 at each Tree::Fill, the branch numbers containing the
 referenced objects are saved to the TBranchRef basket
 When the Tree header is saved (via TTree::Write), the branch
 is saved keeping the information with the pointers to the branches
 having referenced objects.

TBranch* Bronch(const char *name, const char *classname, void *add, Int_t bufsize, Int_t splitlevel)
*-*-*-*-*-*-*-*-*-*-*Create a new TTree BranchElement*-*-*-*-*-*-*-*-*-*-*-*
*-*                  ================================

    WARNING about this new function
    ===============================
    This function is designed to replace the function TTree::Branch above.
    This function is far more powerful than the Branch function.
    It supports the full C++, including STL and has the same behaviour
    in split or non-split mode. classname does not have to derive from TObject.
    The function is based on the new TStreamerInfo.

    Build a TBranchElement for an object of class classname.
    addobj is the address of a pointer to an object of class classname.
    The class dictionary must be available (ClassDef in class header).

    This option requires access to the library where the corresponding class
    is defined. Accessing one single data member in the object implies
    reading the full object.

    By default the branch buffers are stored in the same file as the Tree.
    use TBranch::SetFile to specify a different file

      IMPORTANT NOTE about branch names
    In case two or more master branches contain subbranches with
    identical names, one must add a "." (dot) character at the end
    of the master branch name. This will force the name of the subbranch
    to be master.subbranch instead of simply subbranch.
    This situation happens when the top level object (say event)
    has two or more members referencing the same class.
    For example, if a Tree has two branches B1 and B2 corresponding
    to objects of the same class MyClass, one can do:
       tree.Branch("B1.","MyClass",&b1,8000,1);
       tree.Branch("B2.","MyClass",&b2,8000,1);
    if MyClass has 3 members a,b,c, the two instructions above will generate
    subbranches called B1.a, B1.b ,B1.c, B2.a, B2.b, B2.c

   Use splitlevel < 0 instead of splitlevel=0 when the class
   has a custom Streamer

   Note: if the split level is set to the default (99),  TTree::Branch will
   not issue a warning if the class can not be split.

void Browse(TBrowser *b)

Int_t BuildIndex(const char *majorname, const char *minorname)
 Build a Tree Index (default is TtreeIndex).
 see a description of the parameters and functionality in
  TTreeIndex::TTreeIndex

 The return value is the number of entries in the Index (< 0 indicates failure)

 A TTreeIndex object pointed by fTreeIndex is created.
 This object will be automatically deleted by the TTree destructor
 See also comments in TTree::SetTreeIndex

void SetTreeIndex(TVirtualIndex*index)
 The current TreeIndex is replaced by the new index.
 Note that this function does not delete the previous index.
 This gives the possibility to play with more than one index, eg
 TVirtualIndex *oldIndex = tree.GetTreeIndex();
 tree.SetTreeIndex(newIndex);
 tree.Draw(...);
 tree.SetTreeIndex(oldIndex);
 tree.Draw(); etc

TStreamerInfo* BuildStreamerInfo(TClass *cl, void *pointer)
 Build StreamerInfo for class cl
 pointer is an optional argument that may contain a pointer to an object of cl

TFile* ChangeFile(TFile *file)
 called by TTree::Fill when file has reached its maximum fgMaxTreeSize.
 Create a new file. If the original file is named "myfile.root",
 subsequent files are named "myfile_1.root", "myfile_2.root", etc.

 Return pointer to new file
 Currently, the automatic change of file is restricted
 to the case where the Tree is in the top level directory.
 The file should not contain sub-directories.

 Before switching to a new file, the Tree header is written
 to the current file, then the current file is closed.

 To process the multiple files created by ChangeFile, one must use
 a TChain.

 The new file name has a suffix "_N" where N is equal to fFileNumber+1.
 By default a Root session starts with fFileNumber=0. One can set
 fFileNumber to a different value via TTree::SetFileNumber.
 In case a file named "_N" already exists, the function will try
 a file named "__N", then "___N", etc.

 fgMaxTreeSize can be set via the static function TTree::SetMaxTreeSize.
 The default value of fgMaxTreeSize is 1.9 Gigabytes.

 If the current file contains other objects like TH1 and TTree,
 these objects are automatically moved to the new file.

 IMPORTANT NOTE:
 Be careful when writing the final Tree header to the file!
 Don't do:
  TFile *file = new TFile("myfile.root","recreate");
  TTree *T = new TTree("T","title");
  T->Fill(); //loop
  file->Write();
  file->Close();
 but do the following:
  TFile *file = new TFile("myfile.root","recreate");
  TTree *T = new TTree("T","title");
  T->Fill(); //loop
  file = T->GetCurrentFile(); //to get the pointer to the current file
  file->Write();
  file->Close();

Bool_t CheckBranchAddressType(TBranch *branch, TClass *ptrClass, EDataType datatype, Bool_t ptr)
 Check whether the address described by the last 3 parameters match the
 content of the branch.

TTree* CloneTree(Long64_t nentries, Option_t *)
 Create a clone of this tree and copy nentries
 By default copy all entries
 option is reserved for future use
 Note that only active branches are copied.

 IMPORTANT: The cloned tree stays connected with this tree until this tree
            is deleted.  In particular, any changes in branch addresses
            in this tree are forwarded to the clone trees.  Any changes
            made to the branch addresses of the cloned trees are over-ridden
            anytime this tree changes its branch addresses.
            Once this tree is deleted, all the addresses of the cloned tree
            are reset to their default values.

 For examples of CloneTree, see tutorials
  -copytree:
    Example of Root macro to copy a subset of a Tree to a new Tree
    The input file has been generated by the program in $ROOTSYS/test/Event
    with   Event 1000 1 1 1
  -copytree2:
    Example of Root macro to copy a subset of a Tree to a new Tree
    One branch of the new Tree is written to a separate file
    The input file has been generated by the program in $ROOTSYS/test/Event
    with   Event 1000 1 1 1

void CopyAddresses(TTree *tree)
 Set branch addresses of tree equal to the ones of this tree

Long64_t CopyEntries(TTree *tree, Long64_t nentries)
 Copy nentries from tree to this tree
 By default copy all entries
 Return number of bytes copied to this tree.

TTree* CopyTree(const char *selection, Option_t *option, Long64_t nentries, Long64_t firstentry)
*-*-*-*-*-*-*-*-*copy a Tree with selection*-*-*-*-*-*
*-*              ==========================

 IMPORTANT: The copied tree stays connected with this tree until this tree
            is deleted.  In particular, any changes in branch addresses
            in this tree are forwarded to the clone trees.  Any changes
            made to the branch addresses of the copied trees are over-ridden
            anytime this tree changes its branch addresses.
            Once this tree is deleted, all the addresses of the copied tree
            are reset to their default values.

 For examples of CloneTree, see tutorials
  -copytree:
    Example of Root macro to copy a subset of a Tree to a new Tree
    The input file has been generated by the program in $ROOTSYS/test/Event
    with   Event 1000 1 1 1
  -copytree2:
    Example of Root macro to copy a subset of a Tree to a new Tree
    One branch of the new Tree is written to a separate file
    The input file has been generated by the program in $ROOTSYS/test/Event
    with   Event 1000 1 1 1
  -copytree3:
    Example of Root macro to copy a subset of a Tree to a new Tree
    Only selected entries are copied to the new Tree

 NOTE that only the active branches are copied.

TBasket* CreateBasket(TBranch *branch)
 Create a basket for this implementation of TTree.

void Delete(Option_t *option)
*-*-*-*-*-*-*-*-*Delete this tree from memory or/and disk
*-*              ========================================

  if option == "all" delete Tree object from memory AND from disk
                     all baskets on disk are deleted. All keys with same name
                     are deleted.
  if option =="" only Tree object in memory is deleted.

Long64_t Draw(const char *varexp, const TCut &selection, Option_t *option, Long64_t nentries, Long64_t firstentry)
*-*-*-*-*-*-*-*-*-*-*Draw expression varexp for specified entries-*-*-*-*-*
*-*                  ===========================================

      This function accepts TCut objects as arguments.
      Useful to use the string operator +
         example:
            ntuple.Draw("x",cut1+cut2+cut3);


Long64_t Draw(const char *varexp, const char *selection, Option_t *option,Long64_t nentries, Long64_t firstentry)
*-*-*-*-*-*-*-*-*-*-*Draw expression varexp for specified entries-*-*-*-*-*
*-*                  ===========================================

  varexp is an expression of the general form
   - "e1"           produces a 1-d histogram (TH1F) of expression "e1"
   - "e1:e2"        produces an unbinned 2-d scatter-plot (TGraph) of "e1" versus "e2"
   - "e1:e2:e3"     produces an unbinned 3-d scatter-plot (TPolyMarker3D) of "e1"
                    versus "e2" versus "e3"
   - "e1:e2:e3:e4"  produces an unbinned 3-d scatter-plot (TPolyMarker3D) of "e1"
                    versus "e2" versus "e3" and "e4" mapped on the color number.
  (to create histograms in the 2, 3, and 4 dimesional case, see section "Saving
  the result of Draw to an histogram")

  Example:
     varexp = x     simplest case: 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)
            = px:py:pz:2.5*E  produces a 3-d scatter-plot of px vs py ps pz
              and the color number of each marker will be 2.5*E.
              If the color number is negative it is set to 0.
              If the color number is greater than the current number of colors
                 it is set to the highest color number.
              The default number of colors is 50.
              see TStyle::SetPalette for setting a new color palette.

  Note that the variables e1, e2 or e3 may contain a selection.
  example, if e1= x*(y<0), the value histogrammed will be x if y<0
  and will be 0 otherwise.

  selection is an expression with a combination of the columns.
  In a selection all the C++ operators are authorized.
  The value corresponding to the selection expression is used as a weight
  to fill the histogram.
  If the expression includes only boolean operations, the result
  is 0 or 1. If the result is 0, the histogram is not filled.
  In general, the expression may be of the form:
      value*(boolean expression)
  if boolean expression is true, the histogram is filled with
  a weight = value.
  Examples:
      selection1 = "x<y && sqrt(z)>3.2"
      selection2 = "(x+y)*(sqrt(z)>3.2)"
  selection1 returns a weigth = 0 or 1
  selection2 returns a weight = x+y if sqrt(z)>3.2
             returns a weight = 0 otherwise.

  option is the drawing option.
    - See TH1::Draw for the list of all drawing options.
    - If option COL is specified when varexp has three fields:
            tree.Draw("e1:e2:e3","","col");
      a 2D scatter is produced with e1 vs e2, and e3 is mapped on the color
      table.
    - If option contains the string "goff", no graphics is generated.

  nentries is the number of entries to process (default is all)
  first is the first entry to process (default is 0)

  This function returns the number of selected entries. It returns -1
  if an error occurs.

     Drawing expressions using arrays and array elements
     ===================================================
 Let assumes, a leaf fMatrix, on the branch fEvent, which is a 3 by 3 array,
 or a TClonesArray.
 In a TTree::Draw expression you can now access fMatrix using the following
 syntaxes:

   String passed    What is used for each entry of the tree

   "fMatrix"       the 9 elements of fMatrix
   "fMatrix[][]"   the 9 elements of fMatrix
   "fMatrix[2][2]" only the elements fMatrix[2][2]
   "fMatrix[1]"    the 3 elements fMatrix[1][0], fMatrix[1][1] and fMatrix[1][2]
   "fMatrix[1][]"  the 3 elements fMatrix[1][0], fMatrix[1][1] and fMatrix[1][2]
   "fMatrix[][0]"  the 3 elements fMatrix[0][0], fMatrix[1][0] and fMatrix[2][0]

   "fEvent.fMatrix...." same as "fMatrix..." (unless there is more than one leaf named fMatrix!).

 In summary, if a specific index is not specified for a dimension, TTree::Draw
 will loop through all the indices along this dimension.  Leaving off the
 last (right most) dimension of specifying then with the two characters '[]'
 is equivalent.  For variable size arrays (and TClonesArray) the range
 of the first dimension is recalculated for each entry of the tree.

 TTree::Draw also now properly handling operations involving 2 or more arrays.

 Let assume a second matrix fResults[5][2], here are a sample of some
 of the possible combinations, the number of elements they produce and
 the loop used:

  expression                       element(s)  Loop

  "fMatrix[2][1] - fResults[5][2]"   one     no loop
  "fMatrix[2][]  - fResults[5][2]"   three   on 2nd dim fMatrix
  "fMatrix[2][]  - fResults[5][]"    two     on both 2nd dimensions
  "fMatrix[][2]  - fResults[][1]"    three   on both 1st dimensions
  "fMatrix[][2]  - fResults[][]"     six     on both 1st and 2nd dimensions of
                                             fResults
  "fMatrix[][2]  - fResults[3][]"    two     on 1st dim of fMatrix and 2nd of
                                             fResults (at the same time)
  "fMatrix[][]   - fResults[][]"     six     on 1st dim then on  2nd dim


 In summary, TTree::Draw loops through all un-specified dimensions.  To
 figure out the range of each loop, we match each unspecified dimension
 from left to right (ignoring ALL dimensions for which an index has been
 specified), in the equivalent loop matched dimensions use the same index
 and are restricted to the smallest range (of only the matched dimensions).
 When involving variable arrays, the range can of course be different
 for each entry of the tree.

 So the loop equivalent to "fMatrix[][2] - fResults[3][]" is:

    for (Int_t i0; i < min(3,2); i++) {
       use the value of (fMatrix[i0][2] - fMatrix[3][i0])
    }

 So the loop equivalent to "fMatrix[][2] - fResults[][]" is:

    for (Int_t i0; i < min(3,5); i++) {
       for (Int_t i1; i1 < 2; i1++) {
          use the value of (fMatrix[i0][2] - fMatrix[i0][i1])
       }
    }

 So the loop equivalent to "fMatrix[][] - fResults[][]" is:

    for (Int_t i0; i < min(3,5); i++) {
       for (Int_t i1; i1 < min(3,2); i1++) {
          use the value of (fMatrix[i0][i1] - fMatrix[i0][i1])
       }
    }

     Retrieving the result of Draw
     =============================

  By default the temporary histogram created is called "htemp", but only in
  the one dimensional Draw("e1") it contains the TTree's data points. For
  a two dimensional Draw, the data is filled into a TGraph which is named
  "Graph". They can be retrieved by calling
    TH1F *htemp = (TH1F*)gPad->GetPrimitive("htemp"); // 1D
    TGraph *graph = (TGraph*)gPad->GetPrimitive("graph"); // 2D

  For a three and four dimensional Draw the TPloyMarker3D is unnamed, and
  cannot be retrieved.

  gPad always contains a TH1 derived object called "htemp" which allows to
  access the axes:
    TGraph *graph = (TGraph*)gPad->GetPrimitive("graph"); // 2D
    TH2F   *htemp = (TH2F*)gPad->GetPrimitive("htemp"); // empty, but has axes
    TAxis  *xaxis = htemp->GetXaxis();

     Saving the result of Draw to an histogram
     =========================================

  If varexp0 contains >>hnew (following the variable(s) name(s),
  the new histogram created is called hnew and it is kept in the current
  directory (and also the current pad). This works for all dimensions.
  Example:
    tree.Draw("sqrt(x)>>hsqrt","y>0")
    will draw sqrt(x) and save the histogram as "hsqrt" in the current
    directory.  To retrieve it do:
    TH1F *hsqrt = (TH1F*)gDirectory->Get("hsqrt");

  The binning information is taken from the environment variables

     Hist.Binning.?D.?

  In addition, the name of the histogram can be followed by up to 9
  numbers between '(' and ')', where the numbers describe the
  following:

   1 - bins in x-direction
   2 - lower limit in x-direction
   3 - upper limit in x-direction
   4-6 same for y-direction
   7-9 same for z-direction

   When a new binning is used the new value will become the default.
   Values can be skipped.
  Example:
    tree.Draw("sqrt(x)>>hsqrt(500,10,20)")
          // plot sqrt(x) between 10 and 20 using 500 bins
    tree.Draw("sqrt(x):sin(y)>>hsqrt(100,10,60,50,.1,.5)")
          // plot sqrt(x) against sin(y)
          // 100 bins in x-direction; lower limit on x-axis is 10; upper limit is 60
          //  50 bins in y-direction; lower limit on y-axis is .1; upper limit is .5

  By default, the specified histogram is reset.
  To continue to append data to an existing histogram, use "+" in front
  of the histogram name.
  A '+' in front of the histogram name is ignored, when the name is followed by
  binning information as described in the previous paragraph.
    tree.Draw("sqrt(x)>>+hsqrt","y>0")
      will not reset hsqrt, but will continue filling.
  This works for 1-D, 2-D and 3-D histograms.

     Accessing collection objects
     ============================

  TTree::Draw default's handling of collections is to assume that any
  request on a collection pertain to it content.  For example, if fTracks
  is a collection of Track objects, the following:
      tree->Draw("event.fTracks.fPx");
  will plot the value of fPx for each Track objects inside the collection.
  Also
     tree->Draw("event.fTracks.size()");
  would plot the result of the member function Track::size() for each
  Track object inside the collection.
  To access information about the collection itself, TTree::Draw support
  the '@' notation.  If a variable which points to a collection is prefixed
  or postfixed with '@', the next part of the expression will pertain to
  the collection object.  For example:
     tree->Draw("event.@fTracks.size()");
  will plot the size of the collection refered to by fTracks (i.e the number
  of Track objects).

     Special functions and variables
     ===============================

  Entry$:  A TTree::Draw formula can use the special variable Entry$
  to access the entry number being read.  For example to draw every
  other entry use:
    tree.Draw("myvar","Entry$%2==0");

  Entry$    : return the current entry number (== TTree::GetReadEntry())
  Entries$  : return the total number of entries (== TTree::GetEntries())
  Length$   : return the total number of element of this formula for this
                 entry (==TTreeFormula::GetNdata())
  Iteration$: return the current iteration over this formula for this
                 entry (i.e. varies from 0 to Length$).

  Length$(formula): return the total number of element of the formula given as a
                    parameter.
  Sum$(formula): return the sum of the value of the elements of the formula given
                    as a parameter.  For eaxmple the mean for all the elements in
                    one entry can be calculated with:
                Sum$(formula)/Length$(formula)

  Alt$(primary,alternate) : return the value of "primary" if it is available
                 for the current iteration otherwise return the value of "alternate".
                 For example, with arr1[3] and arr2[2]
    tree->Draw("arr1+Alt$(arr2,0)");
                 will draw arr1[0]+arr2[0] ; arr1[1]+arr2[1] and arr1[2]+0
                 Or with a variable size array arr3
    tree->Draw("Alt$(arr3[0],0)+Alt$(arr3[1],0)+Alt$(arr3[2],0)");
                 will draw the sum arr3 for the index 0 to min(2,actual_size_of_arr3-1)
                 As a comparison
    tree->Draw("arr3[0]+arr3[1]+arr3[2]");
                 will draw the sum arr3 for the index 0 to 2 only if the
                 actual_size_of_arr3 is greater or equal to 3.
                 Note that the array in 'primary' is flatened/linearilized thus using
                 Alt$ with multi-dimensional arrays of different dimensions in unlikely
                 to yield the expected results.  To visualize a bit more what elements
                 would be matched by TTree::Draw, TTree::Scan can be used:
    tree->Scan("arr1:Alt$(arr2,0)");
                 will print on one line the value of arr1 and (arr2,0) that will be
                 matched by
    tree->Draw("arr1-Alt$(arr2,0)");

     Drawing a user function accessing the TTree data directly
     =========================================================

  If the formula contains  a file name, TTree::MakeProxy will be used
  to load and execute this file.   In particular it will draw the
  result of a function with the same name as the file.  The function
  will be executed in a context where the name of the branches can
  be used as a C++ variable.

  For example draw px using the file hsimple.root (generated by the
  hsimple.C tutorial), we need a file named hsimple.cxx:

     double hsimple() {
        return px;
     }

  MakeProxy can then be used indirectly via the TTree::Draw interface
  as follow:
     new TFile("hsimple.root")
     ntuple->Draw("hsimple.cxx");

  A more complete example is available in the tutorials directory:
    h1analysisProxy.cxx , h1analysProxy.h and h1analysisProxyCut.C
  which reimplement the selector found in h1analysis.C

  The main features of this facility are:

    * on-demand loading of branches
    * ability to use the 'branchname' as if it was a data member
    * protection against array out-of-bound
    * ability to use the branch data as object (when the user code is available)

  See TTree::MakeProxy for more details.

     Making a Profile histogram
     ==========================
  In case of a 2-Dim expression, one can generate a TProfile histogram
  instead of a TH2F histogram by specyfying option=prof or option=profs.
  The option=prof is automatically selected in case of y:x>>pf
  where pf is an existing TProfile histogram.

     Saving the result of Draw to a TEventList
     =========================================
  TTree::Draw can be used to fill a TEventList object (list of entry numbers)
  instead of histogramming one variable.
  If varexp0 has the form >>elist , a TEventList object named "elist"
  is created in the current directory. elist will contain the list
  of entry numbers satisfying the current selection.
  Example:
    tree.Draw(">>yplus","y>0")
    will create a TEventList object named "yplus" in the current directory.
    In an interactive session, one can type (after TTree::Draw)
       yplus.Print("all")
    to print the list of entry numbers in the list.

  By default, the specified entry list is reset.
  To continue to append data to an existing list, use "+" in front
  of the list name;
    tree.Draw(">>+yplus","y>0")
      will not reset yplus, but will enter the selected entries at the end
      of the existing list.

      Using a TEventList as Input
      ===========================
  Once a TEventList object has been generated, it can be used as input
  for TTree::Draw. Use TTree::SetEventList to set the current event list
  Example:
     TEventList *elist = (TEventList*)gDirectory->Get("yplus");
     tree->SetEventList(elist);
     tree->Draw("py");

  If arrays are used in the selection critera, the entry entered in the
  list are all the entries that have at least one element of the array that
  satisfy the selection.
  Example:
      tree.Draw(">>pyplus","fTracks.fPy>0");
      tree->SetEventList(pyplus);
      tree->Draw("fTracks.fPy");
  will draw the fPy of ALL tracks in event with at least one track with
  a positive fPy.

  To select only the elements that did match the original selection
  use TEventList::SetReapplyCut.
  Example:
      tree.Draw(">>pyplus","fTracks.fPy>0");
      pyplus->SetReapplyCut(kTRUE);
      tree->SetEventList(pyplus);
      tree->Draw("fTracks.fPy");
  will draw the fPy of only the tracks that have a positive fPy.

  Note: Use tree->SetEventList(0) if you do not want use the list as input.

      How to obtain more info from TTree::Draw
      ========================================

  Once TTree::Draw has been called, it is possible to access useful
  information still stored in the TTree object via the following functions:
    -GetSelectedRows()    // return the number of entries accepted by the
                          //selection expression. In case where no selection
                          //was specified, returns the number of entries processed.
    -GetV1()              //returns a pointer to the double array of V1
    -GetV2()              //returns a pointer to the double array of V2
    -GetV3()              //returns a pointer to the double array of V3
    -GetW()               //returns a pointer to the double array of Weights
                          //where weight equal the result of the selection expression.
   where V1,V2,V3 correspond to the expressions in
   TTree::Draw("V1:V2:V3",selection);

   Example:
    Root > ntuple->Draw("py:px","pz>4");
    Root > TGraph *gr = new TGraph(ntuple->GetSelectedRows(),
                                   ntuple->GetV2(), ntuple->GetV1());
    Root > gr->Draw("ap"); //draw graph in current pad
    creates a TGraph object with a number of points corresponding to the
    number of entries selected by the expression "pz>4", the x points of the graph
    being the px values of the Tree and the y points the py values.

   Important note: By default TTree::Draw creates the arrays obtained
    with GetV1, GetV2, GetV3, GetW with a length corresponding to the
    parameter fEstimate. By default fEstimate=10000 and can be modified
    via TTree::SetEstimate. A possible recipee is to do
       tree->SetEstimate(tree->GetEntries());
    You must call SetEstimate if the expected number of selected rows
    is greater than 10000.

    You can use the option "goff" to turn off the graphics output
    of TTree::Draw in the above example.

           Automatic interface to TTree::Draw via the TTreeViewer
           ======================================================

    A complete graphical interface to this function is implemented
    in the class TTreeViewer.
    To start the TTreeViewer, three possibilities:
       - select TTree context menu item "StartViewer"
       - type the command  "TTreeViewer TV(treeName)"
       - execute statement "tree->StartViewer();"


void DropBaskets()

void DropBuffers(Int_t)
*-*-*-*-*Drop branch buffers to accomodate nbytes below MaxVirtualsize*-*-*-*

Int_t Fill()
*-*-*-*-*Fill all branches of a Tree*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
*-*      ===========================

   This function loops on all the branches of this tree.
   For each branch, it copies to the branch buffer (basket) the current
   values of the leaves data types.
   If a leaf is a simple data type, a simple conversion to a machine
   independent format has to be done.

   The function returns the number of bytes committed to the
   individual branch(es).
   If a write error occurs, the number of bytes returned is -1.
   If no data are written, because e.g. the branch is disabled,
   the number of bytes returned is 0.


TBranch* FindBranch(const char* branchname)
 We already have been visited while recursively looking
 through the friends tree, let return

TLeaf* FindLeaf(const char* searchname)
 We already have been visited while recursively looking
 through the friends tree, let return

Long64_t Fit(const char *funcname ,const char *varexp, const char *selection,Option_t *option ,Option_t *goption,Long64_t nentries, Long64_t firstentry)
*-*-*-*-*-*-*-*-*Fit  a projected item(s) from a Tree*-*-*-*-*-*-*-*-*-*
*-*              ======================================

  funcname is a TF1 function.

  See TTree::Draw for explanations of the other parameters.

  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.

  The function returns the number of selected entries.

  Example:
    tree.Fit(pol4,sqrt(x)>>hsqrt,y>0)
    will fit sqrt(x) and save the histogram as "hsqrt" in the current
    directory.

   See also TTree::UnbinnedFit

const char* GetAlias(const char *aliasName) const
 Returns the expanded value of the alias.  Search in the friends if any

TBranch* GetBranch(const char *name)
 Return pointer to the branch with name in this Tree or the list
 of friends of this tree.

Bool_t GetBranchStatus(const char *branchname) const
 return status of branch with name branchname
 0 if branch is not activated
 1 if branch is activated

Int_t GetBranchStyle()
 static function returning the current branch style
 style = 0 old Branch
 style = 1 new Bronch

TFile* GetCurrentFile() const
*-*-*-*-*-*Return pointer to the current file*-*-*-*-*-*-*-*
*-*        ==================================

Long64_t GetEntriesFriend() const
 return number of entries of this Tree if not zero
 otherwise return the number of entries in the first friend Tree.

Int_t GetEntry(Long64_t entry, Int_t getall)
*-*-*-*-*-*Read all branches of entry and return total number of bytes*-*-*
*-*        ===========================================================
     getall = 0 : get only active branches
     getall = 1 : get all branches

  The function returns the number of bytes read from the input buffer.
  If entry does not exist  the function returns 0.
  If an I/O error occurs,  the function returns -1.

  If the Tree has friends, also read the friends entry

  To activate/deactivate one or more branches, use TBranch::SetBranchStatus
  For example, if you have a Tree with several hundred branches, and you
  are interested only by branches named "u" and "v", do
     mytree.SetBranchStatus("*",0); //disable all branches
     mytree.SetBranchStatus("a",1);
     mytree.SetBranchStatus("b",1);
  when calling mytree.GetEntry(i); only branches "a" and "b" will be read.

  WARNING!!
  If your Tree has been created in split mode with a parent branch "parent",
     mytree.SetBranchStatus("parent",1);
  will not activate the sub-branches of "parent". You should do:
     mytree.SetBranchStatus("parent*",1);

  An alternative is to call directly
     brancha.GetEntry(i)
     branchb.GetEntry(i);

  IMPORTANT NOTE
  ==============
 By default, GetEntry reuses the space allocated by the previous object
 for each branch. You can force the previous object to be automatically
 deleted if you call mybranch.SetAutoDelete(kTRUE) (default is kFALSE).
 Example:
 Consider the example in $ROOTSYS/test/Event.h
 The top level branch in the tree T is declared with:
    Event *event = 0;  //event must be null or point to a valid object
                       //it must be initialized
    T.SetBranchAddress("event",&event);
 When reading the Tree, one can choose one of these 3 options:

   OPTION 1
   --------

    for (Long64_t i=0;i<nentries;i++) {
       T.GetEntry(i);
       // the object event has been filled at this point
    }
   The default (recommended). At the first entry an object of the
   class Event will be created and pointed by event.
   At the following entries, event will be overwritten by the new data.
   All internal members that are TObject* are automatically deleted.
   It is important that these members be in a valid state when GetEntry
   is called. Pointers must be correctly initialized.
   However these internal members will not be deleted if the characters "->"
   are specified as the first characters in the comment field of the data
   member declaration.
   If "->" is specified, the pointer member is read via pointer->Streamer(buf).
   In this case, it is assumed that the pointer is never null (case
   of pointer TClonesArray *fTracks in the Event example).
   If "->" is not specified, the pointer member is read via buf >> pointer.
   In this case the pointer may be null. Note that the option with "->"
   is faster to read or write and it also consumes less space in the file.

   OPTION 2
   --------
  The option AutoDelete is set
   TBranch *branch = T.GetBranch("event");
   branch->SetAddress(&event);
   branch->SetAutoDelete(kTRUE);
    for (Long64_t i=0;i<nentries;i++) {
       T.GetEntry(i);
       // the objrect event has been filled at this point
    }
   In this case, at each iteration, the object event is deleted by GetEntry
   and a new instance of Event is created and filled.

   OPTION 3
   --------
   Same as option 1, but you delete yourself the event.
    for (Long64_t i=0;i<nentries;i++) {
       delete event;
       event = 0;  // EXTREMELY IMPORTANT
       T.GetEntry(i);
       // the objrect event has been filled at this point
    }

  It is strongly recommended to use the default option 1. It has the
  additional advantage that functions like TTree::Draw (internally
  calling TTree::GetEntry) will be functional even when the classes in the
  file are not available.

Long64_t GetEntryNumber(Long64_t entry) const
*-*-*-*-*-*Return entry number corresponding to entry*-*-*
*-*        ==========================================
     if no selection list returns entry
     else returns the entry number corresponding to the list index=entry

Long64_t GetEntryNumberWithBestIndex(Int_t major, Int_t minor) const
 Return entry number corresponding to major and minor number
 Note that this function returns only the entry number, not the data
 To read the data corresponding to an entry number, use TTree::GetEntryWithIndex
 the BuildIndex function has created a table of Double_t* of sorted values
 corresponding to val = major + minor*1e-9;
 The function performs binary search in this sorted table.
 If it finds a pair that maches val, it returns directly the
 index in the table.
 If an entry corresponding to major and minor is not found, the function
 returns the index of the major,minor pair immediatly lower than the
 requested value, ie it will return -1 if the pair is lower than
 the first entry in the index.

 See also GetEntryNumberWithIndex

Long64_t GetEntryNumberWithIndex(Int_t major, Int_t minor) const
 Return entry number corresponding to major and minor number
 Note that this function returns only the entry number, not the data
 To read the data corresponding to an entry number, use TTree::GetEntryWithIndex
 the BuildIndex function has created a table of Double_t* of sorted values
 corresponding to val = major + minor*1e-9;
 The function performs binary search in this sorted table.
 If it finds a pair that maches val, it returns directly the
 index in the table, otherwise it returns -1.

 See also GetEntryNumberWithBestIndex

Int_t GetEntryWithIndex(Int_t major, Int_t minor)
  Read entry  corresponding to major and minor number
  The function returns the total number of bytes read.
  If the Tree has friend trees, the corresponding entry with
  the index values (major,minor) is read. Note that the master Tree
  and its friend may have different entry serial numbers corresponding
  to (major,minor).

const char* GetFriendAlias(TTree *tree) const
 If the the 'tree' is a friend, this method returns its alias name
 This 'alias' is a an alias for the TTree itself.
 It can be used in conjunction with a branch or leaf name in a TTreeFormula.
  Is can also be used in conjunction with an alias created using
  TTree::SetAlias in a TTreeFormula, eg:
     maintree->Draw("treealias.fPx - treealias.myAlias");
  where fPx is a branch of the friend tree aliased as 'treealias' and 'myAlias;
  was created using TTree::SetAlias on the tree aliases as 'treealias'.

  However, note that 'treealias.myAlias' will be expanded literally, without
  'remembering' it comes from the aliased friend and thus the branch
  name might not be disambiguated properly.

TIterator* GetIteratorOnAllLeaves(Bool_t dir)
 Creates a new iterator that will go through all the leaves on the tree
 itself and its friend.

TLeaf* GetLeaf(const char *aname)
 Return pointer to the 1st Leaf named name in any Branch of this Tree
 or any branch in the list of friend trees.

  aname may be of the form branchname/leafname

Double_t GetMaximum(const char *columname)
*-*-*-*-*-*-*-*-*Return maximum of column with name columname*-*-*-*-*-*-*
*-*              ============================================

Long64_t GetMaxTreeSize()
 static function
 return maximum size of a Tree file

Double_t GetMinimum(const char *columname)
*-*-*-*-*-*-*-*-*Return minimum of column with name columname*-*-*-*-*-*-*
*-*              ============================================

const char* GetNameByIndex(TString &varexp, Int_t *index,Int_t colindex) const
*-*-*-*-*-*-*-*-*Return name corresponding to colindex in varexp*-*-*-*-*-*
*-*              ===============================================

   varexp is a string of names separated by :
   index is an array with pointers to the start of name[i] in varexp


TVirtualTreePlayer* GetPlayer()
 Load the TTreePlayer (if not already done)
 Pointer to player is fPlayer

TList* GetUserInfo()
 return a pointer to the list containing user objects associated to this Tree
 The list is automatically created if it does not exist

void KeepCircular()
 keep a maximum of fMaxEntries in memory

Long64_t LoadTree(Long64_t entry)
*-*-*-*-*-*-*-*-*Set current Tree entry
*-*              ======================

 returns -2 if entry does not exist (just as TChain::LoadTree())

Int_t LoadBaskets(Long64_t maxmemory)
 Read in memory all baskets from all branches up to the limit
 of maxmemory bytes.
 If maxmemory is non null and positive SetMaxVirtualSize is called
 with this value. Default for maxmemory is 2000000000 (2 Gigabytes).
 The function returns the total number of baskets read into memory
 if negative an error occured while loading the branches.
 This method may be called to force branch baskets in memory
 when random access to branch entries is required.
 If random access to only a few branches is required, you should
 call directly TBranch::LoadBaskets.

Long64_t LoadTreeFriend(Long64_t entry, TTree *T)
 called by TTree::LoadTree when TTree *T looks for the entry
 number in a friend Tree (this) corresponding to the entry number in T.
 If the friend Tree has no TTreeIndex, entry in the friend and entry
 in T are the same.
 If the friend Tree has an index, one must find the value pair major,minor
 in T to locate the corresponding entry in the friend Tree.

Int_t MakeSelector(const char *selector)
 Generate skeleton selector class for this Tree

 The following files are produced: selector.h and selector.C.
 If selector is 0, the selector will be called "nameoftree".

 The generated code in selector.h includes the following:
    - Identification of the original Tree and Input file name
    - Definition of selector class (data and functions)
    - The following class functions:
       - constructor and destructor
       - void    Begin(TTree *tree)
       - void    Init(TTree *tree)
       - Bool_t  Notify()
       - Bool_t  Process(Long64_t entry)
       - void    Terminate

 The class selector derives from TSelector.
 The generated code in selector.C includes empty functions defined above:

 To use this function:
    - connect your Tree file (eg: TFile f("myfile.root");)
    - T->MakeSelector("myselect");
 where T is the name of the Tree in file myfile.root
 and myselect.h, myselect.C the name of the files created by this function.
 In a ROOT session, you can do:
    root > T->Process("select.C")

Int_t MakeProxy(const char *proxyClassname, const char *macrofilename, const char *cutfilename, const char *option, Int_t maxUnrolling)
 Generate a skeleton analysis class for this Tree using TBranchProxy.
 TBranchProxy is the base of a class hierarchy implementing an
 indirect access to the content of the branches of a TTree.

 "proxyClassname" is expected to be of the form:
    [path/]fileprefix
 The skeleton will then be generated in the file:
    fileprefix.h
 located in the current directory or in 'path/' if it is specified.
 The class generated will be named 'fileprefix'

 "macrofilename" and optionally "cutfilename" are expected to point
 to source file which will be included in by the generated skeletong.
 Method of the same name as the file(minus the extension and path)
 will be called by the generated skeleton's Process method as follow:
    [if (cutfilename())] htemp->Fill(macrofilename());

 "option" can be used select some of the optional features during
 the code generation.  The possible options are:
    nohist : indicates that the generated ProcessFill should not
             fill the histogram.

 'maxUnrolling' controls how deep in the class hierachy does the
 system 'unroll' class that are not split.  'unrolling' a class
 will allow direct access to its data members a class (this
 emulates the behavior of TTreeFormula).

 The main features of this skeleton are:

    * on-demand loading of branches
    * ability to use the 'branchname' as if it was a data member
    * protection against array out-of-bound
    * ability to use the branch data as object (when the user code is available)

 For example with Event.root, if
    Double_t somepx = fTracks.fPx[2];
 is executed by one of the method of the skeleton,
 somepx will updated with the current value of fPx of the 3rd track.

 Both macrofilename and the optional cutfilename are expected to be
 the name of source files which contain at least a free standing
 function with the signature:
     x_t macrofilename(); // i.e function with the same name as the file
 and
     y_t cutfilename();   // i.e function with the same name as the file

 x_t and y_t needs to be types that can convert respectively to a double
 and a bool (because the skeleton uses:
     if (cutfilename()) htemp->Fill(macrofilename());

 This 2 functions are run in a context such that the branch names are
 available as local variables of the correct (read-only) type.

 Note that if you use the same 'variable' twice, it is more efficient
 to 'cache' the value. For example
   Int_t n = fEventNumber; // Read fEventNumber
   if (n<10 || n>10) { ... }
 is more efficient than
   if (fEventNumber<10 || fEventNumber>10)

 Also, optionally, the generated selector will also call methods named
 macrofilename_methodname in each of 6 main selector methods if the method
 macrofilename_methodname exist (Where macrofilename is stripped of its
 extension).

 Concretely, with the script named h1analysisProxy.C,

 The method         calls the method (if it exist)
 Begin           -> h1analysisProxy_Begin
 SlaveBegin      -> h1analysisProxy_SlaveBegin
 Notify          -> h1analysisProxy_Notify
 Process         -> h1analysisProxy_Process
 SlaveTerminate  -> h1analysisProxy_SlaveTerminate
 Terminate       -> h1analysisProxy_Terminate

 If a file name macrofilename.h (or .hh, .hpp, .hxx, .hPP, .hXX) exist
 it is included before the declaration of the proxy class.  This can
 be used in particular to insure that the include files needed by
 the macro file are properly loaded.

 The default histogram is accessible via the variable named 'htemp'.

 If the library of the classes describing the data in the branch is
 loaded, the skeleton will add the needed #include statements and
 give the ability to access the object stored in the branches.

 To draw px using the file hsimple.root (generated by the
 hsimple.C tutorial), we need a file named hsimple.cxx:

     double hsimple() {
        return px;
     }

 MakeProxy can then be used indirectly via the TTree::Draw interface
 as follow:
     new TFile("hsimple.root")
     ntuple->Draw("hsimple.cxx");

 A more complete example is available in the tutorials directory:
   h1analysisProxy.cxx , h1analysProxy.h and h1analysisProxyCut.C
 which reimplement the selector found in h1analysis.C

Int_t MakeClass(const char *classname, Option_t *option)
 Generate skeleton analysis class for this Tree

 The following files are produced: classname.h and classname.C
 If classname is 0, classname will be called "nameoftree.

 The generated code in classname.h includes the following:
    - Identification of the original Tree and Input file name
    - Definition of analysis class (data and functions)
    - the following class functions:
       - constructor (connecting by default the Tree file)
       - GetEntry(Long64_t entry)
       - Init(TTree *tree) to initialize a new TTree
       - Show(Long64_t entry) to read and Dump entry

 The generated code in classname.C includes only the main
 analysis function Loop.

 To use this function:
    - connect your Tree file (eg: TFile f("myfile.root");)
    - T->MakeClass("MyClass");
 where T is the name of the Tree in file myfile.root
 and MyClass.h, MyClass.C the name of the files created by this function.
 In a ROOT session, you can do:
    root > .L MyClass.C
    root > MyClass t
    root > t.GetEntry(12); // Fill t data members with entry number 12
    root > t.Show();       // Show values of entry 12
    root > t.Show(16);     // Read and show values of entry 16
    root > t.Loop();       // Loop on all entries

Int_t MakeCode(const char *filename)
 Generate skeleton function for this Tree

 The function code is written on filename.
 If filename is 0, filename will be called nameoftree.C

 The generated code includes the following:
    - Identification of the original Tree and Input file name
    - Connection of the Tree file
    - Declaration of Tree variables
    - Setting of branches addresses
    - A skeleton for the entry loop

 To use this function:
    - connect your Tree file (eg: TFile f("myfile.root");)
    - T->MakeCode("anal.C");
 where T is the name of the Tree in file myfile.root
 and anal.C the name of the file created by this function.

 NOTE: Since the implementation of this function, a new and better
       function TTree::MakeClass() has been developped.

void MakeIndex(TString &varexp, Int_t *index)
*-*-*-*-*-*-*-*-*Build Index array for names in varexp*-*-*-*-*-*-*-*-*-*-*
*-*              =====================================

Bool_t MemoryFull(Int_t nbytes)
*-*-*-*-*-*Check if adding nbytes to memory we are still below MaxVirtualsize
*-*        ==================================================================

TTree* MergeTrees(TList *li)
static function merging the Trees in the TList into a new Tree.
Trees in the list can be memory or disk-resident trees
The new tree is created in the current directory (memory if gROOT)

Long64_t Merge(TCollection *li)
function merging the Trees in the TList into this Tree.
 return the total number of entries in the merged Tree

Bool_t Notify()
 function called when loading a new class library

TPrincipal* Principal(const char *varexp, const char *selection, Option_t *option, Long64_t nentries, Long64_t firstentry)
*-*-*-*-*-*-*-*-*Interface to the Principal Components Analysis class*-*-*
*-*              ====================================================

   Create an instance of TPrincipal
   Fill it with the selected variables
   if option "n" is specified, the TPrincipal object is filled with
                 normalized variables.
   If option "p" is specified, compute the principal components
   If option "p" and "d" print results of analysis
   If option "p" and "h" generate standard histograms
   If option "p" and "c" generate code of conversion functions
   return a pointer to the TPrincipal object. It is the user responsability
   to delete this object.
   The option default value is "np"

   see TTree::Draw for explanation of the other parameters.

   The created object is  named "principal" and a reference to it
   is added to the list of specials Root objects.
   you can retrieve a pointer to the created object via:
      TPrincipal *principal =
        (TPrincipal*)gROOT->GetListOfSpecials()->FindObject("principal");


void Print(Option_t *option) const
 Print a summary of the Tree contents.
 if option contains "all" friend trees are also printed.
 if option contains "toponly" only the top level branches are printed.

 Wildcarding can be used to print only a subset of the branches
 eg, T.Print("Elec*") will print all branches with name starting with "Elec"

Long64_t Process(const char *filename,Option_t *option,Long64_t nentries, Long64_t firstentry)
*-*-*-*-*-*-*-*-*Process this tree executing the code in filename*-*-*-*-*
*-*              ================================================

   The code in filename is loaded (interpreted or compiled , see below)
   filename must contain a valid class implementation derived from TSelector.
   where TSelector has the following member functions:

     void TSelector::Begin(). This function is called before looping on the
          events in the Tree. The user can create his histograms in this function.

     Bool_t TSelector::ProcessCut(Long64_t entry). This function is called
          before processing entry. It is the user's responsability to read
          the corresponding entry in memory (may be just a partial read).
          The function returns kTRUE if the entry must be processed,
          kFALSE otherwise.
     void TSelector::ProcessFill(Long64_t entry). This function is called for
          all selected events. User fills histograms in this function.
     void TSelector::Terminate(). This function is called at the end of
          the loop on all events.
     void TTreeProcess::Begin(). This function is called before looping on the
          events in the Tree. The user can create his histograms in this function.

   if filename is of the form file.C, the file will be interpreted.
   if filename is of the form file.C++, the file file.C will be compiled
      and dynamically loaded.
   if filename is of the form file.C+, the file file.C will be compiled
      and dynamically loaded. At next call, if file.C is older than file.o
      and file.so, the file.C is not compiled, only file.so is loaded.

   The function returns the number of processed entries. It returns -1
   in case of an error.

  NOTE1
  It may be more interesting to invoke directly the other Process function
  accepting a TSelector* as argument.eg
     MySelector *selector = (MySelector*)TSelector::GetSelector(filename);
     selector->CallSomeFunction(..);
     mytree.Process(selector,..);

  NOTE2
  One should not call this function twice with the same selector file
  in the same script. If this is required, proceed as indicated in NOTE1,
  by getting a pointer to the corresponding TSelector,eg
    workaround 1
    ------------
void stubs1() {
   TSelector *selector = TSelector::GetSelector("h1test.C");
   TFile *f1 = new TFile("stubs_nood_le1.root");
   TTree *h1 = (TTree*)f1->Get("h1");
   h1->Process(selector);
   TFile *f2 = new TFile("stubs_nood_le1_coarse.root");
   TTree *h2 = (TTree*)f2->Get("h1");
   h2->Process(selector);
}
  or use ACLIC to compile the selector
   workaround 2
   ------------
void stubs2() {
   TFile *f1 = new TFile("stubs_nood_le1.root");
   TTree *h1 = (TTree*)f1->Get("h1");
   h1->Process("h1test.C+");
   TFile *f2 = new TFile("stubs_nood_le1_coarse.root");
   TTree *h2 = (TTree*)f2->Get("h1");
   h2->Process("h1test.C+");
}

Long64_t Process(TSelector *selector,Option_t *option, Long64_t nentries, Long64_t firstentry)
*-*-*-*-*-*-*-*-*Process this tree executing the code in selector*-*-*-*-*
*-*              ================================================

   The TSelector class has the following member functions:

     void TSelector::Begin(). This function is called before looping on the
          events in the Tree. The user can create his histograms in this function.

     Bool_t TSelector::ProcessCut(Long64_t entry). This function is called
          before processing entry. It is the user's responsability to read
          the corresponding entry in memory (may be just a partial read).
          The function returns kTRUE if the entry must be processed,
          kFALSE otherwise.
     void TSelector::ProcessFill(Long64_t entry). This function is called for
          all selected events. User fills histograms in this function.
     void TSelector::Terminate(). This function is called at the end of
          the loop on all events.
     void TTreeProcess::Begin(). This function is called before looping on the
          events in the Tree. The user can create his histograms in this function.

Long64_t Project(const char *hname, const char *varexp, const char *selection, Option_t *option,Long64_t nentries, Long64_t firstentry)
*-*-*-*-*-*-*-*-*Make a projection of a Tree using selections*-*-*-*-*-*-*
*-*              =============================================

   Depending on the value of varexp (described in Draw) a 1-D,2-D,etc
   projection of the Tree will be filled in histogram hname.
   Note that the dimension of hname must match with the dimension of varexp.


TSQLResult* Query(const char *varexp, const char *selection, Option_t *option, Long64_t nentries, Long64_t firstentry)
 Loop on Tree & return TSQLResult object containing entries following selection

Long64_t ReadFile(const char *filename, const char *branchDescriptor)
 Create or simply read branches from filename
 if branchDescriptor = "" (default), it is assumed that the Tree descriptor
    is given in the first line of the file with a syntax like
     A/D:Table[2]/F:Ntracks/I:astring/C
  otherwise branchDescriptor must be specified with the above syntax.

 Lines in the input file starting with "#" are ignored.

 A TBranch object is created for each variable in the expression.
 The total number of rows read from the file is returned.

void Refresh()
  Refresh contents of this Tree and his branches from the current
  Tree status on its file
  One can call this function in case the Tree on its file is being
  updated by another process

void RemoveFriend(TTree *oldFriend)
*-*-*-*-*-*-*-*Remove a friend from the list of friend *-*-*
*-*            =============================================

void Reset(Option_t *option)
*-*-*-*-*-*-*-*Reset buffers and entries count in all branches/leaves*-*-*
*-*            ======================================================

void ResetBranchAddresses()
*-*-*-*-*-*-*-*Reset the address of the branches                *-*-*
*-*            ======================================================

Long64_t Scan(const char *varexp, const char *selection, Option_t *option, Long64_t nentries, Long64_t firstentry)
 Loop on Tree and print entries passing selection. If varexp is 0 (or "")
 then print only first 8 columns. If varexp = "*" print all columns.
 Otherwise a columns selection can be made using "var1:var2:var3".
 see TTreePlayer::Scan for more information

Bool_t SetAlias(const char *aliasName, const char *aliasFormula)
  Set a tree variable alias

  Set an alias for an expression/formula based on the tree 'variables'.

  The content of 'aliasName' can be used in TTreeFormula (i.e. TTree::Draw,
  TTree::Scan, TTreeViewer) and will be evaluated as the content of
  'aliasFormula'.
  If the content of 'aliasFormula only contains symbol names, periods and
  array index specification (for example event.fTracks[3]), then
  the content of 'aliasName' can be used as the start of symbol.

  If the alias 'aliasName' already existed, it is replaced by the new
  value.

  When being used, the alias can be preceded by an eventual 'Friend Alias'
  (see TTree::GetFriendAlias)

  Return true if it was added properly.

  For example:
     tree->SetAlias("x1","(tdc1[1]-tdc1[0])/49");
     tree->SetAlias("y1","(tdc1[3]-tdc1[2])/47");
     tree->SetAlias("x2","(tdc2[1]-tdc2[0])/49");
     tree->SetAlias("y2","(tdc2[3]-tdc2[2])/47");
     tree->Draw("y2-y1:x2-x1");

     tree->SetAlias("theGoodTrack","event.fTracks[3]");
     tree->Draw("theGoodTrack.fPx"); // same as "event.fTracks[3].fPx"

void SetBasketSize(const char *bname, Int_t buffsize)
*-*-*-*-*-*-*-*-*Set branc(es) basket size*-*-*-*-*-*-*-*
*-*              =========================

     bname is the name of a branch.
     if bname="*", apply to all branches.
     if bname="xxx*", apply to all branches with name starting with xxx
     see TRegexp for wildcarding options
     buffsize = branc basket size

void SetBranchAddress(const char *bname, void *add)
*-*-*-*-*-*-*-*-*Set branch address*-*-*-*-*-*-*-*
*-*              ==================

      If object is a TTree, this function is only an interface to TBranch::SetAddress
      Function overloaded by TChain.

void SetBranchAddress(const char *bname, void *add, TClass *ptrClass, EDataType datatype, Bool_t ptr)
  Verify the validity of the type of add before calling SetBranchAddress.

void SetBranchStatus(const char *bname, Bool_t status, UInt_t *found)
*-*-*-*-*-*-*-*-*Set branch status Process or DoNotProcess*-*-*-*-*-*-*-*
*-*              =========================================

  When reading a Tree, by default, all branches are read.
  One can speed up considerably the analysis phase by activating
  only the branches that hold variables involved in a query.

     bname is the name of a branch.
     if bname="*", apply to all branches.
     if bname="xxx*", apply to all branches with name starting with xxx
     see TRegexp for wildcarding options
      status = 1  branch will be processed
             = 0  branch will not be processed
    Example:
  Assume a tree T with sub-branches a,b,c,d,e,f,g,etc..
  when doing T.GetEntry(i) all branches are read for entry i.
  to read only the branches c and e, one can do
    T.SetBranchStatus("*",0); //disable all branches
    T.SetBranchStatus("c",1);
    T.setBranchStatus("e",1);
    T.GetEntry(i);

  WARNING!!
  If your Tree has been created in split mode with a parent branch "parent",
     T.SetBranchStatus("parent",1);
  will not activate the sub-branches of "parent". You should do:
     T.SetBranchStatus("parent*",1);

  An alternative to this function is to read directly and only
  the interesting branches. Example:
    TBranch *brc = T.GetBranch("c");
    TBranch *bre = T.GetBranch("e");
    brc->GetEntry(i);
    bre->GetEntry(i);

  If found is not 0, the number of branch(es) found matching the regular
  expression is returned in *found AND the error message 'unknown branch'
  is suppressed.

void SetBranchStyle(Int_t style)
 static function setting the current branch style
 style = 0 old Branch
 style = 1 new Bronch

void SetCircular(Long64_t maxEntries)
 Organize this Tree with circular buffers, keeping in memory
 a maximum of maxEntries

void SetDebug(Int_t level, Long64_t min, Long64_t max)
 Set the debug level and the debug range
 for entries in the debug range, the functions TBranchElement::Fill
 and TBranchElement::GetEntry will print the number of bytes filled
 or read for each branch.

void SetDirectory(TDirectory *dir)
 Remove reference to this tree from current directory and add
 reference to new directory dir. dir can be 0 in which case the tree
 does not belong to any directory.

Long64_t SetEntries(Long64_t n)
 if n >= 0 Set number of entries in the Tree = n.

 if (n < 0) Set number of entries in the Tree to match the
 number of entries in each branch. (default for n is -1)
 This function should be called only when one fills each branch
 independently via TBranch::Fill without calling TTree::Fill
 Calling TTree::SetEntries() make sense only if the number of entries
 in each branch is identical.  A Warning is issued otherwise.
 The function returns the number of entries.

void SetEstimate(Long64_t n)
*-*-*-*-*-*-*-*-*Set number of entries to estimate variable limits*-*-*-*
*-*              ================================================

void SetFileNumber(Int_t number)
 Set fFileNumber to number.
 fFileNumber is used by TTree::Fill to set the file name
 for a new file to be created when the current file exceeds fgTreeMaxSize.
    (see TTree::ChangeFile)
 if fFileNumber=10, the new file name will have a suffix "_11",
 ie, fFileNumber is incremented before setting the file name

void SetMaxTreeSize(Long64_t maxsize)
 static function
 Set the maximum size of a Tree file.
 In TTree::fill, when the file has a size > fgMaxTreeSize,
 the function closes the current file and starts writing into
 a new file with a name of the style "file_1.root" if the original
 requested file name was "file.root"

void SetName(const char *name)
 Change the name of this Tree


void SetObject(const char *name, const char *title)
 Change the name and title of this Tree


void SetWeight(Double_t w, Option_t *)
  Set tree weight.
  The weight is used by TTree::Draw to automatically weight each
  selected entry in the resulting histogram.
  For example the equivalent of
     T.Draw("x","w")
  is
     T.SetWeight(w);
     T.Draw("x");

 This function is redefined by TChain::SetWeight. In case of a TChain,
 an option "global" may be specified to set the same weight
 for all Trees in the TChain instead of the default behaviour
 using the weights of each Tree in the chain. (see TChain::SetWeight)

void Show(Long64_t entry, Int_t lenmax)
*-*-*-*-*-*Print values of all active leaves for entry*-*-*-*-*-*-*-*
*-*        ===========================================
 if entry==-1, print current entry (default)
 if a leaf is an array, a maximum of lenmax elements is printed.


void StartViewer()
*-*-*-*-*-*-*-*-*Start the TTreeViewer on this TTree*-*-*-*-*-*-*-*-*-*
*-*              ===================================

  ww is the width of the canvas in pixels
  wh is the height of the canvas in pixels

void Streamer(TBuffer &b)
*-*-*-*-*-*-*-*-*Stream a class object*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
*-*              =========================================

Long64_t UnbinnedFit(const char *funcname ,const char *varexp, const char *selection,Option_t *option ,Long64_t nentries, Long64_t firstentry)
*-*-*-*-*-*Unbinned fit of one or more variable(s) from a Tree*-*-*-*-*-*
*-*        ===================================================

  funcname is a TF1 function.

  See TTree::Draw for explanations of the other parameters.

   Fit the variable varexp using the function funcname using the
   selection cuts given by selection.

   The list of fit options is given in parameter option.
      option = "Q" Quiet mode (minimum printing)
             = "V" Verbose mode (default is between Q and V)
             = "E" Perform better Errors estimation using Minos technique
             = "M" More. Improve fit results

   You can specify boundary limits for some or all parameters via
        func->SetParLimits(p_number, parmin, parmax);
   if parmin>=parmax, the parameter is fixed
   Note that you are not forced to fix the limits for all parameters.
   For example, if you fit a function with 6 parameters, you can do:
     func->SetParameters(0,3.1,1.e-6,0.1,-8,100);
     func->SetParLimits(4,-10,-4);
     func->SetParLimits(5, 1,1);
   With this setup, parameters 0->3 can vary freely
   Parameter 4 has boundaries [-10,-4] with initial value -8
   Parameter 5 is fixed to 100.

   For the fit to be meaningful, the function must be self-normalized.

   i.e. It must have the same integral regardless of the parameter
   settings.  Otherwise the fit will effectively just maximize the
   area.

   It is mandatory to have a normalization variable
   which is fixed for the fit.  e.g.

     TF1* f1 = new TF1("f1", "gaus(0)/sqrt(2*3.14159)/[2]", 0, 5);
     f1->SetParameters(1, 3.1, 0.01);
     f1->SetParLimits(0, 1, 1); // fix the normalization parameter to 1
     data->UnbinnedFit("f1", "jpsimass", "jpsipt>3.0");
   

   1, 2 and 3 Dimensional fits are supported.
   See also TTree::Fit

void UseCurrentStyle()
*-*-*-*-*-*Replace current attributes by current style*-*-*-*-*
*-*        ===========================================



Inline Functions


                   void AddTotBytes(Int_t tot)
                   void AddZipBytes(Int_t zip)
                  Int_t Debug() const
               Long64_t Draw(const char* varexp, const char* selection, Option_t* option = "", Long64_t nentries = 1000000000, Long64_t firstentry = 0)
            TBranchRef* GetBranchRef() const
               Long64_t GetChainEntryNumber(Long64_t entry) const
               Long64_t GetChainOffset() const
               Long64_t GetDebugMax() const
               Long64_t GetDebugMin() const
            TDirectory* GetDirectory() const
               Long64_t GetEntries() const
               Long64_t GetEntriesFast() const
               Long64_t GetEstimate() const
                  Int_t GetEvent(Long64_t entry = 0, Int_t getall = 0)
            TEventList* GetEventList() const
                  Int_t GetFileNumber() const
                   TH1* GetHistogram()
                 Int_t* GetIndex()
              Double_t* GetIndexValues()
                 TList* GetListOfClones()
             TObjArray* GetListOfBranches()
             TObjArray* GetListOfLeaves()
                 TList* GetListOfFriends() const
        TSeqCollection* GetListOfAliases() const
                  Int_t GetMakeClass() const
               Long64_t GetMaxEntryLoop() const
               Long64_t GetMaxVirtualSize() const
                  Int_t GetNbranches()
               TObject* GetNotify() const
                  Int_t GetPacketSize() const
               Long64_t GetReadEntry() const
               Long64_t GetReadEvent() const
                  Int_t GetScanField() const
          TTreeFormula* GetSelect()
               Long64_t GetSelectedRows()
                  Int_t GetTimerInterval() const
               Long64_t GetTotBytes() const
                 TTree* GetTree() const
         TVirtualIndex* GetTreeIndex() const
                  Int_t GetTreeNumber() const
                  Int_t GetUpdate() const
          TTreeFormula* GetVar1()
          TTreeFormula* GetVar2()
          TTreeFormula* GetVar3()
          TTreeFormula* GetVar4()
              Double_t* GetV1()
              Double_t* GetV2()
              Double_t* GetV3()
              Double_t* GetV4()
              Double_t* GetW()
               Double_t GetWeight() const
               Long64_t GetZipBytes() const
                   void IncrementTotalBuffers(Int_t nbytes)
                 Bool_t IsFolder() const
                   void SetAutoSave(Long64_t autos = 10000000)
                   void SetChainOffset(Int_t offset = 0)
                   void SetEventList(TEventList* list)
                   void SetMakeClass(Int_t make)
                   void SetMaxEntryLoop(Long64_t maxev = 1000000000)
                   void SetMaxVirtualSize(Long64_t size = 0)
                   void SetNotify(TObject* obj)
                   void SetScanField(Int_t n = 50)
                   void SetTimerInterval(Int_t msec = 333)
                   void SetUpdate(Int_t freq = 0)
                TClass* Class()
                TClass* IsA() const
                   void ShowMembers(TMemberInspector& insp, char* parent)
                   void StreamerNVirtual(TBuffer& b)
               TBranch* Branch(const char* name, void*** addobj, Int_t bufsize = 32000, Int_t splitlevel = 99)


Author: Rene Brun 12/01/96
Last update: root/tree:$Name: $:$Id: TTree.cxx,v 1.264 2005/09/08 14:22:17 brun Exp $
Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *


ROOT page - Class index - Class Hierarchy - Top of the page

This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.