ROOT logo

From $ROOTSYS/tutorials/hist/rebin.C

//this tutorial illustrates how to:
//  -create a variable binwidth histogram with a binning such
//   that the population per bin is about the same.
//  -rebin a variable binwidth histogram into another one.
//Author: Rene Brun
   
#include "TH1.h"
#include "TCanvas.h"
void rebin() {
   //create a fix bin histogram
   TH1F *h = new TH1F("h","test rebin",100,-3,3);
   Int_t nentries = 1000;
   h->FillRandom("gaus",nentries);
   Double_t xbins[1001];
   Int_t k=0;
   TAxis *axis = h->GetXaxis();
   for (Int_t i=1;i<=100;i++) {
      Int_t y = (Int_t)h->GetBinContent(i);
      if (y <=0) continue;
      Double_t dx = axis->GetBinWidth(i)/y;
      Double_t xmin = axis->GetBinLowEdge(i);
      for (Int_t j=0;j<y;j++) {
         xbins[k] = xmin +j*dx;
         k++;
      }
   }
   xbins[k] = axis->GetXmax();
   //create a variable binwidth histogram out of fix bin histogram
   //new rebinned histogram should have about 10 entries per bin
   TH1F *hnew = new TH1F("hnew","rebinned",k,xbins);
   hnew->FillRandom("gaus",10*nentries);
   
   //rebin hnew keeping only 50% of the bins
   Double_t xbins2[501];
   Int_t kk=0;
   for (Int_t j=0;j<k;j+=2) {
      xbins2[kk] = xbins[j];
      kk++;
   }
   xbins2[kk] = xbins[k];
   TH1F *hnew2 = (TH1F*)hnew->Rebin(kk,"hnew2",xbins2);

   //draw the 3 histograms
   TCanvas *c1 = new TCanvas("c1","c1",800,1000);
   c1->Divide(1,3);
   c1->cd(1);
   h->Draw();
   c1->cd(2);
   hnew->Draw();
   c1->cd(3);
   hnew2->Draw();
}
 rebin.C:1
 rebin.C:2
 rebin.C:3
 rebin.C:4
 rebin.C:5
 rebin.C:6
 rebin.C:7
 rebin.C:8
 rebin.C:9
 rebin.C:10
 rebin.C:11
 rebin.C:12
 rebin.C:13
 rebin.C:14
 rebin.C:15
 rebin.C:16
 rebin.C:17
 rebin.C:18
 rebin.C:19
 rebin.C:20
 rebin.C:21
 rebin.C:22
 rebin.C:23
 rebin.C:24
 rebin.C:25
 rebin.C:26
 rebin.C:27
 rebin.C:28
 rebin.C:29
 rebin.C:30
 rebin.C:31
 rebin.C:32
 rebin.C:33
 rebin.C:34
 rebin.C:35
 rebin.C:36
 rebin.C:37
 rebin.C:38
 rebin.C:39
 rebin.C:40
 rebin.C:41
 rebin.C:42
 rebin.C:43
 rebin.C:44
 rebin.C:45
 rebin.C:46
 rebin.C:47
 rebin.C:48
 rebin.C:49
 rebin.C:50
 rebin.C:51
 rebin.C:52
 rebin.C:53