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