Skip to content

Commit

Permalink
feat: Experiment cuts to dublets in seed finder (#2326)
Browse files Browse the repository at this point in the history
This PR adds an option to include experimental specific cuts to discard bottom-middle dublets in a certain (r, eta) region of the detector.

This PR is related to the implementation of the ITk "fast tracking" configuration from PR #2166
(I am splitting the changes to make it more organised @CarloVarni @noemina)

For the integration with Athena we will need to implement something like that in Athena:
```
m_finderCfg.experimentCuts.connect(
	[](const void*, const float& bottomRadius, const float& cotTheta) -> bool {

        if (bottomRadius < fastTrackingRMin and
               (cotTheta > fastTrackingCotThetaMax or
                cotTheta < -fastTrackingCotThetaMax)) {
             return false;
        }
        return true;
}
```

where:

```
float fastTrackingRMin = 50. * Acts::UnitConstants::mm;
float fastTrackingCotThetaMax = 1.5;
```
  • Loading branch information
LuisFelipeCoelho authored Nov 24, 2023
1 parent b77e16e commit f32f90f
Show file tree
Hide file tree
Showing 5 changed files with 47 additions and 4 deletions.
5 changes: 5 additions & 0 deletions Core/include/Acts/Seeding/SeedConfirmationRangeConfig.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -15,6 +15,11 @@

namespace Acts {

// defaults experimental cuts to no operation in both seeding algorithms
inline bool noopExperimentCuts(float /*bottomRadius*/, float /*cotTheta*/) {
return true;
}

/// @brief contains parameters for seed confirmation
struct SeedConfirmationRangeConfig {
// z minimum and maximum of middle component of the seed used to define the
Expand Down
28 changes: 24 additions & 4 deletions Core/include/Acts/Seeding/SeedFinder.ipp
Original file line number Diff line number Diff line change
Expand Up @@ -206,7 +206,11 @@ SeedFinder<external_spacepoint_t, platform_t>::getCompatibleDoublets(
const float deltaRMinSP, const float deltaRMaxSP, const float uIP,
const float uIP2, const float cosPhiM, const float sinPhiM) const {
float impactMax = m_config.impactMax;
if constexpr (candidateType == Acts::SpacePointCandidateType::eBottom) {

constexpr bool isBottomCandidate =
candidateType == Acts::SpacePointCandidateType::eBottom;

if constexpr (isBottomCandidate) {
impactMax = -impactMax;
}

Expand Down Expand Up @@ -252,7 +256,7 @@ SeedFinder<external_spacepoint_t, platform_t>::getCompatibleDoublets(
// the iterator so we don't need to look at the other SPs again
for (; min_itr != otherSPs.end(); ++min_itr) {
const auto& otherSP = *min_itr;
if constexpr (candidateType == Acts::SpacePointCandidateType::eBottom) {
if constexpr (isBottomCandidate) {
// if r-distance is too big, try next SP in bin
if ((rM - otherSP->radius()) <= deltaRMaxSP) {
break;
Expand All @@ -272,7 +276,7 @@ SeedFinder<external_spacepoint_t, platform_t>::getCompatibleDoublets(
for (; min_itr != otherSPs.end(); ++min_itr) {
const auto& otherSP = *min_itr;

if constexpr (candidateType == Acts::SpacePointCandidateType::eBottom) {
if constexpr (isBottomCandidate) {
deltaR = (rM - otherSP->radius());

// if r-distance is too small, try next SP in bin
Expand All @@ -288,7 +292,7 @@ SeedFinder<external_spacepoint_t, platform_t>::getCompatibleDoublets(
}
}

if constexpr (candidateType == Acts::SpacePointCandidateType::eBottom) {
if constexpr (isBottomCandidate) {
deltaZ = (zM - otherSP->z());
} else {
deltaZ = (otherSP->z() - zM);
Expand Down Expand Up @@ -379,6 +383,14 @@ SeedFinder<external_spacepoint_t, platform_t>::getCompatibleDoublets(
const float iDeltaR = std::sqrt(iDeltaR2);
const float cotTheta = deltaZ * iDeltaR;

// discard bottom-middle dublets in a certain (r, eta) region according
// to detector specific cuts
if constexpr (isBottomCandidate) {
if (!m_config.experimentCuts(otherSP->radius(), cotTheta)) {
continue;
}
}

const float Er =
((varianceZM + otherSP->varianceZ()) +
(cotTheta * cotTheta) * (varianceRM + otherSP->varianceR())) *
Expand Down Expand Up @@ -420,6 +432,14 @@ SeedFinder<external_spacepoint_t, platform_t>::getCompatibleDoublets(
const float iDeltaR = std::sqrt(iDeltaR2);
const float cotTheta = deltaZ * iDeltaR;

// discard bottom-middle dublets in a certain (r, eta) region according
// to detector specific cuts
if constexpr (isBottomCandidate) {
if (!m_config.experimentCuts(otherSP->radius(), cotTheta)) {
continue;
}
}

const float Er =
((varianceZM + otherSP->varianceZ()) +
(cotTheta * cotTheta) * (varianceRM + otherSP->varianceR())) *
Expand Down
4 changes: 4 additions & 0 deletions Core/include/Acts/Seeding/SeedFinderConfig.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -173,6 +173,10 @@ struct SeedFinderConfig {
// Returns position of the center of the top strip.
Delegate<Acts::Vector3(const SpacePoint&)> getTopStripCenterPosition;

// Delegate to apply experiment specific cuts
Delegate<bool(float /*bottomRadius*/, float /*cotTheta*/)> experimentCuts{
DelegateFuncTag<&noopExperimentCuts>{}};

bool isInInternalUnits = false;

SeedFinderConfig toInternalUnits() const {
Expand Down
9 changes: 9 additions & 0 deletions Core/include/Acts/Seeding/SeedFinderOrthogonal.ipp
Original file line number Diff line number Diff line change
Expand Up @@ -226,6 +226,15 @@ bool SeedFinderOrthogonal<external_spacepoint_t>::validTuple(
return false;
}

/*
* Cut: Ensure that inner-middle dublet is in a certain (r, eta) region of the
* detector according to detector specific cuts.
*/
const float rInner = (isMiddleInverted) ? rH : rL;
if (!m_config.experimentCuts(rInner, cotTheta)) {
return false;
}

return true;
}

Expand Down
5 changes: 5 additions & 0 deletions Core/include/Acts/Seeding/SeedFinderOrthogonalConfig.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -11,6 +11,7 @@
#include "Acts/Definitions/Algebra.hpp"
#include "Acts/Definitions/Units.hpp"
#include "Acts/Seeding/SeedConfirmationRangeConfig.hpp"
#include "Acts/Utilities/Delegate.hpp"

#include <memory>

Expand Down Expand Up @@ -109,6 +110,10 @@ struct SeedFinderOrthogonalConfig {
float highland = 0;
float maxScatteringAngle2 = 0;

// Delegate to apply experiment specific cuts
Delegate<bool(float /*bottomRadius*/, float /*cotTheta*/)> experimentCuts{
DelegateFuncTag<&noopExperimentCuts>{}};

bool isInInternalUnits = false;

SeedFinderOrthogonalConfig calculateDerivedQuantities() const {
Expand Down

0 comments on commit f32f90f

Please sign in to comment.