Skip to content
Draft
Show file tree
Hide file tree
Changes from all 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
9 changes: 6 additions & 3 deletions RecoTracker/LSTCore/interface/alpaka/Common.h
Original file line number Diff line number Diff line change
Expand Up @@ -87,9 +87,12 @@ namespace ALPAKA_ACCELERATOR_NAMESPACE::lst {
HOST_DEVICE_CONSTANT float kEta_norm = 2.5f;
HOST_DEVICE_CONSTANT float kZ_max = 267.2349854f;
HOST_DEVICE_CONSTANT float kR_max = 110.1099396f;
HOST_DEVICE_CONSTANT float kWp[kPtBins][kEtaBins] = {
{0.4493f, 0.4939f, 0.5715f, 0.6488f, 0.5709f, 0.5938f, 0.7164f, 0.7565f, 0.8103f, 0.8593f},
{0.4488f, 0.4448f, 0.5067f, 0.5929f, 0.4836f, 0.4112f, 0.4968f, 0.4403f, 0.5597f, 0.5067f}};
HOST_DEVICE_CONSTANT float kWp_prompt[kPtBins][kEtaBins] = {
{0.398084f, 0.403121f, 0.354281f, 0.358435f, 0.367765f, 0.349150f, 0.378903f, 0.373203f, 0.406608f, 0.374461f},
{0.097341f, 0.106300f, 0.077975f, 0.103842f, 0.124790f, 0.052589f, 0.054577f, 0.112889f, 0.179196f, 0.200328f}};
HOST_DEVICE_CONSTANT float kWp_displaced[kPtBins][kEtaBins] = {
{0.421372f, 0.466420f, 0.611434f, 0.646523f, 0.543573f, 0.552745f, 0.577597f, 0.630934f, 0.620449f, 0.665398f},
{0.674274f, 0.768983f, 0.888747f, 0.847046f, 0.566714f, 0.946426f, 0.912854f, 0.522756f, 0.538444f, 0.473106f}};
} // namespace t5dnn

namespace pt3dnn {
Expand Down
102 changes: 56 additions & 46 deletions RecoTracker/LSTCore/src/alpaka/NeuralNetwork.h
Original file line number Diff line number Diff line change
Expand Up @@ -197,10 +197,17 @@ namespace ALPAKA_ACCELERATOR_NAMESPACE::lst {
const unsigned int mdIndex5,
const float innerRadius,
const float outerRadius,
const float bridgeRadius) {
const float bridgeRadius,
const float fakeScore1,
const float promptScore1,
const float dispScore1,
const float fakeScore2,
const float promptScore2,
const float dispScore2) {
// Constants
constexpr unsigned int kInputFeatures = 23;
constexpr unsigned int kInputFeatures = 30;
constexpr unsigned int kHiddenFeatures = 32;
constexpr unsigned int kOutputFeatures = 3;

float eta1 = alpaka::math::abs(acc, mds.anchorEta()[mdIndex1]); // inner T3 anchor hit 1 eta
float eta2 = alpaka::math::abs(acc, mds.anchorEta()[mdIndex2]); // inner T3 anchor hit 2 eta
Expand All @@ -220,51 +227,52 @@ namespace ALPAKA_ACCELERATOR_NAMESPACE::lst {
float z4 = alpaka::math::abs(acc, mds.anchorZ()[mdIndex4]); // outer T3 anchor hit 4 z
float z5 = alpaka::math::abs(acc, mds.anchorZ()[mdIndex5]); // outer T3 anchor hit 5 z

float r1 = mds.anchorRt()[mdIndex1]; // inner T3 anchor hit 1 r
float r2 = mds.anchorRt()[mdIndex2]; // inner T3 anchor hit 2 r
float r3 = mds.anchorRt()[mdIndex3]; // inner T3 anchor hit 3 r
float r4 = mds.anchorRt()[mdIndex4]; // outer T3 anchor hit 4 r
float r5 = mds.anchorRt()[mdIndex5]; // outer T3 anchor hit 5 r
float r1 = mds.anchorRt()[mdIndex1];
float r2 = mds.anchorRt()[mdIndex2];
float r3 = mds.anchorRt()[mdIndex3];
float r4 = mds.anchorRt()[mdIndex4];
float r5 = mds.anchorRt()[mdIndex5];

// Build the input feature vector using pairwise differences after the first hit
float x[kInputFeatures] = {
eta1 / dnn::t5dnn::kEta_norm, // inner T3: First hit eta normalized
alpaka::math::abs(acc, phi1) / dnn::kPhi_norm, // inner T3: First hit phi normalized
z1 / dnn::t5dnn::kZ_max, // inner T3: First hit z normalized
r1 / dnn::t5dnn::kR_max, // inner T3: First hit r normalized

eta2 - eta1, // inner T3: Difference in eta between hit 2 and 1
cms::alpakatools::deltaPhi(acc, phi2, phi1) /
dnn::kPhi_norm, // inner T3: Difference in phi between hit 2 and 1
(z2 - z1) / dnn::t5dnn::kZ_max, // inner T3: Difference in z between hit 2 and 1 normalized
(r2 - r1) / dnn::t5dnn::kR_max, // inner T3: Difference in r between hit 2 and 1 normalized

eta3 - eta2, // inner T3: Difference in eta between hit 3 and 2
cms::alpakatools::deltaPhi(acc, phi3, phi2) /
dnn::kPhi_norm, // inner T3: Difference in phi between hit 3 and 2
(z3 - z2) / dnn::t5dnn::kZ_max, // inner T3: Difference in z between hit 3 and 2 normalized
(r3 - r2) / dnn::t5dnn::kR_max, // inner T3: Difference in r between hit 3 and 2 normalized

eta4 - eta3, // outer T3: Difference in eta between hit 4 and 3
cms::alpakatools::deltaPhi(acc, phi4, phi3) /
dnn::kPhi_norm, // outer T3: Difference in phi between hit 4 and 3
(z4 - z3) / dnn::t5dnn::kZ_max, // outer T3: Difference in z between hit 4 and 3 normalized
(r4 - r3) / dnn::t5dnn::kR_max, // outer T3: Difference in r between hit 4 and 3 normalized

eta5 - eta4, // outer T3: Difference in eta between hit 5 and 4
cms::alpakatools::deltaPhi(acc, phi5, phi4) /
dnn::kPhi_norm, // outer T3: Difference in phi between hit 5 and 4
(z5 - z4) / dnn::t5dnn::kZ_max, // outer T3: Difference in z between hit 5 and 4 normalized
(r5 - r4) / dnn::t5dnn::kR_max, // outer T3: Difference in r between hit 5 and 4 normalized

alpaka::math::log10(acc, innerRadius), // T5 inner radius
alpaka::math::log10(acc, bridgeRadius), // T5 bridge radius
alpaka::math::log10(acc, outerRadius) // T5 outer radius
};
float x[kInputFeatures] = {eta1 / dnn::t5dnn::kEta_norm,
alpaka::math::cos(acc, phi1),
alpaka::math::sin(acc, phi1),
z1 / dnn::t5dnn::kZ_max,
r1 / dnn::t5dnn::kR_max,

eta2 - eta1,
cms::alpakatools::deltaPhi(acc, phi2, phi1),
(z2 - z1) / dnn::t5dnn::kZ_max,
(r2 - r1) / dnn::t5dnn::kR_max,

eta3 - eta2,
cms::alpakatools::deltaPhi(acc, phi3, phi2),
(z3 - z2) / dnn::t5dnn::kZ_max,
(r3 - r2) / dnn::t5dnn::kR_max,

eta4 - eta3,
cms::alpakatools::deltaPhi(acc, phi4, phi3),
(z4 - z3) / dnn::t5dnn::kZ_max,
(r4 - r3) / dnn::t5dnn::kR_max,

eta5 - eta4,
cms::alpakatools::deltaPhi(acc, phi5, phi4),
(z5 - z4) / dnn::t5dnn::kZ_max,
(r5 - r4) / dnn::t5dnn::kR_max,

1.0f / innerRadius,
1.0f / bridgeRadius,
1.0f / outerRadius,

fakeScore1,
promptScore1,
dispScore1,
(fakeScore2 - fakeScore1),
(promptScore2 - promptScore1),
(dispScore2 - dispScore1)};

float x_1[kHiddenFeatures]; // Layer 1 output
float x_2[kHiddenFeatures]; // Layer 2 output
float x_3[1]; // Layer 3 linear output
float x_3[kOutputFeatures];

// Layer 1: Linear + Relu
linear_layer<kInputFeatures, kHiddenFeatures>(x, x_1, dnn::t5dnn::wgtT_layer1, dnn::t5dnn::bias_layer1);
Expand All @@ -275,8 +283,9 @@ namespace ALPAKA_ACCELERATOR_NAMESPACE::lst {
relu_activation<kHiddenFeatures>(x_2);

// Layer 3: Linear + Sigmoid
linear_layer<kHiddenFeatures, 1>(x_2, x_3, dnn::t5dnn::wgtT_output_layer, dnn::t5dnn::bias_output_layer);
float x_5 = sigmoid_activation(acc, x_3[0]);
linear_layer<kHiddenFeatures, kOutputFeatures>(
x_2, x_3, dnn::t5dnn::wgtT_output_layer, dnn::t5dnn::bias_output_layer);
softmax_activation<dnn::t3dnn::kOutputFeatures>(acc, x_3);

// Get the bin index based on abs(eta) of first hit and t5_pt
float t5_pt = innerRadius * lst::k2Rinv1GeVf * 2;
Expand All @@ -285,7 +294,8 @@ namespace ALPAKA_ACCELERATOR_NAMESPACE::lst {
uint8_t bin_index = (eta1 > 2.5f) ? (dnn::kEtaBins - 1) : static_cast<unsigned int>(eta1 / dnn::kEtaSize);

// Compare x_5 to the cut value for the relevant bin
return x_5 > dnn::t5dnn::kWp[pt_index][bin_index];
return x_3[1] > dnn::t5dnn::kWp_prompt[pt_index][bin_index] ||
x_3[2] > dnn::t5dnn::kWp_displaced[pt_index][bin_index];
}
} // namespace t5dnn

Expand Down
8 changes: 7 additions & 1 deletion RecoTracker/LSTCore/src/alpaka/Quintuplet.h
Original file line number Diff line number Diff line change
Expand Up @@ -1540,7 +1540,13 @@ namespace ALPAKA_ACCELERATOR_NAMESPACE::lst {
fifthMDIndex,
innerRadius,
outerRadius,
bridgeRadius);
bridgeRadius,
triplets.fakeScore()[innerTripletIndex],
triplets.promptScore()[innerTripletIndex],
triplets.displacedScore()[innerTripletIndex],
triplets.fakeScore()[outerTripletIndex],
triplets.promptScore()[outerTripletIndex],
triplets.displacedScore()[outerTripletIndex]);
tightCutFlag = tightCutFlag and inference; // T5-in-TC cut
if (!inference) // T5-building cut
return false;
Expand Down
Loading