Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-06-14 07:48:14

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/RootAthenaDumpReader.hpp"
0010 
0011 #include "Acts/Definitions/Units.hpp"
0012 #include "Acts/EventData/SourceLink.hpp"
0013 #include "Acts/Geometry/GeometryIdentifier.hpp"
0014 #include "Acts/Utilities/Zip.hpp"
0015 #include "ActsExamples/Digitization/MeasurementCreation.hpp"
0016 #include "ActsExamples/EventData/Cluster.hpp"
0017 #include "ActsExamples/EventData/IndexSourceLink.hpp"
0018 #include "ActsExamples/EventData/SimParticle.hpp"
0019 #include "ActsExamples/EventData/SpacePoint.hpp"
0020 
0021 #include <algorithm>
0022 #include <numeric>
0023 
0024 #include <TChain.h>
0025 #include <boost/container/static_vector.hpp>
0026 
0027 using namespace Acts::UnitLiterals;
0028 
0029 namespace ActsExamples {
0030 
0031 namespace {
0032 
0033 /// In cases when there is built up a particle collection in an iterative way it
0034 /// can be way faster to build up a vector and afterwards use a special
0035 /// constructor to speed up the set creation.
0036 inline auto particleVectorToSet(std::vector<SimParticle>& particles) {
0037   using namespace ActsExamples;
0038   auto cmp = [](const auto& a, const auto& b) {
0039     return a.particleId() == b.particleId();
0040   };
0041 
0042   std::ranges::sort(particles, detail::CompareParticleId{});
0043   particles.erase(std::unique(particles.begin(), particles.end(), cmp),
0044                   particles.end());
0045 
0046   return SimParticleContainer(boost::container::ordered_unique_range_t{},
0047                               particles.begin(), particles.end());
0048 }
0049 
0050 }  // namespace
0051 
0052 enum SpacePointType { ePixel = 1, eStrip = 2 };
0053 
0054 RootAthenaDumpReader::RootAthenaDumpReader(
0055     const RootAthenaDumpReader::Config& config, Acts::Logging::Level level)
0056     : IReader(),
0057       m_cfg(config),
0058       m_logger(Acts::getDefaultLogger(name(), level)) {
0059   if (m_cfg.inputfiles.empty()) {
0060     throw std::invalid_argument("Empty input file list");
0061   }
0062   if (m_cfg.treename.empty()) {
0063     throw std::invalid_argument("Missing tree name");
0064   }
0065 
0066   m_inputchain = std::make_shared<TChain>(m_cfg.treename.c_str());
0067 
0068   m_outputPixelSpacePoints.initialize(m_cfg.outputPixelSpacePoints);
0069   m_outputStripSpacePoints.initialize(m_cfg.outputStripSpacePoints);
0070   m_outputSpacePoints.initialize(m_cfg.outputSpacePoints);
0071   if (!m_cfg.onlySpacePoints) {
0072     m_outputMeasurements.initialize(m_cfg.outputMeasurements);
0073     m_outputMeasurementSubset.initialize(m_cfg.outputMeasurementSubset);
0074     m_outputClusters.initialize(m_cfg.outputClusters);
0075     if (!m_cfg.noTruth) {
0076       m_outputParticles.initialize(m_cfg.outputParticles);
0077       m_outputMeasParticleMap.initialize(m_cfg.outputMeasurementParticlesMap);
0078       m_outputParticleMeasMap.initialize(m_cfg.outputParticleMeasurementsMap);
0079     }
0080   }
0081 
0082   if (m_inputchain->GetBranch("SPtopStripDirection") == nullptr) {
0083     ACTS_INFO("Additional SP strip features not available");
0084     m_haveStripFeatures = false;
0085   }
0086 
0087   if (m_inputchain->GetBranch("CLetas") == nullptr) {
0088     ACTS_INFO("Cluster cell data (CLetas/CLphis/CLtots) not available");
0089     m_haveCellData = false;
0090   }
0091 
0092   // Set the branches
0093   m_inputchain->SetBranchAddress("run_number", &run_number);
0094   m_inputchain->SetBranchAddress("event_number", &event_number);
0095   m_inputchain->SetBranchAddress("nSE", &nSE);
0096   m_inputchain->SetBranchAddress("SEID", SEID);
0097 
0098   // Cluster features
0099   m_inputchain->SetBranchAddress("nCL", &nCL);
0100   m_inputchain->SetBranchAddress("CLindex", CLindex);
0101   m_inputchain->SetBranchAddress("CLhardware", &CLhardware.get());
0102   m_inputchain->SetBranchAddress("CLx", CLx);
0103   m_inputchain->SetBranchAddress("CLy", CLy);
0104   m_inputchain->SetBranchAddress("CLz", CLz);
0105   m_inputchain->SetBranchAddress("CLbarrel_endcap", CLbarrel_endcap);
0106   m_inputchain->SetBranchAddress("CLlayer_disk", CLlayer_disk);
0107   m_inputchain->SetBranchAddress("CLeta_module", CLeta_module);
0108   m_inputchain->SetBranchAddress("CLphi_module", CLphi_module);
0109   m_inputchain->SetBranchAddress("CLside", CLside);
0110   m_inputchain->SetBranchAddress("CLmoduleID", CLmoduleID);
0111   if (m_haveCellData) {
0112     m_inputchain->SetBranchAddress("CLphis", &CLphis.get());
0113     m_inputchain->SetBranchAddress("CLetas", &CLetas.get());
0114     m_inputchain->SetBranchAddress("CLtots", &CLtots.get());
0115   }
0116   m_inputchain->SetBranchAddress("CLloc_direction1", CLloc_direction1);
0117   m_inputchain->SetBranchAddress("CLloc_direction2", CLloc_direction2);
0118   m_inputchain->SetBranchAddress("CLloc_direction3", CLloc_direction3);
0119   m_inputchain->SetBranchAddress("CLJan_loc_direction1", CLJan_loc_direction1);
0120   m_inputchain->SetBranchAddress("CLJan_loc_direction2", CLJan_loc_direction2);
0121   m_inputchain->SetBranchAddress("CLJan_loc_direction3", CLJan_loc_direction3);
0122   m_inputchain->SetBranchAddress("CLpixel_count", CLpixel_count);
0123   m_inputchain->SetBranchAddress("CLcharge_count", CLcharge_count);
0124   m_inputchain->SetBranchAddress("CLloc_eta", CLloc_eta);
0125   m_inputchain->SetBranchAddress("CLloc_phi", CLloc_phi);
0126   m_inputchain->SetBranchAddress("CLglob_eta", CLglob_eta);
0127   m_inputchain->SetBranchAddress("CLglob_phi", CLglob_phi);
0128   m_inputchain->SetBranchAddress("CLeta_angle", CLeta_angle);
0129   m_inputchain->SetBranchAddress("CLphi_angle", CLphi_angle);
0130   m_inputchain->SetBranchAddress("CLnorm_x", CLnorm_x);
0131   m_inputchain->SetBranchAddress("CLnorm_y", CLnorm_y);
0132   m_inputchain->SetBranchAddress("CLnorm_z", CLnorm_z);
0133   m_inputchain->SetBranchAddress("CLlocal_cov", &CLlocal_cov.get());
0134   if (!m_cfg.noTruth) {
0135     m_inputchain->SetBranchAddress("CLparticleLink_eventIndex",
0136                                    &CLparticleLink_eventIndex.get());
0137     m_inputchain->SetBranchAddress("CLparticleLink_barcode",
0138                                    &CLparticleLink_barcode.get());
0139     m_inputchain->SetBranchAddress("CLbarcodesLinked", &CLbarcodesLinked.get());
0140     m_inputchain->SetBranchAddress("CLparticle_charge",
0141                                    &CLparticle_charge.get());
0142   }
0143 
0144   // Particle features
0145   if (!m_cfg.noTruth) {
0146     m_inputchain->SetBranchAddress("nPartEVT", &nPartEVT);
0147     m_inputchain->SetBranchAddress("Part_event_number", Part_event_number);
0148     m_inputchain->SetBranchAddress("Part_barcode", Part_barcode);
0149     m_inputchain->SetBranchAddress("Part_px", Part_px);
0150     m_inputchain->SetBranchAddress("Part_py", Part_py);
0151     m_inputchain->SetBranchAddress("Part_pz", Part_pz);
0152     m_inputchain->SetBranchAddress("Part_pt", Part_pt);
0153     m_inputchain->SetBranchAddress("Part_eta", Part_eta);
0154     m_inputchain->SetBranchAddress("Part_vx", Part_vx);
0155     m_inputchain->SetBranchAddress("Part_vy", Part_vy);
0156     m_inputchain->SetBranchAddress("Part_vz", Part_vz);
0157     m_inputchain->SetBranchAddress("Part_radius", Part_radius);
0158     m_inputchain->SetBranchAddress("Part_status", Part_status);
0159     m_inputchain->SetBranchAddress("Part_charge", Part_charge);
0160     m_inputchain->SetBranchAddress("Part_pdg_id", Part_pdg_id);
0161     m_inputchain->SetBranchAddress("Part_passed", Part_passed);
0162     m_inputchain->SetBranchAddress("Part_vProdNin", Part_vProdNin);
0163     m_inputchain->SetBranchAddress("Part_vProdNout", Part_vProdNout);
0164     m_inputchain->SetBranchAddress("Part_vProdStatus", Part_vProdStatus);
0165     m_inputchain->SetBranchAddress("Part_vProdBarcode", Part_vProdBarcode);
0166     m_inputchain->SetBranchAddress("Part_vParentID", &Part_vParentID.get());
0167     m_inputchain->SetBranchAddress("Part_vParentBarcode",
0168                                    &Part_vParentBarcode.get());
0169   }
0170 
0171   // Space point features
0172   m_inputchain->SetBranchAddress("nSP", &nSP);
0173   m_inputchain->SetBranchAddress("SPindex", SPindex);
0174   m_inputchain->SetBranchAddress("SPx", SPx);
0175   m_inputchain->SetBranchAddress("SPy", SPy);
0176   m_inputchain->SetBranchAddress("SPz", SPz);
0177   m_inputchain->SetBranchAddress("SPCL1_index", SPCL1_index);
0178   m_inputchain->SetBranchAddress("SPCL2_index", SPCL2_index);
0179   m_inputchain->SetBranchAddress("SPisOverlap", SPisOverlap);
0180   if (m_haveStripFeatures) {
0181     m_inputchain->SetBranchAddress("SPradius", SPradius);
0182     m_inputchain->SetBranchAddress("SPcovr", SPcovr);
0183     m_inputchain->SetBranchAddress("SPcovz", SPcovz);
0184     m_inputchain->SetBranchAddress("SPhl_topstrip", SPhl_topstrip);
0185     m_inputchain->SetBranchAddress("SPhl_botstrip", SPhl_botstrip);
0186     m_inputchain->SetBranchAddress("SPtopStripDirection",
0187                                    &SPtopStripDirection.get());
0188     m_inputchain->SetBranchAddress("SPbottomStripDirection",
0189                                    &SPbottomStripDirection.get());
0190     m_inputchain->SetBranchAddress("SPstripCenterDistance",
0191                                    &SPstripCenterDistance.get());
0192     m_inputchain->SetBranchAddress("SPtopStripCenterPosition",
0193                                    &SPtopStripCenterPosition.get());
0194   }
0195 
0196   // These quantities are not used currently and thus commented out
0197   // I would like to keep the code, since it is always a pain to write it
0198   /*
0199   m_inputchain->SetBranchAddress("nTRK", &nTRK);
0200   m_inputchain->SetBranchAddress("TRKindex", TRKindex);
0201   m_inputchain->SetBranchAddress("TRKtrack_fitter", TRKtrack_fitter);
0202   m_inputchain->SetBranchAddress("TRKparticle_hypothesis",
0203                                  TRKparticle_hypothesis);
0204   m_inputchain->SetBranchAddress("TRKproperties", &TRKproperties);
0205   m_inputchain->SetBranchAddress("TRKpattern", &TRKpattern);
0206   m_inputchain->SetBranchAddress("TRKndof", TRKndof);
0207   m_inputchain->SetBranchAddress("TRKmot", TRKmot);
0208   m_inputchain->SetBranchAddress("TRKoot", TRKoot);
0209   m_inputchain->SetBranchAddress("TRKchiSq", TRKchiSq);
0210   m_inputchain->SetBranchAddress("TRKmeasurementsOnTrack_pixcl_sctcl_index",
0211                                  &TRKmeasurementsOnTrack_pixcl_sctcl_index);
0212   m_inputchain->SetBranchAddress("TRKoutliersOnTrack_pixcl_sctcl_index",
0213                                  &TRKoutliersOnTrack_pixcl_sctcl_index);
0214   m_inputchain->SetBranchAddress("TRKcharge", TRKcharge);
0215   m_inputchain->SetBranchAddress("TRKperigee_position", &TRKperigee_position);
0216   m_inputchain->SetBranchAddress("TRKperigee_momentum", &TRKperigee_momentum);
0217   m_inputchain->SetBranchAddress("TTCindex", TTCindex);
0218   m_inputchain->SetBranchAddress("TTCevent_index", TTCevent_index);
0219   m_inputchain->SetBranchAddress("TTCparticle_link", TTCparticle_link);
0220   m_inputchain->SetBranchAddress("TTCprobability", TTCprobability);
0221   m_inputchain->SetBranchAddress("nDTT", &nDTT);
0222   m_inputchain->SetBranchAddress("DTTindex", DTTindex);
0223   m_inputchain->SetBranchAddress("DTTsize", DTTsize);
0224   m_inputchain->SetBranchAddress("DTTtrajectory_eventindex",
0225                                  &DTTtrajectory_eventindex);
0226   m_inputchain->SetBranchAddress("DTTtrajectory_barcode",
0227                                  &DTTtrajectory_barcode);
0228   m_inputchain->SetBranchAddress("DTTstTruth_subDetType",
0229                                  &DTTstTruth_subDetType);
0230   m_inputchain->SetBranchAddress("DTTstTrack_subDetType",
0231                                  &DTTstTrack_subDetType);
0232   m_inputchain->SetBranchAddress("DTTstCommon_subDetType",
0233                                  &DTTstCommon_subDetType);
0234   */
0235 
0236   for (const auto& file : m_cfg.inputfiles) {
0237     ACTS_DEBUG("Adding file '" << file << "' to tree " << m_cfg.treename);
0238     m_inputchain->Add(file.c_str());
0239   }
0240 
0241   m_events = m_inputchain->GetEntries();
0242 
0243   ACTS_DEBUG("End of constructor. In total available events=" << m_events);
0244 }  // constructor
0245 
0246 SimParticleContainer RootAthenaDumpReader::readParticles() const {
0247   std::vector<SimParticle> particles;
0248   particles.reserve(nPartEVT);
0249 
0250   for (auto ip = 0; ip < nPartEVT; ++ip) {
0251     if (m_cfg.onlyPassedParticles && !static_cast<bool>(Part_passed[ip])) {
0252       continue;
0253     }
0254 
0255     SimBarcode dummyBarcode =
0256         SimBarcode()
0257             .withVertexPrimary(
0258                 static_cast<SimBarcode::PrimaryVertexId>(Part_event_number[ip]))
0259             .withVertexSecondary(static_cast<SimBarcode::SecondaryVertexId>(
0260                 Part_barcode[ip] < s_maxBarcodeForPrimary ? 0 : 1))
0261             .withParticle(
0262                 static_cast<SimBarcode::ParticleId>(Part_barcode[ip]));
0263     SimParticleState particle(dummyBarcode,
0264                               static_cast<Acts::PdgParticle>(Part_pdg_id[ip]));
0265 
0266     Acts::Vector3 p = Acts::Vector3{Part_px[ip], Part_py[ip], Part_pz[ip]} *
0267                       Acts::UnitConstants::MeV;
0268     particle.setAbsoluteMomentum(p.norm());
0269 
0270     particle.setDirection(p.normalized());
0271 
0272     auto x = Acts::Vector4{Part_vx[ip], Part_vy[ip], Part_vz[ip], 0.0};
0273     particle.setPosition4(x);
0274 
0275     particles.push_back(SimParticle(particle, particle));
0276   }
0277 
0278   ACTS_DEBUG("Created " << particles.size() << " particles");
0279   auto before = particles.size();
0280 
0281   auto particlesSet = particleVectorToSet(particles);
0282 
0283   if (particlesSet.size() < before) {
0284     ACTS_WARNING("Particle IDs not unique for " << before - particles.size()
0285                                                 << " particles!");
0286   }
0287 
0288   return particlesSet;
0289 }
0290 
0291 std::tuple<ClusterContainer, MeasurementContainer, MeasurementParticlesMap,
0292            std::unordered_map<int, std::size_t>>
0293 RootAthenaDumpReader::readMeasurements(
0294     SimParticleContainer& particles, const Acts::GeometryContext& gctx) const {
0295   ClusterContainer clusters;
0296   clusters.reserve(nCL);
0297 
0298   MeasurementContainer measurements;
0299   measurements.reserve(nCL);
0300 
0301   std::size_t nTotalTotZero = 0;
0302 
0303   const auto prevParticlesSize = particles.size();
0304   MeasurementParticlesMap measPartMap;
0305 
0306   // We cannot use im for the index since we might skip measurements
0307   std::unordered_map<int, std::size_t> imIdxMap;
0308   imIdxMap.reserve(nCL);
0309 
0310   for (int im = 0; im < nCL; im++) {
0311     if (!(CLhardware->at(im) == "PIXEL" || CLhardware->at(im) == "STRIP")) {
0312       ACTS_ERROR("hardware is neither 'PIXEL' or 'STRIP', skip particle");
0313       continue;
0314     }
0315     ACTS_VERBOSE("Cluster " << im << ": " << CLhardware->at(im));
0316 
0317     auto type = (CLhardware->at(im) == "PIXEL") ? ePixel : eStrip;
0318 
0319     // Make cluster
0320     // TODO refactor Cluster class so it is not so tedious
0321     Cluster cluster;
0322     if (m_cfg.readCellData && m_haveCellData) {
0323       const auto& etas = CLetas->at(im);
0324       const auto& phis = CLetas->at(im);
0325       const auto& tots = CLtots->at(im);
0326 
0327       const auto totalTot = std::accumulate(tots.begin(), tots.end(), 0);
0328 
0329       const auto [minEta, maxEta] =
0330           std::minmax_element(etas.begin(), etas.end());
0331       const auto [minPhi, maxPhi] =
0332           std::minmax_element(phis.begin(), phis.end());
0333 
0334       cluster.channels.reserve(etas.size());
0335 
0336       cluster.sizeLoc0 = *maxEta - *minEta;
0337       cluster.sizeLoc1 = *maxPhi - *minPhi;
0338 
0339       if (totalTot == 0.0) {
0340         ACTS_VERBOSE(
0341             "total time over threshold is 0, set all activations to 0");
0342         nTotalTotZero++;
0343       }
0344 
0345       for (const auto& [eta, phi, tot] : Acts::zip(etas, phis, tots)) {
0346         // Make best out of what we have:
0347         // Weight the overall collected charge corresponding to the
0348         // time-over-threshold of each cell Use this as activation (does this
0349         // make sense?)
0350         auto activation =
0351             (totalTot != 0.0) ? CLcharge_count[im] * tot / totalTot : 0.0;
0352 
0353         // This bases every cluster at zero, but shouldn't matter right now
0354         ActsFatras::Segmentizer::Bin2D bin;
0355         bin[0] = eta - *minEta;
0356         bin[1] = phi - *minPhi;
0357 
0358         // Of course we have no Segment2D because this is not Fatras
0359         cluster.channels.emplace_back(bin, ActsFatras::Segmentizer::Segment2D{},
0360                                       activation);
0361       }
0362 
0363       ACTS_VERBOSE("- shape: " << cluster.channels.size()
0364                                << "cells, dimensions: " << cluster.sizeLoc0
0365                                << ", " << cluster.sizeLoc1);
0366     }
0367 
0368     cluster.globalPosition = {CLx[im], CLy[im], CLz[im]};
0369     cluster.localDirection = {CLloc_direction1[im], CLloc_direction2[im],
0370                               CLloc_direction3[im]};
0371     cluster.lengthDirection = {CLJan_loc_direction1[im],
0372                                CLJan_loc_direction2[im],
0373                                CLJan_loc_direction3[im]};
0374     cluster.localEta = CLloc_eta[im];
0375     cluster.localPhi = CLloc_phi[im];
0376     cluster.globalEta = CLglob_eta[im];
0377     cluster.globalPhi = CLglob_phi[im];
0378     cluster.etaAngle = CLeta_angle[im];
0379     cluster.phiAngle = CLphi_angle[im];
0380 
0381     // Measurement creation
0382     Acts::GeometryIdentifier geoId;
0383     std::vector<double> localParams;
0384     if (m_cfg.geometryIdMap && m_cfg.trackingGeometry) {
0385       const auto& geoIdMap = m_cfg.geometryIdMap->left;
0386       if (geoIdMap.find(CLmoduleID[im]) == geoIdMap.end()) {
0387         ACTS_WARNING("Missing geo id for " << CLmoduleID[im] << ", skip hit");
0388         continue;
0389       }
0390 
0391       geoId = m_cfg.geometryIdMap->left.at(CLmoduleID[im]);
0392 
0393       auto surface = m_cfg.trackingGeometry->findSurface(geoId);
0394       if (surface == nullptr) {
0395         ACTS_WARNING("Did not find " << geoId
0396                                      << " in tracking geometry, skip hit");
0397         continue;
0398       }
0399 
0400       bool inside =
0401           surface->isOnSurface(gctx, cluster.globalPosition, {},
0402                                Acts::BoundaryTolerance::AbsoluteEuclidean(
0403                                    m_cfg.absBoundaryTolerance),
0404                                std::numeric_limits<double>::max());
0405 
0406       if (!inside) {
0407         const Acts::Vector3 v =
0408             surface->localToGlobalTransform(gctx).inverse() *
0409             cluster.globalPosition;
0410         ACTS_WARNING("Projected position is not in surface bounds for "
0411                      << surface->geometryId() << ", skip hit");
0412         ACTS_WARNING("Position in local coordinates: " << v.transpose());
0413         ACTS_WARNING("Surface details:\n" << surface->toStream(gctx));
0414         continue;
0415       }
0416 
0417       auto loc = surface->globalToLocal(gctx, cluster.globalPosition, {},
0418                                         Acts::s_onSurfaceTolerance);
0419 
0420       if (!loc.ok()) {
0421         const Acts::Vector3 v =
0422             surface->localToGlobalTransform(gctx).inverse() *
0423             cluster.globalPosition;
0424         ACTS_WARNING("Global-to-local fit failed on "
0425                      << geoId << " (z dist: " << v[2]
0426                      << ", projected on surface: " << std::boolalpha << inside
0427                      << ") , skip hit");
0428         continue;
0429       }
0430 
0431       // TODO is this in strip coordinates or in polar coordinates for annulus
0432       // bounds?
0433       localParams = std::vector<double>(loc->begin(), loc->end());
0434     } else {
0435       geoId = Acts::GeometryIdentifier(CLmoduleID[im]);
0436       localParams = {CLloc_direction1[im], CLloc_direction2[im]};
0437     }
0438 
0439     const auto& locCov = CLlocal_cov->at(im);
0440 
0441     DigitizedParameters digiPars;
0442     if (type == ePixel) {
0443       digiPars.indices = {Acts::eBoundLoc0, Acts::eBoundLoc1};
0444       assert(locCov.size() == 4);
0445       digiPars.variances = {locCov[0], locCov[3]};
0446       digiPars.values = localParams;
0447     } else {
0448       assert(!locCov.empty());
0449       // Barrel-endcap index can be -2/2 for endcap or 0 for barrel
0450       // We need to choose the coordinate of local measurement depending on that
0451       const static std::array boundLoc = {Acts::eBoundLoc0, Acts::eBoundLoc1};
0452       auto i = CLbarrel_endcap[im] == 0 ? 0 : 1;
0453       digiPars.variances = {locCov[i]};
0454       digiPars.values = {localParams[i]};
0455       digiPars.indices = {boundLoc[i]};
0456     }
0457 
0458     std::size_t measIndex = measurements.size();
0459     ACTS_VERBOSE("Add measurement with index " << measIndex);
0460     imIdxMap.emplace(im, measIndex);
0461     createMeasurement(measurements, geoId, digiPars);
0462     clusters.push_back(cluster);
0463 
0464     if (!m_cfg.noTruth) {
0465       // Create measurement particles map and particles container
0466       for (const auto& [subevt, barcode] :
0467            Acts::zip(CLparticleLink_eventIndex->at(im),
0468                      CLparticleLink_barcode->at(im))) {
0469         SimBarcode dummyBarcode =
0470             SimBarcode()
0471                 .withVertexPrimary(
0472                     static_cast<SimBarcode::PrimaryVertexId>(subevt))
0473                 .withVertexSecondary(static_cast<SimBarcode::SecondaryVertexId>(
0474                     barcode < s_maxBarcodeForPrimary ? 0 : 1))
0475                 .withParticle(static_cast<SimBarcode::ParticleId>(barcode));
0476         // If we don't find the particle, create one with default values
0477         if (particles.find(dummyBarcode) == particles.end()) {
0478           ACTS_VERBOSE("Particle with subevt "
0479                        << subevt << ", barcode " << barcode
0480                        << "not found, create dummy one");
0481           particles.emplace(dummyBarcode, Acts::PdgParticle::eInvalid);
0482         }
0483         measPartMap.insert(
0484             std::pair<Index, SimBarcode>{measIndex, dummyBarcode});
0485       }
0486     }
0487   }
0488 
0489   if (measurements.size() < static_cast<std::size_t>(nCL)) {
0490     ACTS_WARNING("Could not convert " << nCL - measurements.size() << " / "
0491                                       << nCL << " measurements");
0492   }
0493 
0494   if (particles.size() - prevParticlesSize > 0) {
0495     ACTS_DEBUG("Created " << particles.size() - prevParticlesSize
0496                           << " dummy particles");
0497   }
0498 
0499   if (nTotalTotZero > 0) {
0500     ACTS_DEBUG(nTotalTotZero << " / " << nCL
0501                              << " clusters have zero time-over-threshold");
0502   }
0503 
0504   return {std::move(clusters), std::move(measurements), std::move(measPartMap),
0505           std::move(imIdxMap)};
0506 }
0507 
0508 std::tuple<SpacePointContainer, SpacePointContainer, SpacePointContainer>
0509 RootAthenaDumpReader::readSpacePoints(
0510     const std::optional<std::unordered_map<int, std::size_t>>& imIdxMap) const {
0511   SpacePointContainer pixelSpacePoints(
0512       SpacePointColumns::SourceLinks | SpacePointColumns::X |
0513       SpacePointColumns::Y | SpacePointColumns::Z |
0514       SpacePointColumns::VarianceZ | SpacePointColumns::VarianceR);
0515   pixelSpacePoints.reserve(nSP);
0516 
0517   SpacePointContainer stripSpacePoints(
0518       SpacePointColumns::SourceLinks | SpacePointColumns::X |
0519       SpacePointColumns::Y | SpacePointColumns::Z |
0520       SpacePointColumns::VarianceZ | SpacePointColumns::VarianceR |
0521       SpacePointColumns::Strip);
0522   stripSpacePoints.reserve(nSP);
0523 
0524   SpacePointContainer spacePoints(
0525       SpacePointColumns::SourceLinks | SpacePointColumns::X |
0526       SpacePointColumns::Y | SpacePointColumns::Z |
0527       SpacePointColumns::VarianceZ | SpacePointColumns::VarianceR |
0528       SpacePointColumns::Strip);
0529   spacePoints.reserve(nSP);
0530 
0531   // Loop on space points
0532   std::size_t skippedSpacePoints = 0;
0533   for (int isp = 0; isp < nSP; isp++) {
0534     auto isPhiOverlap = (SPisOverlap[isp] == 2) || (SPisOverlap[isp] == 3);
0535     auto isEtaOverlap = (SPisOverlap[isp] == 1) || (SPisOverlap[isp] == 3);
0536     if (m_cfg.skipOverlapSPsPhi && isPhiOverlap) {
0537       ++skippedSpacePoints;
0538       continue;
0539     }
0540     if (m_cfg.skipOverlapSPsEta && isEtaOverlap) {
0541       ++skippedSpacePoints;
0542       continue;
0543     }
0544 
0545     const Acts::Vector3 globalPos{SPx[isp], SPy[isp], SPz[isp]};
0546     const double spCovr = SPcovr[isp];
0547     const double spCovz = SPcovz[isp];
0548 
0549     // PIX=1  STRIP = 2
0550     auto type = SPCL2_index[isp] == -1 ? ePixel : eStrip;
0551 
0552     ACTS_VERBOSE("SP:: " << type << " [" << globalPos.transpose() << "] "
0553                          << spCovr << " " << spCovz);
0554 
0555     boost::container::static_vector<Acts::SourceLink, 2> sLinks;
0556 
0557     const auto cl1Index = SPCL1_index[isp];
0558     assert(cl1Index >= 0 && cl1Index < nCL);
0559 
0560     auto getGeoId =
0561         [&](auto athenaId) -> std::optional<Acts::GeometryIdentifier> {
0562       if (m_cfg.geometryIdMap == nullptr) {
0563         return Acts::GeometryIdentifier{athenaId};
0564       }
0565       if (m_cfg.geometryIdMap->left.find(athenaId) ==
0566           m_cfg.geometryIdMap->left.end()) {
0567         return std::nullopt;
0568       }
0569       return m_cfg.geometryIdMap->left.at(athenaId);
0570     };
0571 
0572     auto cl1GeoId = getGeoId(CLmoduleID[cl1Index]);
0573     if (!cl1GeoId) {
0574       ACTS_WARNING("Could not find geoId for space point cluster 1");
0575       continue;
0576     }
0577 
0578     if (imIdxMap && !imIdxMap->contains(cl1Index)) {
0579       ACTS_WARNING("Measurement 1 for space point " << isp << " not created");
0580       continue;
0581     }
0582 
0583     IndexSourceLink first(*cl1GeoId,
0584                           imIdxMap ? imIdxMap->at(cl1Index) : cl1Index);
0585     sLinks.emplace_back(first);
0586 
0587     // First create pixel space point here, later maybe overwrite with strip
0588     // space point
0589     auto sp = spacePoints.createSpacePoint();
0590     sp.x() = globalPos.x();
0591     sp.y() = globalPos.y();
0592     sp.z() = globalPos.z();
0593     sp.varianceR() = spCovr;
0594     sp.varianceZ() = spCovz;
0595 
0596     if (type == ePixel) {
0597       auto pixelSp = pixelSpacePoints.createSpacePoint();
0598       pixelSp.assignSourceLinks(sLinks);
0599       pixelSp.x() = globalPos.x();
0600       pixelSp.y() = globalPos.y();
0601       pixelSp.z() = globalPos.z();
0602       pixelSp.varianceR() = spCovr;
0603       pixelSp.varianceZ() = spCovz;
0604     } else {
0605       const auto cl2Index = SPCL2_index[isp];
0606       assert(cl2Index >= 0 && cl2Index < nCL);
0607 
0608       auto cl2GeoId = getGeoId(CLmoduleID[cl2Index]);
0609       if (!cl2GeoId) {
0610         ACTS_WARNING("Could not find geoId for space point cluster 2");
0611         continue;
0612       }
0613 
0614       if (imIdxMap && !imIdxMap->contains(cl2Index)) {
0615         ACTS_WARNING("Measurement 2 for space point " << isp << " not created");
0616         continue;
0617       }
0618 
0619       IndexSourceLink second(*cl2GeoId,
0620                              imIdxMap ? imIdxMap->at(cl2Index) : cl2Index);
0621       sLinks.emplace_back(second);
0622 
0623       Eigen::Vector3f topStripDirection = Eigen::Vector3f::Zero();
0624       Eigen::Vector3f bottomStripDirection = Eigen::Vector3f::Zero();
0625       Eigen::Vector3f stripCenterDistance = Eigen::Vector3f::Zero();
0626       Eigen::Vector3f topStripCenterPosition = Eigen::Vector3f::Zero();
0627 
0628       if (m_haveStripFeatures) {
0629         topStripDirection = {SPtopStripDirection->at(isp).at(0),
0630                              SPtopStripDirection->at(isp).at(1),
0631                              SPtopStripDirection->at(isp).at(2)};
0632         bottomStripDirection = {SPbottomStripDirection->at(isp).at(0),
0633                                 SPbottomStripDirection->at(isp).at(1),
0634                                 SPbottomStripDirection->at(isp).at(2)};
0635         stripCenterDistance = {SPstripCenterDistance->at(isp).at(0),
0636                                SPstripCenterDistance->at(isp).at(1),
0637                                SPstripCenterDistance->at(isp).at(2)};
0638         topStripCenterPosition = {SPtopStripCenterPosition->at(isp).at(0),
0639                                   SPtopStripCenterPosition->at(isp).at(1),
0640                                   SPtopStripCenterPosition->at(isp).at(2)};
0641       }
0642 
0643       auto stripSp = stripSpacePoints.createSpacePoint();
0644       stripSp.assignSourceLinks(sLinks);
0645       stripSp.x() = globalPos.x();
0646       stripSp.y() = globalPos.y();
0647       stripSp.z() = globalPos.z();
0648       stripSp.varianceR() = spCovr;
0649       stripSp.varianceZ() = spCovz;
0650       Eigen::Map<Eigen::Vector3f>(stripSp.topStripVector().data()) =
0651           topStripDirection.cast<float>();
0652       Eigen::Map<Eigen::Vector3f>(stripSp.bottomStripVector().data()) =
0653           bottomStripDirection.cast<float>();
0654       Eigen::Map<Eigen::Vector3f>(stripSp.stripCenterDistance().data()) =
0655           stripCenterDistance.cast<float>();
0656       Eigen::Map<Eigen::Vector3f>(stripSp.topStripCenter().data()) =
0657           topStripCenterPosition.cast<float>();
0658 
0659       sp.topStripVector() = stripSp.topStripVector();
0660       sp.bottomStripVector() = stripSp.bottomStripVector();
0661       sp.stripCenterDistance() = stripSp.stripCenterDistance();
0662       sp.topStripCenter() = stripSp.topStripCenter();
0663     }
0664 
0665     sp.assignSourceLinks(sLinks);
0666   }
0667 
0668   if (m_cfg.skipOverlapSPsEta || m_cfg.skipOverlapSPsPhi) {
0669     ACTS_DEBUG("Skipped " << skippedSpacePoints
0670                           << " because of eta/phi overlaps");
0671   }
0672   if (spacePoints.size() <
0673       (static_cast<std::size_t>(nSP) - skippedSpacePoints)) {
0674     ACTS_WARNING("Could not convert " << nSP - spacePoints.size() << " of "
0675                                       << nSP << " space points");
0676   }
0677 
0678   ACTS_DEBUG("Created " << spacePoints.size() << " overall space points");
0679   ACTS_DEBUG("Created " << pixelSpacePoints.size() << " "
0680                         << " pixel space points");
0681   ACTS_DEBUG("Created " << stripSpacePoints.size() << " "
0682                         << " strip space points");
0683 
0684   return {std::move(spacePoints), std::move(pixelSpacePoints),
0685           std::move(stripSpacePoints)};
0686 }
0687 
0688 std::pair<SimParticleContainer, MeasurementParticlesMap>
0689 RootAthenaDumpReader::reprocessParticles(
0690     const SimParticleContainer& particles,
0691     const MeasurementParticlesMap& measPartMap) const {
0692   std::vector<SimParticle> newParticles;
0693   newParticles.reserve(particles.size());
0694   MeasurementParticlesMap newMeasPartMap;
0695   newMeasPartMap.reserve(measPartMap.size());
0696 
0697   const auto partMeasMap = invertIndexMultimap(measPartMap);
0698 
0699   std::uint16_t primaryCount = 0;
0700   std::uint16_t secondaryCount = 0;
0701 
0702   for (const auto& particle : particles) {
0703     const auto [begin, end] = partMeasMap.equal_range(particle.particleId());
0704 
0705     if (begin == end) {
0706       ACTS_VERBOSE("Particle " << particle.particleId()
0707                                << " has no measurements");
0708       continue;
0709     }
0710 
0711     auto primary = particle.particleId().vertexSecondary() == 0;
0712 
0713     // vertex primary shouldn't be zero for a valid particle
0714     SimBarcode fatrasBarcode = SimBarcode().withVertexPrimary(1);
0715     if (primary) {
0716       fatrasBarcode =
0717           fatrasBarcode.withVertexSecondary(0).withParticle(primaryCount);
0718       assert(primaryCount < std::numeric_limits<std::uint16_t>::max());
0719       primaryCount++;
0720     } else {
0721       fatrasBarcode =
0722           fatrasBarcode.withVertexSecondary(1).withParticle(secondaryCount);
0723       assert(primaryCount < std::numeric_limits<std::uint16_t>::max());
0724       secondaryCount++;
0725     }
0726 
0727     auto newParticle = particle.withParticleId(fatrasBarcode);
0728     newParticle.finalState().setNumberOfHits(std::distance(begin, end));
0729     newParticles.push_back(newParticle);
0730 
0731     for (auto it = begin; it != end; ++it) {
0732       newMeasPartMap.insert(
0733           std::pair<Index, SimBarcode>{it->second, fatrasBarcode});
0734     }
0735   }
0736 
0737   ACTS_DEBUG("After reprocessing particles " << newParticles.size() << " of "
0738                                              << particles.size() << " remain");
0739   return {particleVectorToSet(newParticles), std::move(newMeasPartMap)};
0740 }
0741 
0742 ProcessCode RootAthenaDumpReader::read(const AlgorithmContext& ctx) {
0743   ACTS_DEBUG("Reading event " << ctx.eventNumber);
0744   auto entry = ctx.eventNumber;
0745   if (entry >= m_events) {
0746     ACTS_ERROR("event out of bounds");
0747     return ProcessCode::ABORT;
0748   }
0749 
0750   std::lock_guard<std::mutex> lock(m_read_mutex);
0751 
0752   m_inputchain->GetEntry(entry);
0753 
0754   std::optional<std::unordered_map<int, std::size_t>> optImIdxMap;
0755 
0756   if (!m_cfg.onlySpacePoints) {
0757     SimParticleContainer candidateParticles;
0758 
0759     if (!m_cfg.noTruth) {
0760       candidateParticles = readParticles();
0761     }
0762 
0763     auto [clusters, measurements, candidateMeasPartMap, imIdxMap] =
0764         readMeasurements(candidateParticles, ctx.geoContext);
0765     optImIdxMap.emplace(std::move(imIdxMap));
0766 
0767     m_outputClusters(ctx, std::move(clusters));
0768     const auto& storedMeasurements =
0769         m_outputMeasurements(ctx, std::move(measurements));
0770     std::vector<MeasurementContainer::Index> allIndices(
0771         storedMeasurements.size());
0772     std::iota(allIndices.begin(), allIndices.end(),
0773               MeasurementContainer::Index{0});
0774     m_outputMeasurementSubset(
0775         ctx, MeasurementSubset(storedMeasurements, std::move(allIndices)));
0776 
0777     if (!m_cfg.noTruth) {
0778       auto [particles, measPartMap] =
0779           reprocessParticles(candidateParticles, candidateMeasPartMap);
0780 
0781       m_outputParticles(ctx, std::move(particles));
0782       m_outputParticleMeasMap(ctx, invertIndexMultimap(measPartMap));
0783       m_outputMeasParticleMap(ctx, std::move(measPartMap));
0784     }
0785   }
0786 
0787   auto [spacePoints, pixelSpacePoints, stripSpacePoints] =
0788       readSpacePoints(optImIdxMap);
0789 
0790   m_outputPixelSpacePoints(ctx, std::move(pixelSpacePoints));
0791   m_outputStripSpacePoints(ctx, std::move(stripSpacePoints));
0792   m_outputSpacePoints(ctx, std::move(spacePoints));
0793 
0794   return ProcessCode::SUCCESS;
0795 }
0796 }  // namespace ActsExamples