00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022 #include "TTree.h"
00023 #include "TFile.h"
00024 #include "TRandom.h"
00025 #include "TTree.h"
00026 #include "TTree.h"
00027
00028 Int_t Run, Event;
00029 Float_t x,y,z;
00030
00031 void CreateParentTree() {
00032
00033
00034 TFile *f = new TFile("treeparent.root","recreate");
00035 TTree *T = new TTree("T","test friend trees");
00036 T->Branch("Run",&Run,"Run/I");
00037 T->Branch("Event",&Event,"Event/I");
00038 T->Branch("x",&x,"x/F");
00039 T->Branch("y",&y,"y/F");
00040 T->Branch("z",&z,"z/F");
00041 TRandom r;
00042 for (Int_t i=0;i<10000;i++) {
00043 if (i < 5000) Run = 1;
00044 else Run = 2;
00045 Event = i;
00046 x = r.Gaus(10,1);
00047 y = r.Gaus(20,2);
00048 z = r.Landau(2,1);
00049 T->Fill();
00050 }
00051 T->Print();
00052 T->Write();
00053 delete f;
00054 }
00055 void CreateFriendTree() {
00056
00057
00058
00059
00060
00061
00062 TFile *f = new TFile("treeparent.root");
00063 TTree *T = (TTree*)f->Get("T");
00064 TFile *ff = new TFile("treefriend.root","recreate");
00065 TTree *TF = T->CopyTree("z<10");
00066 TF->SetName("TF");
00067 TF->BuildIndex("Run","Event");
00068 TF->Write();
00069 TF->Print();
00070 delete ff;
00071 }
00072
00073 void CompareTrees() {
00074
00075
00076
00077
00078 TFile *f = new TFile("treeparent.root");
00079 TTree *T = (TTree*)f->Get("T");
00080 TFile *ff = new TFile("treefriend.root");
00081 TTree *TF = (TTree*)ff->Get("TF");
00082 Int_t fRun,fEvent;
00083 Float_t fx,fy,fz;
00084 T->SetBranchAddress("Run",&Run);
00085 T->SetBranchAddress("Event",&Event);
00086 T->SetBranchAddress("x",&x);
00087 T->SetBranchAddress("y",&y);
00088 T->SetBranchAddress("z",&z);
00089 TF->SetBranchAddress("Run",&fRun);
00090 TF->SetBranchAddress("Event",&fEvent);
00091 TF->SetBranchAddress("x",&fx);
00092 TF->SetBranchAddress("y",&fy);
00093 TF->SetBranchAddress("z",&fz);
00094 T->AddFriend(TF);
00095
00096 Int_t nentries = (Int_t)T->GetEntries();
00097 Int_t nok = 0;
00098 for (Int_t i=0;i<nentries;i++) {
00099 T->GetEntry(i);
00100 if (fRun == Run && fEvent==Event && x==fx && y==fy &&z==fz) {
00101 nok++;
00102 } else {
00103 if (TF->GetEntryWithIndex(Run,Event) > 0) {
00104 if (i <100) printf("i=%d, 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);
00105 }
00106 }
00107 }
00108 printf("nok = %d, fentries=%lld\n",nok,TF->GetEntries());
00109
00110 delete f;
00111 delete ff;
00112 }
00113
00114 void DrawFriend() {
00115
00116
00117
00118
00119 TFile *f = TFile::Open("treeparent.root");
00120 TTree *T = (TTree*)f->Get("T");
00121 T->AddFriend("TF","treefriend.root");
00122 T->Draw("x:TF.x");
00123 }
00124
00125 void treefriend() {
00126 CreateParentTree();
00127 CreateFriendTree();
00128 CompareTrees();
00129 DrawFriend();
00130 }