Project
Loading...
Searching...
No Matches
TPCLoopers.cxx
Go to the documentation of this file.
1// Copyright 2024-2025 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
13
16#include "CCDB/CcdbApi.h"
18#include "TF1.h"
19#include <filesystem>
22#include <iostream>
23#include <fstream>
24#include "TDatabasePDG.h"
25
26// Static Ort::Env instance for multiple onnx model loading
27Ort::Env global_env(ORT_LOGGING_LEVEL_WARNING, "GlobalEnv");
28
29// This class is responsible for loading the scaler parameters from a JSON file
30// and applying the inverse transformation to the generated data.
31
32void Scaler::load(const std::string& filename)
33{
34 std::ifstream file(filename);
35 if (!file.is_open()) {
36 throw std::runtime_error("Error: Could not open scaler file!");
37 }
38
39 std::string json_str((std::istreambuf_iterator<char>(file)), std::istreambuf_iterator<char>());
40 file.close();
41
42 rapidjson::Document doc;
43 doc.Parse(json_str.c_str());
44
45 if (doc.HasParseError()) {
46 throw std::runtime_error("Error: JSON parsing failed!");
47 }
48
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"]);
53}
54
55std::vector<double> Scaler::inverse_transform(const std::vector<double>& input)
56{
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]);
61 } else {
62 output.push_back(input[i] * outlier_scale[i - (input.size() - 2)] + outlier_center[i - (input.size() - 2)]);
63 }
64 }
65
66 return output;
67}
68
69std::vector<double> Scaler::jsonArrayToVector(const rapidjson::Value& jsonArray)
70{
71 std::vector<double> vec;
72 for (int i = 0; i < jsonArray.Size(); ++i) {
73 vec.push_back(jsonArray[i].GetDouble());
74 }
75 return vec;
76}
77
78// This class loads the ONNX model and generates samples using it.
79
80ONNXGenerator::ONNXGenerator(Ort::Env& shared_env, const std::string& model_path)
81 : env(shared_env), session(env, model_path.c_str(), Ort::SessionOptions{})
82{
83 // Create session options
84 Ort::SessionOptions session_options;
85 session = Ort::Session(env, model_path.c_str(), session_options);
86}
87
88std::vector<double> ONNXGenerator::generate_sample()
89{
90 Ort::AllocatorWithDefaultOptions allocator;
91
92 // Generate a latent vector (z)
93 std::vector<float> z(100);
94 for (auto& v : z) {
95 v = rand_gen.Gaus(0.0, 1.0);
96 }
97
98 // Prepare input tensor
99 std::vector<int64_t> input_shape = {1, 100};
100 // Get memory information
101 Ort::MemoryInfo memory_info = Ort::MemoryInfo::CreateCpu(OrtArenaAllocator, OrtMemTypeDefault);
102
103 // Create input tensor correctly
104 Ort::Value input_tensor = Ort::Value::CreateTensor<float>(
105 memory_info, z.data(), z.size(), input_shape.data(), input_shape.size());
106 // Run inference
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);
110
111 // Extract output
112 float* output_data = output_tensors.front().GetTensorMutableData<float>();
113 // Get the size of the output tensor
114 auto output_tensor_info = output_tensors.front().GetTensorTypeAndShapeInfo();
115 size_t output_data_size = output_tensor_info.GetElementCount(); // Total number of elements in the tensor
116 std::vector<double> output;
117 for (int i = 0; i < output_data_size; ++i) {
118 output.push_back(output_data[i]);
119 }
120
121 return output;
122}
123
124namespace o2
125{
126namespace eventgen
127{
128
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)
132{
133 // Checking if the model files exist and are not empty
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!";
139 exit(1);
140 }
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!";
143 exit(1);
144 }
145 model_file[0].close();
146 model_file[1].close();
147 // Checking if the scaler files exist and are not empty
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!";
153 exit(1);
154 }
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!";
157 exit(1);
158 }
159 scaler_file[0].close();
160 scaler_file[1].close();
161 // Checking if the poisson file exists and it's not empty
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!";
166 exit(1);
167 } else {
168 poisson_file >> mPoisson[0] >> mPoisson[1] >> mPoisson[2];
169 poisson_file.close();
170 mPoissonSet = true;
171 }
172 }
173 // Checking if the gauss file exists and it's not empty
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!";
178 exit(1);
179 } else {
180 gauss_file >> mGauss[0] >> mGauss[1] >> mGauss[2] >> mGauss[3];
181 gauss_file.close();
182 mGaussSet = true;
183 }
184 }
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);
191}
192
193Bool_t GenTPCLoopers::generateEvent()
194{
195 // Clear the vector of pairs
196 mGenPairs.clear();
197 // Clear the vector of compton electrons
198 mGenElectrons.clear();
199 if (mFlatGas) {
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();
205 // With flat gas the number of loopers are adapted based on time interval widths
206 // The denominator is either the LHC orbit (if mFlatGasOrbit is true) or the mean interaction time record interval
207 nLoopers = mFlatGasOrbit ? (mFlatGasNumber * (mTimeLimit / o2::constants::lhc::LHCOrbitNS)) : (mFlatGasNumber * (mTimeLimit / mIntTimeRecMean));
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);
213 mCurrentEvent++;
214 } else {
215 // Set number of loopers if poissonian params are available
216 if (mPoissonSet) {
217 mNLoopersPairs = static_cast<unsigned int>(std::round(mMultiplier[0] * PoissonPairs()));
218 LOG(debug) << "Generated loopers pairs (Poisson): " << mNLoopersPairs;
219 }
220 if (mGaussSet) {
221 mNLoopersCompton = static_cast<unsigned int>(std::round(mMultiplier[1] * GaussianElectrons()));
222 LOG(debug) << "Generated compton electrons (Gauss): " << mNLoopersCompton;
223 }
224 // Generate pairs
225 for (int i = 0; i < mNLoopersPairs; ++i) {
226 std::vector<double> pair = mONNX_pair->generate_sample();
227 // Apply the inverse transformation using the scaler
228 std::vector<double> transformed_pair = mScaler_pair->inverse_transform(pair);
229 mGenPairs.push_back(transformed_pair);
230 }
231 // Generate compton electrons
232 for (int i = 0; i < mNLoopersCompton; ++i) {
233 std::vector<double> electron = mONNX_compton->generate_sample();
234 // Apply the inverse transformation using the scaler
235 std::vector<double> transformed_electron = mScaler_compton->inverse_transform(electron);
236 mGenElectrons.push_back(transformed_electron);
237 }
238 }
239 return true;
240}
241
242Bool_t GenTPCLoopers::generateEvent(double time_limit)
243{
244 LOG(info) << "Time constraint for loopers: " << time_limit << " ns";
245 // Generate pairs
246 for (int i = 0; i < mNLoopersPairs; ++i) {
247 std::vector<double> pair = mONNX_pair->generate_sample();
248 // Apply the inverse transformation using the scaler
249 std::vector<double> transformed_pair = mScaler_pair->inverse_transform(pair);
250 transformed_pair[9] = gRandom->Uniform(0., time_limit); // Regenerate time, scaling is not needed because time_limit is already in nanoseconds
251 mGenPairs.push_back(transformed_pair);
252 }
253 // Generate compton electrons
254 for (int i = 0; i < mNLoopersCompton; ++i) {
255 std::vector<double> electron = mONNX_compton->generate_sample();
256 // Apply the inverse transformation using the scaler
257 std::vector<double> transformed_electron = mScaler_compton->inverse_transform(electron);
258 transformed_electron[6] = gRandom->Uniform(0., time_limit); // Regenerate time, scaling is not needed because time_limit is already in nanoseconds
259 mGenElectrons.push_back(transformed_electron);
260 }
261 LOG(info) << "Generated Particles with time limit";
262 return true;
263}
264
265std::vector<TParticle> GenTPCLoopers::importParticles()
266{
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();
270 // Get looper pairs from the event
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;
275 px_e = pair[0];
276 py_e = pair[1];
277 pz_e = pair[2];
278 px_p = pair[3];
279 py_p = pair[4];
280 pz_p = pair[5];
281 vx = pair[6];
282 vy = pair[7];
283 vz = pair[8];
284 time = pair[9];
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);
287 // Push the electron
288 TParticle electron(11, 1, -1, -1, -1, -1, px_e, py_e, pz_e, e_etot, vx, vy, vz, time / 1e9);
289 electron.SetStatusCode(o2::mcgenstatus::MCGenStatusEncoding(electron.GetStatusCode(), 0).fullEncoding);
290 electron.SetBit(ParticleStatus::kToBeDone, //
291 o2::mcgenstatus::getHepMCStatusCode(electron.GetStatusCode()) == 1);
292 particles.push_back(electron);
293 // Push the positron
294 TParticle positron(-11, 1, -1, -1, -1, -1, px_p, py_p, pz_p, p_etot, vx, vy, vz, time / 1e9);
295 positron.SetStatusCode(o2::mcgenstatus::MCGenStatusEncoding(positron.GetStatusCode(), 0).fullEncoding);
296 positron.SetBit(ParticleStatus::kToBeDone, //
297 o2::mcgenstatus::getHepMCStatusCode(positron.GetStatusCode()) == 1);
298 particles.push_back(positron);
299 }
300 // Get compton electrons from the event
301 for (auto& compton : mGenElectrons) {
302 double px, py, pz;
303 double vx, vy, vz, time;
304 double etot;
305 px = compton[0];
306 py = compton[1];
307 pz = compton[2];
308 vx = compton[3];
309 vy = compton[4];
310 vz = compton[5];
311 time = compton[6];
312 etot = TMath::Sqrt(px * px + py * py + pz * pz + mass_e * mass_e);
313 // Push the electron
314 TParticle electron(11, 1, -1, -1, -1, -1, px, py, pz, etot, vx, vy, vz, time / 1e9);
315 electron.SetStatusCode(o2::mcgenstatus::MCGenStatusEncoding(electron.GetStatusCode(), 0).fullEncoding);
316 electron.SetBit(ParticleStatus::kToBeDone, //
317 o2::mcgenstatus::getHepMCStatusCode(electron.GetStatusCode()) == 1);
318 particles.push_back(electron);
319 }
320
321 return particles;
322}
323
324unsigned int GenTPCLoopers::PoissonPairs()
325{
326 unsigned int poissonValue;
327 do {
328 // Generate a Poisson-distributed random number with mean mPoisson[0]
329 poissonValue = mRandGen.Poisson(mPoisson[0]);
330 } while (poissonValue < mPoisson[1] || poissonValue > mPoisson[2]); // Regenerate if out of range
331
332 return poissonValue;
333}
334
335unsigned int GenTPCLoopers::GaussianElectrons()
336{
337 unsigned int gaussValue;
338 do {
339 // Generate a Normal-distributed random number with mean mGass[0] and stddev mGauss[1]
340 gaussValue = mRandGen.Gaus(mGauss[0], mGauss[1]);
341 } while (gaussValue < mGauss[2] || gaussValue > mGauss[3]); // Regenerate if out of range
342
343 return gaussValue;
344}
345
346void GenTPCLoopers::SetNLoopers(unsigned int nsig_pair, unsigned int nsig_compton)
347{
348 if (mFlatGas) {
349 mNLoopersPairs = nsig_pair;
350 mNLoopersCompton = nsig_compton;
351 } else {
352 if (mPoissonSet) {
353 LOG(info) << "Poissonian parameters correctly loaded.";
354 } else {
355 mNLoopersPairs = nsig_pair;
356 }
357 if (mGaussSet) {
358 LOG(info) << "Gaussian parameters correctly loaded.";
359 } else {
360 mNLoopersCompton = nsig_compton;
361 }
362 }
363}
364
365void GenTPCLoopers::SetMultiplier(const std::array<float, 2>& mult)
366{
367 // Multipliers will work only if the poissonian and gaussian parameters are set
368 // otherwise they will be ignored
369 if (mult[0] < 0 || mult[1] < 0) {
370 LOG(fatal) << "Error: Multiplier values must be non-negative!";
371 exit(1);
372 } else {
373 LOG(info) << "Multiplier values set to: Pair = " << mult[0] << ", Compton = " << mult[1];
374 mMultiplier[0] = mult[0];
375 mMultiplier[1] = mult[1];
376 }
377}
378
379void GenTPCLoopers::setFlatGas(Bool_t flat, Int_t number, Int_t nloopers_orbit)
380{
381 mFlatGas = flat;
382 if (mFlatGas) {
383 if (nloopers_orbit > 0) {
384 mFlatGasOrbit = true;
385 mFlatGasNumber = nloopers_orbit;
386 LOG(info) << "Flat gas loopers will be generated using orbit reference.";
387 } else {
388 mFlatGasOrbit = false;
389 if (number < 0) {
390 LOG(warn) << "Warning: Number of loopers per event must be non-negative! Switching option off.";
391 mFlatGas = false;
392 mFlatGasNumber = -1;
393 } else {
394 mFlatGasNumber = number;
395 }
396 }
397 if (mFlatGas) {
398 mContextFile = std::filesystem::exists("collisioncontext.root") ? TFile::Open("collisioncontext.root") : nullptr;
399 mCollisionContext = mContextFile ? (o2::steer::DigitizationContext*)mContextFile->Get("DigitizationContext") : 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!";
403 exit(1);
404 } else {
405 LOG(info) << "Interaction Time records has " << mInteractionTimeRecords.size() << " entries.";
406 mCollisionContext->printCollisionSummary();
407 }
408 for (int c = 0; c < mInteractionTimeRecords.size() - 1; c++) {
409 mIntTimeRecMean += mInteractionTimeRecords[c + 1].bc2ns() - mInteractionTimeRecords[c].bc2ns();
410 }
411 mIntTimeRecMean /= (mInteractionTimeRecords.size() - 1); // Average interaction time record used as reference
412 const auto& hbfUtils = o2::raw::HBFUtils::Instance();
413 // Get the start time of the second orbit after the last interaction record
414 const auto& lastIR = mInteractionTimeRecords.back();
415 o2::InteractionRecord finalOrbitIR(0, lastIR.orbit + 2); // Final orbit, BC = 0
416 mTimeEnd = finalOrbitIR.bc2ns();
417 LOG(debug) << "Final orbit start time: " << mTimeEnd << " ns while last interaction record time is " << mInteractionTimeRecords.back().bc2ns() << " ns";
418 }
419 } else {
420 mFlatGasNumber = -1;
421 }
422 LOG(info) << "Flat gas loopers: " << (mFlatGas ? "ON" : "OFF") << ", Reference loopers number per " << (mFlatGasOrbit ? "orbit " : "event ") << mFlatGasNumber;
423}
424
425void GenTPCLoopers::setFractionPairs(float fractionPairs)
426{
427 if (fractionPairs < 0 || fractionPairs > 1) {
428 LOG(fatal) << "Error: Loops fraction for pairs must be in the range [0, 1].";
429 exit(1);
430 }
431 mLoopsFractionPairs = fractionPairs;
432 LOG(info) << "Pairs fraction set to: " << mLoopsFractionPairs;
433}
434
435void GenTPCLoopers::SetRate(const std::string& rateFile, bool isPbPb = true, int intRate)
436{
437 // Checking if the rate file exists and is not empty
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!";
441 exit(1);
442 }
443 const char* fitName = isPbPb ? "fitPbPb" : "fitpp";
444 auto fit = (TF1*)rate_file.Get(fitName);
445 if (!fit) {
446 LOG(fatal) << "Error: Could not find fit function '" << fitName << "' in rate file!";
447 exit(1);
448 }
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!";
454 exit(1);
455 }
456 mCollisionContext = (o2::steer::DigitizationContext*)mContextFile->Get("DigitizationContext");
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!";
461 exit(1);
462 }
463 }
464 auto ref = static_cast<int>(std::floor(fit->Eval(mInteractionRate / 1000.))); // fit expects rate in kHz
465 rate_file.Close();
466 if (ref <= 0) {
467 LOG(fatal) << "Computed flat gas number reference per orbit is <=0";
468 exit(1);
469 } else {
470 LOG(info) << "Set flat gas number to " << ref << " loopers per orbit using " << fitName << " from " << mInteractionRate << " Hz interaction rate.";
471 auto flat = true;
472 setFlatGas(flat, -1, ref);
473 }
474}
475
476void GenTPCLoopers::SetAdjust(float adjust)
477{
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;
482 }
483}
484
485} // namespace eventgen
486} // namespace o2
std::ostringstream debug
int16_t time
Definition RawEventData.h:4
int32_t i
void output(const std::map< std::string, ChannelStat > &channels)
Definition rawdump.cxx:197
@ kToBeDone
uint32_t c
Definition RawData.h:2
Ort::Env global_env(ORT_LOGGING_LEVEL_WARNING, "GlobalEnv")
const GLdouble * v
Definition glcorearb.h:832
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
constexpr double LHCOrbitNS
TFitResultPtr fit(const size_t nBins, const T *arr, const T xMin, const T xMax, TF1 &func, std::string_view option="")
Definition fit.h:60
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)
Definition common.h:52
std::string filename()
std::vector< o2::ctf::BufferType > vec
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"