The script calculates the bandwidth for filling and retrieving bin contents (in million entries per second) for these two histogramming techniques, where "seconds" is CPU and real time.
The first line of the plots contains the bandwidth based on the CPU time (THnSpase, TH1/2/3/nF*, ratio), the second line shows the plots for real time, and the third line shows the fraction of filled bins and memory used by THnSparse vs. TH1/2/3/nF.
The timing depends on the distribution and the amount of entries in the histograms; here, a Gaussian distribution (center is contained in the histograms) is used to fill each histogram with 1000 entries. The filling and reading is repeated until enough statistics have been collected.
This macro should be run in compiled mode due to the many nested loops that force CINT to disable its optimization. If run interpreted one would not benchmark THnSparse but CINT.
#ifndef INT_MAX
#define INT_MAX std::numeric_limits<int>::max()
#endif
class TTimeHists {
public:
enum EHist { kHist, kSparse, kNumHist };
enum ETime {
kReal, kCPU, kNumTime };
fValue(0), fDim(dim), fBins(bins), fNum(num),
fSparse(0), fHist(0), fHn(0) {}
~TTimeHists();
bool Run();
Double_t GetTime(EHist hist, ETime time)
const {
if (time ==
kReal)
return fTime[hist][0];
return fTime[hist][1]; }
static void SetDebug(
Int_t lvl) { fgDebug = lvl; }
THnSparse* GetSparse()
const {
return fSparse; }
protected:
void Fill(EHist hist);
void SetupHist(EHist hist);
void NextValues();
void SetupValues();
private:
};
Int_t TTimeHists::fgDebug = 0;
TTimeHists::~TTimeHists()
{
delete fSparse;
delete fHist;
delete fHn;
}
bool TTimeHists::Run()
{
for (
int h = 0; h < 2; ++
h) {
SetupValues();
try {
SetupHist((EHist) h);
do {
Fill((EHist) h);
check[
h] = Check((EHist) h);
|| (h && rep[0] > 0 && rep[1] < rep[0]));
fTime[
h][0] = (1.* fNum * rep[
h]) / w.
RealTime() / 1E6;
fTime[
h][1] = (1.* fNum * rep[
h]) / w.
CpuTime() / 1E6;
if (h == 1 && (fTime[h][0] > 1E20 || fTime[h][1] > 1E20)) {
do {
Fill((EHist) h);
Check((EHist) h);
fTime[
h][0] = (1.* fNum * rep[
h]) / w.
RealTime() / 1E6;
fTime[
h][1] = (1.* fNum * rep[
h]) / w.
CpuTime() / 1E6;
}
if (fTime[h][0] > 1E20) fTime[
h][0] = 1E20;
if (fTime[h][1] > 1E20) fTime[
h][1] = 1E20;
}
catch (std::exception&) {
fTime[
h][0] = fTime[
h][1] = -1.;
}
}
if (check[0] != check[1])
if (check[0] != -1.)
printf(
"ERROR: mismatch of histogram (%g) and sparse histogram (%g) for dim=%d, bins=%d!\n",
check[0], check[1], fDim, fBins);
return (check[0] == check[1]);
}
void TTimeHists::NextValues()
{
for (
Int_t d = 0; d < fDim; ++
d)
}
void TTimeHists::SetupValues()
{
}
void TTimeHists::Fill(EHist hist)
{
for (
Long_t n = 0; n < fNum; ++
n) {
NextValues();
if (fgDebug > 1) {
printf(
"%ld: fill %s", n, hist == kHist? (fDim < 4 ?
"hist" :
"arr") :
"sparse");
for (
Int_t d = 0; d < fDim; ++
d)
}
if (hist == kHist) {
switch (fDim) {
case 1: fHist->Fill(
fValue[0]);
break;
default: fHn->Fill(
fValue);
break;
}
} else {
}
}
}
void TTimeHists::SetupHist(EHist hist)
{
if (hist == kHist) {
switch (fDim) {
case 1: fHist =
new TH1F(
"h1",
"h1", fBins, -1., 1.);
break;
case 2: fHist =
new TH2F(
"h2",
"h2", fBins, -1., 1., fBins, -1., 1.);
break;
case 3: fHist =
new TH3F(
"h3",
"h3", fBins, -1., 1., fBins, -1., 1., fBins, -1., 1.);
break;
default:
{
for (
Int_t d = 0; d < fDim; ++
d) {
throw std::bad_alloc();
size *= (fBins + 2);
}
throw std::bad_alloc();
for (
Int_t d = 0; d < fDim; ++
d) {
}
fHn =
new THnF(
"hn",
"hn", fDim, bins, xmin, xmax);
}
}
} else {
for (
Int_t d = 0; d < fDim; ++
d) {
}
fSparse =
new THnSparseF(
"hs",
"hs", fDim, bins, xmin, xmax);
}
}
{
memset(x, 0,
sizeof(
Int_t) * fDim);
if (hist == kHist) {
for (
Int_t d = 0; d < fDim; ++
d)
size *= (fBins + 2);
while (x[0] <= fBins + 1) {
if (fDim < 4) {
if (fDim == 2) histidx = fHist->GetBin(x[0], x[1]);
else if (fDim == 3) histidx = fHist->GetBin(x[0], x[1], x[2]);
v = fHist->GetBinContent(histidx);
}
else v = fHn->GetBinContent(x);
if (v)
for (
Int_t d = 0; d < fDim; ++
d)
checkx += x[d];
if (fgDebug > 2 || (fgDebug > 1 && v)) {
printf(
"%s%d", fDim < 4 ?
"hist" :
"arr", fDim);
for (
Int_t d = 0; d < fDim; ++
d)
}
++x[fDim - 1];
for (
Int_t d = fDim - 1; d > 0; --
d) {
if (x[d] > fBins + 1) {
++x[d - 1];
}
}
++idx;
}
} else {
for (
Long64_t i = 0; i < fSparse->GetNbins(); ++i) {
Double_t v = fSparse->GetBinContent(i, x);
for (
Int_t d = 0; d < fDim; ++
d)
checkx += x[d];
if (fgDebug > 1) {
for (
Int_t d = 0; d < fDim; ++
d)
}
}
}
check /= fNum;
if (fgDebug > 0)
printf(
"check %s%d = %g\n", hist == kHist ? (fDim < 4 ?
"hist" :
"arr") :
"sparse", fDim, check);
return check;
}
void sparsehist() {
#if defined (__CLING__)
printf(
"Please run this script in compiled mode by running \".x sparsehist.C+\"\n");
return;
#endif
TH2F* htime[TTimeHists::kNumHist][TTimeHists::kNumTime];
for (
int h = 0; h < TTimeHists::kNumHist; ++
h)
for (
int t = 0; t < TTimeHists::kNumTime; ++
t) {
if (h == 0) name += "arr";
else name += "sp";
if (t == 0) name += "_r";
title.
Form(
"Throughput (fill,get) %s (%s, 1M entries/sec);dim;bins;1M entries/sec", h == 0 ?
"TH1/2/3/nF" :
"THnSparseF", t == 0 ?
"real" :
"CPU");
htime[
h][
t] =
new TH2F(name, title, 6, 0.5, 6.5, 10, 5, 105);
}
TH2F* hsparse_mem =
new TH2F(
"hsparse_mem",
"Fractional memory usage;dim;bins;fraction of memory used", 6, 0.5, 6.5, 10, 5, 105);
TH2F* hsparse_bins =
new TH2F(
"hsparse_bins",
"Fractional number of used bins;dim;bins;fraction of filled bins", 6, 0.5, 6.5, 10, 5, 105);
for (
Int_t dim = 1; dim < 7; ++dim) {
printf(
"Processing dimension %d", dim);
timer.Run();
for (
int h = 0; h < TTimeHists::kNumHist; ++
h) {
for (
int t = 0; t < TTimeHists::kNumTime; ++
t) {
if (time >= 0.)
}
}
hsparse_mem->
Fill(dim,
bins, timer.GetSparse()->GetSparseFractionMem());
hsparse_bins->
Fill(dim,
bins, timer.GetSparse()->GetSparseFractionBins());
fflush(stdout);
}
}
TH2F* htime_ratio[TTimeHists::kNumTime];
for (
int t = 0; t < TTimeHists::kNumTime; ++
t) {
const char* name = t ? "htime_ratio" : "htime_ratio_r";
htime_ratio[
t] = (
TH2F*) htime[TTimeHists::kSparse][t]->Clone(name);
title.
Form(
"Relative speed improvement (%s, 1M entries/sec): sparse/hist;dim;bins;#Delta 1M entries/sec", t == 0 ?
"real" :
"CPU");
htime_ratio[
t]->
Divide(htime[TTimeHists::kHist][t]);
}
TFile* f =
new TFile(
"sparsehist.root",
"RECREATE");
const char* opt = "TEXT COL";
for (
int t = 0; t < TTimeHists::kNumTime; ++
t) {
for (
int h = 0; h < TTimeHists::kNumHist; ++
h) {
}
}
canv->
cd(7); hsparse_mem->
Draw(opt);
canv->
cd(8); hsparse_bins->
Draw(opt);
}