File indexing completed on 2025-08-05 08:09:51
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "ActsExamples/Vertexing/IterativeVertexFinderAlgorithm.hpp"
0010
0011 #include "Acts/Definitions/Algebra.hpp"
0012 #include "Acts/Propagator/EigenStepper.hpp"
0013 #include "Acts/Propagator/VoidNavigator.hpp"
0014 #include "Acts/Utilities/Logger.hpp"
0015 #include "Acts/Utilities/Result.hpp"
0016 #include "Acts/Vertexing/IterativeVertexFinder.hpp"
0017 #include "Acts/Vertexing/TrackAtVertex.hpp"
0018 #include "Acts/Vertexing/Vertex.hpp"
0019 #include "ActsExamples/EventData/ProtoVertex.hpp"
0020 #include "ActsExamples/Framework/AlgorithmContext.hpp"
0021
0022 #include <chrono>
0023 #include <ostream>
0024 #include <stdexcept>
0025 #include <system_error>
0026
0027 #include "VertexingHelpers.hpp"
0028
0029 ActsExamples::IterativeVertexFinderAlgorithm::IterativeVertexFinderAlgorithm(
0030 const Config& config, Acts::Logging::Level level)
0031 : ActsExamples::IAlgorithm("IterativeVertexFinder", level), m_cfg(config) {
0032 if (m_cfg.inputTrackParameters.empty()) {
0033 throw std::invalid_argument("Missing input track parameter collection");
0034 }
0035 if (m_cfg.outputProtoVertices.empty()) {
0036 throw std::invalid_argument("Missing output proto vertices collection");
0037 }
0038 if (m_cfg.outputVertices.empty()) {
0039 throw std::invalid_argument("Missing output vertices collection");
0040 }
0041
0042 m_inputTrackParameters.initialize(m_cfg.inputTrackParameters);
0043 m_outputProtoVertices.initialize(m_cfg.outputProtoVertices);
0044 m_outputVertices.initialize(m_cfg.outputVertices);
0045 }
0046
0047 ActsExamples::ProcessCode ActsExamples::IterativeVertexFinderAlgorithm::execute(
0048 const ActsExamples::AlgorithmContext& ctx) const {
0049
0050
0051 const auto& inputTrackParameters = m_inputTrackParameters(ctx);
0052
0053 auto inputTracks = makeInputTracks(inputTrackParameters);
0054
0055 if (inputTrackParameters.size() != inputTracks.size()) {
0056 ACTS_ERROR("Input track containers do not align: "
0057 << inputTrackParameters.size() << " != " << inputTracks.size());
0058 }
0059
0060 for (const auto& trk : inputTrackParameters) {
0061 if (trk.covariance() && trk.covariance()->determinant() <= 0) {
0062
0063
0064 ACTS_WARNING("input track " << trk << " has det(cov) = "
0065 << trk.covariance()->determinant());
0066 }
0067 }
0068
0069
0070 Acts::EigenStepper<> stepper(m_cfg.bField);
0071
0072
0073 auto propagator = std::make_shared<Propagator>(
0074 stepper, Acts::VoidNavigator{}, logger().cloneWithSuffix("Propagator"));
0075
0076 Fitter::Config vertexFitterCfg;
0077 vertexFitterCfg.extractParameters
0078 .connect<&Acts::InputTrack::extractParameters>();
0079
0080 Linearizer::Config linearizerCfg;
0081 linearizerCfg.bField = m_cfg.bField;
0082 linearizerCfg.propagator = propagator;
0083 Linearizer linearizer(linearizerCfg,
0084 logger().cloneWithSuffix("HelicalTrackLinearizer"));
0085
0086 vertexFitterCfg.trackLinearizer.connect<&Linearizer::linearizeTrack>(
0087 &linearizer);
0088 Fitter vertexFitter(vertexFitterCfg,
0089 logger().cloneWithSuffix("FullBilloirVertexFitter"));
0090
0091
0092 Acts::ImpactPointEstimator::Config ipEstCfg(m_cfg.bField, propagator);
0093 Acts::ImpactPointEstimator ipEst(
0094 ipEstCfg, logger().cloneWithSuffix("ImpactPointEstimator"));
0095
0096 Acts::GaussianTrackDensity::Config densityCfg;
0097 densityCfg.extractParameters.connect<&Acts::InputTrack::extractParameters>();
0098 auto seeder = std::make_shared<Seeder>(Seeder::Config{{densityCfg}});
0099
0100 Finder::Config finderCfg(std::move(vertexFitter), seeder, ipEst);
0101 finderCfg.trackLinearizer.connect<&Linearizer::linearizeTrack>(&linearizer);
0102
0103 finderCfg.maxVertices = 200;
0104 finderCfg.reassignTracksAfterFirstFit = false;
0105 finderCfg.extractParameters.connect<&Acts::InputTrack::extractParameters>();
0106 finderCfg.field = m_cfg.bField;
0107 Finder finder(std::move(finderCfg), logger().clone());
0108 Acts::IVertexFinder::State state{std::in_place_type<Finder::State>,
0109 *m_cfg.bField, ctx.magFieldContext};
0110 Options finderOpts(ctx.geoContext, ctx.magFieldContext);
0111
0112
0113 auto result = finder.find(inputTracks, finderOpts, state);
0114
0115 VertexCollection vertices;
0116 if (result.ok()) {
0117 vertices = std::move(result.value());
0118 } else {
0119 ACTS_ERROR("Error in vertex finder: " << result.error().message());
0120 }
0121
0122
0123 ACTS_INFO("Found " << vertices.size() << " vertices in event");
0124 for (const auto& vtx : vertices) {
0125 ACTS_INFO("Found vertex at " << vtx.fullPosition().transpose() << " with "
0126 << vtx.tracks().size() << " tracks.");
0127 }
0128
0129
0130 m_outputProtoVertices(ctx, makeProtoVertices(inputTracks, vertices));
0131
0132
0133 m_outputVertices(ctx, std::move(vertices));
0134
0135 return ActsExamples::ProcessCode::SUCCESS;
0136 }