24#include "TDatabasePDG.h"
27Ort::Env
global_env(ORT_LOGGING_LEVEL_WARNING,
"GlobalEnv");
32void Scaler::load(
const std::string&
filename)
35 if (!
file.is_open()) {
36 throw std::runtime_error(
"Error: Could not open scaler file!");
39 std::string json_str((std::istreambuf_iterator<char>(
file)), std::istreambuf_iterator<char>());
42 rapidjson::Document doc;
43 doc.Parse(json_str.c_str());
45 if (doc.HasParseError()) {
46 throw std::runtime_error(
"Error: JSON parsing failed!");
49 normal_min = jsonArrayToVector(doc[
"normal"][
"min"]);
50 normal_max = jsonArrayToVector(doc[
"normal"][
"max"]);
51 outlier_center = jsonArrayToVector(doc[
"outlier"][
"center"]);
52 outlier_scale = jsonArrayToVector(doc[
"outlier"][
"scale"]);
55std::vector<double> Scaler::inverse_transform(
const std::vector<double>& input)
57 std::vector<double>
output;
58 for (
int i = 0;
i < input.size(); ++
i) {
59 if (
i < input.size() - 2) {
60 output.push_back(input[
i] * (normal_max[
i] - normal_min[
i]) + normal_min[
i]);
62 output.push_back(input[
i] * outlier_scale[
i - (input.size() - 2)] + outlier_center[
i - (input.size() - 2)]);
69std::vector<double> Scaler::jsonArrayToVector(
const rapidjson::Value& jsonArray)
71 std::vector<double>
vec;
72 for (
int i = 0;
i < jsonArray.Size(); ++
i) {
73 vec.push_back(jsonArray[
i].GetDouble());
80ONNXGenerator::ONNXGenerator(Ort::Env& shared_env,
const std::string& model_path)
81 : env(shared_env), session(env, model_path.c_str(),
Ort::SessionOptions{})
84 Ort::SessionOptions session_options;
85 session = Ort::Session(env, model_path.c_str(), session_options);
88std::vector<double> ONNXGenerator::generate_sample()
90 Ort::AllocatorWithDefaultOptions allocator;
93 std::vector<float>
z(100);
95 v = rand_gen.Gaus(0.0, 1.0);
99 std::vector<int64_t> input_shape = {1, 100};
101 Ort::MemoryInfo memory_info = Ort::MemoryInfo::CreateCpu(OrtArenaAllocator, OrtMemTypeDefault);
104 Ort::Value input_tensor = Ort::Value::CreateTensor<float>(
105 memory_info,
z.data(),
z.size(), input_shape.data(), input_shape.size());
107 const char* input_names[] = {
"z"};
108 const char* output_names[] = {
"output"};
109 auto output_tensors = session.Run(Ort::RunOptions{
nullptr}, input_names, &input_tensor, 1, output_names, 1);
112 float* output_data = output_tensors.front().GetTensorMutableData<
float>();
114 auto output_tensor_info = output_tensors.front().GetTensorTypeAndShapeInfo();
115 size_t output_data_size = output_tensor_info.GetElementCount();
116 std::vector<double>
output;
117 for (
int i = 0;
i < output_data_size; ++
i) {
118 output.push_back(output_data[
i]);
129GenTPCLoopers::GenTPCLoopers(std::string model_pairs, std::string model_compton,
130 std::string poisson, std::string gauss, std::string scaler_pair,
131 std::string scaler_compton)
134 std::ifstream model_file[2];
135 model_file[0].open(model_pairs);
136 model_file[1].open(model_compton);
137 if (!model_file[0].is_open() || model_file[0].peek() == std::ifstream::traits_type::eof()) {
138 LOG(fatal) <<
"Error: Pairs model file is empty or does not exist!";
141 if (!model_file[1].is_open() || model_file[1].peek() == std::ifstream::traits_type::eof()) {
142 LOG(fatal) <<
"Error: Compton model file is empty or does not exist!";
145 model_file[0].close();
146 model_file[1].close();
148 std::ifstream scaler_file[2];
149 scaler_file[0].open(scaler_pair);
150 scaler_file[1].open(scaler_compton);
151 if (!scaler_file[0].is_open() || scaler_file[0].peek() == std::ifstream::traits_type::eof()) {
152 LOG(fatal) <<
"Error: Pairs scaler file is empty or does not exist!";
155 if (!scaler_file[1].is_open() || scaler_file[1].peek() == std::ifstream::traits_type::eof()) {
156 LOG(fatal) <<
"Error: Compton scaler file is empty or does not exist!";
159 scaler_file[0].close();
160 scaler_file[1].close();
162 if (poisson !=
"" && poisson !=
"None" && poisson !=
"none") {
163 std::ifstream poisson_file(poisson);
164 if (!poisson_file.is_open() || poisson_file.peek() == std::ifstream::traits_type::eof()) {
165 LOG(fatal) <<
"Error: Poisson file is empty or does not exist!";
168 poisson_file >> mPoisson[0] >> mPoisson[1] >> mPoisson[2];
169 poisson_file.close();
174 if (gauss !=
"" && gauss !=
"None" && gauss !=
"none") {
175 std::ifstream gauss_file(gauss);
176 if (!gauss_file.is_open() || gauss_file.peek() == std::ifstream::traits_type::eof()) {
177 LOG(fatal) <<
"Error: Gauss file is empty or does not exist!";
180 gauss_file >> mGauss[0] >> mGauss[1] >> mGauss[2] >> mGauss[3];
185 mONNX_pair = std::make_unique<ONNXGenerator>(
global_env, model_pairs);
186 mScaler_pair = std::make_unique<Scaler>();
187 mScaler_pair->load(scaler_pair);
188 mONNX_compton = std::make_unique<ONNXGenerator>(
global_env, model_compton);
189 mScaler_compton = std::make_unique<Scaler>();
190 mScaler_compton->load(scaler_compton);
193Bool_t GenTPCLoopers::generateEvent()
198 mGenElectrons.clear();
200 unsigned int nLoopers, nLoopersPairs, nLoopersCompton;
201 LOG(
debug) <<
"mCurrentEvent is " << mCurrentEvent;
202 LOG(
debug) <<
"Current event time: " << ((mCurrentEvent < mInteractionTimeRecords.size() - 1) ?
std::to_string(mInteractionTimeRecords[mCurrentEvent + 1].bc2ns() - mInteractionTimeRecords[mCurrentEvent].bc2ns()) :
std::
to_string(mTimeEnd - mInteractionTimeRecords[mCurrentEvent].bc2ns())) <<
" ns";
203 LOG(
debug) <<
"Current time offset wrt BC: " << mInteractionTimeRecords[mCurrentEvent].getTimeOffsetWrtBC() <<
" ns";
204 mTimeLimit = (mCurrentEvent < mInteractionTimeRecords.size() - 1) ? mInteractionTimeRecords[mCurrentEvent + 1].bc2ns() - mInteractionTimeRecords[mCurrentEvent].bc2ns() : mTimeEnd - mInteractionTimeRecords[mCurrentEvent].bc2ns();
208 nLoopersPairs =
static_cast<unsigned int>(std::round(nLoopers * mLoopsFractionPairs));
209 nLoopersCompton = nLoopers - nLoopersPairs;
210 SetNLoopers(nLoopersPairs, nLoopersCompton);
211 LOG(info) <<
"Flat gas loopers: " << nLoopers <<
" (pairs: " << nLoopersPairs <<
", compton: " << nLoopersCompton <<
")";
212 generateEvent(mTimeLimit);
217 mNLoopersPairs =
static_cast<unsigned int>(std::round(mMultiplier[0] * PoissonPairs()));
218 LOG(
debug) <<
"Generated loopers pairs (Poisson): " << mNLoopersPairs;
221 mNLoopersCompton =
static_cast<unsigned int>(std::round(mMultiplier[1] * GaussianElectrons()));
222 LOG(
debug) <<
"Generated compton electrons (Gauss): " << mNLoopersCompton;
225 for (
int i = 0;
i < mNLoopersPairs; ++
i) {
226 std::vector<double> pair = mONNX_pair->generate_sample();
228 std::vector<double> transformed_pair = mScaler_pair->inverse_transform(pair);
229 mGenPairs.push_back(transformed_pair);
232 for (
int i = 0;
i < mNLoopersCompton; ++
i) {
233 std::vector<double> electron = mONNX_compton->generate_sample();
235 std::vector<double> transformed_electron = mScaler_compton->inverse_transform(electron);
236 mGenElectrons.push_back(transformed_electron);
242Bool_t GenTPCLoopers::generateEvent(
double time_limit)
244 LOG(info) <<
"Time constraint for loopers: " << time_limit <<
" ns";
246 for (
int i = 0;
i < mNLoopersPairs; ++
i) {
247 std::vector<double> pair = mONNX_pair->generate_sample();
249 std::vector<double> transformed_pair = mScaler_pair->inverse_transform(pair);
250 transformed_pair[9] = gRandom->Uniform(0., time_limit);
251 mGenPairs.push_back(transformed_pair);
254 for (
int i = 0;
i < mNLoopersCompton; ++
i) {
255 std::vector<double> electron = mONNX_compton->generate_sample();
257 std::vector<double> transformed_electron = mScaler_compton->inverse_transform(electron);
258 transformed_electron[6] = gRandom->Uniform(0., time_limit);
259 mGenElectrons.push_back(transformed_electron);
261 LOG(info) <<
"Generated Particles with time limit";
265std::vector<TParticle> GenTPCLoopers::importParticles()
267 std::vector<TParticle> particles;
268 const double mass_e = TDatabasePDG::Instance()->GetParticle(11)->Mass();
269 const double mass_p = TDatabasePDG::Instance()->GetParticle(-11)->Mass();
271 for (
auto& pair : mGenPairs) {
272 double px_e, py_e, pz_e, px_p, py_p, pz_p;
273 double vx, vy, vz,
time;
274 double e_etot, p_etot;
285 e_etot = TMath::Sqrt(px_e * px_e + py_e * py_e + pz_e * pz_e + mass_e * mass_e);
286 p_etot = TMath::Sqrt(px_p * px_p + py_p * py_p + pz_p * pz_p + mass_p * mass_p);
288 TParticle electron(11, 1, -1, -1, -1, -1, px_e, py_e, pz_e, e_etot, vx, vy, vz,
time / 1e9);
292 particles.push_back(electron);
294 TParticle positron(-11, 1, -1, -1, -1, -1, px_p, py_p, pz_p, p_etot, vx, vy, vz,
time / 1e9);
298 particles.push_back(positron);
301 for (
auto& compton : mGenElectrons) {
303 double vx, vy, vz,
time;
312 etot = TMath::Sqrt(px * px + py * py + pz * pz + mass_e * mass_e);
314 TParticle electron(11, 1, -1, -1, -1, -1, px, py, pz, etot, vx, vy, vz,
time / 1e9);
318 particles.push_back(electron);
324unsigned int GenTPCLoopers::PoissonPairs()
326 unsigned int poissonValue;
329 poissonValue = mRandGen.Poisson(mPoisson[0]);
330 }
while (poissonValue < mPoisson[1] || poissonValue > mPoisson[2]);
335unsigned int GenTPCLoopers::GaussianElectrons()
337 unsigned int gaussValue;
340 gaussValue = mRandGen.Gaus(mGauss[0], mGauss[1]);
341 }
while (gaussValue < mGauss[2] || gaussValue > mGauss[3]);
346void GenTPCLoopers::SetNLoopers(
unsigned int nsig_pair,
unsigned int nsig_compton)
349 mNLoopersPairs = nsig_pair;
350 mNLoopersCompton = nsig_compton;
353 LOG(info) <<
"Poissonian parameters correctly loaded.";
355 mNLoopersPairs = nsig_pair;
358 LOG(info) <<
"Gaussian parameters correctly loaded.";
360 mNLoopersCompton = nsig_compton;
365void GenTPCLoopers::SetMultiplier(
const std::array<float, 2>& mult)
369 if (mult[0] < 0 || mult[1] < 0) {
370 LOG(fatal) <<
"Error: Multiplier values must be non-negative!";
373 LOG(info) <<
"Multiplier values set to: Pair = " << mult[0] <<
", Compton = " << mult[1];
374 mMultiplier[0] = mult[0];
375 mMultiplier[1] = mult[1];
379void GenTPCLoopers::setFlatGas(Bool_t flat, Int_t number, Int_t nloopers_orbit)
383 if (nloopers_orbit > 0) {
384 mFlatGasOrbit =
true;
385 mFlatGasNumber = nloopers_orbit;
386 LOG(info) <<
"Flat gas loopers will be generated using orbit reference.";
388 mFlatGasOrbit =
false;
390 LOG(warn) <<
"Warning: Number of loopers per event must be non-negative! Switching option off.";
394 mFlatGasNumber = number;
398 mContextFile = std::filesystem::exists(
"collisioncontext.root") ? TFile::Open(
"collisioncontext.root") : nullptr;
400 mInteractionTimeRecords = mCollisionContext ? mCollisionContext->getEventRecords() : std::vector<o2::InteractionTimeRecord>{};
401 if (mInteractionTimeRecords.empty()) {
402 LOG(error) <<
"Error: No interaction time records found in the collision context!";
405 LOG(info) <<
"Interaction Time records has " << mInteractionTimeRecords.size() <<
" entries.";
406 mCollisionContext->printCollisionSummary();
408 for (
int c = 0;
c < mInteractionTimeRecords.size() - 1;
c++) {
409 mIntTimeRecMean += mInteractionTimeRecords[
c + 1].bc2ns() - mInteractionTimeRecords[
c].bc2ns();
411 mIntTimeRecMean /= (mInteractionTimeRecords.size() - 1);
414 const auto& lastIR = mInteractionTimeRecords.back();
416 mTimeEnd = finalOrbitIR.bc2ns();
417 LOG(
debug) <<
"Final orbit start time: " << mTimeEnd <<
" ns while last interaction record time is " << mInteractionTimeRecords.back().bc2ns() <<
" ns";
422 LOG(info) <<
"Flat gas loopers: " << (mFlatGas ?
"ON" :
"OFF") <<
", Reference loopers number per " << (mFlatGasOrbit ?
"orbit " :
"event ") << mFlatGasNumber;
425void GenTPCLoopers::setFractionPairs(
float fractionPairs)
427 if (fractionPairs < 0 || fractionPairs > 1) {
428 LOG(fatal) <<
"Error: Loops fraction for pairs must be in the range [0, 1].";
431 mLoopsFractionPairs = fractionPairs;
432 LOG(info) <<
"Pairs fraction set to: " << mLoopsFractionPairs;
435void GenTPCLoopers::SetRate(
const std::string& rateFile,
bool isPbPb =
true,
int intRate)
438 TFile rate_file(rateFile.c_str(),
"READ");
439 if (!rate_file.IsOpen() || rate_file.IsZombie()) {
440 LOG(fatal) <<
"Error: Rate file is empty or does not exist!";
443 const char* fitName = isPbPb ?
"fitPbPb" :
"fitpp";
444 auto fit = (
TF1*)rate_file.Get(fitName);
446 LOG(fatal) <<
"Error: Could not find fit function '" << fitName <<
"' in rate file!";
449 mInteractionRate = intRate;
450 if (mInteractionRate < 0) {
451 mContextFile = std::filesystem::exists(
"collisioncontext.root") ? TFile::Open(
"collisioncontext.root") : nullptr;
452 if (!mContextFile || mContextFile->IsZombie()) {
453 LOG(fatal) <<
"Error: Interaction rate not provided and collision context file not found!";
457 mInteractionRate = std::floor(mCollisionContext->getDigitizerInteractionRate());
458 LOG(info) <<
"Interaction rate retrieved from collision context: " << mInteractionRate <<
" Hz";
459 if (mInteractionRate < 0) {
460 LOG(fatal) <<
"Error: Invalid interaction rate retrieved from collision context!";
464 auto ref =
static_cast<int>(std::floor(
fit->Eval(mInteractionRate / 1000.)));
467 LOG(fatal) <<
"Computed flat gas number reference per orbit is <=0";
470 LOG(info) <<
"Set flat gas number to " <<
ref <<
" loopers per orbit using " << fitName <<
" from " << mInteractionRate <<
" Hz interaction rate.";
472 setFlatGas(flat, -1,
ref);
476void GenTPCLoopers::SetAdjust(
float adjust)
478 if (mFlatGas && mFlatGasOrbit && adjust >= -1.f && adjust != 0.f) {
479 LOG(info) <<
"Adjusting flat gas number per orbit by " << adjust * 100.f <<
"%";
480 mFlatGasNumber =
static_cast<int>(std::round(mFlatGasNumber * (1.f + adjust)));
481 LOG(info) <<
"New flat gas number per orbit: " << mFlatGasNumber;
Ort::Env global_env(ORT_LOGGING_LEVEL_WARNING, "GlobalEnv")
static const HBFUtils & Instance()
GLdouble GLdouble GLdouble z
constexpr double LHCOrbitNS
TFitResultPtr fit(const size_t nBins, const T *arr, const T xMin, const T xMax, TF1 &func, std::string_view option="")
int getHepMCStatusCode(MCGenStatusEncoding enc)
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
Defining DataPointCompositeObject explicitly as copiable.
std::string to_string(gsl::span< T, Size > span)
std::vector< o2::ctf::BufferType > vec
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"