ATLAS Offline Software
Loading...
Searching...
No Matches
FPGATrackSimGNNGraphConstructionTool.h
Go to the documentation of this file.
1// Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
2
3#ifndef FPGATRACKSIMGNNGRAPHCONSTRUCTIONTOOL_H
4#define FPGATRACKSIMGNNGRAPHCONSTRUCTIONTOOL_H
5
17
19
23
26#include <onnxruntime_cxx_api.h>
27
28#include <vector>
29#include <string>
30#include <memory>
31#include <array>
32#include <unordered_map>
33
35{
36 public:
37
39 // AthAlgTool
40
41 FPGATrackSimGNNGraphConstructionTool(const std::string&, const std::string&, const IInterface*);
42
43 virtual StatusCode initialize() override;
44
46 // Functions
47
48 virtual StatusCode getEdges(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits,
49 std::vector<std::shared_ptr<FPGATrackSimGNNEdge>> & edges);
50
51 private:
52
54 // Handles
55
56 ToolHandle<AthOnnx::IOnnxRuntimeInferenceTool> m_MLInferenceTool {this, "MLInferenceTool", "AthOnnx::OnnxRuntimeInferenceTool"};
58
60 // Properties
61
62 Gaudi::Property<std::string> m_graphTool { this, "graphTool", "", "Tool for graph construction" };
63 Gaudi::Property<std::string> m_moduleMapType { this, "moduleMapType", "", "Type for Module Map for graph construction" };
64 Gaudi::Property<std::string> m_moduleMapFunc { this, "moduleMapFunc", "", "Function for Module Map for graph construction" };
65 Gaudi::Property<float> m_moduleMapTol { this, "moduleMapTol", 0.0, "Tolerance value for Module Map cut calculations" };
66 Gaudi::Property<float> m_moduleMapRMSThresholdFactor { this, "moduleMapRMSThresholdFactor", 0.0, "RMS Threshold value for Module Map cut calculations" };
67 Gaudi::Property<float> m_metricLearningR { this, "metricLearningR", 0.0, "Clustering radius for Metric Learning"};
68 Gaudi::Property<int> m_metricLearningMaxN { this, "metricLearningMaxN", 1, "Max number of neighbours for Metric Learning"};
69
71 // Convenience
72
73 // Module Map Information
74 std::string m_moduleMapPath;
75
77 unsigned mid1 = 0;
78 unsigned mid2 = 0;
79 unsigned mid3 = 0; // Only exists for triplet map
80 unsigned occurence = 0; // Only exists for triplet map
81
82 struct FeatureCuts {
83 float min = 0.0f;
84 float max = 0.0f;
85 float sum = 0.0f;
86 float sumSq = 0.0f;
87 float mean = 0.0f;
88 float rms = 0.0f;
89 };
90
97
102
103 std::array<DoubletCuts, 2> doubletCuts = {};
105 };
106
107 struct TripletKey {
108 unsigned mid1;
109 unsigned mid2;
110 unsigned mid3;
111
112 bool operator==(const TripletKey& other) const {
113 return mid1 == other.mid1 && mid2 == other.mid2 && mid3 == other.mid3;
114 }
115 };
116
118 std::size_t operator()(const TripletKey& k) const {
119 return std::hash<unsigned>()(k.mid1) ^
120 (std::hash<unsigned>()(k.mid2) << 1) ^
121 (std::hash<unsigned>()(k.mid3) << 2);
122 }
123 };
124
125 std::vector<ModuleMapConfig> m_cfgs;
126 std::unordered_map<TripletKey, const ModuleMapConfig*, TripletKeyHash> m_tripletMap;
127
129 // Helpers
130
133 void doModuleMap(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits,
134 std::vector<std::shared_ptr<FPGATrackSimGNNEdge>> & edges);
135 void getTripletEdges(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits,
136 std::vector<std::shared_ptr<FPGATrackSimGNNEdge>> & edges);
137 void getDoubletEdges(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits,
138 std::vector<std::shared_ptr<FPGATrackSimGNNEdge>> & edges,
139 int cutIndex);
140 bool applyDoubletCuts(const std::shared_ptr<FPGATrackSimGNNHit> & hit1,
141 const std::shared_ptr<FPGATrackSimGNNHit> & hit2,
142 const ModuleMapConfig::DoubletCuts& cuts);
143 bool applyTripletCuts(const std::shared_ptr<FPGATrackSimGNNHit> & hit1,
144 const std::shared_ptr<FPGATrackSimGNNHit> & hit2,
145 const std::shared_ptr<FPGATrackSimGNNHit> & hit3,
146 const ModuleMapConfig::TripletCuts& cuts);
147 bool doMask(float val, const ModuleMapConfig::FeatureCuts& cuts);
148 bool doMinMaxMask(float val, const ModuleMapConfig::FeatureCuts& cuts);
149 bool doMeanRMSMask(float val, const ModuleMapConfig::FeatureCuts& cuts);
150 float featureSign(float feature);
151 void doMetricLearning(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits, std::vector<std::shared_ptr<FPGATrackSimGNNEdge>> & edges);
152 std::vector<float> getNodeFeatures(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits);
153 std::vector<float> embed(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits);
154 void doClustering(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits, std::vector<std::shared_ptr<FPGATrackSimGNNEdge>> & edges, std::vector<float> & gEmbedded);
155 // Metric Learning Properties
156 StringArrayProperty m_MLFeatureNamesVec{
157 this, "MLFeatureNames",
158 {"r", "phi", "z"},
159 "Feature names for the Metric Learning model"};
160 FloatArrayProperty m_MLFeatureScalesVec{
161 this, "MLFeatureScales",
162 {1000.0, 3.14159265359, 1000.0},
163 "Feature scales for the Metric Learning model"};
164};
165
166#endif // FPGATRACKSIMGNNGRAPHCONSTRUCTIONTOOL_H
FPGATrackSim-specific class to represent an edge as a connection between two hits in the detector use...
FPGATrackSim-specific class to represent an hit in the detector used for GNN pattern recognition.
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
bool doMask(float val, const ModuleMapConfig::FeatureCuts &cuts)
bool doMinMaxMask(float val, const ModuleMapConfig::FeatureCuts &cuts)
ToolHandle< AthOnnx::IOnnxRuntimeInferenceTool > m_MLInferenceTool
bool applyTripletCuts(const std::shared_ptr< FPGATrackSimGNNHit > &hit1, const std::shared_ptr< FPGATrackSimGNNHit > &hit2, const std::shared_ptr< FPGATrackSimGNNHit > &hit3, const ModuleMapConfig::TripletCuts &cuts)
virtual StatusCode getEdges(const std::vector< std::shared_ptr< FPGATrackSimGNNHit > > &hits, std::vector< std::shared_ptr< FPGATrackSimGNNEdge > > &edges)
void getDoubletEdges(const std::vector< std::shared_ptr< FPGATrackSimGNNHit > > &hits, std::vector< std::shared_ptr< FPGATrackSimGNNEdge > > &edges, int cutIndex)
std::vector< float > getNodeFeatures(const std::vector< std::shared_ptr< FPGATrackSimGNNHit > > &hits)
void getTripletEdges(const std::vector< std::shared_ptr< FPGATrackSimGNNHit > > &hits, std::vector< std::shared_ptr< FPGATrackSimGNNEdge > > &edges)
std::vector< float > embed(const std::vector< std::shared_ptr< FPGATrackSimGNNHit > > &hits)
void doClustering(const std::vector< std::shared_ptr< FPGATrackSimGNNHit > > &hits, std::vector< std::shared_ptr< FPGATrackSimGNNEdge > > &edges, std::vector< float > &gEmbedded)
void doMetricLearning(const std::vector< std::shared_ptr< FPGATrackSimGNNHit > > &hits, std::vector< std::shared_ptr< FPGATrackSimGNNEdge > > &edges)
FPGATrackSimGNNGraphConstructionTool(const std::string &, const std::string &, const IInterface *)
void doModuleMap(const std::vector< std::shared_ptr< FPGATrackSimGNNHit > > &hits, std::vector< std::shared_ptr< FPGATrackSimGNNEdge > > &edges)
bool doMeanRMSMask(float val, const ModuleMapConfig::FeatureCuts &cuts)
bool applyDoubletCuts(const std::shared_ptr< FPGATrackSimGNNHit > &hit1, const std::shared_ptr< FPGATrackSimGNNHit > &hit2, const ModuleMapConfig::DoubletCuts &cuts)
ServiceHandle< IFPGATrackSimMappingSvc > m_FPGATrackSimMapping
std::unordered_map< TripletKey, const ModuleMapConfig *, TripletKeyHash > m_tripletMap