swh:1:snp:af87cd67498ef4fe47c76ed3e7caffe5b61facaf
Raw File
Tip revision: c9efe024ca039cdac6f9fe015dd231f95f7c00c2 authored by Unknown Author on 18 January 2001, 11:26:51 UTC
This commit was manufactured by cvs2svn to create tag 'v3-00-02'.
Tip revision: c9efe02
tcl.C
// Example of macro illustrating how to write a TClonesArray to a TTree
// the following tests can be run
//    Interactive tests
// Root > .x tcl.C        //no-split interpreted
// Root > .x tcl.C(1)     //split    interpreted
// Root > .x tcl.C++      //no-split compiled
// Root > .x tcl.C++(1)   //split    compiled
//    batch tests: same as above but with no graphics
// root -b -q tcl.C
// root -b -q tcl.C++
// root -b -q "tcl.C(1)"
// root -b -q "tcl.C++(1)"

#include "TFile.h"
#include "TClonesArray.h"
#include "TH2.h"
#include "TLine.h"
#include "TTree.h"
#include "TBenchmark.h"
#include "TRandom.h"

void tclwrite(Int_t split)
{
// Generate a Tree with a TClonesArray
// The array can be split or not
   TFile f("tcl.root","recreate");
   f.SetCompressionLevel(1); //try level 2 also
   TTree T("T","test tcl");
   TClonesArray *arr = new TClonesArray("TLine");
   TClonesArray &ar = *arr;
   T.Branch("tcl",&arr,256000,split);
   for (Int_t ev=0;ev<10000;ev++) {
      ar.Clear();
      Int_t nlines = Int_t(gRandom->Gaus(50,10));
      if(nlines < 0) nlines = 1;
      for (Int_t i=0;i<nlines;i++) {
         Float_t x1 = gRandom->Rndm();
         Float_t y1 = gRandom->Rndm();
         Float_t x2 = gRandom->Rndm();
         Float_t y2 = gRandom->Rndm();
         new(ar[i]) TLine(x1,y1,x2,y2);
      }
      T.Fill();
   }
   T.Print();
   T.Write();
}

void tclread()
{
// read file generated by tclwrite
// loop on all entries.
// histogram center of lines
   TFile *f = new TFile("tcl.root");
   TTree *T = (TTree*)f->Get("T");
   TH2F *h2 = new TH2F("h2","center of lines",40,0,1,40,0,1);

   TClonesArray *arr = 0;
   T->GetBranch("tcl")->SetAutoDelete(kFALSE);
   T->SetBranchAddress("tcl",&arr);
   Int_t nentries = (Int_t)(T->GetEntries());
   for (Int_t ev=0;ev<nentries;ev++) {
      arr->Clear();
      T->GetEntry(ev);
      Int_t nlines = arr->GetEntriesFast();
      for (Int_t i=0;i<nlines;i++) {
         TLine *line = (TLine*)arr->At(i);
         h2->Fill(0.5*(line->GetX1()+line->GetX2()), 0.5*(line->GetY1()+line->GetY2()));
      }
   }
   h2->Draw("lego");
}

void tcl(Int_t split=0)
{
   gBenchmark->Start("tcl");
   tclwrite(split);
   tclread();
   gBenchmark->Show("tcl");
}
back to top