Skip to content

Commit bcba1a3

Browse files
authored
Merge pull request #49164 from kk428/threshold
LST: Introducing a condition before checking cuts in the T5 counting kernel
2 parents b152705 + 788d8b7 commit bcba1a3

File tree

3 files changed

+155
-3
lines changed

3 files changed

+155
-3
lines changed

RecoTracker/LSTCore/interface/alpaka/Common.h

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -42,6 +42,7 @@ namespace ALPAKA_ACCELERATOR_NAMESPACE::lst {
4242
HOST_DEVICE_CONSTANT float kWidth2S = 0.009;
4343
HOST_DEVICE_CONSTANT float kWidthPS = 0.01;
4444
HOST_DEVICE_CONSTANT float kPt_betaMax = 7.0;
45+
HOST_DEVICE_CONSTANT int kNTripletThreshold = 1000;
4546
// To be updated with std::numeric_limits<float>::infinity() in the code and data files
4647
HOST_DEVICE_CONSTANT float kVerticalModuleSlope = 123456789.0;
4748

RecoTracker/LSTCore/src/alpaka/LSTEvent.dev.cc

Lines changed: 3 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -746,10 +746,12 @@ void LSTEvent::createQuintuplets() {
746746
countConn_workDiv,
747747
CountTripletConnections{},
748748
modules_.const_view<ModulesSoA>(),
749+
miniDoubletsDC_->const_view<MiniDoubletsSoA>(),
749750
segmentsDC_->const_view<SegmentsSoA>(),
750751
tripletsDC_->view<TripletsSoA>(),
751752
tripletsDC_->const_view<TripletsOccupancySoA>(),
752-
rangesDC_->const_view());
753+
rangesDC_->const_view(),
754+
ptCut_);
753755

754756
auto const createEligibleModulesListForQuintuplets_workDiv = cms::alpakatools::make_workdiv<Acc1D>(1, 1024);
755757

RecoTracker/LSTCore/src/alpaka/Quintuplet.h

Lines changed: 151 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -1747,6 +1747,110 @@ namespace ALPAKA_ACCELERATOR_NAMESPACE::lst {
17471747
if (miniDoublet3Index != outerInnerInnerMiniDoubletIndex)
17481748
continue;
17491749

1750+
// If densely connected, do not attempt parallel processing to avoid truncation
1751+
if (nInnerTriplets >= kNTripletThreshold || nOuterTriplets >= kNTripletThreshold) {
1752+
uint16_t lowerModule2 = triplets.lowerModuleIndices()[innerTripletIndex][1];
1753+
uint16_t lowerModule4 = triplets.lowerModuleIndices()[outerTripletIndex][1];
1754+
uint16_t lowerModule5 = triplets.lowerModuleIndices()[outerTripletIndex][2];
1755+
1756+
float innerRadius, outerRadius, bridgeRadius, regressionCenterX, regressionCenterY, regressionRadius,
1757+
rzChiSquared, chiSquared, nonAnchorChiSquared, dBeta1, dBeta2; //required for making distributions
1758+
1759+
float t5Embed[Params_T5::kEmbed] = {0.f};
1760+
1761+
bool tightCutFlag = false;
1762+
1763+
bool success = runQuintupletDefaultAlgo(acc,
1764+
modules,
1765+
mds,
1766+
segments,
1767+
triplets,
1768+
lowerModule1,
1769+
lowerModule2,
1770+
lowerModule3,
1771+
lowerModule4,
1772+
lowerModule5,
1773+
innerTripletIndex,
1774+
outerTripletIndex,
1775+
innerRadius,
1776+
outerRadius,
1777+
bridgeRadius,
1778+
regressionCenterX,
1779+
regressionCenterY,
1780+
regressionRadius,
1781+
rzChiSquared,
1782+
chiSquared,
1783+
nonAnchorChiSquared,
1784+
dBeta1,
1785+
dBeta2,
1786+
tightCutFlag,
1787+
t5Embed,
1788+
ptCut);
1789+
if (success) {
1790+
int totOccupancyQuintuplets =
1791+
alpaka::atomicAdd(acc,
1792+
&quintupletsOccupancy.totOccupancyQuintuplets()[lowerModule1],
1793+
1u,
1794+
alpaka::hierarchy::Threads{});
1795+
if (totOccupancyQuintuplets >= ranges.quintupletModuleOccupancy()[lowerModule1]) {
1796+
#ifdef WARNINGS
1797+
printf("Quintuplet excess alert! Module index = %d, Occupancy = %d\n",
1798+
lowerModule1,
1799+
totOccupancyQuintuplets);
1800+
#endif
1801+
} else {
1802+
int quintupletModuleIndex = alpaka::atomicAdd(
1803+
acc, &quintupletsOccupancy.nQuintuplets()[lowerModule1], 1u, alpaka::hierarchy::Threads{});
1804+
if (ranges.quintupletModuleIndices()[lowerModule1] == -1) {
1805+
#ifdef WARNINGS
1806+
printf("Quintuplets : no memory for module at module index = %d\n", lowerModule1);
1807+
#endif
1808+
} else {
1809+
unsigned int quintupletIndex =
1810+
ranges.quintupletModuleIndices()[lowerModule1] + quintupletModuleIndex;
1811+
float phi = mds.anchorPhi()[segments.mdIndices()[triplets.segmentIndices()[innerTripletIndex][0]]
1812+
[layer2_adjustment]];
1813+
float eta = mds.anchorEta()[segments.mdIndices()[triplets.segmentIndices()[innerTripletIndex][0]]
1814+
[layer2_adjustment]];
1815+
float pt = (innerRadius + outerRadius) * k2Rinv1GeVf;
1816+
float scores = chiSquared + nonAnchorChiSquared;
1817+
addQuintupletToMemory(triplets,
1818+
quintuplets,
1819+
innerTripletIndex,
1820+
outerTripletIndex,
1821+
lowerModule1,
1822+
lowerModule2,
1823+
lowerModule3,
1824+
lowerModule4,
1825+
lowerModule5,
1826+
innerRadius,
1827+
bridgeRadius,
1828+
outerRadius,
1829+
regressionCenterX,
1830+
regressionCenterY,
1831+
regressionRadius,
1832+
rzChiSquared,
1833+
chiSquared,
1834+
nonAnchorChiSquared,
1835+
dBeta1,
1836+
dBeta2,
1837+
pt,
1838+
eta,
1839+
phi,
1840+
scores,
1841+
layer,
1842+
quintupletIndex,
1843+
t5Embed,
1844+
tightCutFlag);
1845+
1846+
triplets.partOfT5()[quintuplets.tripletIndices()[quintupletIndex][0]] = true;
1847+
triplets.partOfT5()[quintuplets.tripletIndices()[quintupletIndex][1]] = true;
1848+
}
1849+
}
1850+
}
1851+
continue;
1852+
}
1853+
17501854
// Match inner Sg and Outer Sg
17511855
int mIdx = alpaka::atomicAdd(acc, &matchCount, 1, alpaka::hierarchy::Threads{});
17521856
unsigned int quintupletIndex = ranges.quintupletModuleIndices()[lowerModule1] + mIdx;
@@ -1892,10 +1996,12 @@ namespace ALPAKA_ACCELERATOR_NAMESPACE::lst {
18921996
struct CountTripletConnections {
18931997
ALPAKA_FN_ACC void operator()(Acc3D const& acc,
18941998
ModulesConst modules,
1999+
MiniDoubletsConst mds,
18952000
SegmentsConst segments,
18962001
Triplets triplets,
18972002
TripletsOccupancyConst tripletsOcc,
1898-
ObjectRangesConst ranges) const {
2003+
ObjectRangesConst ranges,
2004+
const float ptCut) const {
18992005
// The atomicAdd below with hierarchy::Threads{} requires one block in x, y dimensions.
19002006
ALPAKA_ASSERT_ACC((alpaka::getWorkDiv<alpaka::Grid, alpaka::Blocks>(acc)[1] == 1) &&
19012007
(alpaka::getWorkDiv<alpaka::Grid, alpaka::Blocks>(acc)[2] == 1));
@@ -1929,7 +2035,50 @@ namespace ALPAKA_ACCELERATOR_NAMESPACE::lst {
19292035
const unsigned int thirdMDInner = mdIndices[thirdSegIdx][0];
19302036

19312037
if (secondMDOuter == thirdMDInner) {
1932-
alpaka::atomicAdd(acc, &triplets.connectedMax()[innerTripletIndex], 1u, alpaka::hierarchy::Threads{});
2038+
// Will only perform runQuintupletDefaultAlgorithm() checks if densely connected
2039+
if (nInnerTriplets < kNTripletThreshold && nOuterTriplets < kNTripletThreshold) {
2040+
alpaka::atomicAdd(acc, &triplets.connectedMax()[innerTripletIndex], 1u, alpaka::hierarchy::Threads{});
2041+
} else {
2042+
const uint16_t lowerModule2 = lmIdx[innerTripletIndex][1];
2043+
const uint16_t lowerModule4 = lmIdx[outerTripletIndex][1];
2044+
const uint16_t lowerModule5 = lmIdx[outerTripletIndex][2];
2045+
2046+
float innerRadius, outerRadius, bridgeRadius;
2047+
float regCx, regCy, regR;
2048+
float rzChi2, chi2, nonAnchorChi2, dBeta1, dBeta2;
2049+
float t5Embed[Params_T5::kEmbed] = {0.f};
2050+
bool tightFlag = false;
2051+
2052+
const bool ok = runQuintupletDefaultAlgo(acc,
2053+
modules,
2054+
mds,
2055+
segments,
2056+
triplets,
2057+
lowerModule1,
2058+
lowerModule2,
2059+
lowerModule3,
2060+
lowerModule4,
2061+
lowerModule5,
2062+
innerTripletIndex,
2063+
outerTripletIndex,
2064+
innerRadius,
2065+
outerRadius,
2066+
bridgeRadius,
2067+
regCx,
2068+
regCy,
2069+
regR,
2070+
rzChi2,
2071+
chi2,
2072+
nonAnchorChi2,
2073+
dBeta1,
2074+
dBeta2,
2075+
tightFlag,
2076+
t5Embed,
2077+
ptCut);
2078+
if (ok) {
2079+
alpaka::atomicAdd(acc, &triplets.connectedMax()[innerTripletIndex], 1u, alpaka::hierarchy::Threads{});
2080+
}
2081+
}
19332082
}
19342083
}
19352084
}

0 commit comments

Comments
 (0)