htest.C

Go to the documentation of this file.
00001 //Example illustrating how to save histograms in Tree branches.
00002 //To run this example, do
00003 // root > .L htest.C
00004 // root > htw()
00005 // root > htr1()
00006 // root > htr2()
00007 // root > htr3()
00008 //
00009 //Author: Rene Brun
00010       
00011 void htw() {
00012    //create a Tree with a few branches of type histogram
00013    //25000 entries are filled in the Tree
00014    //For each entry, the copy of 3 histograms is written
00015    //The data base will contain 75000 histograms.
00016    gBenchmark->Start("hsimple");
00017    TFile f("ht.root","recreate");
00018    TTree *T     = new TTree("T","test");
00019    TH1F *hpx    = new TH1F("hpx","This is the px distribution",100,-4,4);
00020    TH2F *hpxpy  = new TH2F("hpxpy","py vs px",40,-4,4,40,-4,4);
00021    TProfile *hprof  = new TProfile("hprof","Profile of pz versus px",100,-4,4,0,20);
00022    T->Branch("hpx","TH1F",&hpx,32000,0);
00023    T->Branch("hpxpy","TH2F",&hpxpy,32000,0);
00024    T->Branch("hprof","TProfile",&hprof,32000,0);
00025    Float_t px, py, pz;
00026    for (Int_t i = 0; i < 25000; i++) {
00027       if (i%1000 == 0) printf("at entry: %d\n",i);
00028       gRandom->Rannor(px,py);
00029       pz = px*px + py*py;
00030       hpx->Fill(px);
00031       hpxpy->Fill(px,py);
00032       hprof->Fill(px,pz);
00033       T->Fill();
00034    }
00035    T->Print();
00036    f.Write();
00037    gBenchmark->Show("hsimple");
00038 }
00039 void htr1() {
00040    //connect Tree generated by htw and show histograms for entry 12345
00041    TFile *f = new TFile("ht.root");
00042    TTree *T = (TTree*)f->Get("T");
00043    TH1F *hpx = 0;
00044    TH2F *hpxpy = 0;
00045    TProfile *hprof = 0;
00046    T->SetBranchAddress("hpx",&hpx);
00047    T->SetBranchAddress("hpxpy",&hpxpy);
00048    T->SetBranchAddress("hprof",&hprof);
00049    T->GetEntry(12345);
00050    TCanvas *c1 = new TCanvas("c1","test",10,10,600,1000);
00051    c1->Divide(1,3);
00052    c1->cd(1);
00053    hpx->Draw();
00054    c1->cd(2);
00055    hpxpy->Draw();
00056    c1->cd(3);
00057    hprof->Draw();
00058 }
00059 void htr2() {
00060    //connect Tree generated by htw and show histograms for entry 12345
00061    // a variant of htr1
00062    TFile *f = new TFile("ht.root");
00063    TTree *T = (TTree*)f->Get("T");
00064    TCanvas *c1 = new TCanvas("c1","test",10,10,600,1000);
00065    c1->Divide(1,3);
00066    c1->cd(1);
00067    T->Draw("hpx.Draw()","","goff",1,12345);
00068    c1->cd(2);
00069    T->Draw("hpxpy.Draw()","","goff",1,12345);
00070    c1->cd(3);
00071    T->Draw("hprof.Draw()","","goff",1,12345);
00072 }
00073 void htr3() {
00074    //connect Tree generated by htw
00075    //read all histograms and plot the RMS of hpx versus the Mean of hprof
00076    //for each of the 25000 entries
00077    TFile *f = new TFile("ht.root");
00078    TTree *T = (TTree*)f->Get("T");
00079    T->Draw("hpx.GetRMS():hprof.GetMean()");
00080 }
00081 void htest() {
00082    htw();
00083    htr1();
00084    htr2();
00085    htr3();
00086 }

Generated on Tue Jul 5 15:45:11 2011 for ROOT_528-00b_version by  doxygen 1.5.1