Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

feat: Validate track parameters #3756

Merged
merged 20 commits into from
Oct 29, 2024
Merged
Show file tree
Hide file tree
Changes from 17 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
7 changes: 5 additions & 2 deletions Core/include/Acts/EventData/GenericBoundTrackParameters.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -9,6 +9,7 @@
#pragma once

#include "Acts/Definitions/Tolerance.hpp"
#include "Acts/EventData/TrackParameterHelpers.hpp"
#include "Acts/EventData/TransformationHelpers.hpp"
#include "Acts/EventData/detail/PrintParameters.hpp"
#include "Acts/Surfaces/Surface.hpp"
Expand All @@ -18,7 +19,6 @@
#include <cassert>
#include <cmath>
#include <memory>
#include <type_traits>

namespace Acts {

Expand Down Expand Up @@ -61,7 +61,10 @@ class GenericBoundTrackParameters {
m_cov(std::move(cov)),
m_surface(std::move(surface)),
m_particleHypothesis(std::move(particleHypothesis)) {
assert(m_surface);
// TODO set `validateAngleRange` to `true` after fixing caller code
assert(isBoundVectorValid(m_params, false) &&
"Invalid bound parameters vector");
assert(m_surface != nullptr && "Reference surface must not be null");
normalizePhiTheta();
}

Expand Down
11 changes: 8 additions & 3 deletions Core/include/Acts/EventData/GenericFreeTrackParameters.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -11,17 +11,16 @@
#include "Acts/Definitions/Algebra.hpp"
#include "Acts/Definitions/Common.hpp"
#include "Acts/Definitions/TrackParametrization.hpp"
#include "Acts/EventData/TrackParameterHelpers.hpp"
#include "Acts/EventData/TrackParametersConcept.hpp"
#include "Acts/EventData/TransformationHelpers.hpp"
#include "Acts/EventData/detail/PrintParameters.hpp"
#include "Acts/Utilities/MathHelpers.hpp"
#include "Acts/Utilities/UnitVectors.hpp"
#include "Acts/Utilities/VectorHelpers.hpp"

#include <cassert>
#include <cmath>
#include <optional>
#include <type_traits>

namespace Acts {

Expand Down Expand Up @@ -55,7 +54,9 @@ class GenericFreeTrackParameters {
ParticleHypothesis particleHypothesis)
: m_params(params),
m_cov(std::move(cov)),
m_particleHypothesis(std::move(particleHypothesis)) {}
m_particleHypothesis(std::move(particleHypothesis)) {
assert(isFreeVectorValid(m_params) && "Invalid free parameters vector");
}

/// Construct from four-position, direction, absolute momentum, and charge.
///
Expand All @@ -78,6 +79,8 @@ class GenericFreeTrackParameters {
m_params[eFreeDir1] = dir[eMom1];
m_params[eFreeDir2] = dir[eMom2];
m_params[eFreeQOverP] = qOverP;

assert(isFreeVectorValid(m_params) && "Invalid free parameters vector");
}

/// Construct from four-position, angles, absolute momentum, and charge.
Expand All @@ -103,6 +106,8 @@ class GenericFreeTrackParameters {
m_params[eFreeDir1] = dir[eMom1];
m_params[eFreeDir2] = dir[eMom2];
m_params[eFreeQOverP] = qOverP;

assert(isFreeVectorValid(m_params) && "Invalid free parameters vector");
}

/// Converts a free track parameter with a different hypothesis.
Expand Down
30 changes: 30 additions & 0 deletions Core/include/Acts/EventData/TrackParameterHelpers.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -11,8 +11,38 @@
#include "Acts/Definitions/TrackParametrization.hpp"
#include "Acts/Utilities/detail/periodic.hpp"

#include <limits>

namespace Acts {

/// Check if a bound vector is valid. This checks the following:
/// - All values are finite
/// - (optionally) The phi value is in the range [-pi, pi)
/// - (optionally) The theta value is in the range [0, pi]
///
/// @param v The bound vector to check
/// @param validateAngleRange If true, the phi and theta values are range checked
/// @param epsilon The epsilon to use for the checks
/// @param maxAbsEta The maximum allowed eta value
///
/// @return True if the bound vector is valid
bool isBoundVectorValid(
const BoundVector& v, bool validateAngleRange, double epsilon = 1e-6,
double maxAbsEta = std::numeric_limits<double>::infinity());

/// Check if a free vector is valid. This checks the following:
/// - All values are finite
/// - Direction is normalized
///
/// @param v The free vector to check
/// @param epsilon The epsilon to use for the checks
/// @param maxAbsEta The maximum allowed eta value
///
/// @return True if the free vector is valid
bool isFreeVectorValid(
const FreeVector& v, double epsilon = 1e-6,
double maxAbsEta = std::numeric_limits<double>::infinity());

/// Normalize the bound parameter angles
///
/// @param boundParams The bound parameters to normalize
Expand Down
3 changes: 1 addition & 2 deletions Core/include/Acts/EventData/TrackParametersConcept.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -12,9 +12,7 @@
#include "Acts/Definitions/TrackParametrization.hpp"
#include "Acts/Geometry/GeometryContext.hpp"

#include <concepts>
#include <optional>
#include <type_traits>

namespace Acts {

Expand Down Expand Up @@ -68,4 +66,5 @@ concept BoundTrackParametersConcept =
{ p.position(c) } -> std::same_as<Vector3>;
};
};

} // namespace Acts
8 changes: 8 additions & 0 deletions Core/include/Acts/EventData/detail/GenerateParameters.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -195,6 +195,10 @@ struct GenerateBoundParametersOptions {
GenerateQoverPOptions qOverP;
};

inline BoundVector someBoundParametersA() {
return {0.0, 0.0, 0.0, std::numbers::pi / 2, 1.0, 0.0};
}

template <typename generator_t>
inline BoundVector generateBoundParameters(
generator_t& rng, const GenerateBoundParametersOptions& options) {
Expand Down Expand Up @@ -242,6 +246,10 @@ struct GenerateFreeParametersOptions {
GenerateQoverPOptions qOverP;
};

inline FreeVector someFreeParametersA() {
return {0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 1.0};
}

template <typename generator_t>
inline FreeVector generateFreeParameters(
generator_t& rng, const GenerateFreeParametersOptions& options) {
Expand Down
6 changes: 3 additions & 3 deletions Core/include/Acts/EventData/detail/TestTrackState.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -38,11 +38,11 @@ struct TestTrackState {
: surface(
CurvilinearSurface(Vector3::Zero(), Vector3::UnitZ()).surface()),
// set bogus parameters first since they are not default-constructible
predicted(surface, BoundVector::Zero(), std::nullopt,
predicted(surface, someBoundParametersA(), std::nullopt,
ParticleHypothesis::pion()),
filtered(surface, BoundVector::Zero(), std::nullopt,
filtered(surface, someBoundParametersA(), std::nullopt,
ParticleHypothesis::pion()),
smoothed(surface, BoundVector::Zero(), std::nullopt,
smoothed(surface, someBoundParametersA(), std::nullopt,
ParticleHypothesis::pion()),
jacobian(BoundMatrix::Identity()),
chi2(std::chi_squared_distribution<double>(measdim)(rng)),
Expand Down
42 changes: 42 additions & 0 deletions Core/include/Acts/Utilities/AngleHelpers.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -9,16 +9,51 @@
#pragma once

#include <cmath>
#include <numbers>

namespace Acts::AngleHelpers {

template <typename Scalar>
struct EtaThetaConversionTraits {
static constexpr Scalar minTheta = 0;
static constexpr Scalar maxTheta = std::numbers::pi_v<Scalar>;

static constexpr Scalar minEta = -std::numeric_limits<Scalar>::infinity();
static constexpr Scalar maxEta = std::numeric_limits<Scalar>::infinity();
};

template <>
struct EtaThetaConversionTraits<float> {
static constexpr float minTheta = 0;
static constexpr float maxTheta = std::numbers::pi_v<float> - 1e-6f;
andiwand marked this conversation as resolved.
Show resolved Hide resolved

static constexpr float minEta = -std::numeric_limits<float>::infinity();
static constexpr float maxEta = std::numeric_limits<float>::infinity();
};

template <>
struct EtaThetaConversionTraits<double> {
static constexpr double minTheta = 0;
static constexpr double maxTheta = std::numbers::pi;

static constexpr double minEta = -std::numeric_limits<double>::infinity();
static constexpr double maxEta = std::numeric_limits<double>::infinity();
};

/// Calculate the pseudorapidity from the polar angle theta.
///
/// @param theta is the polar angle in radian towards the z-axis.
///
/// @return the pseudorapidity towards the z-axis.
template <typename Scalar>
Scalar etaFromTheta(Scalar theta) {
if (theta < EtaThetaConversionTraits<Scalar>::minTheta) {
return std::numeric_limits<Scalar>::infinity();
}
if (theta > EtaThetaConversionTraits<Scalar>::maxTheta) {
return -std::numeric_limits<Scalar>::infinity();
}

return -std::log(std::tan(0.5 * theta));
}

Expand All @@ -29,6 +64,13 @@ Scalar etaFromTheta(Scalar theta) {
/// @return the polar angle in radian towards the z-axis.
template <typename Scalar>
Scalar thetaFromEta(Scalar eta) {
if (eta < EtaThetaConversionTraits<Scalar>::minEta) {
return std::numbers::pi_v<Scalar>;
}
if (eta > EtaThetaConversionTraits<Scalar>::maxEta) {
return 0;
}

andiwand marked this conversation as resolved.
Show resolved Hide resolved
return 2 * std::atan(std::exp(-eta));
}

Expand Down
1 change: 1 addition & 0 deletions Core/src/EventData/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -8,4 +8,5 @@ target_sources(
TrackStatePropMask.cpp
VectorMultiTrajectory.cpp
VectorTrackContainer.cpp
TrackParameterHelpers.cpp
)
64 changes: 64 additions & 0 deletions Core/src/EventData/TrackParameterHelpers.cpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,64 @@
// This file is part of the ACTS project.
//
// Copyright (C) 2016 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 https://mozilla.org/MPL/2.0/.

#include "Acts/EventData/TrackParameterHelpers.hpp"

#include "Acts/Utilities/AngleHelpers.hpp"
#include "Acts/Utilities/VectorHelpers.hpp"

#include <numbers>

bool Acts::isBoundVectorValid(const BoundVector& v, bool validateAngleRange,
double epsilon, double maxAbsEta) {
constexpr auto pi = std::numbers::pi_v<double>;

bool loc0Valid = std::isfinite(v[eBoundLoc0]);
bool loc1Valid = std::isfinite(v[eBoundLoc1]);
bool phiValid = std::isfinite(v[eBoundPhi]);
bool thetaValid = std::isfinite(v[eBoundTheta]);
bool qOverPValid = std::isfinite(v[eBoundQOverP]);
bool timeValid = std::isfinite(v[eBoundTime]);

if (validateAngleRange) {
phiValid = phiValid && (v[eBoundPhi] + epsilon >= -pi) &&
(v[eBoundPhi] - epsilon < pi);
thetaValid = thetaValid && (v[eBoundTheta] + epsilon >= 0) &&
(v[eBoundTheta] - epsilon <= pi);
}

bool etaValid = true;
if (std::isfinite(maxAbsEta)) {
double eta = AngleHelpers::etaFromTheta(v[eBoundTheta]);
etaValid = std::isfinite(eta) && (std::abs(eta) - epsilon <= maxAbsEta);
}

return loc0Valid && loc1Valid && phiValid && thetaValid && qOverPValid &&
timeValid && etaValid;
}

bool Acts::isFreeVectorValid(const FreeVector& v, double epsilon,
double maxAbsEta) {
bool pos0Valid = std::isfinite(v[eFreePos0]);
bool pos1Valid = std::isfinite(v[eFreePos1]);
bool pos2Valid = std::isfinite(v[eFreePos2]);
bool dir0Valid = std::isfinite(v[eFreeDir0]);
bool dir1Valid = std::isfinite(v[eFreeDir1]);
bool dir2Valid = std::isfinite(v[eFreeDir2]);
bool dirValid = (std::abs(v.segment<3>(eFreeDir0).norm() - 1) - epsilon <= 0);
bool qOverPValid = std::isfinite(v[eFreeQOverP]);
bool timeValid = std::isfinite(v[eFreeTime]);

bool etaValid = true;
if (std::isfinite(maxAbsEta)) {
double eta = VectorHelpers::eta(v.segment<3>(eFreeDir0));
etaValid = std::isfinite(eta) && (std::abs(eta) - epsilon <= maxAbsEta);
}

return pos0Valid && pos1Valid && pos2Valid && dir0Valid && dir1Valid &&
dir2Valid && dirValid && qOverPValid && timeValid && etaValid;
}
10 changes: 10 additions & 0 deletions Tests/UnitTests/Core/EventData/TrackParameterHelpersTests.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -14,6 +14,16 @@

BOOST_AUTO_TEST_SUITE(TrackParameterHelpers)

BOOST_AUTO_TEST_CASE(isBoundVectorValid) {
BOOST_CHECK(!Acts::isBoundVectorValid({1, 2, 3, 4, 5, 6}, true));
BOOST_CHECK(Acts::isBoundVectorValid({1, 2, 1, 1, 5, 6}, true));
}

BOOST_AUTO_TEST_CASE(isFreeVectorValid) {
BOOST_CHECK(!Acts::isFreeVectorValid({1, 2, 3, 4, 5, 6, 7, 8}));
BOOST_CHECK(Acts::isFreeVectorValid({1, 2, 3, 4, 1, 0, 0, 8}));
}

BOOST_AUTO_TEST_CASE(normalizeBoundParameters) {
CHECK_CLOSE_OR_SMALL(Acts::normalizeBoundParameters({1, 2, 3, 4, 5, 6}),
Acts::BoundVector(1, 2, -0.141593, 2.28319, 5, 6), 1e-3,
Expand Down
Loading
Loading