Skip to content

Commit

Permalink
refactor: avoid unnecessary calculations in triplet compatibility (#1929
Browse files Browse the repository at this point in the history
)

This PR moves the calculation of `iSinTheta2 * options.sigmapT2perRadius` to outside the loop over the top triplet component to avoid recalculating this variable for every iteration.
It also avoids the calculation of the seed pT if not necessary.
This is not supposed to change the physics performance.
  • Loading branch information
LuisFelipeCoelho committed Mar 15, 2023
1 parent 49e2c9e commit 447f0a2
Show file tree
Hide file tree
Showing 3 changed files with 28 additions and 32 deletions.
2 changes: 0 additions & 2 deletions Core/include/Acts/Seeding/SeedFinder.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -51,8 +51,6 @@ class SeedFinder {
std::vector<InternalSpacePoint<external_spacepoint_t>*> topSpVec;
std::vector<float> curvatures;
std::vector<float> impactParameters;
std::vector<float> etaVec;
std::vector<float> ptVec;

// managing seed candidates for SpM
CandidatesForMiddleSp<InternalSpacePoint<external_spacepoint_t>>
Expand Down
32 changes: 14 additions & 18 deletions Core/include/Acts/Seeding/SeedFinder.ipp
Original file line number Diff line number Diff line change
Expand Up @@ -342,6 +342,7 @@ void SeedFinder<external_spacepoint_t, platform_t>::filterCandidates(

// 1+(cot^2(theta)) = 1/sin^2(theta)
float iSinTheta2 = (1. + cotThetaB * cotThetaB);
float sigmaSquaredSPtDependent = iSinTheta2 * options.sigmapT2perRadius;
// calculate max scattering for min momentum at the seed's theta angle
// scaling scatteringAngle^2 by sin^2(theta) to convert pT^2 to p^2
// accurate would be taking 1/atan(thetaBottom)-1/atan(thetaTop) <
Expand Down Expand Up @@ -553,19 +554,21 @@ void SeedFinder<external_spacepoint_t, platform_t>::filterCandidates(
// the two seed segments using a scattering term scaled by the actual
// measured pT (p2scatterSigma)
float iHelixDiameter2 = B2 / S2;
// calculate scattering for p(T) calculated from seed curvature
float pT2scatterSigma = iHelixDiameter2 * options.sigmapT2perRadius;
// if pT > maxPtScattering, calculate allowed scattering angle using
// maxPtScattering instead of pt.
float pT = options.pTPerHelixRadius * std::sqrt(S2 / B2) / 2.;
if (pT > m_config.maxPtScattering) {
float pTscatterSigma = (m_config.highland / m_config.maxPtScattering) *
m_config.sigmaScattering;
pT2scatterSigma = pTscatterSigma * pTscatterSigma;
}
// convert p(T) to p scaling by sin^2(theta) AND scale by 1/sin^4(theta)
// from rad to deltaCotTheta
float p2scatterSigma = pT2scatterSigma * iSinTheta2;
float p2scatterSigma = iHelixDiameter2 * sigmaSquaredSPtDependent;
if (!std::isinf(m_config.maxPtScattering)) {
// if pT > maxPtScattering, calculate allowed scattering angle using
// maxPtScattering instead of pt.
float pT = options.pTPerHelixRadius * std::sqrt(S2 / B2) / 2.;
if (pT > m_config.maxPtScattering) {
float pTscatterSigma =
(m_config.highland / m_config.maxPtScattering) *
m_config.sigmaScattering;
p2scatterSigma = pTscatterSigma * pTscatterSigma * iSinTheta2;
}
}

// if deltaTheta larger than allowed scattering for calculated pT, skip
if (deltaCotTheta2 > (error2 + p2scatterSigma)) {
if (not m_config.skipPreviousTopSP) {
Expand Down Expand Up @@ -593,13 +596,6 @@ void SeedFinder<external_spacepoint_t, platform_t>::filterCandidates(
// positive/negative in phi
state.curvatures.push_back(B / std::sqrt(S2));
state.impactParameters.push_back(Im);

// evaluate eta and pT of the seed
float cotThetaAvg = std::sqrt(cotThetaAvg2);
float theta = std::atan(1. / cotThetaAvg);
float eta = -std::log(std::tan(0.5 * theta));
state.etaVec.push_back(eta);
state.ptVec.push_back(pT);
} // loop on tops

// continue if number of top SPs is smaller than minimum required for filter
Expand Down
26 changes: 14 additions & 12 deletions Core/include/Acts/Seeding/SeedFinderOrthogonal.ipp
Original file line number Diff line number Diff line change
Expand Up @@ -321,6 +321,7 @@ void SeedFinderOrthogonal<external_spacepoint_t>::filterCandidates(

// 1+(cot^2(theta)) = 1/sin^2(theta)
float iSinTheta2 = (1. + cotThetaB * cotThetaB);
float sigmaSquaredSPtDependent = iSinTheta2 * options.sigmapT2perRadius;
// calculate max scattering for min momentum at the seed's theta angle
// scaling scatteringAngle^2 by sin^2(theta) to convert pT^2 to p^2
// accurate would be taking 1/atan(thetaBottom)-1/atan(thetaTop) <
Expand Down Expand Up @@ -408,21 +409,22 @@ void SeedFinderOrthogonal<external_spacepoint_t>::filterCandidates(

// 1/helixradius: (B/sqrt(S2))*2 (we leave everything squared)
float iHelixDiameter2 = B2 / S2;
// calculate scattering for p(T) calculated from seed curvature
float pT2scatter = 4 * iHelixDiameter2 * options.pT2perRadius;
// if pT > maxPtScattering, calculate allowed scattering angle using
// maxPtScattering instead of pt.
float pT = options.pTPerHelixRadius * std::sqrt(S2 / B2) / 2.;
if (pT > m_config.maxPtScattering) {
float pTscatter = m_config.highland / m_config.maxPtScattering;
pT2scatter = pTscatter * pTscatter;
}
// convert p(T) to p scaling by sin^2(theta) AND scale by 1/sin^4(theta)
// from rad to deltaCotTheta
float p2scatterSigma = iHelixDiameter2 * sigmaSquaredSPtDependent;
if (!std::isinf(m_config.maxPtScattering)) {
// if pT > maxPtScattering, calculate allowed scattering angle using
// maxPtScattering instead of pt.
float pT = options.pTPerHelixRadius * std::sqrt(S2 / B2) / 2.;
if (pT > m_config.maxPtScattering) {
float pTscatterSigma =
(m_config.highland / m_config.maxPtScattering) *
m_config.sigmaScattering;
p2scatterSigma = pTscatterSigma * pTscatterSigma * iSinTheta2;
}
}
// if deltaTheta larger than allowed scattering for calculated pT, skip
if (deltaCotTheta2 >
(error2 + (pT2scatter * iSinTheta2 * m_config.sigmaScattering *
m_config.sigmaScattering))) {
if (deltaCotTheta2 > (error2 + p2scatterSigma)) {
if (m_config.skipPreviousTopSP) {
if (cotThetaB - cotThetaT < 0) {
break;
Expand Down

0 comments on commit 447f0a2

Please sign in to comment.