Logo ROOT  
Reference Guide
htest.C File Reference

Detailed Description

View in nbviewer Open in SWAN Save histograms in Tree branches

To run this example, do

root > .L htest.C
root > htw()
root > htr1()
root > htr2()
root > htr3()
void htw() {
// Create a Tree with a few branches of type histogram
// 25000 entries are filled in the Tree
// For each entry, the copy of 3 histograms is written
// The data base will contain 75000 histograms.
gBenchmark->Start("hsimple");
TFile f("ht.root","recreate");
auto T = new TTree("T","test");
auto hpx = new TH1F("hpx","This is the px distribution",100,-4,4);
auto hpxpy = new TH2F("hpxpy","py vs px",40,-4,4,40,-4,4);
auto hprof = new TProfile("hprof","Profile of pz versus px",100,-4,4,0,20);
T->Branch("hpx","TH1F",&hpx,32000,0);
T->Branch("hpxpy","TH2F",&hpxpy,32000,0);
T->Branch("hprof","TProfile",&hprof,32000,0);
Float_t px, py, pz;
for (Int_t i = 0; i < 25000; i++) {
if (i%1000 == 0) printf("at entry: %d\n",i);
gRandom->Rannor(px,py);
pz = px*px + py*py;
hpx->Fill(px);
hpxpy->Fill(px,py);
hprof->Fill(px,pz);
T->Fill();
}
T->Print();
f.Write();
gBenchmark->Show("hsimple");
}
void htr1() {
// Connect Tree generated by htw and show histograms for entry 12345
auto f = new TFile("ht.root");
auto T = (TTree*)f->Get("T");
TH1F *hpx = nullptr;
TH2F *hpxpy = nullptr;
TProfile *hprof = nullptr;
T->SetBranchAddress("hpx",&hpx);
T->SetBranchAddress("hpxpy",&hpxpy);
T->SetBranchAddress("hprof",&hprof);
T->GetEntry(12345);
auto c1 = new TCanvas("c1","test",10,10,600,1000);
c1->Divide(1,3);
c1->cd(1);
hpx->Draw();
c1->cd(2);
hpxpy->Draw();
c1->cd(3);
hprof->Draw();
c1->Print("htr1.png");
}
void htr2() {
// Connect Tree generated by htw and show histograms for entry 12345
// a variant of htr1
auto f = new TFile("ht.root");
auto T = (TTree*)f->Get("T");
auto c1 = new TCanvas("c1","test",10,10,600,1000);
c1->Divide(1,3);
c1->cd(1);
T->Draw("hpx.Draw()","","goff",1,12345);
c1->cd(2);
T->Draw("hpxpy.Draw()","","goff",1,12345);
c1->cd(3);
T->Draw("hprof.Draw()","","goff",1,12345);
c1->Print("htr2.png");
}
void htr3() {
// Connect Tree generated by htw
// read all histograms and plot the RMS of hpx versus the Mean of hprof
// for each of the 25000 entries
auto f = new TFile("ht.root");
auto T = (TTree*)f->Get("T");
auto c1 = new TCanvas("c1","test",10,10,600,400);
T->Draw("hpx.GetRMS():hprof.GetMean()");
c1->Print("htr3.png");
}
void htest() {
htw();
htr1();
htr2();
htr3();
}
Author
Rene Brun

Definition in file htest.C.

f
#define f(i)
Definition: RSha256.hxx:122
TBenchmark::Start
virtual void Start(const char *name)
Starts Benchmark with the specified name.
Definition: TBenchmark.cxx:172
TH2F
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:251
TTree
Definition: TTree.h:79
Float_t
float Float_t
Definition: RtypesCore.h:57
Int_t
int Int_t
Definition: RtypesCore.h:45
TRandom::Rannor
virtual void Rannor(Float_t &a, Float_t &b)
Return 2 numbers distributed following a gaussian with mean=0 and sigma=1.
Definition: TRandom.cxx:489
gBenchmark
R__EXTERN TBenchmark * gBenchmark
Definition: TBenchmark.h:59
TProfile::Fill
Int_t Fill(const Double_t *v)
Definition: TProfile.h:54
TH1::Fill
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3274
gRandom
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
TFile
Definition: TFile.h:54
TH2::Fill
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:294
TProfile
Definition: TProfile.h:32
TCanvas
Definition: TCanvas.h:23
TH1F
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:572
ROOT::Math::Chebyshev::T
double T(double x)
Definition: ChebyshevPol.h:52
TBenchmark::Show
virtual void Show(const char *name)
Stops Benchmark name and Prints results.
Definition: TBenchmark.cxx:155
TH1::Draw
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2997
c1
return c1
Definition: legend1.C:41