Evaluate the performance of THnSparse vs TH1/2/3/nF for different numbers of dimensions and bins per dimension.
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.
tutorials/io/tree/tree143_drawsparse.C shows an example for visualizing a THnSparse. It creates a TTree which is then drawn using TParallelCoord.
This macro should be run in compiled mode due to the many nested loops that force CLING to disable its optimization. If run interpreted one would not benchmark THnSparse but CLING.
#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(nullptr), fDim(
dim), fBins(bins), fNum(num), fSparse(nullptr), fHist(nullptr), fHn(nullptr)
{
}
~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 SetupHist(EHist hist);
void NextValues();
void SetupValues();
private:
THnSparse *fSparse;
TH1 *fHist;
THn *fHn;
};
Int_t TTimeHists::fgDebug = 0;
TTimeHists::~TTimeHists()
{
delete[] fValue;
delete fSparse;
delete fHist;
delete fHn;
}
bool TTimeHists::Run()
{
for (
int h = 0;
h < 2; ++
h) {
SetupValues();
try {
w.Start();
w.Stop();
do {
check[
h] = Check((EHist)
h);
w.Stop();
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 {
w.Stop();
} while (w.RealTime() < 0.1);
fTime[
h][0] = (1. * fNum * rep[
h]) / w.RealTime() / 1E6;
fTime[
h][1] = (1. * fNum * rep[
h]) / w.CpuTime() / 1E6;
}
} 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()
{
}
void TTimeHists::SetupValues()
{
if (!fValue)
}
void TTimeHists::Fill(EHist hist)
{
NextValues();
if (fgDebug > 1) {
printf(
"%ld: fill %s",
n, hist == kHist ? (fDim < 4 ?
"hist" :
"arr") :
"sparse");
printf(
"[%g]", fValue[
d]);
printf("\n");
}
if (hist == kHist) {
switch (fDim) {
case 1: fHist->Fill(fValue[0]); break;
case 2: ((
TH2F *)fHist)->Fill(fValue[0], fValue[1]);
break;
case 3: ((
TH3F *)fHist)->Fill(fValue[0], fValue[1], fValue[2]);
break;
default: fHn->Fill(fValue); break;
}
} else {
fSparse->Fill(fValue);
}
}
}
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: {
throw std::bad_alloc();
}
throw std::bad_alloc();
}
}
}
} else {
}
}
}
{
memset(
x, 0,
sizeof(
Int_t) * fDim);
if (hist == kHist) {
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 (fgDebug > 2 || (fgDebug > 1 &&
v)) {
printf("%s%d", fDim < 4 ? "hist" : "arr", fDim);
}
}
}
++idx;
}
} else {
for (
Long64_t i = 0; i < fSparse->GetNbins(); ++i) {
if (fgDebug > 1) {
printf("sparse%d", fDim);
}
}
}
check /= fNum;
if (fgDebug > 0)
printf("check %s%d = %g\n", hist == kHist ? (fDim < 4 ? "hist" : "arr") : "sparse", fDim, check);
return check;
}
void hist103_THnSparse_hist()
{
#if defined(__CLING__)
printf("Please run this script in compiled mode by running \".x hist103_THnSparse_hist.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) {
else
if (t == 0)
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);
}
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);
printf("Processing dimension %d", dim);
for (
Int_t bins = 10; bins <= 100; bins += 10) {
TTimeHists timer(dim, bins, 1000);
timer.Run();
for (
int h = 0;
h < TTimeHists::kNumHist; ++
h) {
for (int t = 0; t < TTimeHists::kNumTime; ++t) {
Double_t time = timer.GetTime((TTimeHists::EHist)
h, (TTimeHists::ETime)t);
if (time >= 0.)
htime[
h][t]->
Fill(dim, bins, time);
}
}
hsparse_mem->
Fill(dim, bins, timer.GetSparse()->GetSparseFractionMem());
hsparse_bins->
Fill(dim, bins, timer.GetSparse()->GetSparseFractionBins());
if (max < timer.GetTime(TTimeHists::kSparse, TTimeHists::kReal))
max = timer.GetTime(TTimeHists::kSparse, TTimeHists::kReal);
printf(".");
fflush(stdout);
}
printf(" done\n");
}
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]);
}
gStyle->SetPalette(8,
nullptr);
gStyle->SetPaintTextFormat(
".2g");
const char *opt = "TEXT COL";
for (int t = 0; t < TTimeHists::kNumTime; ++t) {
for (
int h = 0;
h < TTimeHists::kNumHist; ++
h) {
}
htime_ratio[t]->
Draw(opt);
}
}
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
int Int_t
Signed integer 4 bytes (int).
long Long_t
Signed long integer 4 bytes (long). Size depends on architecture.
double Double_t
Double 8 bytes.
long long Long64_t
Portable signed long integer 8 bytes.
float Float_t
Float 4 bytes (float).
THnSparseT< TArrayF > THnSparseF
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
TVirtualPad * cd(Int_t subpadnumber=0) override
Set current canvas & pad.
A file, usually with extension .root, that stores data and code in the form of serialized objects in ...
1-D histogram with a float per channel (see TH1 documentation)
void SetTitle(const char *title) override
Change/set the title.
virtual void SetMaximum(Double_t maximum=-1111)
void Draw(Option_t *option="") override
Draw this histogram with options.
virtual void SetMinimum(Double_t minimum=-1111)
virtual Bool_t Divide(TF1 *f1, Double_t c1=1)
Performs the operation: this = this/(c1*f1) if errors are defined (see TH1::Sumw2),...
2-D histogram with a float per channel (see TH1 documentation)
Int_t Fill(Double_t) override
Invalid Fill method.
3-D histogram with a float per channel (see TH1 documentation)
virtual Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0) override
Automatic pad generation by division.
__device__ AFloat max(AFloat x, AFloat y)
void Fill(float *output, float value, int size)