Loading [MathJax]/extensions/tex2jax.js
Logo ROOT  
Reference Guide
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
tcl.C File Reference

Detailed Description

View in nbviewer Open in SWAN 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
const char * Root
Definition: TXMLSetup.cxx:46

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);
//By default a TClonesArray is created with its BypassStreamer bit set.
//However, because TLine has a custom Streamer, this bit was reset
//by TTree::Branch above. We set again this bit because the current
//version of TLine uses the automatic Streamer.
//BypassingStreamer saves space and time.
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 y1 = 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 = new TClonesArray("TLine");
T->GetBranch("tcl")->SetAutoDelete(kFALSE);
T->SetBranchAddress("tcl",&arr);
Long64_t nentries = T->GetEntries();
for (Long64_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");
}
#define f(i)
Definition: RSha256.hxx:104
static const double x2[5]
static const double x1[5]
int Int_t
Definition: RtypesCore.h:43
const Bool_t kFALSE
Definition: RtypesCore.h:90
long long Long64_t
Definition: RtypesCore.h:71
float Float_t
Definition: RtypesCore.h:55
R__EXTERN TBenchmark * gBenchmark
Definition: TBenchmark.h:59
int nentries
Definition: THbookFile.cxx:89
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
virtual void Start(const char *name)
Starts Benchmark with the specified name.
Definition: TBenchmark.cxx:174
virtual void Show(const char *name)
Stops Benchmark name and Prints results.
Definition: TBenchmark.cxx:157
An array of clone (identical) objects.
Definition: TClonesArray.h:32
void BypassStreamer(Bool_t bypass=kTRUE)
When the kBypassStreamer bit is set, the automatically generated Streamer can call directly TClass::W...
virtual void Clear(Option_t *option="")
Clear the clones array.
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition: TFile.h:53
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2998
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:251
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:294
A simple line.
Definition: TLine.h:23
Double_t GetY1() const
Definition: TLine.h:53
Double_t GetX2() const
Definition: TLine.h:52
Double_t GetX1() const
Definition: TLine.h:51
Double_t GetY2() const
Definition: TLine.h:54
Int_t GetEntriesFast() const
Definition: TObjArray.h:64
TObject * At(Int_t idx) const
Definition: TObjArray.h:166
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
Definition: TRandom.cxx:263
virtual Double_t Rndm()
Machine independent random number generator.
Definition: TRandom.cxx:541
A TTree represents a columnar dataset.
Definition: TTree.h:78
TLine * line
double T(double x)
Definition: ChebyshevPol.h:34
Author
Rene Brun

Definition in file tcl.C.