This macro shows the usage of TMPIFile to simulate event reconstruction and merging them in parallel.
#include "TMPIFile.h"
#ifdef TMPI_SECOND_RUN
#include <chrono>
#include <sstream>
void test_tmpi()
{
Int_t events_per_rank = 6;
std::string treename = "test_tmpi";
std::string branchname = "event";
std::stringstream smpifname;
smpifname << "/tmp/merged_output_" << getpid() << ".root";
TMPIFile *newfile = new TMPIFile(smpifname.str().c_str(), "RECREATE", N_collectors);
if (newfile->GetMPIGlobalRank() == 0) {
Info(
"test_tmpi",
" running with parallel ranks: %d", newfile->GetMPIGlobalSize());
Info(
"test_tmpi",
" running with collecting ranks: %d", N_collectors);
Info(
"test_tmpi",
" running with working ranks: %d", (newfile->GetMPIGlobalSize() - N_collectors));
Info(
"test_tmpi",
" running with sync rate: %d", sync_rate);
Info(
"test_tmpi",
" running with events per rank: %d", events_per_rank);
Info(
"test_tmpi",
" running with sleep mean: %d", sleep_mean);
Info(
"test_tmpi",
" running with sleep sigma: %d", sleep_sigma);
}
if (newfile->IsCollector()) {
Info(
"Collector",
"[%d]\troot output filename = %s", newfile->GetMPIGlobalRank(), smpifname.str().c_str());
}
if (newfile->IsCollector()) {
newfile->RunCollector();
} else {
tree->SetAutoFlush(sync_rate);
tree->Branch(branchname.c_str(),
"JetEvent", &event, 8000, 2);
auto sync_start = std::chrono::high_resolution_clock::now();
for (int i = 0; i < events_per_rank; i++) {
auto start = std::chrono::high_resolution_clock::now();
event->Build(jetm, trackm, hitam, hitbm);
auto evt_built = std::chrono::high_resolution_clock::now();
double build_time = std::chrono::duration_cast<std::chrono::duration<double>>(evt_built - start).count();
Info(
"Rank",
"[%d] [%d]\tevt = %d;\tbuild_time = %f", newfile->GetMPIColor(), newfile->GetMPILocalRank(), i,
build_time);
auto adjusted_sleep = (
int)(sleep_mean - build_time);
std::this_thread::sleep_for(std::chrono::seconds(
int(
sleep)));
if ((i + 1) % sync_rate == 0) {
newfile->Sync();
auto end = std::chrono::high_resolution_clock::now();
double sync_time = std::chrono::duration_cast<std::chrono::duration<double>>(end - sync_start).count();
Info(
"Rank",
"[%d] [%d]\tevent collection time: %f", newfile->GetMPIColor(), newfile->GetMPILocalRank(),
sync_time);
sync_start = std::chrono::high_resolution_clock::now();
}
}
if (events_per_rank % sync_rate != 0) {
newfile->Sync();
}
}
Info(
"Rank",
"[%d] [%d]\tclosing file", newfile->GetMPIColor(), newfile->GetMPILocalRank());
newfile->Close();
if (newfile->GetMPILocalRank() == 0) {
TString filename = newfile->GetMPIFilename();
Info(
"Rank",
"[%d] [%d]\topening file: %s", newfile->GetMPIColor(), newfile->GetMPILocalRank(), filename.
Data());
Info(
"Rank",
"[%d] [%d]\tfile should have %d events and has %lld", newfile->GetMPIColor(),
newfile->GetMPILocalRank(), (newfile->GetMPILocalSize() - 1) * events_per_rank,
tree->GetEntries());
}
}
}
void testTMPIFile(
Bool_t secRun)
{
auto start = std::chrono::high_resolution_clock::now();
test_tmpi();
auto end = std::chrono::high_resolution_clock::now();
double time = std::chrono::duration_cast<std::chrono::duration<double>>(end - start).count();
std::string msg = "Total elapsed time: ";
msg += std::to_string(time);
Info(
"testTMPIFile",
"%s", msg.c_str());
Info(
"testTMPIFile",
"exiting");
}
#else
void testTMPIFile()
{
MPI_Initialized(&flag);
if (!flag) {
MPI_Init(NULL, NULL);
}
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
MPI_Comm_size(MPI_COMM_WORLD, &size);
if (rank == 0) {
gROOT->ProcessLine(
".L JetEvent.cxx+");
}
MPI_Barrier(MPI_COMM_WORLD);
gROOT->ProcessLine(
"#define TMPI_SECOND_RUN yes");
gROOT->ProcessLine(
"#include \"" __FILE__
"\"");
gROOT->ProcessLine(
"testTMPIFile(true)");
MPI_Finalized(&finalized);
if (!finalized) {
MPI_Finalize();
}
}
#endif
#define R__LOAD_LIBRARY(LIBRARY)
void Info(const char *location, const char *msgfmt,...)
Use this function for informational messages.
R__EXTERN C unsigned int sleep(unsigned int seconds)
R__EXTERN TRandom * gRandom
R__EXTERN TSystem * gSystem
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
virtual void SetSeed(ULong_t seed=0)
Set the random generator seed.
virtual UInt_t GetSeed() const
Get the random generator seed.
const char * Data() const
virtual Int_t Exec(const char *shellcmd)
Execute a command.
A TTree represents a columnar dataset.