File indexing completed on 2025-01-18 09:11:43
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/SympyStepper.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 namespace ActsExamples {
0030
0031 IterativeVertexFinderAlgorithm::IterativeVertexFinderAlgorithm(
0032 const Config& config, Acts::Logging::Level level)
0033 : IAlgorithm("IterativeVertexFinder", level), m_cfg(config) {
0034 if (m_cfg.inputTrackParameters.empty()) {
0035 throw std::invalid_argument("Missing input track parameter collection");
0036 }
0037 if (m_cfg.outputProtoVertices.empty()) {
0038 throw std::invalid_argument("Missing output proto vertices collection");
0039 }
0040 if (m_cfg.outputVertices.empty()) {
0041 throw std::invalid_argument("Missing output vertices collection");
0042 }
0043
0044 m_inputTrackParameters.initialize(m_cfg.inputTrackParameters);
0045 m_outputProtoVertices.initialize(m_cfg.outputProtoVertices);
0046 m_outputVertices.initialize(m_cfg.outputVertices);
0047 }
0048
0049 ProcessCode IterativeVertexFinderAlgorithm::execute(
0050 const AlgorithmContext& ctx) const {
0051
0052
0053 const auto& inputTrackParameters = m_inputTrackParameters(ctx);
0054
0055 auto inputTracks = makeInputTracks(inputTrackParameters);
0056
0057 if (inputTrackParameters.size() != inputTracks.size()) {
0058 ACTS_ERROR("Input track containers do not align: "
0059 << inputTrackParameters.size() << " != " << inputTracks.size());
0060 }
0061
0062 for (const auto& trk : inputTrackParameters) {
0063 if (trk.covariance() && trk.covariance()->determinant() <= 0) {
0064
0065
0066 ACTS_WARNING("input track " << trk << " has det(cov) = "
0067 << trk.covariance()->determinant());
0068 }
0069 }
0070
0071
0072 Acts::SympyStepper stepper(m_cfg.bField);
0073
0074
0075 auto propagator = std::make_shared<Propagator>(
0076 stepper, Acts::VoidNavigator{}, logger().cloneWithSuffix("Propagator"));
0077
0078 Fitter::Config vertexFitterCfg;
0079 vertexFitterCfg.extractParameters
0080 .connect<&Acts::InputTrack::extractParameters>();
0081
0082 Linearizer::Config linearizerCfg;
0083 linearizerCfg.bField = m_cfg.bField;
0084 linearizerCfg.propagator = propagator;
0085 Linearizer linearizer(linearizerCfg,
0086 logger().cloneWithSuffix("HelicalTrackLinearizer"));
0087
0088 vertexFitterCfg.trackLinearizer.connect<&Linearizer::linearizeTrack>(
0089 &linearizer);
0090 Fitter vertexFitter(vertexFitterCfg,
0091 logger().cloneWithSuffix("FullBilloirVertexFitter"));
0092
0093
0094 Acts::ImpactPointEstimator::Config ipEstCfg(m_cfg.bField, propagator);
0095 Acts::ImpactPointEstimator ipEst(
0096 ipEstCfg, logger().cloneWithSuffix("ImpactPointEstimator"));
0097
0098 Acts::GaussianTrackDensity::Config densityCfg;
0099 densityCfg.extractParameters.connect<&Acts::InputTrack::extractParameters>();
0100 auto seeder = std::make_shared<Seeder>(Seeder::Config{{densityCfg}});
0101
0102 Finder::Config finderCfg(std::move(vertexFitter), seeder, ipEst);
0103 finderCfg.trackLinearizer.connect<&Linearizer::linearizeTrack>(&linearizer);
0104
0105 finderCfg.maxVertices = m_cfg.maxIterations;
0106 finderCfg.reassignTracksAfterFirstFit = false;
0107 finderCfg.extractParameters.connect<&Acts::InputTrack::extractParameters>();
0108 finderCfg.field = m_cfg.bField;
0109 Finder finder(std::move(finderCfg), logger().clone());
0110 Acts::IVertexFinder::State state{std::in_place_type<Finder::State>,
0111 *m_cfg.bField, ctx.magFieldContext};
0112 Options finderOpts(ctx.geoContext, ctx.magFieldContext);
0113
0114
0115 auto result = finder.find(inputTracks, finderOpts, state);
0116
0117 VertexCollection vertices;
0118 if (result.ok()) {
0119 vertices = std::move(result.value());
0120 } else {
0121 ACTS_ERROR("Error in vertex finder: " << result.error().message());
0122 }
0123
0124
0125 ACTS_INFO("Found " << vertices.size() << " vertices in event");
0126 for (const auto& vtx : vertices) {
0127 ACTS_DEBUG("Found vertex at " << vtx.fullPosition().transpose() << " with "
0128 << vtx.tracks().size() << " tracks.");
0129 }
0130
0131
0132 m_outputProtoVertices(ctx, makeProtoVertices(inputTracks, vertices));
0133
0134
0135 m_outputVertices(ctx, std::move(vertices));
0136
0137 return ProcessCode::SUCCESS;
0138 }
0139
0140 }