Skip to content
Snippets Groups Projects
Commit ed13fb44 authored by Gerardo Ganis's avatar Gerardo Ganis Committed by Danilo Piparo
Browse files

tutorials/multicore: two new tutorials to illustrate autobin in power-of-2 mode

parent 3d9152ed
No related branches found
No related tags found
No related merge requests found
/// \file
/// \ingroup tutorial_multicore
/// Fill histograms in parallel with automatic binning.
/// Illustrates use of power-of-two autobin algorithm
///
/// \macro_code
///
/// \author Gerardo Ganis
/// \date November 2017
#include "TH1D.h"
#include "TH2D.h"
#include "TH3D.h"
#include "TList.h"
#include "TRandom3.h"
#include "TDirectory.h"
#include "TROOT.h"
#include "TCanvas.h"
#include "TString.h"
#include "TStyle.h"
#include "ROOT/TSeq.hxx"
#include "ROOT/TThreadedObject.hxx"
#include <thread>
#include <iostream>
// The number of workers
const UInt_t nWorkers = 8U;
// Reference boundaries
const Double_t xmiref = -1.;
const Double_t xmaref = 7.;
Int_t mt304_fillHistos(UInt_t nNumbers = 1001)
{
// The first, fundamental operation to be performed in order to make ROOT
// thread-aware.
ROOT::EnableThreadSafety();
// Histograms to be filled in parallel
ROOT::TThreadedObject<TH1D> h1d("h1d", "1D test histogram", 64, 0., -1.);
ROOT::TThreadedObject<TH1D> h1dr("h1dr", "1D test histogram w/ ref boundaries", 64, xmiref, xmaref);
// We define our work item
auto workItem = [&](UInt_t workerID) {
// One generator, file and ntuple per worker
TRandom3 workerRndm(workerID); // Change the seed
auto wh1d = h1d.Get();
wh1d->SetBit(TH1::kAutoBinPTwo);
auto wh1dr = h1dr.Get();
Double_t x;
for (UInt_t i = 0; i < nNumbers; ++i) {
x = workerRndm.Gaus(3.);
wh1d->Fill(x);
wh1dr->Fill(x);
}
};
// Create the collection which will hold the threads, our "pool"
std::vector<std::thread> workers;
// Fill the "pool" with workers
for (auto workerID : ROOT::TSeqI(nWorkers)) {
workers.emplace_back(workItem, workerID);
}
// Now join them
for (auto && worker : workers) worker.join();
// Merge
auto fh1d = h1d.Merge();
auto fh1dr = h1dr.Merge();
// Make the canvas
TCanvas *c = new TCanvas("c", "c", 800, 800);
c->Divide(1,2);
gStyle->SetOptStat(111110);
c->cd(1);
fh1d->DrawCopy();
c->cd(2);
fh1dr->DrawCopy();
c->Update();
gROOTMutex = 0;
return 0;
}
/// \file
/// \ingroup tutorial_multicore
/// Fill multiple histograms with different functions and automatic binning.
/// Illustrates merging with the power-of-two autobin algorithm
///
/// \macro_code
///
/// \author Gerardo Ganis
/// \date November 2017
#include "TF1.h"
#include "TH1D.h"
#include "TMath.h"
#include "TF1.h"
#include "TCanvas.h"
#include "TRandom3.h"
#include "TStatistic.h"
#include "TFile.h"
#include "TStyle.h"
TF1 *gam = new TF1("gam", "1/(1+0.1*x*0.1*x)", -100., 100.);
TF1 *gam1 = new TF1("gam", "1/(1+0.1*x*0.1*x)", -1., .25);
TF1 *iga = new TF1("inv gam", "1.-1/(1+0.1*x*0.1*x)", -100., 100.);
TF1 *iga1 = new TF1("inv gam", "1.-1/(1+0.1*x*0.1*x)", -.5, 1.);
void mt305_fillHistosAutobin(unsigned opt = 1, unsigned n = 1001)
{
UInt_t nh = 10;
UInt_t bsize = 1000;
TRandom3 rndm((Long64_t)time(0));
// Standard autobinning reference
auto href = new TH1D("myhref", "current", 50, 0., -1.);
href->SetBuffer(bsize);
// New autobinning 1-histo reference
auto href2 = new TH1D("myhref", "Auto P2, sequential", 50, 0., -1.);
href2->SetBit(TH1::kAutoBinPTwo);
href2->SetBuffer(bsize);
TList *hlist = new TList;
Int_t nbins = 50;
TStatistic x("min"), y("max"), d("dif"), a("mean"), r("rms");
for (UInt_t j = 0; j < nh; ++j) {
Double_t xmi = 1e15, xma = -1e15;
TStatistic xw("work");
TString hname = TString::Format("myh%d", j);
auto hw = new TH1D(hname.Data(), "Auto P2, merged", nbins, 0., -1.);
hw->SetBit(TH1::kAutoBinPTwo);
hw->SetBuffer(bsize);
Double_t xhma, xhmi, ovf, unf;
Bool_t emptied = kFALSE, tofill = kTRUE;
Bool_t buffering = kTRUE;
for (UInt_t i = 0; i < n; ++i) {
Double_t xx;
switch (opt) {
case 1: xx = rndm.Gaus(3, 1); break;
case 2: xx = rndm.Rndm() * 100. - 50.; break;
case 3: xx = gam->GetRandom(); break;
case 4: xx = gam1->GetRandom(); break;
case 5: xx = iga->GetRandom(); break;
case 6: xx = iga1->GetRandom(); break;
default: xx = rndm.Gaus(0, 1);
}
if (buffering) {
if (xx > xma)
xma = xx;
if (xx < xmi)
xmi = xx;
xw.Fill(xx);
}
hw->Fill(xx);
href->Fill(xx);
href2->Fill(xx);
if (!hw->GetBuffer()) {
// Not buffering anymore
buffering = kFALSE;
}
}
x.Fill(xmi);
y.Fill(xma);
d.Fill(xma - xmi);
a.Fill(xw.GetMean());
r.Fill(xw.GetRMS());
hlist->Add(hw);
}
x.Print();
y.Print();
d.Print();
a.Print();
r.Print();
TH1D *h0 = (TH1D *)hlist->First();
hlist->Remove(h0);
if (!h0->Merge(hlist))
return;
gStyle->SetOptStat(111110);
if (gROOT->GetListOfCanvases()->FindObject("c3"))
delete gROOT->GetListOfCanvases()->FindObject("c3");
TCanvas *c3 = new TCanvas("c3", "c3", 800, 800);
c3->Divide(1, 3);
c3->cd(1);
h0->StatOverflows();
h0->DrawClone("HIST");
c3->cd(2);
href2->StatOverflows();
href2->DrawClone();
c3->cd(3);
href->StatOverflows();
href->DrawClone();
c3->Update();
std::cout << " ent: " << h0->GetEntries() << "\n";
h0->Print();
href->Print();
hlist->SetOwner(kTRUE);
delete hlist;
delete href;
delete href2;
delete h0;
}
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment