ATLAS Offline Software
Loading...
Searching...
No Matches
FPGATrackSimNNPathfinderExtensionTool.h
Go to the documentation of this file.
1// Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
2
3#ifndef FPGATrackPATHFINDEREXTENSION_H
4#define FPGATrackPATHFINDEREXTENSION_H
5
12
13#include "GaudiKernel/ServiceHandle.h"
15
16// add public header directory to algorithms for this?
18
30#include "GaudiKernel/ITHistSvc.h"
31
32#include "GaudiKernel/IChronoStatSvc.h"
33
34#include <vector>
35#include <memory>
36
37 // internal object for book-keeping during the tree branching, basically just a vector of hits with helper functions - NOTHING else
38 struct miniRoad {
39
41 miniRoad(unsigned nLayers) {m_hits.resize(nLayers);}
42
43 // only one objects that we storek a vector of the hits
44 std::vector<std::shared_ptr<const FPGATrackSimHit>> m_hits;
45
46 // getter functions
47 std::shared_ptr<const FPGATrackSimHit> getHit(size_t layer) const {
48 if (layer < m_hits.size()) return m_hits[layer];
49 const FPGATrackSimHit dummyHit;
50 return std::make_shared<const FPGATrackSimHit>(dummyHit);
51 }
52
53 std::vector<std::shared_ptr<const FPGATrackSimHit>>& getHits() {
54 return m_hits;
55 }
56
57 std::vector<std::vector<std::shared_ptr<const FPGATrackSimHit>>> getVecHits() const {
58 std::vector<std::vector<std::shared_ptr<const FPGATrackSimHit>>> vecHits;
59 vecHits.resize(m_hits.size());
60 for (unsigned layer = 0; layer < m_hits.size(); layer++) {
61 std::vector<std::shared_ptr<const FPGATrackSimHit>> thislayerVec;
62 thislayerVec.push_back(m_hits[layer]);
63 vecHits[layer] = std::move(thislayerVec);
64 }
65 return vecHits;
66 }
67
68 unsigned getNHits() const {
69 return m_hits.size();
70 }
71
72 size_t getNLayers() const { return m_hits.size(); }
73 size_t getNHitLayers() const {
74 return getNLayers() - getNWCLayers();
75 }
76 size_t getNWCLayers() const {
77 size_t nwc = 0;
78 for (auto hit : m_hits) {
79 if (!hit->isReal()) nwc++;
80 }
81 return nwc;
82 }
83
85 layer_bitmask_t wcLayers = 0;
86 for (unsigned layer = 0; layer < m_hits.size(); layer++)
87 if (!((*m_hits[layer]).isReal()))
88 wcLayers |= (0x1 << layer);
89 return wcLayers;
90 }
91
93 layer_bitmask_t hitLayers = 0;
94 for (unsigned layer = 0; layer < m_hits.size(); layer++)
95 if ((*m_hits[layer]).isReal())
96 hitLayers |= (0x1 << layer);
97 return hitLayers;
98 }
99
100
101 // setters
102 void setNLayers(unsigned layers) { m_hits.resize(layers); }
103 void setHits(std::vector<std::shared_ptr<const FPGATrackSimHit>> &hits) {
104 m_hits = hits;
105 }
106 void addHits(const std::vector<std::shared_ptr<const FPGATrackSimHit>> &hits) {
107 m_hits.insert(m_hits.end(), hits.begin(), hits.end());
108 }
109
110 void setHit(unsigned layer, const std::shared_ptr<const FPGATrackSimHit> & hit) {
111 if (layer >= m_hits.size()) m_hits.resize(layer + 1);
112 m_hits[layer] = hit;
113 }
114
115 void addHit(const std::shared_ptr<const FPGATrackSimHit> & hit) {
116 m_hits.push_back(hit);
117 }
118
119
121 // Utility
122
123
124 };
125
126
127
128class FPGATrackSimNNPathfinderExtensionTool : public extends <AthAlgTool, IFPGATrackSimTrackExtensionTool>
129{
130 public:
132 using base_class::base_class;
133
134 virtual StatusCode initialize() override;
135
136 virtual StatusCode extendTracks(const std::vector<std::shared_ptr<const FPGATrackSimHit>> & hits,
137 const FPGATrackSimTrackCollection & tracks,
138 std::vector<FPGATrackSimRoad> & roads) override;
139
140 // We don't have a "union" tool that sits in front of the extension tool, so this is needed here.
141 virtual StatusCode setupSlices(FPGATrackSimLogicalEventInputHeader *slicedHitHeader) override {
142 m_slicedHitHeader = slicedHitHeader;
143 return StatusCode::SUCCESS;
144 };
145
146
147 private:
148 ServiceHandle<IFPGATrackSimMappingSvc> m_FPGATrackSimMapping {this, "FPGATrackSimMappingSvc", "FPGATrackSimMappingSvc"};
149 ServiceHandle<ITHistSvc> m_tHistSvc {this, "THistSvc", "THistSvc"};
150 ServiceHandle<IChronoStatSvc> m_chronoSvc{this, "ChronoStatSvc", "ChronoStatSvc"};
151
152 Gaudi::Property<unsigned> m_maxMiss { this, "threshold", 2, "Maximum number of missing hits to reject a road"};
153 Gaudi::Property <std::string> m_region { this, "OutputRegion", "", "region ID"};
154
155 // Options only needed for sector assignment.
156 // The eta pattern option here should probably be dropped, because we're not using it
157 // and supporting it requires having two sets of eta patterns (one for the first stage, one for the second)
158 // and then running the eta pattern filter a second time.
159 Gaudi::Property<std::vector<float>> m_windowR { this, "windowR", {20.0}, "Window Size to search in for r, either pass one value for all layers or use the number of layers"};
160 Gaudi::Property<std::vector<float>> m_windowZ { this, "windowZ", {20.0}, "Window Size to search in for z, either pass one value for all layers or use the number of layers"};
161 Gaudi::Property<std::vector<float>> m_windowPhi { this, "windowPhi", {0.4}, "Window Size to search in for phi, either pass one value for all layers or use the number of layers"};
162 Gaudi::Property<std::vector<int>> m_windowFineID { this, "windowFineID", {0}, "Fine ID indexing for windows"};
163 Gaudi::Property <float> m_lowPtValueForWindowRScaling { this, "lowPtValueWindowR", -1, "Value in MeV below which we scale the r window size"};
164 Gaudi::Property <float> m_lowPtWindowRScaling {this, "lowPtRScaling", 1.0, "Scaling factor for low pt in R"};
165 Gaudi::Property <float> m_lowPtValueForWindowZScaling { this, "lowPtValueWindowZ", -1, "Value in MeV below which we scale the r window size"};
166 Gaudi::Property <float> m_lowPtWindowZScaling {this, "lowPtZScaling", 1.0, "Scaling factor for low pt in Z"};
167 Gaudi::Property <float> m_lowPtValueForWindowPhiScaling { this, "lowPtValueWindowPhi", -1, "Value in MeV below which we scale the phi window size"};
168 Gaudi::Property <float> m_lowPtWindowPhiScaling {this, "lowPtPhiScaling", 1.0, "Scaling factor for low pt in Phi"};
169
170 Gaudi::Property <float> m_missedHitRScaling {this, "missedHitRScaling", -1, "Amount to scale R window if previous hit was missed. Negative means this is disabled"};
171 Gaudi::Property <float> m_missedHitZScaling {this, "missedHitZScaling", -1, "Amount to scale Z window if previous hit was missed. Negative means this is disabled"};
172 Gaudi::Property <float> m_missedHitPhiScaling {this, "missedHitPhiScaling", -1, "Amount to scale Phi window if previous hit was missed. Negative means this is disabled"};
173 Gaudi::Property <int> m_maxBranches { this, "maxBranches", -1, "Max number of branches before we stop, if negative this is disabled"};
174 Gaudi::Property <bool> m_doOutsideIn { this, "doOutsideIn", true, "Setup the tool so it's doing outside in extrap"};
175 Gaudi::Property <int> m_predictionWindowLength { this, "predictionWindowLength", 3, "Length of hits needed for prediction"};
176 Gaudi::Property <bool> m_useCartesian { this, "useCartesian", true, "If true, NNs use Cartestian coordinates. If false,they use cylindrical coordiantes"};
177 Gaudi::Property <unsigned int> m_batchSize { this, "batchSize", 1, "Batch size for NN inference (1 = sequential, >1 for true batching on GPU)"};
178
179
180 std::vector<FPGATrackSimRoad> m_roads;
181 unsigned m_nLayers_1stStage = 0;
182 unsigned m_nLayers_2ndStage = 0;
183
184 static float getXScale() { return 1015.;};
185 static float getYScale() { return 1015.;};
186 static float getZScale() { return 3000.;};
187 static float getRScale() {return 1015.;};
188 static float getPhiScale() {return 3.15;};
189
190 bool m_debugEvent = false;
191
192 // Internal storage for the sliced hits (implemented as a LogicalEventInputHeader,
193 // so we can easily copy to the output ROOT file).
195
198
199 StatusCode fillInputTensorForNN(miniRoad& thisRoad, std::vector<float>& inputTensorValues);
200 StatusCode getPredictedHit(std::vector<float>& inputTensorValues, std::vector<float>& outputTensorValues, long& fineID);
201 StatusCode getPredictedHitBatched(const std::vector<std::vector<float>>& batchInputTensors,
202 std::vector<std::vector<float>>& batchOutputTensors,
203 std::vector<long>& batchFineIDs);
204 StatusCode addHitToRoad(miniRoad& newroad, miniRoad& currentRoad, const std::vector<std::shared_ptr<const FPGATrackSimHit>>& hits);
205 StatusCode getFakeHit(miniRoad& currentRoad, std::vector<float>& predhit, const long& fineID, std::vector<std::shared_ptr<const FPGATrackSimHit>>& hits);
206 StatusCode getLastLayer(miniRoad& currentRoad, unsigned& lastHitLayer, std::shared_ptr<const FPGATrackSimHit>& lastHit);
207
208 StatusCode findHitinNextStripLayer(std::shared_ptr<const FPGATrackSimHit> hit, const std::vector<std::shared_ptr<const FPGATrackSimHit>>& hitList, std::vector<std::shared_ptr<const FPGATrackSimHit>>& hits);
209
210 void printRoad(miniRoad& currentRoad);
211
212
213};
214
215#endif
: FPGATrackSim-specific class to represent an hit in the detector.
Utilize NN score to build track candidates.
Maps physical layers to logical layers.
Maps ITK module indices to FPGATrackSim regions.
Defines a class for roads.
std::vector< FPGATrackSimTrack > FPGATrackSimTrackCollection
uint32_t layer_bitmask_t
Interface declaration for track extension tools (tracks + hits -> roads)
ServiceHandle< IFPGATrackSimMappingSvc > m_FPGATrackSimMapping
StatusCode getPredictedHit(std::vector< float > &inputTensorValues, std::vector< float > &outputTensorValues, long &fineID)
StatusCode getPredictedHitBatched(const std::vector< std::vector< float > > &batchInputTensors, std::vector< std::vector< float > > &batchOutputTensors, std::vector< long > &batchFineIDs)
StatusCode findHitinNextStripLayer(std::shared_ptr< const FPGATrackSimHit > hit, const std::vector< std::shared_ptr< const FPGATrackSimHit > > &hitList, std::vector< std::shared_ptr< const FPGATrackSimHit > > &hits)
StatusCode getFakeHit(miniRoad &currentRoad, std::vector< float > &predhit, const long &fineID, std::vector< std::shared_ptr< const FPGATrackSimHit > > &hits)
Gaudi::Property< std::vector< float > > m_windowPhi
virtual StatusCode extendTracks(const std::vector< std::shared_ptr< const FPGATrackSimHit > > &hits, const FPGATrackSimTrackCollection &tracks, std::vector< FPGATrackSimRoad > &roads) override
StatusCode getLastLayer(miniRoad &currentRoad, unsigned &lastHitLayer, std::shared_ptr< const FPGATrackSimHit > &lastHit)
StatusCode addHitToRoad(miniRoad &newroad, miniRoad &currentRoad, const std::vector< std::shared_ptr< const FPGATrackSimHit > > &hits)
FPGATrackSimLogicalEventInputHeader * m_slicedHitHeader
StatusCode fillInputTensorForNN(miniRoad &thisRoad, std::vector< float > &inputTensorValues)
virtual StatusCode setupSlices(FPGATrackSimLogicalEventInputHeader *slicedHitHeader) override
void setHit(unsigned layer, const std::shared_ptr< const FPGATrackSimHit > &hit)
layer_bitmask_t getHitLayers() const
void addHits(const std::vector< std::shared_ptr< const FPGATrackSimHit > > &hits)
void addHit(const std::shared_ptr< const FPGATrackSimHit > &hit)
std::vector< std::shared_ptr< const FPGATrackSimHit > > & getHits()
std::vector< std::shared_ptr< const FPGATrackSimHit > > m_hits
std::shared_ptr< const FPGATrackSimHit > getHit(size_t layer) const
std::vector< std::vector< std::shared_ptr< const FPGATrackSimHit > > > getVecHits() const
void setNLayers(unsigned layers)
layer_bitmask_t getWCLayers() const
void setHits(std::vector< std::shared_ptr< const FPGATrackSimHit > > &hits)