ATLAS Offline Software
FPGATrackSimGNNGraphConstructionTool.cxx
Go to the documentation of this file.
1 // Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
2 
4 
5 #include <TFile.h>
6 #include <TTree.h>
8 
10 // AthAlgTool
11 
12 FPGATrackSimGNNGraphConstructionTool::FPGATrackSimGNNGraphConstructionTool(const std::string& algname, const std::string &name, const IInterface *ifc)
13  : AthAlgTool(algname, name, ifc) {}
14 
16 {
17  if(m_graphTool == "ModuleMap") {
18  if(m_moduleMapPath == "") { // Require a path provided for the Module Map
19  ATH_MSG_FATAL("ERROR! No Module Map provided. Please provide a valid path to a ROOT file.");
20  return StatusCode::FAILURE;
21  }
22  if(m_moduleMapType == "doublet") {
23  loadDoubletModuleMap(); // Load the doublet module map and store entry branches in vectors
24  }
25  }
26  else if(m_graphTool == "MetricLearning") {
27  ATH_CHECK( m_MLInferenceTool.retrieve() );
28  m_MLInferenceTool->printModelInfo();
29  assert(m_MLFeatureNamesVec.size() == m_MLFeatureScalesVec.size());
30  }
31 
32  return StatusCode::SUCCESS;
33 }
34 
36 // Functions
37 
38 StatusCode FPGATrackSimGNNGraphConstructionTool::getEdges(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits, std::vector<std::shared_ptr<FPGATrackSimGNNEdge>> & edges)
39 {
40  if(m_graphTool == "ModuleMap") {
41  doModuleMap(hits, edges);
42  }
43  else if(m_graphTool == "MetricLearning") {
44  doMetricLearning(hits, edges);
45  }
46 
47  return StatusCode::SUCCESS;
48 }
49 
51 {
52  std::unique_ptr<TFile> file(TFile::Open(m_moduleMapPath.value().c_str()));
53  std::unique_ptr<TTree> tree(static_cast<TTree*>(file->Get("TreeModuleDoublet")));
54 
55  unsigned int mid1_value = 0;
56  unsigned int mid2_value = 0;
57  float z0min_12_value = 0.0;
58  float dphimin_12_value = 0.0;
59  float phislopemin_12_value = 0.0;
60  float detamin_12_value = 0.0;
61  float z0max_12_value = 0.0;
62  float dphimax_12_value = 0.0;
63  float phislopemax_12_value = 0.0;
64  float detamax_12_value = 0.0;
65 
66  tree->SetBranchAddress("Module1", &mid1_value);
67  tree->SetBranchAddress("Module2", &mid2_value);
68  tree->SetBranchAddress("z0min_12", &z0min_12_value);
69  tree->SetBranchAddress("dphimin_12", &dphimin_12_value);
70  tree->SetBranchAddress("phiSlopemin_12", &phislopemin_12_value);
71  tree->SetBranchAddress("detamin_12", &detamin_12_value);
72  tree->SetBranchAddress("z0max_12", &z0max_12_value);
73  tree->SetBranchAddress("dphimax_12", &dphimax_12_value);
74  tree->SetBranchAddress("phiSlopemax_12", &phislopemax_12_value);
75  tree->SetBranchAddress("detamax_12", &detamax_12_value);
76 
77  int64_t nEntries = tree->GetEntries();
78  for (int64_t i = 0; i < nEntries; ++i) {
79  tree->GetEntry(i);
80  m_mid1.emplace_back(mid1_value);
81  m_mid2.emplace_back(mid2_value);
82  m_z0min_12.emplace_back(z0min_12_value);
83  m_dphimin_12.emplace_back(dphimin_12_value);
84  m_phislopemin_12.emplace_back(phislopemin_12_value);
85  m_detamin_12.emplace_back(detamin_12_value);
86  m_z0max_12.emplace_back(z0max_12_value);
87  m_dphimax_12.emplace_back(dphimax_12_value);
88  m_phislopemax_12.emplace_back(phislopemax_12_value);
89  m_detamax_12.emplace_back(detamax_12_value);
90  }
91 }
92 
93 void FPGATrackSimGNNGraphConstructionTool::doModuleMap(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits, std::vector<std::shared_ptr<FPGATrackSimGNNEdge>> & edges)
94 {
95  // Use Module Map method for edge building
96  // Two types of module maps: Doublet and Triplet
97  // For each type of module map there is three functions: minmax, meanrms, and hybrid
98  // Use the proper configuration set by the input script and passed as Gaudi::Property variables
99  // Currently only Doublet Module Map with minmax cuts exist, but others can be implemented later on as desired
100 
101  if(m_moduleMapType == "doublet") {
102  getDoubletEdges(hits, edges);
103  }
104 }
105 
106 void FPGATrackSimGNNGraphConstructionTool::getDoubletEdges(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits, std::vector<std::shared_ptr<FPGATrackSimGNNEdge>> & edges)
107 {
108  // Take the list of hits and use the doublet module map to generate all the edges between hits that pass the doublet cuts
109 
110  for (size_t i = 0; i < m_mid2.size(); i++) {
111  std::vector<std::shared_ptr<FPGATrackSimGNNHit>> hit1_matches;
112  std::vector<std::shared_ptr<FPGATrackSimGNNHit>> hit2_matches;
113 
114  std::vector<int> hit1_indices;
115  std::vector<int> hit2_indices;
116 
117  for (size_t j = 0; j < hits.size(); j++) {
118  if (hits[j]->getIdentifier() == m_mid1[i]) {
119  hit1_matches.emplace_back(hits[j]);
120  hit1_indices.emplace_back(j);
121  }
122  if (hits[j]->getIdentifier() == m_mid2[i]) {
123  hit2_matches.emplace_back(hits[j]);
124  hit2_indices.emplace_back(j);
125  }
126  }
127 
128  for (size_t h1 = 0; h1 < hit1_matches.size(); h1++) {
129  for (size_t h2 = 0; h2 < hit2_matches.size(); h2++) {
130  applyDoubletCuts(hit1_matches[h1], hit2_matches[h2], edges, hit1_indices[h1], hit2_indices[h2], i);
131  }
132  }
133  }
134 }
135 
136 void FPGATrackSimGNNGraphConstructionTool::applyDoubletCuts(const std::shared_ptr<FPGATrackSimGNNHit> & hit1, const std::shared_ptr<FPGATrackSimGNNHit> & hit2, std::vector<std::shared_ptr<FPGATrackSimGNNEdge>> & edges, int hit1_index, int hit2_index, unsigned int modulemap_id)
137 {
138  // Four types of doublet cuts (dEta, z0, dPhi, phiSlope)
139  // If an edge passes all four, then it is a valid edge and can be stored
140 
141  // delta_eta cuts
142  float deta = hit1->getEta() - hit2->getEta();
143  if(!doMask(deta, m_detamin_12[modulemap_id], m_detamax_12[modulemap_id])) return;
144 
145  // z0 cuts
146  float dz = hit2->getZ() - hit1->getZ();
147  float dr = hit2->getR() - hit1->getR();
148  float z0 = dr==0. ? 0. : hit1->getZ() - (hit1->getR() * dz / dr);
149  if(!doMask(z0, m_z0min_12[modulemap_id], m_z0max_12[modulemap_id])) return;
150 
151  // delta_phi cuts
152  float dphi = P4Helpers::deltaPhi(hit2->getPhi(),hit1->getPhi());
153  if(!doMask(dphi, m_dphimin_12[modulemap_id], m_dphimax_12[modulemap_id])) return;
154 
155  // phislope cuts
156  float phislope = dr==0. ? 0. : dphi / dr;
157  if(!doMask(phislope, m_phislopemin_12[modulemap_id], m_phislopemax_12[modulemap_id])) return;
158 
159  // if pass all doublet cuts, then record the edge information
160  std::shared_ptr<FPGATrackSimGNNEdge> edge = std::make_shared<FPGATrackSimGNNEdge>();
161  edge->setEdgeIndex1(hit1_index);
162  edge->setEdgeIndex2(hit2_index);
163  edges.emplace_back(edge);
164 }
165 
167 {
168  bool mask = false;
169  if(m_moduleMapFunc == "minmax") {
170  mask = doMinMaxMask(val, min, max);
171  }
172 
173  return mask;
174 }
175 
177 {
178  bool mask = false;
179 
180  if((val <= max * (1.0 + featureSign(max) * m_moduleMapTol)) &&
181  (val >= min * (1.0 - featureSign(min) * m_moduleMapTol))) {
182  mask = true;
183  }
184 
185  return mask;
186 }
187 
189 {
190  if(feature < 0.0) { return -1.0; }
191  else if(feature > 0.0) { return 1.0; }
192  else { return 0.0; }
193 }
194 
195 void FPGATrackSimGNNGraphConstructionTool::doMetricLearning(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits, std::vector<std::shared_ptr<FPGATrackSimGNNEdge>> & edges)
196 {
197  // Use Metric Learning for edge construction
198  // Clustering properties can be set in the input scripta as Gaudi::Property variables
199  std::vector<float> gNodeFeatures = getNodeFeatures(hits);
200  std::vector<float> gEmbedded = embed(hits);
201  doClustering(hits, edges, gEmbedded);
202 }
203 
204 std::vector<float> FPGATrackSimGNNGraphConstructionTool::getNodeFeatures(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits)
205 {
206  std::vector<float> gNodeFeatures;
207 
208  for(auto hit : hits) {
209  std::map<std::string, float> features;
210  features["r"] = hit->getR();
211  features["phi"] = hit->getPhi();
212  features["z"] = hit->getZ();
213 
214  for(size_t i = 0; i < m_MLFeatureNamesVec.size(); i++){
215  gNodeFeatures.push_back(
217  }
218  }
219  return gNodeFeatures;
220 }
221 
222 std::vector<float> FPGATrackSimGNNGraphConstructionTool::embed(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits)
223 {
224  // Use the ML network to embed the hits in a 12-dim latent space
225  std::vector<float> gNodeFeatures = getNodeFeatures(hits);
226  std::vector<float> gEmbedded;
227 
228  std::vector<Ort::Value> gInputTensor;
229  StatusCode s = m_MLInferenceTool->addInput(gInputTensor, gNodeFeatures, 0, hits.size());
230  std::vector<Ort::Value> gOutputTensor;
231  s = m_MLInferenceTool->addOutput(gOutputTensor, gEmbedded, 0, hits.size());
232  s = m_MLInferenceTool->inference(gInputTensor, gOutputTensor);
233 
234  return gEmbedded;
235 }
236 
237 void FPGATrackSimGNNGraphConstructionTool::doClustering(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits, std::vector<std::shared_ptr<FPGATrackSimGNNEdge>> & edges,
238  std::vector<float> & gEmbedded)
239 {
240  // Create graph edges based on the hits distance in the latent space
241  // Creates a directed graph
242  int n_dim = 12;
243  int size = hits.size();
244  float r_squared = m_metricLearningR*m_metricLearningR;
245  int index1 = 0;
246  int index2 = 0;
247  int count = 0;
248  float distance = 0.;
249  std::vector<float> start(n_dim);
250 
251  // Loop over all hits
252  for(int k = 0; k < size; ++k){
253  count = 0;
254  // Setup current hit
255  for(int j = 0; j < n_dim; ++j){
256  start[j] = gEmbedded[k*n_dim + j];
257  }
258  // Loop over the hits not yet checked
259  for (int i = k + 1; i < size; ++i){
260  distance = 0.;
261  for(int d = 0; d < n_dim; ++d){
262  distance += (start[d] - gEmbedded[i*n_dim + d]) * (start[d] - gEmbedded[i*n_dim + d]);
263  }
264  // Store edge if the distance between the hits meets is below the limit
265  if(distance < r_squared){
266  std::shared_ptr<FPGATrackSimGNNEdge> edge = std::make_shared<FPGATrackSimGNNEdge>();
267  // Set order of edge indices to make a directed graph
268  float d_i_sq = (hits[i]->getR() * hits[i]->getR()) + (hits[i]->getZ() * hits[i]->getZ());
269  float d_k_sq = (hits[k]->getR() * hits[k]->getR()) + (hits[k]->getZ() * hits[k]->getZ());
270  if (d_i_sq < d_k_sq){
271  index1 = i;
272  index2 = k;
273  } else {
274  index1 = k;
275  index2 = i;
276  }
277 
278  edge->setEdgeIndex1(index1);
279  edge->setEdgeIndex2(index2);
280  edges.emplace_back(edge);
281  ++count;
282  }
283  // Upper limit for connections of the same hit
285  break;
286  }
287  }
288 
289  }
290 }
getMenu.algname
algname
Definition: getMenu.py:54
ATH_MSG_FATAL
#define ATH_MSG_FATAL(x)
Definition: AthMsgStreamMacros.h:34
FPGATrackSimGNNGraphConstructionTool::getNodeFeatures
std::vector< float > getNodeFeatures(const std::vector< std::shared_ptr< FPGATrackSimGNNHit >> &hits)
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:204
TRTCalib_Extractor.hits
hits
Definition: TRTCalib_Extractor.py:35
FPGATrackSimGNNGraphConstructionTool::m_MLFeatureScalesVec
FloatArrayProperty m_MLFeatureScalesVec
Definition: FPGATrackSimGNNGraphConstructionTool.h:100
FPGATrackSimGNNGraphConstructionTool::m_metricLearningMaxN
Gaudi::Property< int > m_metricLearningMaxN
Definition: FPGATrackSimGNNGraphConstructionTool.h:60
hist_file_dump.d
d
Definition: hist_file_dump.py:142
max
constexpr double max()
Definition: ap_fixedTest.cxx:33
FPGATrackSimGNNGraphConstructionTool::embed
std::vector< float > embed(const std::vector< std::shared_ptr< FPGATrackSimGNNHit >> &hits)
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:222
min
constexpr double min()
Definition: ap_fixedTest.cxx:26
FPGATrackSimGNNGraphConstructionTool::m_dphimin_12
std::vector< float > m_dphimin_12
Definition: FPGATrackSimGNNGraphConstructionTool.h:69
mergePhysValFiles.start
start
Definition: DataQuality/DataQualityUtils/scripts/mergePhysValFiles.py:13
tree
TChain * tree
Definition: tile_monitor.h:30
FPGATrackSimGNNGraphConstructionTool::m_moduleMapFunc
Gaudi::Property< std::string > m_moduleMapFunc
Definition: FPGATrackSimGNNGraphConstructionTool.h:56
FPGATrackSimGNNGraphConstructionTool::m_moduleMapTol
Gaudi::Property< float > m_moduleMapTol
Definition: FPGATrackSimGNNGraphConstructionTool.h:57
python.TurnDataReader.dr
dr
Definition: TurnDataReader.py:111
FPGATrackSimGNNEdge::setEdgeIndex1
void setEdgeIndex1(int v)
Definition: FPGATrackSimGNNEdge.h:27
FPGATrackSimGNNGraphConstructionTool::FPGATrackSimGNNGraphConstructionTool
FPGATrackSimGNNGraphConstructionTool(const std::string &, const std::string &, const IInterface *)
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:12
read_hist_ntuple.h1
h1
Definition: read_hist_ntuple.py:21
FPGATrackSimGNNGraphConstructionTool::m_MLFeatureNamesVec
StringArrayProperty m_MLFeatureNamesVec
Definition: FPGATrackSimGNNGraphConstructionTool.h:96
FPGATrackSimGNNEdge::setEdgeIndex2
void setEdgeIndex2(int v)
Definition: FPGATrackSimGNNEdge.h:28
FPGATrackSimGNNGraphConstructionTool.h
Implements graph construction tool to build edges (connections) between hits.
XMLtoHeader.count
count
Definition: XMLtoHeader.py:84
python.utils.AtlRunQueryLookup.mask
string mask
Definition: AtlRunQueryLookup.py:459
FPGATrackSimGNNGraphConstructionTool::applyDoubletCuts
void applyDoubletCuts(const std::shared_ptr< FPGATrackSimGNNHit > &hit1, const std::shared_ptr< FPGATrackSimGNNHit > &hit2, std::vector< std::shared_ptr< FPGATrackSimGNNEdge >> &edges, int hit1_index, int hit2_index, unsigned int modulemap_id)
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:136
FPGATrackSimGNNHit::getR
float getR() const
Definition: FPGATrackSimGNNHit.h:44
FPGATrackSimGNNGraphConstructionTool::m_z0min_12
std::vector< float > m_z0min_12
Definition: FPGATrackSimGNNGraphConstructionTool.h:68
python.setupRTTAlg.size
int size
Definition: setupRTTAlg.py:39
P4Helpers::deltaPhi
double deltaPhi(double phiA, double phiB)
delta Phi in range [-pi,pi[
Definition: P4Helpers.h:34
Trk::index1
@ index1
Definition: BoundarySurfaceFace.h:48
lumiFormat.i
int i
Definition: lumiFormat.py:85
FPGATrackSimGNNGraphConstructionTool::doMetricLearning
void doMetricLearning(const std::vector< std::shared_ptr< FPGATrackSimGNNHit >> &hits, std::vector< std::shared_ptr< FPGATrackSimGNNEdge >> &edges)
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:195
FPGATrackSimGNNGraphConstructionTool::doMinMaxMask
bool doMinMaxMask(float val, float min, float max)
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:176
FPGATrackSimGNNHit::getEta
float getEta() const
Definition: FPGATrackSimGNNHit.h:46
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
FPGATrackSimGNNGraphConstructionTool::initialize
virtual StatusCode initialize() override
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:15
FPGATrackSimGNNGraphConstructionTool::m_phislopemin_12
std::vector< float > m_phislopemin_12
Definition: FPGATrackSimGNNGraphConstructionTool.h:70
FPGATrackSimGNNGraphConstructionTool::m_MLInferenceTool
ToolHandle< AthOnnx::IOnnxRuntimeInferenceTool > m_MLInferenceTool
Definition: FPGATrackSimGNNGraphConstructionTool.h:49
FPGATrackSimGNNGraphConstructionTool::m_detamin_12
std::vector< float > m_detamin_12
Definition: FPGATrackSimGNNGraphConstructionTool.h:71
file
TFile * file
Definition: tile_monitor.h:29
FPGATrackSimGNNGraphConstructionTool::getDoubletEdges
void getDoubletEdges(const std::vector< std::shared_ptr< FPGATrackSimGNNHit >> &hits, std::vector< std::shared_ptr< FPGATrackSimGNNEdge >> &edges)
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:106
FPGATrackSimGNNHit::getZ
float getZ() const
Definition: FPGATrackSimGNNHit.h:43
FPGATrackSimGNNGraphConstructionTool::m_phislopemax_12
std::vector< float > m_phislopemax_12
Definition: FPGATrackSimGNNGraphConstructionTool.h:74
ATH_CHECK
#define ATH_CHECK
Definition: AthCheckMacros.h:40
FPGATrackSimGNNGraphConstructionTool::m_moduleMapType
Gaudi::Property< std::string > m_moduleMapType
Definition: FPGATrackSimGNNGraphConstructionTool.h:55
FPGATrackSimGNNGraphConstructionTool::m_z0max_12
std::vector< float > m_z0max_12
Definition: FPGATrackSimGNNGraphConstructionTool.h:72
TRT::Track::z0
@ z0
Definition: InnerDetector/InDetCalibEvent/TRT_CalibData/TRT_CalibData/TrackInfo.h:63
FPGATrackSimGNNGraphConstructionTool::m_mid2
std::vector< unsigned int > m_mid2
Definition: FPGATrackSimGNNGraphConstructionTool.h:67
P4Helpers.h
Trk::index2
@ index2
Definition: BoundarySurfaceFace.h:49
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:240
FPGATrackSimGNNGraphConstructionTool::doModuleMap
void doModuleMap(const std::vector< std::shared_ptr< FPGATrackSimGNNHit >> &hits, std::vector< std::shared_ptr< FPGATrackSimGNNEdge >> &edges)
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:93
FPGATrackSimGNNGraphConstructionTool::m_detamax_12
std::vector< float > m_detamax_12
Definition: FPGATrackSimGNNGraphConstructionTool.h:75
FPGATrackSimGNNGraphConstructionTool::getEdges
virtual StatusCode getEdges(const std::vector< std::shared_ptr< FPGATrackSimGNNHit >> &hits, std::vector< std::shared_ptr< FPGATrackSimGNNEdge >> &edges)
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:38
FPGATrackSimGNNGraphConstructionTool::m_mid1
std::vector< unsigned int > m_mid1
Definition: FPGATrackSimGNNGraphConstructionTool.h:66
FPGATrackSimGNNGraphConstructionTool::doClustering
void doClustering(const std::vector< std::shared_ptr< FPGATrackSimGNNHit >> &hits, std::vector< std::shared_ptr< FPGATrackSimGNNEdge >> &edges, std::vector< float > &gEmbedded)
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:237
FPGATrackSimGNNGraphConstructionTool::doMask
bool doMask(float val, float min, float max)
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:166
Pythia8_RapidityOrderMPI.val
val
Definition: Pythia8_RapidityOrderMPI.py:14
FPGATrackSimGNNGraphConstructionTool::m_metricLearningR
Gaudi::Property< float > m_metricLearningR
Definition: FPGATrackSimGNNGraphConstructionTool.h:59
python.SystemOfUnits.s
float s
Definition: SystemOfUnits.py:147
FPGATrackSimGNNGraphConstructionTool::m_moduleMapPath
Gaudi::Property< std::string > m_moduleMapPath
Definition: FPGATrackSimGNNGraphConstructionTool.h:58
FPGATrackSimGNNHit::getPhi
float getPhi() const
Definition: FPGATrackSimGNNHit.h:45
FPGATrackSimGNNGraphConstructionTool::featureSign
float featureSign(float feature)
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:188
AthAlgTool
Definition: AthAlgTool.h:26
dqBeamSpot.nEntries
int nEntries
Definition: dqBeamSpot.py:72
Amg::distance
float distance(const Amg::Vector3D &p1, const Amg::Vector3D &p2)
calculates the distance between two point in 3D space
Definition: GeoPrimitivesHelpers.h:54
FPGATrackSimGNNGraphConstructionTool::m_graphTool
Gaudi::Property< std::string > m_graphTool
Definition: FPGATrackSimGNNGraphConstructionTool.h:54
FPGATrackSimGNNGraphConstructionTool::loadDoubletModuleMap
void loadDoubletModuleMap()
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:50
FPGATrackSimGNNGraphConstructionTool::m_dphimax_12
std::vector< float > m_dphimax_12
Definition: FPGATrackSimGNNGraphConstructionTool.h:73
fitman.k
k
Definition: fitman.py:528