Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-10-13 08:38:07

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 <sys/resource.h>
0020 
0021 #include <HepMC3/ReaderFactory.h>
0022 #include <HepMC3/WriterAscii.h>
0023 #include "HepMC3/WriterRootTree.h"
0024 #include "HepMC3/GenRunInfo.h"
0025 #include "HepMC3/GenEvent.h"
0026 #include "HepMC3/Print.h"
0027 
0028 #include "argparse/argparse.hpp"
0029 
0030 using std::cout;
0031 using std::cerr;
0032 using std::endl;
0033 using std::string;
0034 
0035 
0036 // =============================================================
0037 /**
0038     Combine signal and up to four background HEPMC files.
0039     
0040     Typical usage:
0041     ./SignalBackgroundMerger --signalFile dis.hepmc3.tree.root --signalFreq 0 \
0042             --bgFile hgas.hepmc3.tree.root 2000 0 2000 \
0043         --bgFile egastouschk.hepmc3.tree.root 20 0 3000 \
0044             --bgFile egascouloumb.hepmc3.tree.root 20 0 4000 \
0045         --bgFile egasbrems.hepmc3.tree.root 20 0 5000 \
0046             --bgFile synrad.hepmc3.tree.root 25 0 6000
0047 **/    
0048 
0049 struct BackgroundConfig {
0050     std::string file;
0051     double frequency=0;
0052     int skip=0;
0053     int status=0;
0054 } ;
0055 
0056 class SignalBackgroundMerger {
0057 
0058 private:
0059   // more private data at the end; pulling these more complicated objects up for readability
0060   std::shared_ptr<HepMC3::Reader> sigAdapter;
0061   double sigFreq = 0;
0062   int sigStatus = 0;
0063   std::map<std::string, std::shared_ptr<HepMC3::Reader>> freqAdapters;
0064   std::map<std::string, double> freqs;
0065   std::map<std::string, int> baseStatuses;
0066 
0067   std::map<std::string,
0068       std::tuple<std::vector<HepMC3::GenEvent>,
0069               std::piecewise_constant_distribution<>,
0070               double>
0071        > weightDict;
0072 
0073   // just keep count of some numbers, could be more sophisticated
0074   typedef struct{
0075     long eventCount;
0076     long particleCount;
0077   } stats;
0078   std::map<std::string, stats > infoDict;
0079 
0080 public:
0081 
0082   SignalBackgroundMerger(int argc, char* argv[]) {
0083     auto t0 = std::chrono::high_resolution_clock::now();
0084 
0085     // Parse arguments, print banner, open files, initialize rng
0086     digestArgs(argc, argv);
0087     rng.seed( rngSeed );
0088     banner();
0089     if (outputFile != "" ) {
0090       outputFileName = outputFile;
0091     } else {
0092       outputFileName = nameGen();
0093     }
0094     std::cout << "\n==================================================================\n";
0095     cout << "Writing to " << outputFileName << endl;
0096 
0097     PrepData ( signalFile, signalFreq, signalSkip, signalStatus, true );
0098     for (const auto& bg : backgroundFiles) {
0099     PrepData ( bg.file, bg.frequency, bg.skip, bg.status, false );
0100     }
0101     
0102     
0103     auto t1 = std::chrono::high_resolution_clock::now();
0104     std::cout << "Initiation time: " << std::round(std::chrono::duration<double, std::chrono::seconds::period>(t1 - t0).count()) << " sec" << std::endl;
0105     std::cout << "\n==================================================================\n" << std::endl;
0106 
0107   }
0108 
0109   // Helper to parse raw strings into BackgroundConfig structs
0110   std::vector<BackgroundConfig>
0111   parse_backgrounds(const std::vector<std::string> &raw_args_list) {
0112     std::vector<BackgroundConfig> backgrounds;
0113     auto is_pure_integer = [](const std::string &str) {
0114       if (str.empty()) return false;
0115       for (char c : str) {
0116         if (!std::isdigit(c)) return false;
0117       }
0118       return true;
0119     };
0120       
0121     // Group strings into sets of 2-4 arguments per background
0122     for (size_t i = 0; i < raw_args_list.size();) {
0123       // Determine how many arguments this background has
0124       size_t args_count = 2; // minimum
0125 
0126       // Look ahead to see if next strings can be parsed as numbers (skip/status)
0127       if (i + 2 < raw_args_list.size() && is_pure_integer(raw_args_list[i + 2])) {
0128           args_count = 3;
0129 
0130           if (i + 3 < raw_args_list.size() && is_pure_integer(raw_args_list[i + 3])) {
0131               args_count = 4;
0132           }
0133       }
0134 
0135       // Ensure we don't go beyond the vector bounds
0136       if (i + args_count > raw_args_list.size()) {
0137         args_count = raw_args_list.size() - i;
0138       }
0139 
0140       if (args_count < 2) {
0141         throw std::runtime_error("Background file " +
0142                                  std::to_string(backgrounds.size()) +
0143                                  " must have at least 2 arguments");
0144       }
0145 
0146       try {
0147         BackgroundConfig bg;
0148         bg.file = raw_args_list[i];
0149         bg.frequency = std::stod(raw_args_list[i + 1]);
0150         bg.skip = (args_count > 2) ? std::stoi(raw_args_list[i + 2]) : 0;
0151         bg.status = (args_count > 3) ? std::stoi(raw_args_list[i + 3]) : 0;
0152         backgrounds.push_back(bg);
0153       } catch (const std::exception &e) {
0154         throw std::runtime_error("Error parsing background file " +
0155                                  std::to_string(backgrounds.size()) + ": " +
0156                                  e.what());
0157       }
0158 
0159       i += args_count;
0160     }
0161 
0162     return backgrounds;
0163   }
0164 
0165   void merge(){
0166     auto t1 = std::chrono::high_resolution_clock::now();
0167 
0168     // Open output file
0169     std::shared_ptr<HepMC3::Writer> f;
0170     if (rootFormat)
0171       f = std::make_shared<HepMC3::WriterRootTree>(outputFileName);
0172     else
0173       f = std::make_shared<HepMC3::WriterAscii>(outputFileName);
0174     
0175     // Slice loop
0176     int i = 0;
0177     for (i = 0; i<nSlices; ++i ) {
0178       if (i % 1000 == 0 || verbose ) squawk(i);
0179       auto hepSlice = mergeSlice(i);
0180       if (!hepSlice) {
0181     std::cout << "Exhausted signal source." << std::endl;
0182     break;
0183       }
0184       hepSlice->set_event_number(i);
0185       f->write_event(*hepSlice);
0186     }
0187     std::cout << "Finished all requested slices." << std::endl;
0188 
0189     int slicesDone = i;
0190     auto t2 = std::chrono::high_resolution_clock::now();
0191 
0192     std::cout << "Slice loop time: " << std::round(std::chrono::duration<double, std::chrono::minutes::period>(t2 - t1).count()) << " min" << std::endl;
0193     std::cout << " -- " << std::round(std::chrono::duration<double, std::chrono::microseconds::period>(t2 - t1).count() / i) << " us / slice" << std::endl;
0194 
0195     for (auto info : infoDict) {
0196       std::cout << "From " << info.first << std::endl;
0197       std::cout << "  placed " << info.second.eventCount << " events" << std::endl; 
0198       std::cout << "  --> on average " << std::setprecision(3) << info.second.eventCount / float(nSlices) << std::endl;
0199       std::cout << "  placed " << info.second.particleCount << " final state particles" << std::endl;
0200       std::cout << "  --> on average " << std::setprecision(3) << info.second.particleCount / float(nSlices) << std::endl;
0201       
0202     }
0203 
0204     struct rusage r_usage;
0205     getrusage(RUSAGE_SELF, &r_usage);
0206 
0207     // NOTE: Reported in kB on Linux, bytes in Mac/Darwin
0208     // Could try to explicitly catch __linux__ as well
0209     // Unclear in BSD, I've seen conflicting reports
0210 #ifdef __MACH__
0211     float mbsize = 1024 * 1024;
0212 #else // Linux
0213     float mbsize = 1024;
0214 #endif
0215   
0216 
0217     std::cout << endl << "Maximum Resident Memory " << r_usage.ru_maxrss / mbsize << " MB" << std::endl;
0218     // clean up, close all files
0219     sigAdapter->close();
0220     for (auto& it : freqAdapters) {
0221       it.second->close();
0222     }
0223     f->close();
0224     
0225   }
0226   
0227   // ---------------------------------------------------------------------------
0228   void digestArgs(int argc, char* argv[]) {
0229     // Handle the command line tedium
0230     // ArgumentParser is meant to be used in a single function.
0231     // ArgumentParser internally uses std::string_views,
0232     // references, iterators, etc.
0233     // Many of these elements become invalidated after a copy or move.
0234     argparse::ArgumentParser args ("Merge signal events with up to four background sources.");
0235     
0236     args.add_argument("-i", "--signalFile")
0237       .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"))
0238       .help("Name of the HEPMC file with the signal events");
0239     
0240     args.add_argument("-sf", "--signalFreq")
0241       .default_value(0.0)
0242       .scan<'g', double>()
0243       .help("Signal frequency in kHz. Default is 0 to have exactly one signal event per slice. Set to the estimated DIS rate to randomize.");
0244 
0245     args.add_argument("-S", "--signalSkip")
0246       .default_value(0)
0247     .scan<'i', int>()
0248     .help("Number of signals events to skip. Default is 0.");
0249 
0250     args.add_argument("-St", "--signalStatus")
0251       .default_value(0)
0252     .scan<'i', int>()
0253     .help("Apply shift on particle generatorStatus code for signal. Default is 0. ");
0254 
0255     args.add_argument("-b","--bgFile")
0256       .nargs(2,4)
0257       .append()
0258       .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).");
0259     
0260     args.add_argument("-o", "--outputFile")
0261       .default_value(std::string("bgmerged.hepmc3.tree.root"))
0262       .help("Specify the output file name. By default bgmerged.hepmc3.tree.root is used");
0263 
0264     args.add_argument("-r", "--rootFormat")
0265       .default_value(true)
0266       .implicit_value(true)
0267       .help("Use hepmc.root output format, default is true.");
0268     
0269     args.add_argument("-w", "--intWindow")
0270       .default_value(2000.0)
0271       .scan<'g', double>()
0272       .help("Length of the integration window in nanoseconds. Default is 2000.");
0273     
0274     args.add_argument("-N", "--nSlices")
0275       .default_value(10000)
0276       .scan<'i', int>()
0277       .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.");
0278     
0279     args.add_argument("--squashTime")
0280       .default_value(false)
0281       .implicit_value(true)
0282       .help("Integration is performed but no time information is associated to vertices.");
0283     
0284     args.add_argument("--rngSeed")
0285       .default_value(0)
0286       .action([](const std::string& value) { return std::stoi(value); })
0287       .help("Random seed, default is None");
0288     
0289     args.add_argument("-v", "--verbose")
0290       .default_value(false)
0291       .implicit_value(true)
0292       .help("Display details for every slice.");
0293     
0294     try {
0295       args.parse_args(argc, argv);
0296     }
0297     catch (const std::runtime_error& err) {
0298       std::cout << err.what() << std::endl;
0299       std::cout << args;
0300       exit(EXIT_FAILURE);
0301     }
0302     // Access arguments using args.get method
0303     signalFile = args.get<std::string>("--signalFile");
0304     signalFreq = args.get<double>("--signalFreq");
0305     signalSkip = args.get<int>("--signalSkip");
0306     signalStatus = args.get<int>("--signalStatus");
0307     backgroundFiles = parse_backgrounds(args.get<std::vector<std::string>>("--bgFile"));
0308     outputFile = args.get<std::string>("--outputFile");
0309     rootFormat = args.get<bool>("--rootFormat");
0310     intWindow  = args.get<double>("--intWindow");
0311     nSlices    = args.get<int>("--nSlices");
0312     squashTime = args.get<bool>("--squashTime");
0313     rngSeed    = args.get<int>("--rngSeed");
0314     verbose    = args.get<bool>("--verbose");
0315 
0316     
0317   }
0318   
0319   // ---------------------------------------------------------------------------
0320   void banner() {
0321     // Print banner
0322     std::cout << "==================================================================" << std::endl;
0323     std::cout << "=== EPIC HEPMC MERGER ===" << std::endl;
0324     std::cout << "authors: Benjamen Sterwerf* (bsterwerf@berkeley.edu), Kolja Kauder** (kkauder@bnl.gov), Reynier Cruz-Torres***" << std::endl;
0325     std::cout << "* University of California, Berkeley" << std::endl;
0326     std::cout << "** Brookhaven National Laboratory" << std::endl;
0327     std::cout << "*** formerly Lawrence Berkeley National Laboratory" << std::endl;
0328     std::cout << "\nFor more information, run \n./signal_background_merger --help" << std::endl;
0329 
0330     std::string statusMessage = "Shifting all particle status codes from this source by ";
0331     std::vector<int> statusList_stable, statusList_decay;
0332 
0333     std::cout << "Number of Slices:" << nSlices << endl;
0334     std::string freqTerm = signalFreq > 0 ? std::to_string(signalFreq) + " kHz" : "(one event per time slice)";
0335     std::string statusTerm = signalStatus > 0 ? statusMessage + std::to_string(signalStatus): "";
0336     if (signalStatus>0){
0337       statusList_stable.push_back(signalStatus+1);
0338       statusList_decay.push_back(signalStatus+2);
0339     }
0340     std::cout << "Signal events file and frequency:\n";
0341     std::cout << "\t- " << signalFile << "\t" << freqTerm << "\n" << statusTerm << "\n";
0342     
0343     std::cout << "\nBackground files and their respective frequencies:\n";
0344 
0345     for (const auto& bg : backgroundFiles) {
0346       if (!bg.file.empty()) {
0347         freqTerm = bg.frequency > 0 ? std::to_string(bg.frequency) + " kHz" : "(from weights)";
0348         statusTerm = bg.status > 0 ? statusMessage + std::to_string(bg.status) : "";
0349         std::cout << "\t- " << bg.file << "\t" << freqTerm << "\n" << statusTerm << "\n";
0350         if (bg.status>0){
0351           statusList_stable.push_back(bg.status+1);
0352           statusList_decay.push_back(bg.status+2);
0353         }
0354       }
0355     }
0356     
0357     auto join = [](const std::vector<int>& vec) {
0358         return std::accumulate(vec.begin(), vec.end(), std::string(),
0359             [](const std::string& a, int b) {
0360                 return a.empty() ? std::to_string(b) : a + " " + std::to_string(b);
0361             });
0362     };
0363     std::string stableStatuses = join(statusList_stable);
0364     std::string decayStatuses = join(statusList_decay);
0365     std::string message = "\n!!!Attention!!!\n To proceed the shifted particles statuses in DD4hep, please add the following options to ddsim:\n"
0366                           "--physics.alternativeStableStatuses=\"" + stableStatuses + 
0367                           "\"  --physics.alternativeDecayStatuses=\"" + decayStatuses + "\"\n";
0368     std::cout << message<<std::endl;           
0369   }
0370   
0371   // ---------------------------------------------------------------------------  
0372   void PrepData(const std::string& fileName, double freq, int skip=0, int baseStatus=0, bool signal=false) {
0373     if (fileName.empty()) return;
0374 
0375     cout << "Prepping " << fileName << endl;
0376     std::shared_ptr<HepMC3::Reader> adapter;
0377     try {
0378       adapter = HepMC3::deduce_reader(fileName);
0379       if (!adapter) {
0380         throw std::runtime_error("Failed to open file");
0381       }
0382     } catch (const std::runtime_error& e) {
0383       std::cerr << "Opening " << fileName << " failed: " << e.what() << std::endl;
0384       exit(EXIT_FAILURE);
0385     }
0386     
0387     infoDict[fileName] = {0,0};
0388 
0389     if (signal) {
0390       sigAdapter = adapter;
0391       sigFreq = freq;
0392       sigStatus = baseStatus;
0393       sigAdapter->skip(skip);
0394       return;
0395     }
0396 
0397     // Now catch the weighted case
0398     if (freq <= 0) {
0399       std::cout << "Reading in all events from " << fileName << std::endl;
0400       std::vector<HepMC3::GenEvent> events;
0401       std::vector<double> weights;
0402 
0403       while(!adapter->failed()) {
0404         HepMC3::GenEvent evt(HepMC3::Units::GEV,HepMC3::Units::MM);
0405         adapter->read_event(evt);
0406 
0407         // remove events with 0 weight - note that this does change avgRate = <weight> (by a little)
0408         if (double w=evt.weight() > 0){
0409           events.push_back(evt);
0410           weights.push_back(evt.weight());
0411           }
0412       }
0413       adapter->close();
0414       
0415       double avgRate = 0.0;
0416       for ( auto w : weights ){ avgRate += w;}
0417       avgRate /= weights.size();
0418       avgRate *= 1e-9; // convert to 1/ns == GHz
0419       std::cout << "Average rate is " << avgRate << " GHz" << std::endl;
0420 
0421       std::vector<int> indices (weights.size());
0422       std::iota (std::begin(indices), std::end(indices), 0); // [ 0 , ... , N ] <- draw randomly from this
0423       
0424       // Replacing python's compact toPlace = self.rng.choice( a=events, size=nEvents, p=probs, replace=False )
0425       // is tricky. Possibly more elegant or faster versions exist,
0426       // https://stackoverflow.com/questions/42926209/equivalent-function-to-numpy-random-choice-in-c
0427       // we'll do it rather bluntly, since the need for this code should go away soon with new SR infrastructure
0428       // https://stackoverflow.com/questions/1761626/weighted-random-numbers
0429       // Normalizing is not necessary for this method 
0430       // for ( auto& w : weights ) {
0431       //    w /= avgRate;
0432       // }
0433       std::piecewise_constant_distribution<> weightedDist(std::begin(indices),std::end(indices),
0434                               std::begin(weights));
0435       weightDict[fileName] = { std::make_tuple(events, weightedDist, avgRate) };
0436 
0437       return;
0438     }
0439 
0440     // Not signal and not weighted --> prepare frequency backgrounds
0441     adapter->skip(skip);
0442     freqAdapters[fileName] = adapter;
0443     freqs[fileName] = freq;
0444     baseStatuses[fileName] = baseStatus;
0445   }
0446 
0447    // ---------------------------------------------------------------------------
0448   bool hasEnding (std::string const &fullString, std::string const &ending) {
0449     if (fullString.length() >= ending.length()) {
0450       return (0 == fullString.compare (fullString.length() - ending.length(), ending.length(), ending));
0451     } else {
0452       return false;
0453     }
0454   }
0455   
0456   // ---------------------------------------------------------------------------
0457   std::string nameGen() {
0458     // Generate a name for the output file
0459     // It's too simplistic for input with directories
0460     std::string name = signalFile;
0461     if (nSlices > 0) {
0462         size_t pos = name.find(".hepmc");
0463         if (pos != std::string::npos) {
0464       name.replace(pos, 6, "_n_" + std::to_string(nSlices) + ".hepmc");
0465         }
0466     }
0467 
0468     if ( rootFormat && !hasEnding(name,".root")){
0469       name.append(".root");
0470     }
0471     name = "bgmerged_" + name;
0472 
0473     return name;
0474   }
0475 
0476   // ---------------------------------------------------------------------------
0477   void squawk(int i) {
0478 
0479     // More fine-grained info about current usage
0480 #ifdef __MACH__
0481     task_basic_info_data_t info;
0482     mach_msg_type_number_t size = sizeof(info);
0483     kern_return_t kerr = task_info(mach_task_self(),
0484                                    TASK_BASIC_INFO,
0485                                    (task_info_t)&info,
0486                                    &size);
0487 
0488     long memory_usage = -1;
0489     if (kerr == KERN_SUCCESS) {
0490       memory_usage = info.resident_size  / 1024 / 1024;
0491     }
0492 #else // Linux
0493     std::ifstream statm("/proc/self/statm");
0494     long size, resident, share, text, lib, data, dt;
0495     statm >> size >> resident >> share >> text >> lib >> data >> dt;
0496     statm.close();
0497 
0498     long page_size = sysconf(_SC_PAGESIZE);  // in case x86-64 is configured to use 2MB pages
0499     long memory_usage = resident * page_size  / 1024 / 1024 ;    
0500 #endif
0501   
0502     
0503     std::cout << "Working on slice " << i + 1 << std::endl;
0504     std::cout << "Current memory usage: " << memory_usage << " MB" << std::endl;
0505 
0506   }
0507   // ---------------------------------------------------------------------------
0508 
0509   std::unique_ptr<HepMC3::GenEvent> mergeSlice(int i) {
0510     auto hepSlice = std::make_unique<HepMC3::GenEvent>(HepMC3::Units::GEV, HepMC3::Units::MM);
0511     
0512     addFreqEvents(signalFile, sigAdapter, sigFreq, hepSlice, signalStatus, true);
0513     
0514     for (const auto& freqBgs : freqAdapters) {
0515       auto fileName=freqBgs.first;
0516       addFreqEvents(fileName, freqAdapters[fileName], freqs[fileName], hepSlice, baseStatuses[fileName], false);
0517     }
0518     
0519     for (const auto& fileName : weightDict) {
0520       addWeightedEvents(fileName.first, hepSlice, baseStatuses[fileName.first]);
0521     }
0522 
0523     return hepSlice;
0524   };
0525 
0526   // ---------------------------------------------------------------------------
0527 
0528   void addFreqEvents(std::string fileName, std::shared_ptr<HepMC3::Reader>& adapter, const double freq,
0529              std::unique_ptr<HepMC3::GenEvent>& hepSlice, int baseStatus = 0, bool signal = false) {
0530 
0531     // First, create a timeline
0532     // Signals can be different
0533     std::vector<double> timeline;
0534 
0535     std::uniform_real_distribution<> uni(0, intWindow);
0536     if (freq == 0){
0537       if (!signal) {
0538         std::cerr << "frequency can't be 0 for background files" << std::endl;
0539         exit(EXIT_FAILURE);
0540       }
0541       // exactly one signal event, at an arbitrary point
0542       timeline.push_back(uni(rng));
0543     } else {
0544       // Generate poisson-distributed times to place events
0545       timeline = poissonTimes(freq, intWindow);
0546     }
0547     
0548     if ( verbose) std::cout << "Placing " << timeline.size() << " events from " << fileName << std::endl;
0549 
0550     if (timeline.empty()) return;
0551     long particleCount = 0;
0552 
0553     // Insert events at all specified locations
0554     for (double time : timeline) {
0555       if(adapter->failed()) {
0556       try{
0557         if (signal) { // Exhausted signal events
0558             throw std::ifstream::failure("EOF");
0559           } else { // background file reached its end, reset to the start
0560             std::cout << "Cycling back to the start of " << fileName << std::endl;
0561             adapter->close();
0562             adapter = HepMC3::deduce_reader(fileName);
0563           }
0564         } catch (std::ifstream::failure& e) {
0565             continue; // just need to suppress the error
0566         }
0567       }
0568 
0569       HepMC3::GenEvent inevt;
0570       adapter->read_event(inevt);
0571       if (signal && (signalFreq == 0.0)){
0572         hepSlice->weights() = inevt.weights();
0573       }
0574 
0575       if (squashTime) time = 0;
0576       particleCount += insertHepmcEvent( inevt, hepSlice, time, baseStatus, signal);
0577     }
0578 
0579     infoDict[fileName].eventCount += timeline.size();
0580     infoDict[fileName].particleCount += particleCount;
0581 
0582 
0583     return;
0584   }
0585 
0586   // ---------------------------------------------------------------------------
0587 
0588   void addWeightedEvents(std::string fileName, std::unique_ptr<HepMC3::GenEvent>& hepSlice, int baseStatus=0, bool signal = false) {
0589     auto& [events, weightedDist, avgRate ] = weightDict[fileName];
0590 
0591     // How many events? Assume Poisson distribution
0592     int nEvents;
0593     std::poisson_distribution<> d( intWindow * avgRate );
0594 
0595     // Small SR files may not have enough photons (example or test files). Could use them all or reroll
0596     // Choosing the latter, neither is physical
0597     while (true) {
0598       nEvents = d(rng);
0599       if (nEvents > events.size()) {
0600           std::cout << "WARNING: Trying to place " << nEvents << " events from " << fileName
0601                   << " but the file doesn't have enough. Rerolling, but this is not physical." << std::endl;
0602         continue;
0603       }
0604       break;
0605     }
0606 
0607     if (verbose) std::cout << "Placing " << nEvents << " events from " << fileName << std::endl;
0608     
0609     // Get randomized event indices
0610     // Note: Could change to drawing without replacing ( if ( not in toPLace) ...) , not worth the effort
0611     std::vector<HepMC3::GenEvent> toPlace(nEvents);
0612     for ( auto& e : toPlace ){
0613       auto i = static_cast<int> (weightedDist(rng));
0614       e = events.at(i);
0615     }
0616     
0617     // Place at random times
0618     std::vector<double> timeline;
0619     std::uniform_real_distribution<> uni(0, intWindow);
0620     long particleCount = 0;
0621     if (!squashTime) {
0622       for ( auto& e : toPlace ){
0623           double time = squashTime ? 0 : uni(rng);
0624         particleCount += insertHepmcEvent( e, hepSlice, time, baseStatus, signal);
0625       }
0626     }
0627 
0628     infoDict[fileName].eventCount += nEvents;
0629     infoDict[fileName].particleCount += particleCount;
0630 
0631     return;
0632 }
0633 
0634   // ---------------------------------------------------------------------------
0635   long insertHepmcEvent( const HepMC3::GenEvent& inevt,
0636              std::unique_ptr<HepMC3::GenEvent>& hepSlice, double time=0, int baseStatus=0, bool signal = false) {
0637     // Unit conversion
0638     double timeHepmc = c_light * time;
0639     
0640     std::vector<HepMC3::GenParticlePtr> particles;
0641     std::vector<HepMC3::GenVertexPtr> vertices;
0642 
0643     // Stores the vertices of the event inside a vertex container. These vertices are in increasing order
0644     // so we can index them with [abs(vertex_id)-1]
0645     for (auto& vertex : inevt.vertices()) {
0646       HepMC3::FourVector position = vertex->position();
0647       position.set_t(position.t() + timeHepmc);
0648       auto v1 = std::make_shared<HepMC3::GenVertex>(position);
0649       vertices.push_back(v1);
0650     }
0651       
0652     // copies the particles and attaches them to their corresponding vertices
0653     long finalParticleCount = 0;
0654     for (auto& particle : inevt.particles()) {
0655       HepMC3::FourVector momentum = particle->momentum();
0656       int status = particle->status();
0657       if (status == 1 ) finalParticleCount++;
0658       int pid = particle->pid();
0659       status += baseStatus;
0660       auto p1 = std::make_shared<HepMC3::GenParticle> (momentum, pid, status);
0661       p1->set_generated_mass(particle->generated_mass());
0662       particles.push_back(p1);
0663       // since the beam particles do not have a production vertex they cannot be attached to a production vertex
0664       if (particle->production_vertex()->id() < 0) {
0665           int production_vertex = particle->production_vertex()->id();
0666           vertices[abs(production_vertex) - 1]->add_particle_out(p1);
0667           hepSlice->add_particle(p1);
0668       }
0669     
0670       // Adds particles with an end vertex to their end vertices
0671       if (particle->end_vertex()) {
0672           int end_vertex = particle->end_vertex()->id();
0673           vertices.at(abs(end_vertex) - 1)->add_particle_in(p1);    
0674       }
0675     }
0676 
0677     // Adds the vertices with the attached particles to the event
0678     for (auto& vertex : vertices) {
0679       hepSlice->add_vertex(vertex);
0680     }
0681     
0682     return finalParticleCount;
0683   }
0684 
0685   // ---------------------------------------------------------------------------
0686 
0687   std::vector<double> poissonTimes(double mu, double endTime) {
0688     std::exponential_distribution<> exp(mu);
0689     
0690     double t = 0;
0691     std::vector<double> ret;
0692     while (true) {
0693       double delt = exp(rng)*1e6;
0694       // cout << delt <<endl;
0695       t += delt;
0696       if (t >= endTime) {
0697     break;
0698       }
0699       ret.push_back(t);
0700     }
0701     return ret;
0702 }
0703   // ---------------------------------------------------------------------------
0704 
0705   
0706   // private:
0707   std::mt19937 rng;
0708   string signalFile;
0709   double signalFreq;
0710   int signalSkip;
0711   int signalStatus;
0712   std::vector<BackgroundConfig> backgroundFiles;  
0713   string outputFile;
0714   string outputFileName;
0715   bool rootFormat;
0716   double intWindow;
0717   int nSlices; // should be long, but argparse cannot read that
0718   bool squashTime;
0719   int rngSeed;  // should be unsigned, but argparse cannot read that
0720   bool verbose;
0721   
0722   const double c_light = 299.792458; // speed of light = 299.792458 mm/ns to get mm  
0723 };
0724 
0725 // =============================================================
0726 int main(int argc, char* argv[]) {
0727 
0728   auto t0 = std::chrono::high_resolution_clock::now();
0729   // Create an instance of SignalBackgroundMerger
0730   SignalBackgroundMerger sbm (argc, argv);
0731 
0732 
0733   sbm.merge();
0734 
0735   std::cout << "\n==================================================================\n";
0736   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;
0737   
0738   return EXIT_SUCCESS;
0739 }