File indexing completed on 2026-05-12 08:15:19
0001 #ifdef __MACH__
0002 #include <mach/mach.h>
0003 #endif
0004
0005 #include <iostream>
0006 #include <fstream>
0007 #include <string>
0008 #include <random>
0009 #include <ctime>
0010 #include <cstdlib>
0011 #include <vector>
0012 #include <algorithm>
0013 #include <numeric>
0014 #include <stdexcept>
0015 #include <chrono>
0016 #include <cmath>
0017 #include <random>
0018 #include <tuple>
0019 #include <unistd.h>
0020 #include <sys/resource.h>
0021
0022 #include <HepMC3/ReaderFactory.h>
0023 #include <HepMC3/WriterAscii.h>
0024 #include "HepMC3/WriterRootTree.h"
0025 #include "HepMC3/GenRunInfo.h"
0026 #include "HepMC3/GenEvent.h"
0027 #include "HepMC3/Print.h"
0028
0029 #include "argparse/argparse.hpp"
0030
0031 using std::cout;
0032 using std::cerr;
0033 using std::endl;
0034 using std::string;
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050 struct BackgroundConfig {
0051 std::string file;
0052 double frequency=0;
0053 int skip=0;
0054 int status=0;
0055 } ;
0056
0057 class SignalBackgroundMerger {
0058
0059 private:
0060
0061 std::shared_ptr<HepMC3::Reader> sigAdapter;
0062 double sigFreq = 0;
0063 int sigStatus = 0;
0064 std::map<std::string, std::shared_ptr<HepMC3::Reader>> freqAdapters;
0065 std::map<std::string, double> freqs;
0066 std::map<std::string, int> baseStatuses;
0067
0068 std::map<std::string,
0069 std::tuple<std::vector<HepMC3::GenEvent>,
0070 std::piecewise_constant_distribution<>,
0071 double>
0072 > weightDict;
0073
0074
0075 typedef struct{
0076 long eventCount;
0077 long particleCount;
0078 } stats;
0079 std::map<std::string, stats > infoDict;
0080
0081 public:
0082
0083 SignalBackgroundMerger(int argc, char* argv[]) {
0084 auto t0 = std::chrono::high_resolution_clock::now();
0085
0086
0087 digestArgs(argc, argv);
0088 rng.seed( rngSeed );
0089 banner();
0090 if (outputFile != "" ) {
0091 outputFileName = outputFile;
0092 } else {
0093 outputFileName = nameGen();
0094 }
0095 std::cout << "\n==================================================================\n";
0096 cout << "Writing to " << outputFileName << endl;
0097
0098 PrepData ( signalFile, signalFreq, signalSkip, signalStatus, true );
0099 for (const auto& bg : backgroundFiles) {
0100 PrepData ( bg.file, bg.frequency, bg.skip, bg.status, false );
0101 }
0102
0103
0104 auto t1 = std::chrono::high_resolution_clock::now();
0105 std::cout << "Initiation time: " << std::round(std::chrono::duration<double, std::chrono::seconds::period>(t1 - t0).count()) << " sec" << std::endl;
0106 std::cout << "\n==================================================================\n" << std::endl;
0107
0108 }
0109
0110
0111 std::vector<BackgroundConfig>
0112 parse_backgrounds(const std::vector<std::string> &raw_args_list) {
0113 std::vector<BackgroundConfig> backgrounds;
0114 auto is_pure_integer = [](const std::string &str) {
0115 if (str.empty()) return false;
0116 for (char c : str) {
0117 if (!std::isdigit(c)) return false;
0118 }
0119 return true;
0120 };
0121
0122
0123 for (size_t i = 0; i < raw_args_list.size();) {
0124
0125 size_t args_count = 2;
0126
0127
0128 if (i + 2 < raw_args_list.size() && is_pure_integer(raw_args_list[i + 2])) {
0129 args_count = 3;
0130
0131 if (i + 3 < raw_args_list.size() && is_pure_integer(raw_args_list[i + 3])) {
0132 args_count = 4;
0133 }
0134 }
0135
0136
0137 if (i + args_count > raw_args_list.size()) {
0138 args_count = raw_args_list.size() - i;
0139 }
0140
0141 if (args_count < 2) {
0142 throw std::runtime_error("Background file " +
0143 std::to_string(backgrounds.size()) +
0144 " must have at least 2 arguments");
0145 }
0146
0147 try {
0148 BackgroundConfig bg;
0149 bg.file = raw_args_list[i];
0150 bg.frequency = std::stod(raw_args_list[i + 1]);
0151 bg.skip = (args_count > 2) ? std::stoi(raw_args_list[i + 2]) : 0;
0152 bg.status = (args_count > 3) ? std::stoi(raw_args_list[i + 3]) : 0;
0153 backgrounds.push_back(bg);
0154 } catch (const std::exception &e) {
0155 throw std::runtime_error("Error parsing background file " +
0156 std::to_string(backgrounds.size()) + ": " +
0157 e.what());
0158 }
0159
0160 i += args_count;
0161 }
0162
0163 return backgrounds;
0164 }
0165
0166 void merge(){
0167 auto t1 = std::chrono::high_resolution_clock::now();
0168
0169
0170 auto runInfo = std::make_shared<HepMC3::GenRunInfo>();
0171
0172 runInfo->add_attribute("hepmc_merger_signal_file",
0173 std::make_shared<HepMC3::StringAttribute>(signalFile));
0174 runInfo->add_attribute("hepmc_merger_signal_frequency_kHz",
0175 std::make_shared<HepMC3::DoubleAttribute>(signalFreq));
0176
0177 std::string bgFiles, bgFreqs, bgAvgRates;
0178 for (size_t bi = 0; bi < backgroundFiles.size(); ++bi) {
0179 if (bi) { bgFiles += ";"; bgFreqs += ";"; bgAvgRates += ";"; }
0180 bgFiles += backgroundFiles[bi].file;
0181 bgFreqs += std::to_string(backgroundFiles[bi].frequency);
0182 double avgRate = 0.0;
0183 if (backgroundFiles[bi].frequency <= 0.0) {
0184 auto it = weightDict.find(backgroundFiles[bi].file);
0185 if (it != weightDict.end())
0186 avgRate = std::get<2>(it->second) * 1e6;
0187 }
0188 bgAvgRates += std::to_string(avgRate);
0189 }
0190 runInfo->add_attribute("hepmc_merger_background_files",
0191 std::make_shared<HepMC3::StringAttribute>(bgFiles));
0192 runInfo->add_attribute("hepmc_merger_background_frequencies_kHz",
0193 std::make_shared<HepMC3::StringAttribute>(bgFreqs));
0194 runInfo->add_attribute("hepmc_merger_background_avg_rates_kHz",
0195 std::make_shared<HepMC3::StringAttribute>(bgAvgRates));
0196
0197 runInfo->add_attribute("hepmc_merger_integration_window_ns",
0198 std::make_shared<HepMC3::DoubleAttribute>(intWindow));
0199 runInfo->add_attribute("hepmc_merger_n_slices",
0200 std::make_shared<HepMC3::IntAttribute>(nSlices));
0201
0202
0203 std::shared_ptr<HepMC3::Writer> f;
0204 if (rootFormat)
0205 f = std::make_shared<HepMC3::WriterRootTree>(outputFileName, runInfo);
0206 else
0207 f = std::make_shared<HepMC3::WriterAscii>(outputFileName, runInfo);
0208
0209
0210 int i = 0;
0211 for (i = 0; i<nSlices; ++i ) {
0212 if (i % 1000 == 0 || verbose ) squawk(i);
0213 auto hepSlice = mergeSlice(i);
0214 if (!hepSlice) {
0215 std::cout << "Exhausted signal source." << std::endl;
0216 break;
0217 }
0218 hepSlice->set_event_number(i);
0219 f->write_event(*hepSlice);
0220 }
0221 std::cout << "Finished all requested slices." << std::endl;
0222
0223 int slicesDone = i;
0224 auto t2 = std::chrono::high_resolution_clock::now();
0225
0226 std::cout << "Slice loop time: " << std::round(std::chrono::duration<double, std::chrono::minutes::period>(t2 - t1).count()) << " min" << std::endl;
0227 std::cout << " -- " << std::round(std::chrono::duration<double, std::chrono::microseconds::period>(t2 - t1).count() / i) << " us / slice" << std::endl;
0228
0229 for (auto info : infoDict) {
0230 std::cout << "From " << info.first << std::endl;
0231 std::cout << " placed " << info.second.eventCount << " events" << std::endl;
0232 std::cout << " --> on average " << std::setprecision(3) << info.second.eventCount / float(nSlices) << std::endl;
0233 std::cout << " placed " << info.second.particleCount << " final state particles" << std::endl;
0234 std::cout << " --> on average " << std::setprecision(3) << info.second.particleCount / float(nSlices) << std::endl;
0235
0236 }
0237
0238 struct rusage r_usage;
0239 getrusage(RUSAGE_SELF, &r_usage);
0240
0241
0242
0243
0244 #ifdef __MACH__
0245 float mbsize = 1024 * 1024;
0246 #else
0247 float mbsize = 1024;
0248 #endif
0249
0250
0251 std::cout << endl << "Maximum Resident Memory " << r_usage.ru_maxrss / mbsize << " MB" << std::endl;
0252
0253 sigAdapter->close();
0254 for (auto& it : freqAdapters) {
0255 it.second->close();
0256 }
0257 f->close();
0258
0259 }
0260
0261
0262 void digestArgs(int argc, char* argv[]) {
0263
0264
0265
0266
0267
0268 argparse::ArgumentParser args ("Merge signal events with up to four background sources.");
0269
0270 args.add_argument("-i", "--signalFile")
0271 .default_value(std::string("root://dtn-eic.jlab.org//volatile/eic/EPIC/EVGEN/SIDIS/pythia6-eic/1.0.0/10x100/q2_0to1/pythia_ep_noradcor_10x100_q2_0.000000001_1.0_run1.ab.hepmc3.tree.root"))
0272 .help("Name of the HEPMC file with the signal events");
0273
0274 args.add_argument("-sf", "--signalFreq")
0275 .default_value(0.0)
0276 .scan<'g', double>()
0277 .help("Signal frequency in kHz. Default is 0 to have exactly one signal event per slice. Set to the estimated DIS rate to randomize.");
0278
0279 args.add_argument("-S", "--signalSkip")
0280 .default_value(0)
0281 .scan<'i', int>()
0282 .help("Number of signals events to skip. Default is 0.");
0283
0284 args.add_argument("-St", "--signalStatus")
0285 .default_value(0)
0286 .scan<'i', int>()
0287 .help("Apply shift on particle generatorStatus code for signal. Default is 0. ");
0288
0289 args.add_argument("-b","--bgFile")
0290 .nargs(2,4)
0291 .append()
0292 .help("Tuple with name of the HEPMC file with background events, background frequency in kHz, number of background events to skip (default 0), shift on particle generatorStatus code (default 0).");
0293
0294 args.add_argument("-o", "--outputFile")
0295 .default_value(std::string("bgmerged.hepmc3.tree.root"))
0296 .help("Specify the output file name. By default bgmerged.hepmc3.tree.root is used");
0297
0298 args.add_argument("-r", "--rootFormat")
0299 .default_value(true)
0300 .implicit_value(true)
0301 .help("Use hepmc.root output format, default is true.");
0302
0303 args.add_argument("-w", "--intWindow")
0304 .default_value(2000.0)
0305 .scan<'g', double>()
0306 .help("Length of the integration window in nanoseconds. Default is 2000.");
0307
0308 args.add_argument("-N", "--nSlices")
0309 .default_value(10000)
0310 .scan<'i', int>()
0311 .help("Number of sampled time slices ('events'). Default is 10000. If set to -1, all events in the signal file will be used and background files cycled as needed.");
0312
0313 args.add_argument("--squashTime")
0314 .default_value(false)
0315 .implicit_value(true)
0316 .help("Integration is performed but no time information is associated to vertices.");
0317
0318 args.add_argument("--rngSeed")
0319 .default_value(0)
0320 .action([](const std::string& value) { return std::stoi(value); })
0321 .help("Random seed, default is None");
0322
0323 args.add_argument("-v", "--verbose")
0324 .default_value(false)
0325 .implicit_value(true)
0326 .help("Display details for every slice.");
0327
0328 try {
0329 args.parse_args(argc, argv);
0330 }
0331 catch (const std::runtime_error& err) {
0332 std::cout << err.what() << std::endl;
0333 std::cout << args;
0334 exit(EXIT_FAILURE);
0335 }
0336
0337 signalFile = args.get<std::string>("--signalFile");
0338 signalFreq = args.get<double>("--signalFreq");
0339 signalSkip = args.get<int>("--signalSkip");
0340 signalStatus = args.get<int>("--signalStatus");
0341 backgroundFiles = parse_backgrounds(args.get<std::vector<std::string>>("--bgFile"));
0342 outputFile = args.get<std::string>("--outputFile");
0343 rootFormat = args.get<bool>("--rootFormat");
0344 intWindow = args.get<double>("--intWindow");
0345 nSlices = args.get<int>("--nSlices");
0346 squashTime = args.get<bool>("--squashTime");
0347 rngSeed = args.get<int>("--rngSeed");
0348 verbose = args.get<bool>("--verbose");
0349
0350
0351 }
0352
0353
0354 void banner() {
0355
0356 std::cout << "==================================================================" << std::endl;
0357 std::cout << "=== EPIC HEPMC MERGER ===" << std::endl;
0358 std::cout << "authors: Benjamen Sterwerf* (bsterwerf@berkeley.edu), Kolja Kauder** (kkauder@bnl.gov), Reynier Cruz-Torres***" << std::endl;
0359 std::cout << "* University of California, Berkeley" << std::endl;
0360 std::cout << "** Brookhaven National Laboratory" << std::endl;
0361 std::cout << "*** formerly Lawrence Berkeley National Laboratory" << std::endl;
0362 std::cout << "\nFor more information, run \n./signal_background_merger --help" << std::endl;
0363
0364 std::string statusMessage = "Shifting all particle status codes from this source by ";
0365 std::vector<int> statusList_stable, statusList_decay;
0366
0367 std::cout << "Number of Slices:" << nSlices << endl;
0368 std::string freqTerm = signalFreq > 0 ? std::to_string(signalFreq) + " kHz" : "(one event per time slice)";
0369 std::string statusTerm = signalStatus > 0 ? statusMessage + std::to_string(signalStatus): "";
0370 if (signalStatus>0){
0371 statusList_stable.push_back(signalStatus+1);
0372 statusList_decay.push_back(signalStatus+2);
0373 }
0374 std::cout << "Signal events file and frequency:\n";
0375 std::cout << "\t- " << signalFile << "\t" << freqTerm << "\n" << statusTerm << "\n";
0376
0377 std::cout << "\nBackground files and their respective frequencies:\n";
0378
0379 for (const auto& bg : backgroundFiles) {
0380 if (!bg.file.empty()) {
0381 freqTerm = bg.frequency > 0 ? std::to_string(bg.frequency) + " kHz" : "(from weights)";
0382 statusTerm = bg.status > 0 ? statusMessage + std::to_string(bg.status) : "";
0383 std::cout << "\t- " << bg.file << "\t" << freqTerm << "\n" << statusTerm << "\n";
0384 if (bg.status>0){
0385 statusList_stable.push_back(bg.status+1);
0386 statusList_decay.push_back(bg.status+2);
0387 }
0388 }
0389 }
0390
0391 auto join = [](const std::vector<int>& vec) {
0392 return std::accumulate(vec.begin(), vec.end(), std::string(),
0393 [](const std::string& a, int b) {
0394 return a.empty() ? std::to_string(b) : a + " " + std::to_string(b);
0395 });
0396 };
0397 std::string stableStatuses = join(statusList_stable);
0398 std::string decayStatuses = join(statusList_decay);
0399 std::string message = "\n!!!Attention!!!\n To proceed the shifted particles statuses in DD4hep, please add the following options to ddsim:\n"
0400 "--physics.alternativeStableStatuses=\"" + stableStatuses +
0401 "\" --physics.alternativeDecayStatuses=\"" + decayStatuses + "\"\n";
0402 std::cout << message<<std::endl;
0403 }
0404
0405
0406 void PrepData(const std::string& fileName, double freq, int skip=0, int baseStatus=0, bool signal=false) {
0407 if (fileName.empty()) return;
0408
0409 cout << "Prepping " << fileName << endl;
0410 std::shared_ptr<HepMC3::Reader> adapter;
0411 try {
0412 adapter = HepMC3::deduce_reader(fileName);
0413 if (!adapter) {
0414 throw std::runtime_error("Failed to open file");
0415 }
0416 } catch (const std::runtime_error& e) {
0417 std::cerr << "Opening " << fileName << " failed: " << e.what() << std::endl;
0418 exit(EXIT_FAILURE);
0419 }
0420
0421 infoDict[fileName] = {0,0};
0422
0423 if (signal) {
0424 sigAdapter = adapter;
0425 sigFreq = freq;
0426 sigStatus = baseStatus;
0427 sigAdapter->skip(skip);
0428 return;
0429 }
0430
0431
0432 if (freq <= 0) {
0433 std::cout << "Reading in all events from " << fileName << std::endl;
0434 std::vector<HepMC3::GenEvent> events;
0435 std::vector<double> weights;
0436
0437 while(!adapter->failed()) {
0438 HepMC3::GenEvent evt(HepMC3::Units::GEV,HepMC3::Units::MM);
0439 adapter->read_event(evt);
0440
0441
0442 if (double w=evt.weight() > 0){
0443 events.push_back(evt);
0444 weights.push_back(evt.weight());
0445 }
0446 }
0447 adapter->close();
0448
0449 double avgRate = 0.0;
0450 for ( auto w : weights ){ avgRate += w;}
0451 avgRate /= weights.size();
0452 avgRate *= 1e-9;
0453 std::cout << "Average rate is " << avgRate << " GHz" << std::endl;
0454
0455 std::vector<int> indices (weights.size());
0456 std::iota (std::begin(indices), std::end(indices), 0);
0457
0458
0459
0460
0461
0462
0463
0464
0465
0466
0467 std::piecewise_constant_distribution<> weightedDist(std::begin(indices),std::end(indices),
0468 std::begin(weights));
0469 weightDict[fileName] = { std::make_tuple(events, weightedDist, avgRate) };
0470
0471 return;
0472 }
0473
0474
0475 adapter->skip(skip);
0476 freqAdapters[fileName] = adapter;
0477 freqs[fileName] = freq;
0478 baseStatuses[fileName] = baseStatus;
0479 }
0480
0481
0482 bool hasEnding (std::string const &fullString, std::string const &ending) {
0483 if (fullString.length() >= ending.length()) {
0484 return (0 == fullString.compare (fullString.length() - ending.length(), ending.length(), ending));
0485 } else {
0486 return false;
0487 }
0488 }
0489
0490
0491 std::string nameGen() {
0492
0493
0494 std::string name = signalFile;
0495 if (nSlices > 0) {
0496 size_t pos = name.find(".hepmc");
0497 if (pos != std::string::npos) {
0498 name.replace(pos, 6, "_n_" + std::to_string(nSlices) + ".hepmc");
0499 }
0500 }
0501
0502 if ( rootFormat && !hasEnding(name,".root")){
0503 name.append(".root");
0504 }
0505 name = "bgmerged_" + name;
0506
0507 return name;
0508 }
0509
0510
0511 void squawk(int i) {
0512
0513
0514 #ifdef __MACH__
0515 task_basic_info_data_t info;
0516 mach_msg_type_number_t size = sizeof(info);
0517 kern_return_t kerr = task_info(mach_task_self(),
0518 TASK_BASIC_INFO,
0519 (task_info_t)&info,
0520 &size);
0521
0522 long memory_usage = -1;
0523 if (kerr == KERN_SUCCESS) {
0524 memory_usage = info.resident_size / 1024 / 1024;
0525 }
0526 #else
0527 std::ifstream statm("/proc/self/statm");
0528 long size, resident, share, text, lib, data, dt;
0529 statm >> size >> resident >> share >> text >> lib >> data >> dt;
0530 statm.close();
0531
0532 long page_size = sysconf(_SC_PAGESIZE);
0533 long memory_usage = resident * page_size / 1024 / 1024 ;
0534 #endif
0535
0536
0537 std::cout << "Working on slice " << i + 1 << std::endl;
0538 std::cout << "Current memory usage: " << memory_usage << " MB" << std::endl;
0539
0540 }
0541
0542
0543 std::unique_ptr<HepMC3::GenEvent> mergeSlice(int i) {
0544 auto hepSlice = std::make_unique<HepMC3::GenEvent>(HepMC3::Units::GEV, HepMC3::Units::MM);
0545
0546 addFreqEvents(signalFile, sigAdapter, sigFreq, hepSlice, signalStatus, true);
0547
0548 for (const auto& freqBgs : freqAdapters) {
0549 auto fileName=freqBgs.first;
0550 addFreqEvents(fileName, freqAdapters[fileName], freqs[fileName], hepSlice, baseStatuses[fileName], false);
0551 }
0552
0553 for (const auto& fileName : weightDict) {
0554 addWeightedEvents(fileName.first, hepSlice, baseStatuses[fileName.first]);
0555 }
0556
0557 return hepSlice;
0558 };
0559
0560
0561
0562 void addFreqEvents(std::string fileName, std::shared_ptr<HepMC3::Reader>& adapter, const double freq,
0563 std::unique_ptr<HepMC3::GenEvent>& hepSlice, int baseStatus = 0, bool signal = false) {
0564
0565
0566
0567 std::vector<double> timeline;
0568
0569 std::uniform_real_distribution<> uni(0, intWindow);
0570 if (freq == 0){
0571 if (!signal) {
0572 std::cerr << "frequency can't be 0 for background files" << std::endl;
0573 exit(EXIT_FAILURE);
0574 }
0575
0576 timeline.push_back(uni(rng));
0577 } else {
0578
0579 timeline = poissonTimes(freq, intWindow);
0580 }
0581
0582 if ( verbose) std::cout << "Placing " << timeline.size() << " events from " << fileName << std::endl;
0583
0584 if (timeline.empty()) return;
0585 long particleCount = 0;
0586
0587
0588 for (double time : timeline) {
0589 if(adapter->failed()) {
0590 try{
0591 if (signal) {
0592 throw std::ifstream::failure("EOF");
0593 } else {
0594 std::cout << "Cycling back to the start of " << fileName << std::endl;
0595 adapter->close();
0596 adapter = HepMC3::deduce_reader(fileName);
0597 }
0598 } catch (std::ifstream::failure& e) {
0599 continue;
0600 }
0601 }
0602
0603 HepMC3::GenEvent inevt;
0604 adapter->read_event(inevt);
0605 if (signal && (signalFreq == 0.0)){
0606 hepSlice->weights() = inevt.weights();
0607 }
0608
0609 if (squashTime) time = 0;
0610 particleCount += insertHepmcEvent( inevt, hepSlice, time, baseStatus, signal);
0611 }
0612
0613 infoDict[fileName].eventCount += timeline.size();
0614 infoDict[fileName].particleCount += particleCount;
0615
0616
0617 return;
0618 }
0619
0620
0621
0622 void addWeightedEvents(std::string fileName, std::unique_ptr<HepMC3::GenEvent>& hepSlice, int baseStatus=0, bool signal = false) {
0623 auto& [events, weightedDist, avgRate ] = weightDict[fileName];
0624
0625
0626 int nEvents;
0627 std::poisson_distribution<> d( intWindow * avgRate );
0628
0629
0630
0631 while (true) {
0632 nEvents = d(rng);
0633 if (nEvents > events.size()) {
0634 std::cout << "WARNING: Trying to place " << nEvents << " events from " << fileName
0635 << " but the file doesn't have enough. Rerolling, but this is not physical." << std::endl;
0636 continue;
0637 }
0638 break;
0639 }
0640
0641 if (verbose) std::cout << "Placing " << nEvents << " events from " << fileName << std::endl;
0642
0643
0644
0645 std::vector<HepMC3::GenEvent> toPlace(nEvents);
0646 for ( auto& e : toPlace ){
0647 auto i = static_cast<int> (weightedDist(rng));
0648 e = events.at(i);
0649 }
0650
0651
0652 std::vector<double> timeline;
0653 std::uniform_real_distribution<> uni(0, intWindow);
0654 long particleCount = 0;
0655 if (!squashTime) {
0656 for ( auto& e : toPlace ){
0657 double time = squashTime ? 0 : uni(rng);
0658 particleCount += insertHepmcEvent( e, hepSlice, time, baseStatus, signal);
0659 }
0660 }
0661
0662 infoDict[fileName].eventCount += nEvents;
0663 infoDict[fileName].particleCount += particleCount;
0664
0665 return;
0666 }
0667
0668
0669 long insertHepmcEvent( const HepMC3::GenEvent& inevt,
0670 std::unique_ptr<HepMC3::GenEvent>& hepSlice, double time=0, int baseStatus=0, bool signal = false) {
0671
0672 double timeHepmc = c_light * time;
0673
0674 std::vector<HepMC3::GenParticlePtr> particles;
0675 std::vector<HepMC3::GenVertexPtr> vertices;
0676
0677
0678
0679 for (auto& vertex : inevt.vertices()) {
0680 HepMC3::FourVector position = vertex->position();
0681 position.set_t(position.t() + timeHepmc);
0682 auto v1 = std::make_shared<HepMC3::GenVertex>(position);
0683 vertices.push_back(v1);
0684 }
0685
0686
0687 long finalParticleCount = 0;
0688 for (auto& particle : inevt.particles()) {
0689 HepMC3::FourVector momentum = particle->momentum();
0690 int status = particle->status();
0691 if (status == 1 ) finalParticleCount++;
0692 int pid = particle->pid();
0693 status += baseStatus;
0694 auto p1 = std::make_shared<HepMC3::GenParticle> (momentum, pid, status);
0695 p1->set_generated_mass(particle->generated_mass());
0696 particles.push_back(p1);
0697
0698 if (particle->production_vertex()->id() < 0) {
0699 int production_vertex = particle->production_vertex()->id();
0700 vertices[abs(production_vertex) - 1]->add_particle_out(p1);
0701 hepSlice->add_particle(p1);
0702 }
0703
0704
0705 if (particle->end_vertex()) {
0706 int end_vertex = particle->end_vertex()->id();
0707 vertices.at(abs(end_vertex) - 1)->add_particle_in(p1);
0708 }
0709 }
0710
0711
0712 for (auto& vertex : vertices) {
0713 hepSlice->add_vertex(vertex);
0714 }
0715
0716 return finalParticleCount;
0717 }
0718
0719
0720
0721 std::vector<double> poissonTimes(double mu, double endTime) {
0722 std::exponential_distribution<> exp(mu);
0723
0724 double t = 0;
0725 std::vector<double> ret;
0726 while (true) {
0727 double delt = exp(rng)*1e6;
0728
0729 t += delt;
0730 if (t >= endTime) {
0731 break;
0732 }
0733 ret.push_back(t);
0734 }
0735 return ret;
0736 }
0737
0738
0739
0740
0741 std::mt19937 rng;
0742 string signalFile;
0743 double signalFreq;
0744 int signalSkip;
0745 int signalStatus;
0746 std::vector<BackgroundConfig> backgroundFiles;
0747 string outputFile;
0748 string outputFileName;
0749 bool rootFormat;
0750 double intWindow;
0751 int nSlices;
0752 bool squashTime;
0753 int rngSeed;
0754 bool verbose;
0755
0756 const double c_light = 299.792458;
0757 };
0758
0759
0760 int main(int argc, char* argv[]) {
0761
0762 auto t0 = std::chrono::high_resolution_clock::now();
0763
0764 SignalBackgroundMerger sbm (argc, argv);
0765
0766
0767 sbm.merge();
0768
0769 std::cout << "\n==================================================================\n";
0770 std::cout << "Overall running time: " << std::round(std::chrono::duration<double, std::chrono::minutes::period>(std::chrono::high_resolution_clock::now() - t0).count()) << " min" << std::endl;
0771
0772 return EXIT_SUCCESS;
0773 }