11#include <benchmark/benchmark.h>
23#include <TRandomGen.h>
25#include <boost/histogram.hpp>
27namespace bh = boost::histogram;
31#define MAX_SIZE_COLLECTION (1 << 10)
32#define BENCHMARK_RANGE_COLLECTIONS Arg(1)->Arg(1 << 2)->Arg(1 << 4)->Arg(1 << 6)->Arg(1 << 8)->Arg(1 << 10)
34static void BM_mergingCollectionsTH1I(benchmark::State&
state)
38 const size_t bins = 62500;
40 for (
auto _ :
state) {
41 std::vector<std::unique_ptr<TCollection>> collections;
42 for (
size_t ci = 0; ci < numberOfCollections; ci++) {
43 std::unique_ptr<TCollection> collection = std::make_unique<TObjArray>();
44 collection->SetOwner(
true);
45 TF1* uni =
new TF1(
"uni",
"1", 0, 1000000);
48 h->FillRandom(
"uni", 50000);
51 collections.push_back(std::move(collection));
54 auto m = std::make_unique<TH1I>(
"merged",
"merged",
bins, 0, 1000000);
56 for (
size_t i = 0;
i <
bins;
i++) {
57 m->SetBinContent(
i, 1);
59 auto start = std::chrono::high_resolution_clock::now();
60 for (
const auto& collection : collections) {
61 m->Merge(collection.get(),
"-NOCHECK");
63 auto end = std::chrono::high_resolution_clock::now();
65 auto elapsed_seconds = std::chrono::duration_cast<std::chrono::duration<double>>(
end -
start);
66 state.SetIterationTime(elapsed_seconds.count());
70static void BM_mergingCollectionsTH2I(benchmark::State&
state)
75 const size_t bins = 250;
77 TF2* uni =
new TF2(
"uni",
"1", 0, 1000000, 0, 1000000);
79 for (
auto _ :
state) {
81 std::vector<std::unique_ptr<TCollection>> collections;
83 for (
size_t ci = 0; ci < numberOfCollections; ci++) {
84 std::unique_ptr<TCollection> collection = std::make_unique<TObjArray>();
85 collection->SetOwner(
true);
90 h->FillRandom(
"uni", 50000);
93 collections.push_back(std::move(collection));
95 auto m = std::make_unique<TH2I>(
"merged",
"merged",
bins, 0, 1000000,
bins, 0, 1000000);
97 for (
size_t i = 0;
i <
bins;
i++) {
98 m->SetBinContent(
i, 1);
101 auto start = std::chrono::high_resolution_clock::now();
102 for (
const auto& collection : collections) {
103 m->Merge(collection.get(),
"-NOCHECK");
105 auto end = std::chrono::high_resolution_clock::now();
107 auto elapsed_seconds = std::chrono::duration_cast<std::chrono::duration<double>>(
end -
start);
108 state.SetIterationTime(elapsed_seconds.count());
113static void BM_mergingCollectionsTH3I(benchmark::State&
state)
118 const size_t bins = 40;
120 TF3* uni =
new TF3(
"uni",
"1", 0, 1000000, 0, 1000000, 0, 1000000);
122 for (
auto _ :
state) {
124 std::vector<std::unique_ptr<TCollection>> collections;
126 for (
size_t ci = 0; ci < numberOfCollections; ci++) {
127 std::unique_ptr<TCollection> collection = std::make_unique<TObjArray>();
128 collection->SetOwner(
true);
134 h->FillRandom(
"uni", 50000);
137 collections.push_back(std::move(collection));
139 auto m = std::make_unique<TH3I>(
"merged",
"merged",
bins, 0, 1000000,
bins, 0, 1000000,
bins, 0, 1000000);
141 for (
size_t i = 0;
i <
bins;
i++) {
142 m->SetBinContent(
i, 1);
145 auto start = std::chrono::high_resolution_clock::now();
146 for (
const auto& collection : collections) {
147 m->Merge(collection.get(),
"-NOCHECK");
149 auto end = std::chrono::high_resolution_clock::now();
151 auto elapsed_seconds = std::chrono::duration_cast<std::chrono::duration<double>>(
end -
start);
152 state.SetIterationTime(elapsed_seconds.count());
158static void BM_mergingCollectionsTHNSparse(benchmark::State&
state)
164 const Double_t
min = 0.0;
165 const Double_t
max = 1000000.0;
166 const size_t dim = 10;
167 const Int_t
bins = 250;
173 gen.SetSeed(std::random_device()());
174 Double_t randomArray[dim];
176 for (
auto _ :
state) {
178 std::vector<std::unique_ptr<TCollection>> collections;
180 for (
size_t ci = 0; ci < numberOfCollections; ci++) {
181 std::unique_ptr<TCollection> collection = std::make_unique<TObjArray>();
182 collection->SetOwner(
true);
186 gen.RndmArray(dim, randomArray);
187 for (
double&
r : randomArray) {
190 h->Fill(randomArray);
194 collections.push_back(std::move(collection));
196 auto m = std::make_unique<THnSparseI>(
"merged",
"merged", dim, binsDims, mins, maxs);
198 auto start = std::chrono::high_resolution_clock::now();
199 for (
const auto& collection : collections) {
200 m->Merge(collection.get());
202 auto end = std::chrono::high_resolution_clock::now();
204 auto elapsed_seconds = std::chrono::duration_cast<std::chrono::duration<double>>(
end -
start);
205 state.SetIterationTime(elapsed_seconds.count());
209static void BM_mergingCollectionsTTree(benchmark::State&
state)
223 auto createTree = [&](std::string
name) -> TTree* {
224 TTree*
tree =
new TTree();
226 tree->Branch(
"b1", &branch1,
"a/I:b/L:c/F:d/D");
227 tree->Branch(
"b2", &branch2);
232 gen.SetSeed(std::random_device()());
233 Double_t randomArray[5];
235 for (
auto _ :
state) {
236 std::vector<std::unique_ptr<TCollection>> collections;
238 for (
size_t ci = 0; ci < numberOfCollections; ci++) {
239 std::unique_ptr<TCollection> collection = std::make_unique<TObjArray>();
240 collection->SetOwner(
true);
244 gen.RndmArray(5, randomArray);
245 branch1 = {
static_cast<Int_t
>(randomArray[0]),
static_cast<Long64_t
>(randomArray[1]),
static_cast<Float_t>(randomArray[2]), randomArray[3]};
246 branch2 = randomArray[4];
251 collections.emplace_back(std::move(collection));
253 TTree* merged = createTree(
"merged");
255 auto start = std::chrono::high_resolution_clock::now();
256 for (
const auto& collection : collections) {
257 merged->Merge(collection.get());
259 auto end = std::chrono::high_resolution_clock::now();
261 auto elapsed_seconds = std::chrono::duration_cast<std::chrono::duration<double>>(
end -
start);
262 state.SetIterationTime(elapsed_seconds.count());
268static void BM_mergingPODCollections(benchmark::State&
state)
273 const size_t bins = 62500;
275 using PODHistoType = std::vector<float>;
276 using PODHistoTypePtr = std::unique_ptr<PODHistoType>;
277 TF1* uni =
new TF1(
"uni",
"1", 0, 1000000);
279 const size_t randoms = 50000;
281 gen.SetSeed(std::random_device()());
282 Double_t randomArray[randoms];
284 for (
auto _ :
state) {
285 std::vector<std::unique_ptr<std::vector<PODHistoTypePtr>>> collections;
287 for (
size_t ci = 0; ci < numberOfCollections; ci++) {
288 auto collection = std::make_unique<std::vector<PODHistoTypePtr>>();
290 auto h = PODHistoTypePtr(
new PODHistoType(
bins, 0));
291 gen.RndmArray(randoms, randomArray);
292 for (
double r : randomArray) {
298 collection->emplace_back(std::move(
h));
300 collections.emplace_back(std::move(collection));
302 auto m = PODHistoTypePtr(
new PODHistoType(
bins, 0));
304 for (
size_t i = 0;
i <
bins;
i++) {
308 auto merge = [&](
const std::vector<PODHistoTypePtr>& collection,
size_t i) {
309 auto h = collection[
i].get();
310 for (
size_t b = 0;
b <
bins;
b++) {
315 auto start = std::chrono::high_resolution_clock::now();
316 for (
const auto& collection : collections) {
321 auto end = std::chrono::high_resolution_clock::now();
323 auto elapsed_seconds = std::chrono::duration_cast<std::chrono::duration<double>>(
end -
start);
324 state.SetIterationTime(elapsed_seconds.count());
329static void BM_mergingBoostRegular1DCollections(benchmark::State&
state)
331 const double min = 0.0;
332 const double max = 1000000.0;
336 const size_t bins = 62500;
338 auto merged = bh::make_histogram(bh::axis::regular<>(
bins,
min,
max,
"x"));
340 using HistoType =
decltype(merged);
341 using CollectionHistoType = std::vector<HistoType>;
343 TF1* uni =
new TF1(
"uni",
"1", 0, 1000000);
345 const size_t randoms = 50000;
347 gen.SetSeed(std::random_device()());
348 Double_t randomArray[randoms];
350 for (
auto _ :
state) {
352 std::vector<CollectionHistoType> collections;
354 for (
size_t ci = 0; ci < numberOfCollections; ci++) {
355 CollectionHistoType collection;
357 collection.emplace_back(std::move(bh::make_histogram(bh::axis::regular<>(
bins,
min,
max,
"x"))));
359 auto&
h = collection.back();
360 static_assert(std::is_reference<
decltype(
h)>
::value);
362 gen.RndmArray(randoms, randomArray);
363 for (
double r : randomArray) {
367 collections.emplace_back(std::move(collection));
370 auto merge = [&](
const CollectionHistoType& collection) {
372 merged += collection[
i];
376 auto start = std::chrono::high_resolution_clock::now();
377 for (
size_t ci = 0; ci < numberOfCollections; ci++) {
378 merge(collections[ci]);
380 auto end = std::chrono::high_resolution_clock::now();
382 auto elapsed_seconds = std::chrono::duration_cast<std::chrono::duration<double>>(
end -
start);
383 state.SetIterationTime(elapsed_seconds.count());
389static void BM_mergingBoostRegular2DCollections(benchmark::State&
state)
391 const double min = 0.0;
392 const double max = 1000000.0;
396 const size_t bins = 250;
398 auto merged = bh::make_histogram(bh::axis::regular<>(
bins,
min,
max,
"x"), bh::axis::regular<>(
bins,
min,
max,
"y"));
400 using HistoType =
decltype(merged);
401 using CollectionHistoType = std::vector<HistoType>;
403 TF1* uni =
new TF1(
"uni",
"1", 0, 1000000);
405 const size_t randoms = 50000;
407 gen.SetSeed(std::random_device()());
408 Double_t randomArrayX[randoms];
409 Double_t randomArrayY[randoms];
411 for (
auto _ :
state) {
414 std::vector<CollectionHistoType> collections;
416 for (
size_t ci = 0; ci < numberOfCollections; ci++) {
417 CollectionHistoType collection;
419 collection.emplace_back(std::move(bh::make_histogram(bh::axis::regular<>(
bins,
min,
max,
"x"), bh::axis::regular<>(
bins,
min,
max,
"y"))));
421 auto&
h = collection.back();
422 static_assert(std::is_reference<
decltype(
h)>
::value);
424 gen.RndmArray(randoms, randomArrayX);
425 gen.RndmArray(randoms, randomArrayY);
426 for (
size_t r = 0;
r < randoms;
r++) {
427 h(randomArrayX[
r] *
max, randomArrayY[
r] *
max);
430 collections.emplace_back(std::move(collection));
433 auto merge = [&](
const CollectionHistoType& collection) {
435 merged += collection[
i];
439 auto start = std::chrono::high_resolution_clock::now();
440 for (
size_t ci = 0; ci < numberOfCollections; ci++) {
441 merge(collections[ci]);
443 auto end = std::chrono::high_resolution_clock::now();
445 auto elapsed_seconds = std::chrono::duration_cast<std::chrono::duration<double>>(
end -
start);
446 state.SetIterationTime(elapsed_seconds.count());
453BENCHMARK(BM_mergingCollectionsTH1I)->Arg(1)->UseManualTime();
454BENCHMARK(BM_mergingCollectionsTH1I)->Arg(1)->UseManualTime();
455BENCHMARK(BM_mergingCollectionsTH2I)->Arg(1)->UseManualTime();
456BENCHMARK(BM_mergingCollectionsTH3I)->Arg(1)->UseManualTime();
457BENCHMARK(BM_mergingCollectionsTHNSparse)->Arg(1)->UseManualTime();
458BENCHMARK(BM_mergingPODCollections)->Arg(1)->UseManualTime();
459BENCHMARK(BM_mergingBoostRegular1DCollections)->Arg(1)->UseManualTime();
460BENCHMARK(BM_mergingBoostRegular2DCollections)->Arg(1)->UseManualTime();
461BENCHMARK(BM_mergingCollectionsTTree)->Arg(1)->UseManualTime();
465BENCHMARK(BM_mergingCollectionsTH1I)->BENCHMARK_RANGE_COLLECTIONS->UseManualTime();
466BENCHMARK(BM_mergingCollectionsTH2I)->BENCHMARK_RANGE_COLLECTIONS->UseManualTime();
467BENCHMARK(BM_mergingCollectionsTH3I)->BENCHMARK_RANGE_COLLECTIONS->UseManualTime();
468BENCHMARK(BM_mergingCollectionsTHNSparse)->BENCHMARK_RANGE_COLLECTIONS->UseManualTime();
469BENCHMARK(BM_mergingPODCollections)->BENCHMARK_RANGE_COLLECTIONS->UseManualTime();
470BENCHMARK(BM_mergingBoostRegular1DCollections)->BENCHMARK_RANGE_COLLECTIONS->UseManualTime();
471BENCHMARK(BM_mergingBoostRegular2DCollections)->BENCHMARK_RANGE_COLLECTIONS->UseManualTime();
472BENCHMARK(BM_mergingCollectionsTTree)->BENCHMARK_RANGE_COLLECTIONS->UseManualTime();
default_random_engine gen(dev())
const size_t collectionSize
BENCHMARK(BM_mergingCollectionsTH1I) -> Arg(1) ->UseManualTime()
#define MAX_SIZE_COLLECTION
Class for time synchronization of RawReader instances.
GLuint const GLchar * name
GLboolean GLboolean GLboolean b
GLsizei const GLfloat * value
GLboolean GLboolean GLboolean GLboolean a
std::string to_string(gsl::span< T, Size > span)
std::unique_ptr< TTree > tree((TTree *) flIn.Get(std::string(o2::base::NameConf::CTFTREENAME).c_str()))