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