swh:1:snp:af87cd67498ef4fe47c76ed3e7caffe5b61facaf
Raw File
Tip revision: 59fe6578eda8d4c55e125a80e956ef4431ab5322 authored by Pere Mato on 04 March 2016, 14:47:28 UTC
Update ROOT version files to v6.07/04.
Tip revision: 59fe657
cumulative.C
/// \file
/// \ingroup tutorial_hist
/// Illustrate use of the TH1::GetCumulative method.
///
/// \macro_image
/// \macro_code
///
/// \author M. Schiller

#include <cassert>
#include <cmath>

#include "TH1.h"
#include "TH1D.h"
#include "TCanvas.h"
#include "TRandom.h"

TCanvas* cumulative()
{
   TH1* h = new TH1D("h", "h", 100, -5., 5.);
   gRandom->SetSeed();
   h->FillRandom("gaus", 1u << 16);
   // get the cumulative of h
   TH1* hc = h->GetCumulative();
   // check that c has the "right" contents
   Double_t* integral = h->GetIntegral();
   for (Int_t i = 1; i <= hc->GetNbinsX(); ++i) {
      assert(std::abs(integral[i] * h->GetEntries() - hc->GetBinContent(i)) < 1e-7);
   }
   // draw histogram together with its cumulative distribution
   TCanvas* c = new TCanvas;
   c->Divide(1,2);
   c->cd(1);
   h->Draw();
   c->cd(2);
   hc->Draw();
   c->Update();

   return c;
}
back to top