ROOT logo

From $ROOTSYS/tutorials/tree/treefriend.C

// Illustrates how to use Tree friends:
//   - create a simple TTree
//   - Copy a subset of this TTree to a new TTree
//   - Create a Tree Index
//   - Make a friend TTree
//   - compare two TTrees
//   - Draw a variable from the first tree versus a variable
//     in the friend Tree
//
// You can run this tutorial with:
//  root > .x treefriend.C  (interpreted via CINT)
//  root > .x treefriend.C+ (executed via ACLIC & the native compiler)
// or, variants like:
//  root > .L treefriend.C+
//  root > CreateParentTree();
//  root > CreateFriendTree();
//  root > CompareTrees();
//  root > DrawFriend();
//
//  Author: Rene Brun
      
#include "TTree.h"
#include "TFile.h"
#include "TRandom.h"
#include "TTree.h"
#include "TTree.h"
   
Int_t Run, Event;
Float_t x,y,z;
   
void CreateParentTree() {
   // create a simple TTree with 5 branches
   // Two branches ("Run" and "Event") will be used to index the Tree
   TFile *f = new TFile("treeparent.root","recreate");
   TTree *T = new TTree("T","test friend trees");
   T->Branch("Run",&Run,"Run/I");
   T->Branch("Event",&Event,"Event/I");
   T->Branch("x",&x,"x/F");
   T->Branch("y",&y,"y/F");
   T->Branch("z",&z,"z/F");
   TRandom r;
   for (Int_t i=0;i<10000;i++) {
      if (i < 5000) Run = 1;
      else          Run = 2;
      Event = i;
      x = r.Gaus(10,1);
      y = r.Gaus(20,2);
      z = r.Landau(2,1);
      T->Fill();
   }
   T->Print();
   T->Write();
   delete f;
}
void CreateFriendTree() {
   // Open the file created by CreateParentTree
   // Copy a subset of the TTree into a new TTree
   //   (see also tutorials copytree.C, copytree2.C and copytree3.C)
   // Create an index on the new TTree ("Run","Event")
   // Write the new TTree (including its index)
   
   TFile *f = new TFile("treeparent.root");
   TTree *T = (TTree*)f->Get("T");
   TFile *ff = new TFile("treefriend.root","recreate");
   TTree *TF = T->CopyTree("z<10");
   TF->SetName("TF");
   TF->BuildIndex("Run","Event");
   TF->Write();
   TF->Print();
   delete ff;
}

void CompareTrees() {
   // The two TTrees created above are compared.
   // The subset of entries in the small TTree must be identical
   // to the entries in the original TTree.
   
   TFile *f = new TFile("treeparent.root");
   TTree *T  = (TTree*)f->Get("T");
   TFile *ff = new TFile("treefriend.root");
   TTree *TF = (TTree*)ff->Get("TF");
   Int_t fRun,fEvent;
   Float_t fx,fy,fz;
   T->SetBranchAddress("Run",&Run);
   T->SetBranchAddress("Event",&Event);
   T->SetBranchAddress("x",&x);
   T->SetBranchAddress("y",&y);
   T->SetBranchAddress("z",&z);
   TF->SetBranchAddress("Run",&fRun);
   TF->SetBranchAddress("Event",&fEvent);
   TF->SetBranchAddress("x",&fx);
   TF->SetBranchAddress("y",&fy);
   TF->SetBranchAddress("z",&fz);
   T->AddFriend(TF);
   
   Long64_t nentries = T->GetEntries();
   Int_t nok = 0;
   for (Long64_t i=0;i<nentries;i++) {
      T->GetEntry(i);
      if (fRun == Run && fEvent==Event && x==fx && y==fy &&z==fz) {
         nok++;
      } else {
         if (TF->GetEntryWithIndex(Run,Event) > 0) {
            if (i <100) printf("i=%lld, Run=%d, Event=%d, x=%g, y=%g, z=%g,  : fRun=%d, fEvent=%d, fx=%g, fy=%g, fz=%g\n",i,Run,Event,x,y,z,fRun,fEvent,fx,fy,fz);
         } 
      } 
   }
   printf("nok = %d, fentries=%lld\n",nok,TF->GetEntries());

   delete f;
   delete ff;
}

void DrawFriend() {
  // Draw a scatter plot of variable x in the parent TTree versus
  // the same variable in the subtree.
  // This should produce points along a straight line.
   
   TFile *f  = TFile::Open("treeparent.root");
   TTree *T  = (TTree*)f->Get("T");
   T->AddFriend("TF","treefriend.root");
   T->Draw("x:TF.x");
}

void treefriend() {
   CreateParentTree();
   CreateFriendTree();
   CompareTrees();
   DrawFriend();
}
 treefriend.C:1
 treefriend.C:2
 treefriend.C:3
 treefriend.C:4
 treefriend.C:5
 treefriend.C:6
 treefriend.C:7
 treefriend.C:8
 treefriend.C:9
 treefriend.C:10
 treefriend.C:11
 treefriend.C:12
 treefriend.C:13
 treefriend.C:14
 treefriend.C:15
 treefriend.C:16
 treefriend.C:17
 treefriend.C:18
 treefriend.C:19
 treefriend.C:20
 treefriend.C:21
 treefriend.C:22
 treefriend.C:23
 treefriend.C:24
 treefriend.C:25
 treefriend.C:26
 treefriend.C:27
 treefriend.C:28
 treefriend.C:29
 treefriend.C:30
 treefriend.C:31
 treefriend.C:32
 treefriend.C:33
 treefriend.C:34
 treefriend.C:35
 treefriend.C:36
 treefriend.C:37
 treefriend.C:38
 treefriend.C:39
 treefriend.C:40
 treefriend.C:41
 treefriend.C:42
 treefriend.C:43
 treefriend.C:44
 treefriend.C:45
 treefriend.C:46
 treefriend.C:47
 treefriend.C:48
 treefriend.C:49
 treefriend.C:50
 treefriend.C:51
 treefriend.C:52
 treefriend.C:53
 treefriend.C:54
 treefriend.C:55
 treefriend.C:56
 treefriend.C:57
 treefriend.C:58
 treefriend.C:59
 treefriend.C:60
 treefriend.C:61
 treefriend.C:62
 treefriend.C:63
 treefriend.C:64
 treefriend.C:65
 treefriend.C:66
 treefriend.C:67
 treefriend.C:68
 treefriend.C:69
 treefriend.C:70
 treefriend.C:71
 treefriend.C:72
 treefriend.C:73
 treefriend.C:74
 treefriend.C:75
 treefriend.C:76
 treefriend.C:77
 treefriend.C:78
 treefriend.C:79
 treefriend.C:80
 treefriend.C:81
 treefriend.C:82
 treefriend.C:83
 treefriend.C:84
 treefriend.C:85
 treefriend.C:86
 treefriend.C:87
 treefriend.C:88
 treefriend.C:89
 treefriend.C:90
 treefriend.C:91
 treefriend.C:92
 treefriend.C:93
 treefriend.C:94
 treefriend.C:95
 treefriend.C:96
 treefriend.C:97
 treefriend.C:98
 treefriend.C:99
 treefriend.C:100
 treefriend.C:101
 treefriend.C:102
 treefriend.C:103
 treefriend.C:104
 treefriend.C:105
 treefriend.C:106
 treefriend.C:107
 treefriend.C:108
 treefriend.C:109
 treefriend.C:110
 treefriend.C:111
 treefriend.C:112
 treefriend.C:113
 treefriend.C:114
 treefriend.C:115
 treefriend.C:116
 treefriend.C:117
 treefriend.C:118
 treefriend.C:119
 treefriend.C:120
 treefriend.C:121
 treefriend.C:122
 treefriend.C:123
 treefriend.C:124
 treefriend.C:125
 treefriend.C:126
 treefriend.C:127
 treefriend.C:128
 treefriend.C:129
 treefriend.C:130
 treefriend.C:131