Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:11:40

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/TrackFinding/TrackParamsLookupEstimation.hpp"
0010 
0011 #include "Acts/Surfaces/RectangleBounds.hpp"
0012 #include "ActsExamples/Framework/ProcessCode.hpp"
0013 
0014 ActsExamples::TrackParamsLookupEstimation::TrackParamsLookupEstimation(
0015     const Config& config, Acts::Logging::Level level)
0016     : IAlgorithm("TrackParamsLookupEstimation", level), m_cfg(config) {
0017   // Iterate over the reference layers and create
0018   // track parameter accumulators
0019   for (const auto& [geoId, refSurface] : m_cfg.refLayers) {
0020     // Get bounds to construct the accumulator grid
0021     auto bounds =
0022         dynamic_cast<const Acts::RectangleBounds*>(&refSurface->bounds());
0023 
0024     if (bounds == nullptr) {
0025       throw std::invalid_argument("Only rectangle bounds supported");
0026     }
0027     if (refSurface->type() != Acts::Surface::SurfaceType::Plane) {
0028       throw std::invalid_argument("Only plane surfaces supported");
0029     }
0030 
0031     // Initialize the accumulator grid
0032     auto halfX = bounds->halfLengthX();
0033     auto halfY = bounds->halfLengthY();
0034 
0035     TrackParamsLookupAxisGen axisGen{
0036         {-halfX, halfX}, m_cfg.bins.first, {-halfY, halfY}, m_cfg.bins.second};
0037 
0038     // Each reference layer has its own accumulator
0039     m_accumulators[geoId] = std::make_unique<TrackParamsLookupAccumulator>(
0040         TrackParamsLookupGrid(axisGen()));
0041   }
0042 
0043   m_inputParticles.initialize(m_cfg.inputParticles);
0044   m_inputSimHits.initialize(m_cfg.inputHits);
0045 }
0046 
0047 ActsExamples::ProcessCode
0048 ActsExamples::TrackParamsLookupEstimation::finalize() {
0049   // Finiliaze the lookup tables and write them
0050   ActsExamples::TrackParamsLookup lookup;
0051   for (auto& [id, acc] : m_accumulators) {
0052     lookup.insert({id, acc->finalizeLookup()});
0053   }
0054   for (const auto& writer : m_cfg.trackLookupGridWriters) {
0055     writer->writeLookup(lookup);
0056   }
0057 
0058   return ActsExamples::ProcessCode::SUCCESS;
0059 };
0060 
0061 ActsExamples::ProcessCode ActsExamples::TrackParamsLookupEstimation::execute(
0062     const ActsExamples::AlgorithmContext& ctx) const {
0063   // Get the particles and hits
0064   const auto& particles = m_inputParticles(ctx);
0065   const auto& hits = m_inputSimHits(ctx);
0066 
0067   // Iterate over the reference layer hits and
0068   // accumulate the track parameters
0069   for (const auto& [geoId, refSurface] : m_cfg.refLayers) {
0070     // Get reference layer hits
0071     auto refLayerHits = hits.equal_range(geoId);
0072 
0073     for (auto hit = refLayerHits.first; hit != refLayerHits.second; ++hit) {
0074       // Get the corresponding particle
0075       const auto& id = hit->particleId();
0076       const auto& particle = particles.find(id);
0077 
0078       if (particle == particles.end()) {
0079         throw std::invalid_argument("Particle not found");
0080       }
0081 
0082       // Hit stores the reference layer parameters
0083       auto refLayerPars = Acts::CurvilinearTrackParameters(
0084           hit->fourPosition(), hit->direction(), particle->qOverP(),
0085           std::nullopt, particle->hypothesis());
0086 
0087       // Particle stores the IP parameters
0088       auto ipPars = Acts::CurvilinearTrackParameters(
0089           particle->fourPosition(), particle->direction(), particle->qOverP(),
0090           std::nullopt, particle->hypothesis());
0091 
0092       // Get the local position of the hit
0093       auto localPos = refSurface
0094                           ->globalToLocal(ctx.geoContext, hit->position(),
0095                                           Acts::Vector3{0, 1, 0})
0096                           .value();
0097 
0098       // Add the track parameters to the accumulator grid
0099       m_accumulators.at(geoId)->addTrack(ipPars, refLayerPars, localPos);
0100     }
0101   }
0102 
0103   return ActsExamples::ProcessCode::SUCCESS;
0104 }