Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-03-28 07:46:03

0001 // This file is part of the ACTS project.
0002 //
0003 // Copyright (C) 2016 CERN for the benefit of the ACTS project
0004 //
0005 // This Source Code Form is subject to the terms of the Mozilla Public
0006 // License, v. 2.0. If a copy of the MPL was not distributed with this
0007 // file, You can obtain one at https://mozilla.org/MPL/2.0/.
0008 
0009 #include "ActsExamples/Io/Root/detail/NuclearInteractionParametrisation.hpp"
0010 
0011 #include "Acts/Definitions/Common.hpp"
0012 #include "Acts/Utilities/MathHelpers.hpp"
0013 
0014 #include <cmath>
0015 #include <iterator>
0016 #include <limits>
0017 
0018 #include <Eigen/Eigenvalues>
0019 #include <TMath.h>
0020 #include <TVectorF.h>
0021 #include <TVectorT.h>
0022 
0023 namespace ActsExamples::detail::NuclearInteractionParametrisation {
0024 
0025 namespace {
0026 
0027 /// @brief Evaluate the location in a standard normal distribution for a value
0028 /// from a probability distribution
0029 ///
0030 /// @param [in] histo The probability distribution
0031 /// @param [in] mom The abscissa value in @p histo
0032 ///
0033 /// @return The location in a standard normal distribution
0034 float gaussianValue(TH1F const* histo, const float mom) {
0035   // Get the cumulative probability distribution
0036 
0037   // DrawNormalized and GetCumulative return pointers to histograms where ROOT
0038   // transfers ownership to the caller.
0039   std::unique_ptr<TH1F> normalised(
0040       dynamic_cast<TH1F*>(histo->DrawNormalized()));
0041   std::unique_ptr<TH1F> cumulative(
0042       dynamic_cast<TH1F*>(normalised->GetCumulative()));
0043   assert(cumulative);
0044   assert(normalised);
0045   // Find the cumulative probability
0046   const float binContent = cumulative->GetBinContent(cumulative->FindBin(mom));
0047   // Transform the probability to an entry in a standard normal distribution
0048   const float value = TMath::ErfInverse(2. * binContent - 1.);
0049 
0050   return value;
0051 }
0052 
0053 /// @brief Evaluate the invariant mass of two four vectors
0054 ///
0055 /// @param [in] fourVector1 The one four vector
0056 /// @param [in] fourVector2 The other four vector
0057 ///
0058 /// @return The invariant mass
0059 float invariantMass(const Acts::Vector4& fourVector1,
0060                     const Acts::Vector4& fourVector2) {
0061   const Acts::Vector4 sum = fourVector1 + fourVector2;
0062   const double energy = sum[Acts::eEnergy];
0063   const double momentum = sum.template segment<3>(Acts::eMom0).norm();
0064   return Acts::fastCathetus(energy, momentum);
0065 }
0066 
0067 }  // namespace
0068 
0069 std::pair<Vector, Matrix> calculateMeanAndCovariance(
0070     unsigned int multiplicity, const EventProperties& events) {
0071   // Calculate the mean
0072   Vector mean = Vector::Zero(multiplicity);
0073   for (const std::vector<float>& event : events) {
0074     for (unsigned int j = 0; j < multiplicity; j++) {
0075       mean[j] += event[j];
0076     }
0077   }
0078   mean /= events.size();
0079 
0080   // Calculate the covariance matrix
0081   Matrix covariance = Matrix::Zero(multiplicity, multiplicity);
0082   for (unsigned int i = 0; i < multiplicity; i++) {
0083     for (unsigned int j = 0; j < multiplicity; j++) {
0084       for (unsigned int k = 0; k < events.size(); k++) {
0085         covariance(i, j) += (events[k][i] - mean[i]) * (events[k][j] - mean[j]);
0086       }
0087     }
0088   }
0089   covariance /= events.size();
0090 
0091   return {mean, covariance};
0092 }
0093 
0094 EigenspaceComponents calculateEigenspace(const Vector& mean,
0095                                          const Matrix& covariance) {
0096   // Calculate eigenvalues and eigenvectors
0097   Eigen::EigenSolver<Matrix> es(covariance);
0098   Vector eigenvalues = es.eigenvalues().real();
0099   Matrix eigenvectors = es.eigenvectors().real();
0100   // Transform the mean vector into eigenspace
0101   Vector meanEigenspace = eigenvectors * mean;
0102 
0103   return {eigenvalues, eigenvectors, meanEigenspace};
0104 }
0105 
0106 Parametrisation buildMomentumParameters(const EventCollection& events,
0107                                         unsigned int multiplicity, bool soft,
0108                                         unsigned int nBins) {
0109   // Strip off data
0110   auto momenta = prepareMomenta(events, multiplicity, soft);
0111   if (momenta.empty()) {
0112     return Parametrisation();
0113   }
0114 
0115   // Build histos
0116   ProbabilityDistributions histos = buildMomPerMult(momenta, nBins);
0117 
0118   // Build normal distribution
0119   auto momentaGaussian = convertEventToGaussian(histos, momenta);
0120   auto meanAndCovariance =
0121       calculateMeanAndCovariance(multiplicity + 1, momentaGaussian);
0122   // Calculate the transformation into the eigenspace of the covariance matrix
0123   EigenspaceComponents eigenspaceElements =
0124       calculateEigenspace(meanAndCovariance.first, meanAndCovariance.second);
0125   // Calculate the cumulative distributions
0126   return {eigenspaceElements, histos};
0127 }
0128 
0129 EventProperties prepareMomenta(const EventCollection& events,
0130                                unsigned int multiplicity,
0131                                bool soft)  // TODO: build enum instead of bool
0132 {
0133   EventProperties result;
0134   // Loop over all events
0135   for (const EventFraction& event : events) {
0136     // Test the multiplicity and type of the event
0137     if (event.multiplicity == multiplicity && event.soft == soft) {
0138       const float initialMomentum = event.initialParticle.absoluteMomentum();
0139       float sum = 0.;
0140       std::vector<float> momenta;
0141       momenta.reserve(multiplicity + 1);
0142       // Fill the vector with the scaled momenta
0143       for (const SimParticle& p : event.finalParticles) {
0144         sum += p.absoluteMomentum();
0145         momenta.push_back(p.absoluteMomentum() / initialMomentum);
0146       }
0147       // Add the scaled sum of momenta
0148       momenta.push_back(sum / initialMomentum);
0149       result.push_back(std::move(momenta));
0150     }
0151   }
0152   return result;
0153 }
0154 
0155 ProbabilityDistributions buildMomPerMult(const EventProperties& events,
0156                                          unsigned int nBins) {
0157   // Fast exit
0158   if (events.empty()) {
0159     return {};
0160   }
0161   const unsigned int multMax = events[0].size();
0162 
0163   // Find the range of each histogram
0164   std::vector<float> min(multMax, std::numeric_limits<float>::max());
0165   std::vector<float> max(multMax, 0);
0166   for (const std::vector<float>& event : events) {
0167     for (unsigned int i = 0; i < multMax; i++) {
0168       min[i] = std::min(event[i], min[i]);
0169       max[i] = std::max(event[i], max[i]);
0170     }
0171   }
0172 
0173   // Evaluate the range of the histograms
0174   // This is used to avoid entries in over-/underflow bins
0175   std::vector<float> diff(multMax);
0176   for (unsigned int i = 0; i < multMax; i++) {
0177     diff[i] = (max[i] - min[i]) * 0.1;
0178   }
0179 
0180   // Build the histograms
0181   ProbabilityDistributions histos(multMax);
0182   for (unsigned int i = 0; i < multMax; i++) {
0183     histos[i] = new TH1F("", "", nBins, min[i] - diff[i], max[i] + diff[i]);
0184   }
0185 
0186   // Fill the histograms
0187   for (const std::vector<float>& event : events) {
0188     for (unsigned int i = 0; i < multMax; i++) {
0189       histos[i]->Fill(event[i]);
0190     }
0191   }
0192   return histos;
0193 }
0194 
0195 EventProperties convertEventToGaussian(const ProbabilityDistributions& histos,
0196                                        const EventProperties& events) {
0197   // Fast exit
0198   if (events.empty()) {
0199     return {};
0200   }
0201   const unsigned int multMax = events[0].size();
0202 
0203   // Loop over the events
0204   EventProperties gaussianEvents;
0205   for (const std::vector<float>& event : events) {
0206     // Transform the properties in the events
0207     std::vector<float> gaussianEvent;
0208     for (unsigned int i = 0; i < multMax; i++) {
0209       gaussianEvent.push_back(gaussianValue(histos[i], event[i]));
0210     }
0211     // Store the transformed event
0212     gaussianEvents.push_back(gaussianEvent);
0213   }
0214   return gaussianEvents;
0215 }
0216 
0217 EventProperties prepareInvariantMasses(const EventCollection& events,
0218                                        unsigned int multiplicity, bool soft) {
0219   EventProperties result;
0220   // Loop over all events
0221   for (const EventFraction& event : events) {
0222     // Test the multiplicity and type of the event
0223     if (event.multiplicity == multiplicity && event.soft == soft) {
0224       const auto fourVectorBefore = event.interactingParticle.fourMomentum();
0225       std::vector<float> invariantMasses;
0226       invariantMasses.reserve(multiplicity);
0227       // Fill the vector with the invariant masses
0228       for (const SimParticle& p : event.finalParticles) {
0229         const auto fourVector = p.fourMomentum();
0230         invariantMasses.push_back(invariantMass(fourVectorBefore, fourVector));
0231       }
0232       result.push_back(invariantMasses);
0233     }
0234   }
0235   return result;
0236 }
0237 
0238 Parametrisation buildInvariantMassParameters(const EventCollection& events,
0239                                              unsigned int multiplicity,
0240                                              bool soft, unsigned int nBins) {
0241   // Strip off data
0242   auto invariantMasses = prepareInvariantMasses(events, multiplicity, soft);
0243   if (invariantMasses.empty()) {
0244     return Parametrisation();
0245   }
0246 
0247   // Build histos
0248   ProbabilityDistributions histos = buildMomPerMult(invariantMasses, nBins);
0249 
0250   // Build normal distribution
0251   auto invariantMassesGaussian =
0252       convertEventToGaussian(histos, invariantMasses);
0253   auto meanAndCovariance =
0254       calculateMeanAndCovariance(multiplicity, invariantMassesGaussian);
0255   // Calculate the transformation into the eigenspace of the covariance matrix
0256   EigenspaceComponents eigenspaceElements =
0257       calculateEigenspace(meanAndCovariance.first, meanAndCovariance.second);
0258   // Calculate the cumulative distributions
0259   return {eigenspaceElements, histos};
0260 }
0261 
0262 std::unordered_map<int, std::unordered_map<int, float>>
0263 cumulativePDGprobability(const EventCollection& events) {
0264   std::unordered_map<int, std::unordered_map<int, float>> counter;
0265 
0266   // Count how many and which particles were created by which particle
0267   for (const EventFraction& event : events) {
0268     if (event.finalParticles.empty()) {
0269       continue;
0270     }
0271     if (!event.soft) {
0272       counter[event.initialParticle.pdg()][event.finalParticles[0].pdg()]++;
0273     }
0274     for (unsigned int i = 1; i < event.multiplicity; i++) {
0275       counter[event.finalParticles[i - 1].pdg()]
0276              [event.finalParticles[i].pdg()]++;
0277     }
0278   }
0279 
0280   // Build a cumulative distribution
0281   for (const auto& element : counter) {
0282     float sum = 0;
0283     auto prevIt = counter[element.first].begin();
0284     for (auto it1 = counter[element.first].begin();
0285          it1 != counter[element.first].end(); it1++) {
0286       float binEntry = 0;
0287       if (it1 == counter[element.first].begin()) {
0288         binEntry = it1->second;
0289         prevIt = it1;
0290       } else {
0291         binEntry = it1->second - prevIt->second;
0292         prevIt = it1;
0293       }
0294       // Add content to next bins
0295       for (auto it2 = std::next(it1, 1); it2 != counter[element.first].end();
0296            it2++) {
0297         it2->second += binEntry;
0298         sum = it2->second;
0299       }
0300     }
0301     // Normalise the entry
0302     for (auto it1 = counter[element.first].begin();
0303          it1 != counter[element.first].end(); it1++) {
0304       it1->second /= sum;
0305     }
0306   }
0307   return counter;
0308 }
0309 
0310 std::pair<CumulativeDistribution, CumulativeDistribution>
0311 cumulativeMultiplicityProbability(const EventCollection& events,
0312                                   unsigned int multiplicityMax) {
0313   // Find the range of both histogram
0314   unsigned int minSoft = std::numeric_limits<unsigned int>::max();
0315   unsigned int maxSoft = 0;
0316   unsigned int minHard = std::numeric_limits<unsigned int>::max();
0317   unsigned int maxHard = 0;
0318   for (const EventFraction& event : events) {
0319     if (event.soft) {
0320       minSoft = std::min(event.multiplicity, minSoft);
0321       maxSoft = std::max(event.multiplicity, maxSoft);
0322     } else {
0323       minHard = std::min(event.multiplicity, minHard);
0324       maxHard = std::max(event.multiplicity, maxHard);
0325     }
0326   }
0327 
0328   // Build and fill the histograms
0329   TH1F* softHisto =
0330       new TH1F("", "", std::min(maxSoft, multiplicityMax) + 1 - minSoft,
0331                minSoft, std::min(maxSoft, multiplicityMax) + 1);
0332   TH1F* hardHisto =
0333       new TH1F("", "", std::min(maxHard, multiplicityMax) + 1 - minHard,
0334                minHard, std::min(maxHard, multiplicityMax) + 1);
0335   for (const EventFraction& event : events) {
0336     if (event.multiplicity <= multiplicityMax) {
0337       if (event.soft) {
0338         softHisto->Fill(event.multiplicity);
0339       } else {
0340         hardHisto->Fill(event.multiplicity);
0341       }
0342     }
0343   }
0344 
0345   return {softHisto, hardHisto};
0346 }
0347 
0348 TVectorF softProbability(const EventCollection& events) {
0349   std::size_t countSoft = 0;
0350   // Count the soft events
0351   for (const EventFraction& event : events) {
0352     if (event.soft) {
0353       countSoft++;
0354     }
0355   }
0356 
0357   TVectorF result(1);
0358   result[0] = static_cast<float>(countSoft) / events.size();
0359   return result;
0360 }
0361 
0362 CumulativeDistribution cumulativeNuclearInteractionProbability(
0363     const EventCollection& events, unsigned int interactionProbabilityBins) {
0364   // Find the limits of the histogram
0365   float min = std::numeric_limits<float>::max();
0366   float max = 0.;
0367   for (const EventFraction& event : events) {
0368     min =
0369         std::min(static_cast<float>(event.interactingParticle.pathInL0()), min);
0370     max =
0371         std::max(static_cast<float>(event.interactingParticle.pathInL0()), max);
0372   }
0373 
0374   // Fill the histogram
0375   TH1F* histo = new TH1F("", "", interactionProbabilityBins, min, max);
0376   for (const EventFraction& event : events) {
0377     histo->Fill(event.interactingParticle.pathInL0());
0378   }
0379 
0380   // Build the distributions
0381   return histo;  // TODO: in this case the normalisation is not taking into
0382                  // account
0383 }
0384 
0385 }  // namespace ActsExamples::detail::NuclearInteractionParametrisation