Skip to content

Commit

Permalink
refactor: Factorize track selector tool from algorithm (#2267)
Browse files Browse the repository at this point in the history
This separates out the track selection logic itself from the `TrackSelectorAlgorithm` (formerly called just `TrackSelector`). This is in preparation for #2206.
  • Loading branch information
paulgessinger committed Jul 5, 2023
1 parent f118922 commit c0b9fcc
Show file tree
Hide file tree
Showing 9 changed files with 255 additions and 74 deletions.
104 changes: 104 additions & 0 deletions Core/include/Acts/TrackFinding/TrackSelector.hpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,104 @@
// This file is part of the Acts project.
//
// Copyright (C) 2023 CERN for the benefit of the Acts project
//
// This Source Code Form is subject to the terms of the Mozilla Public
// License, v. 2.0. If a copy of the MPL was not distributed with this
// file, You can obtain one at http://mozilla.org/MPL/2.0/.

#pragma once

#include <cmath>
#include <limits>

namespace Acts {

/// Class which performs filtering of tracks. It accepts an input and an output
/// track container and uses the built-in copy facility to copy tracks into the
/// output container.
class TrackSelector {
public:
struct Config {
// Minimum/maximum local positions.
double loc0Min = -std::numeric_limits<double>::infinity();
double loc0Max = std::numeric_limits<double>::infinity();
double loc1Min = -std::numeric_limits<double>::infinity();
double loc1Max = std::numeric_limits<double>::infinity();
// Minimum/maximum track time.
double timeMin = -std::numeric_limits<double>::infinity();
double timeMax = std::numeric_limits<double>::infinity();
// Direction cuts.
double phiMin = -std::numeric_limits<double>::infinity();
double phiMax = std::numeric_limits<double>::infinity();
double etaMin = -std::numeric_limits<double>::infinity();
double etaMax = std::numeric_limits<double>::infinity();
double absEtaMin = 0.0;
double absEtaMax = std::numeric_limits<double>::infinity();
// Momentum cuts.
double ptMin = 0.0;
double ptMax = std::numeric_limits<double>::infinity();

std::size_t minMeasurements = 0;
};

/// Constructor from a config object
/// @param config is the configuration object
TrackSelector(const Config& config) : m_cfg(config) {}

/// Select tracks from an input container and copy them into an output
/// container
/// @tparam input_tracks_t is the type of the input track container
/// @tparam output_tracks_t is the type of the output track container
/// @param inputTracks is the input track container
/// @param outputTracks is the output track container
template <typename input_tracks_t, typename output_tracks_t>
void selectTracks(const input_tracks_t& inputTracks,
output_tracks_t& outputTracks) const;

/// Helper function to check if a track is valid
/// @tparam track_proxy_t is the type of the track proxy
/// @param track is the track proxy
/// @return true if the track is valid
template <typename track_proxy_t>
bool isValidTrack(const track_proxy_t& track) const;

/// Get readonly access to the config parameters
/// @return the config object
const Config& config() const { return m_cfg; }

private:
Config m_cfg;
};

template <typename input_tracks_t, typename output_tracks_t>
void TrackSelector::selectTracks(const input_tracks_t& inputTracks,
output_tracks_t& outputTracks) const {
for (auto track : inputTracks) {
if (!isValidTrack(track)) {
continue;
}
auto destProxy = outputTracks.getTrack(outputTracks.addTrack());
destProxy.copyFrom(track, false);
destProxy.tipIndex() = track.tipIndex();
}
}

template <typename track_proxy_t>
bool TrackSelector::isValidTrack(const track_proxy_t& track) const {
auto checkMin = [](auto x, auto min) { return min <= x; };
auto within = [](double x, double min, double max) {
return (min <= x) and (x < max);
};
const auto theta = track.theta();
const auto eta = -std::log(std::tan(theta / 2));
return within(track.transverseMomentum(), m_cfg.ptMin, m_cfg.ptMax) and
within(std::abs(eta), m_cfg.absEtaMin, m_cfg.absEtaMax) and
within(eta, m_cfg.etaMin, m_cfg.etaMax) and
within(track.phi(), m_cfg.phiMin, m_cfg.phiMax) and
within(track.loc0(), m_cfg.loc0Min, m_cfg.loc0Max) and
within(track.loc1(), m_cfg.loc1Min, m_cfg.loc1Max) and
within(track.time(), m_cfg.timeMin, m_cfg.timeMax) and
checkMin(track.nMeasurements(), m_cfg.minMeasurements);
}

} // namespace Acts
1 change: 0 additions & 1 deletion Examples/Algorithms/TruthTracking/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -2,7 +2,6 @@ add_library(
ActsExamplesTruthTracking SHARED
ActsExamples/TruthTracking/ParticleSelector.cpp
ActsExamples/TruthTracking/ParticleSmearing.cpp
ActsExamples/TruthTracking/TrackSelector.cpp
ActsExamples/TruthTracking/TrackParameterSelector.cpp
ActsExamples/TruthTracking/TrackModifier.cpp
ActsExamples/TruthTracking/TruthSeedSelector.cpp
Expand Down
1 change: 1 addition & 0 deletions Examples/Algorithms/Utilities/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -3,6 +3,7 @@ add_library(
src/PrototracksToSeeds.cpp
src/SeedsToPrototracks.cpp
src/TrajectoriesToPrototracks.cpp
src/TrackSelectorAlgorithm.cpp
src/TracksToTrajectories.cpp)
target_include_directories(
ActsExamplesUtilities
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -8,6 +8,7 @@

#pragma once

#include "Acts/TrackFinding/TrackSelector.hpp"
#include "Acts/Utilities/Logger.hpp"
#include "ActsExamples/EventData/Track.hpp"
#include "ActsExamples/Framework/DataHandle.hpp"
Expand All @@ -22,37 +23,18 @@ namespace ActsExamples {
struct AlgorithmContext;

/// Select tracks by applying some selection cuts.
class TrackSelector final : public IAlgorithm {
class TrackSelectorAlgorithm final : public IAlgorithm {
public:
struct Config {
/// Input track collection.
std::string inputTracks;
/// Output track collection
std::string outputTracks;

// Minimum/maximum local positions.
double loc0Min = -std::numeric_limits<double>::infinity();
double loc0Max = std::numeric_limits<double>::infinity();
double loc1Min = -std::numeric_limits<double>::infinity();
double loc1Max = std::numeric_limits<double>::infinity();
// Minimum/maximum track time.
double timeMin = -std::numeric_limits<double>::infinity();
double timeMax = std::numeric_limits<double>::infinity();
// Direction cuts.
double phiMin = -std::numeric_limits<double>::infinity();
double phiMax = std::numeric_limits<double>::infinity();
double etaMin = -std::numeric_limits<double>::infinity();
double etaMax = std::numeric_limits<double>::infinity();
double absEtaMin = 0.0;
double absEtaMax = std::numeric_limits<double>::infinity();
// Momentum cuts.
double ptMin = 0.0;
double ptMax = std::numeric_limits<double>::infinity();

std::size_t minMeasurements = 0;
Acts::TrackSelector::Config selectorConfig;
};

TrackSelector(const Config& config, Acts::Logging::Level level);
TrackSelectorAlgorithm(const Config& config, Acts::Logging::Level level);

ProcessCode execute(const AlgorithmContext& ctx) const final;

Expand All @@ -62,6 +44,8 @@ class TrackSelector final : public IAlgorithm {
private:
Config m_cfg;

Acts::TrackSelector m_selector;

ReadDataHandle<ConstTrackContainer> m_inputTrackContainer{this,
"InputTracks"};
WriteDataHandle<ConstTrackContainer> m_outputTrackContainer{this,
Expand Down
76 changes: 76 additions & 0 deletions Examples/Algorithms/Utilities/src/TrackSelectorAlgorithm.cpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,76 @@
// This file is part of the Acts project.
//
// Copyright (C) 2019-2023 CERN for the benefit of the Acts project
//
// This Source Code Form is subject to the terms of the Mozilla Public
// License, v. 2.0. If a copy of the MPL was not distributed with this
// file, You can obtain one at http://mozilla.org/MPL/2.0/.

#include "ActsExamples/Utilities/TrackSelectorAlgorithm.hpp"

#include "Acts/EventData/TrackContainer.hpp"
#include "Acts/EventData/TrackProxy.hpp"
#include "Acts/EventData/VectorMultiTrajectory.hpp"
#include "Acts/EventData/VectorTrackContainer.hpp"
#include "ActsExamples/EventData/Track.hpp"

#include <cmath>
#include <memory>
#include <stdexcept>
#include <utility>

namespace ActsExamples {
struct AlgorithmContext;
} // namespace ActsExamples

ActsExamples::TrackSelectorAlgorithm::TrackSelectorAlgorithm(
const Config& config, Acts::Logging::Level level)
: IAlgorithm("TrackSelector", level),
m_cfg(config),
m_selector(config.selectorConfig) {
if (m_cfg.inputTracks.empty()) {
throw std::invalid_argument("Input track collection is empty");
}

if (m_cfg.outputTracks.empty()) {
throw std::invalid_argument("Output track collection is empty");
}

m_inputTrackContainer.initialize(m_cfg.inputTracks);
m_outputTrackContainer.initialize(m_cfg.outputTracks);
}

ActsExamples::ProcessCode ActsExamples::TrackSelectorAlgorithm::execute(
const ActsExamples::AlgorithmContext& ctx) const {
ACTS_VERBOSE("Reading tracks from: " << m_cfg.inputTracks);

const auto& inputTracks = m_inputTrackContainer(ctx);

std::shared_ptr<Acts::ConstVectorMultiTrajectory> trackStateContainer =
inputTracks.trackStateContainerHolder();

auto trackContainer = std::make_shared<Acts::VectorTrackContainer>();

// temporary empty track state container: we don't change the original one,
// but we need one for filtering
auto tempTrackStateContainer =
std::make_shared<Acts::VectorMultiTrajectory>();

TrackContainer filteredTracks{trackContainer, tempTrackStateContainer};
filteredTracks.ensureDynamicColumns(inputTracks);

ACTS_DEBUG("Track container size before filtering: " << inputTracks.size());

m_selector.selectTracks(inputTracks, filteredTracks);

ACTS_DEBUG("Track container size after filtering: " << filteredTracks.size());

ConstTrackContainer outputTracks{
std::make_shared<Acts::ConstVectorTrackContainer>(
std::move(*trackContainer)),
trackStateContainer};

m_outputTrackContainer(ctx, std::move(outputTracks));

return ProcessCode::SUCCESS;
}
16 changes: 10 additions & 6 deletions Examples/Python/python/acts/examples/reconstruction.py
Original file line number Diff line number Diff line change
Expand Up @@ -1120,13 +1120,10 @@ def addTrackSelection(
inputTracks: str,
outputTracks: str,
logLevel: Optional[acts.logging.Level] = None,
) -> acts.examples.TrackSelector:
) -> acts.examples.TrackSelectorAlgorithm:
customLogLevel = acts.examples.defaultLogging(s, logLevel)

trackSelector = acts.examples.TrackSelector(
level=customLogLevel(),
inputTracks=inputTracks,
outputTracks=outputTracks,
selectorConfig = acts.TrackSelector.Config(
**acts.examples.defaultKWArgs(
loc0Min=trackSelectorConfig.loc0[0],
loc0Max=trackSelectorConfig.loc0[1],
Expand All @@ -1143,7 +1140,14 @@ def addTrackSelection(
ptMin=trackSelectorConfig.pt[0],
ptMax=trackSelectorConfig.pt[1],
minMeasurements=trackSelectorConfig.nMeasurementsMin,
),
)
)

trackSelector = acts.examples.TrackSelectorAlgorithm(
level=customLogLevel(),
inputTracks=inputTracks,
outputTracks=outputTracks,
selectorConfig=selectorConfig,
)

s.addAlgorithm(trackSelector)
Expand Down
57 changes: 56 additions & 1 deletion Examples/Python/src/ExampleAlgorithms.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -8,12 +8,14 @@

#include "Acts/Definitions/Algebra.hpp"
#include "Acts/Plugins/Python/Utilities.hpp"
#include "Acts/TrackFinding/TrackSelector.hpp"
#include "ActsExamples/Fatras/FatrasSimulation.hpp"
#include "ActsExamples/Io/Json/JsonGeometryList.hpp"
#include "ActsExamples/Printers/HitsPrinter.hpp"
#include "ActsExamples/Printers/ParticlesPrinter.hpp"
#include "ActsExamples/Printers/TrackParametersPrinter.hpp"
#include "ActsExamples/Utilities/Range.hpp"
#include "ActsExamples/Utilities/TrackSelectorAlgorithm.hpp"

#include <vector>

Expand All @@ -28,7 +30,7 @@ using namespace Acts;
namespace Acts::Python {

void addExampleAlgorithms(Context& ctx) {
auto mex = ctx.get("examples");
auto [m, mex] = ctx.get("main", "examples");

mex.def("readJsonGeometryList", ActsExamples::readJsonGeometryList);

Expand All @@ -50,6 +52,59 @@ void addExampleAlgorithms(Context& ctx) {

ACTS_PYTHON_DECLARE_ALGORITHM(ActsExamples::TrackParametersPrinter, mex,
"TrackParametersPrinter", inputTrackParameters);

{
using Alg = ActsExamples::TrackSelectorAlgorithm;
using Config = Alg::Config;

auto alg = py::class_<Alg, IAlgorithm, std::shared_ptr<Alg>>(
mex, "TrackSelectorAlgorithm")
.def(py::init<const Alg::Config&, Acts::Logging::Level>(),
py::arg("config"), py::arg("level"))
.def_property_readonly("config", &Alg::config);

auto c = py::class_<Config>(alg, "Config").def(py::init<>());

ACTS_PYTHON_STRUCT_BEGIN(c, Config);
ACTS_PYTHON_MEMBER(inputTracks);
ACTS_PYTHON_MEMBER(outputTracks);
ACTS_PYTHON_MEMBER(selectorConfig);
ACTS_PYTHON_STRUCT_END();
}

{
using Config = Acts::TrackSelector::Config;
auto tool = py::class_<Acts::TrackSelector>(m, "TrackSelector")
.def(py::init<const Config&>(), py::arg("config"));

auto c = py::class_<Config>(tool, "Config").def(py::init<>());

ACTS_PYTHON_STRUCT_BEGIN(c, Config);
ACTS_PYTHON_MEMBER(loc0Min);
ACTS_PYTHON_MEMBER(loc0Max);
ACTS_PYTHON_MEMBER(loc1Min);
ACTS_PYTHON_MEMBER(loc1Max);
ACTS_PYTHON_MEMBER(timeMin);
ACTS_PYTHON_MEMBER(timeMax);
ACTS_PYTHON_MEMBER(phiMin);
ACTS_PYTHON_MEMBER(phiMax);
ACTS_PYTHON_MEMBER(etaMin);
ACTS_PYTHON_MEMBER(etaMax);
ACTS_PYTHON_MEMBER(absEtaMin);
ACTS_PYTHON_MEMBER(absEtaMax);
ACTS_PYTHON_MEMBER(ptMin);
ACTS_PYTHON_MEMBER(ptMax);
ACTS_PYTHON_MEMBER(minMeasurements);
ACTS_PYTHON_STRUCT_END();

pythonRangeProperty(c, "loc0", &Config::loc0Min, &Config::loc0Max);
pythonRangeProperty(c, "loc1", &Config::loc1Min, &Config::loc1Max);
pythonRangeProperty(c, "time", &Config::timeMin, &Config::timeMax);
pythonRangeProperty(c, "phi", &Config::phiMin, &Config::phiMax);
pythonRangeProperty(c, "eta", &Config::etaMin, &Config::etaMax);
pythonRangeProperty(c, "absEta", &Config::absEtaMin, &Config::absEtaMax);
pythonRangeProperty(c, "pt", &Config::ptMin, &Config::ptMax);
}
}

} // namespace Acts::Python

0 comments on commit c0b9fcc

Please sign in to comment.