Loading [MathJax]/extensions/tex2jax.js
ATLAS Offline Software
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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 
27  return StatusCode::SUCCESS;
28 }
29 
31 // Functions
32 
33 StatusCode FPGATrackSimGNNGraphConstructionTool::getEdges(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits, std::vector<std::shared_ptr<FPGATrackSimGNNEdge>> & edges)
34 {
35  if(m_graphTool == "ModuleMap") {
36  doModuleMap(hits, edges);
37  }
38 
39  return StatusCode::SUCCESS;
40 }
41 
43 {
44  std::unique_ptr<TFile> file(TFile::Open(m_moduleMapPath.value().c_str()));
45  std::unique_ptr<TTree> tree(static_cast<TTree*>(file->Get("TreeModuleDoublet")));
46 
47  unsigned int mid1_value = 0;
48  unsigned int mid2_value = 0;
49  float z0min_12_value = 0.0;
50  float dphimin_12_value = 0.0;
51  float phislopemin_12_value = 0.0;
52  float detamin_12_value = 0.0;
53  float z0max_12_value = 0.0;
54  float dphimax_12_value = 0.0;
55  float phislopemax_12_value = 0.0;
56  float detamax_12_value = 0.0;
57 
58  tree->SetBranchAddress("Module1", &mid1_value);
59  tree->SetBranchAddress("Module2", &mid2_value);
60  tree->SetBranchAddress("z0min_12", &z0min_12_value);
61  tree->SetBranchAddress("dphimin_12", &dphimin_12_value);
62  tree->SetBranchAddress("phiSlopemin_12", &phislopemin_12_value);
63  tree->SetBranchAddress("detamin_12", &detamin_12_value);
64  tree->SetBranchAddress("z0max_12", &z0max_12_value);
65  tree->SetBranchAddress("dphimax_12", &dphimax_12_value);
66  tree->SetBranchAddress("phiSlopemax_12", &phislopemax_12_value);
67  tree->SetBranchAddress("detamax_12", &detamax_12_value);
68 
69  int64_t nEntries = tree->GetEntries();
70  for (int64_t i = 0; i < nEntries; ++i) {
71  tree->GetEntry(i);
72  m_mid1.emplace_back(mid1_value);
73  m_mid2.emplace_back(mid2_value);
74  m_z0min_12.emplace_back(z0min_12_value);
75  m_dphimin_12.emplace_back(dphimin_12_value);
76  m_phislopemin_12.emplace_back(phislopemin_12_value);
77  m_detamin_12.emplace_back(detamin_12_value);
78  m_z0max_12.emplace_back(z0max_12_value);
79  m_dphimax_12.emplace_back(dphimax_12_value);
80  m_phislopemax_12.emplace_back(phislopemax_12_value);
81  m_detamax_12.emplace_back(detamax_12_value);
82  }
83 }
84 
85 void FPGATrackSimGNNGraphConstructionTool::doModuleMap(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits, std::vector<std::shared_ptr<FPGATrackSimGNNEdge>> & edges)
86 {
87  // Use Module Map method for edge building
88  // Two types of module maps: Doublet and Triplet
89  // For each type of module map there is three functions: minmax, meanrms, and hybrid
90  // Use the proper configuration set by the input script and passed as Gaudi::Property variables
91  // Currently only Doublet Module Map with minmax cuts exist, but others can be implemented later on as desired
92 
93  if(m_moduleMapType == "doublet") {
94  getDoubletEdges(hits, edges);
95  }
96 }
97 
98 void FPGATrackSimGNNGraphConstructionTool::getDoubletEdges(const std::vector<std::shared_ptr<FPGATrackSimGNNHit>> & hits, std::vector<std::shared_ptr<FPGATrackSimGNNEdge>> & edges)
99 {
100  // Take the list of hits and use the doublet module map to generate all the edges between hits that pass the doublet cuts
101 
102  for (size_t i = 0; i < m_mid2.size(); i++) {
103  std::vector<std::shared_ptr<FPGATrackSimGNNHit>> hit1_matches;
104  std::vector<std::shared_ptr<FPGATrackSimGNNHit>> hit2_matches;
105 
106  std::vector<int> hit1_indices;
107  std::vector<int> hit2_indices;
108 
109  for (size_t j = 0; j < hits.size(); j++) {
110  if (hits[j]->getIdentifier() == m_mid1[i]) {
111  hit1_matches.emplace_back(hits[j]);
112  hit1_indices.emplace_back(j);
113  }
114  if (hits[j]->getIdentifier() == m_mid2[i]) {
115  hit2_matches.emplace_back(hits[j]);
116  hit2_indices.emplace_back(j);
117  }
118  }
119 
120  for (size_t h1 = 0; h1 < hit1_matches.size(); h1++) {
121  for (size_t h2 = 0; h2 < hit2_matches.size(); h2++) {
122  applyDoubletCuts(hit1_matches[h1], hit2_matches[h2], edges, hit1_indices[h1], hit2_indices[h2], i);
123  }
124  }
125  }
126 }
127 
128 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)
129 {
130  // Four types of doublet cuts (dEta, z0, dPhi, phiSlope)
131  // If an edge passes all four, then it is a valid edge and can be stored
132 
133  // delta_eta cuts
134  float deta = hit1->getEta() - hit2->getEta();
135  if(!doMask(deta, m_detamin_12[modulemap_id], m_detamax_12[modulemap_id])) return;
136 
137  // z0 cuts
138  float dz = hit2->getZ() - hit1->getZ();
139  float dr = hit2->getR() - hit1->getR();
140  float z0 = dr==0. ? 0. : hit1->getZ() - (hit1->getR() * dz / dr);
141  if(!doMask(z0, m_z0min_12[modulemap_id], m_z0max_12[modulemap_id])) return;
142 
143  // delta_phi cuts
144  float dphi = P4Helpers::deltaPhi(hit2->getPhi(),hit1->getPhi());
145  if(!doMask(dphi, m_dphimin_12[modulemap_id], m_dphimax_12[modulemap_id])) return;
146 
147  // phislope cuts
148  float phislope = dr==0. ? 0. : dphi / dr;
149  if(!doMask(phislope, m_phislopemin_12[modulemap_id], m_phislopemax_12[modulemap_id])) return;
150 
151  // if pass all doublet cuts, then record the edge information
152  std::shared_ptr<FPGATrackSimGNNEdge> edge = std::make_shared<FPGATrackSimGNNEdge>();
153  edge->setEdgeIndex1(hit1_index);
154  edge->setEdgeIndex2(hit2_index);
155  edge->setEdgeDR(dr);
156  edge->setEdgeDPhi(dphi);
157  edge->setEdgeDZ(dz);
158  edge->setEdgeDEta(deta);
159  edge->setEdgePhiSlope(phislope);
160  edge->setEdgeRPhiSlope(0.5 * (hit2->getR() + hit1->getR()) * phislope);
161  edges.emplace_back(edge);
162 }
163 
165 {
166  bool mask = false;
167  if(m_moduleMapFunc == "minmax") {
168  mask = doMinMaxMask(val, min, max);
169  }
170 
171  return mask;
172 }
173 
175 {
176  bool mask = false;
177 
178  if((val <= max * (1.0 + featureSign(max) * m_moduleMapTol)) &&
179  (val >= min * (1.0 - featureSign(min) * m_moduleMapTol))) {
180  mask = true;
181  }
182 
183  return mask;
184 }
185 
187 {
188  if(feature < 0.0) { return -1.0; }
189  else if(feature > 0.0) { return 1.0; }
190  else { return 0.0; }
191 }
getMenu.algname
algname
Definition: getMenu.py:54
ATH_MSG_FATAL
#define ATH_MSG_FATAL(x)
Definition: AthMsgStreamMacros.h:34
TRTCalib_Extractor.hits
hits
Definition: TRTCalib_Extractor.py:35
max
constexpr double max()
Definition: ap_fixedTest.cxx:33
min
constexpr double min()
Definition: ap_fixedTest.cxx:26
FPGATrackSimGNNGraphConstructionTool::m_dphimin_12
std::vector< float > m_dphimin_12
Definition: FPGATrackSimGNNGraphConstructionTool.h:58
tree
TChain * tree
Definition: tile_monitor.h:30
FPGATrackSimGNNGraphConstructionTool::m_moduleMapFunc
Gaudi::Property< std::string > m_moduleMapFunc
Definition: FPGATrackSimGNNGraphConstructionTool.h:47
FPGATrackSimGNNGraphConstructionTool::m_moduleMapTol
Gaudi::Property< float > m_moduleMapTol
Definition: FPGATrackSimGNNGraphConstructionTool.h:48
python.TurnDataReader.dr
dr
Definition: TurnDataReader.py:112
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
FPGATrackSimGNNEdge::setEdgePhiSlope
void setEdgePhiSlope(float v)
Definition: FPGATrackSimGNNEdge.h:33
FPGATrackSimGNNEdge::setEdgeIndex2
void setEdgeIndex2(int v)
Definition: FPGATrackSimGNNEdge.h:28
FPGATrackSimGNNGraphConstructionTool.h
Implements graph construction tool to build edges (connections) between hits.
python.utils.AtlRunQueryLookup.mask
string mask
Definition: AtlRunQueryLookup.py:460
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:128
FPGATrackSimGNNHit::getR
float getR() const
Definition: FPGATrackSimGNNHit.h:44
FPGATrackSimGNNGraphConstructionTool::m_z0min_12
std::vector< float > m_z0min_12
Definition: FPGATrackSimGNNGraphConstructionTool.h:57
FPGATrackSimGNNEdge::setEdgeDPhi
void setEdgeDPhi(float v)
Definition: FPGATrackSimGNNEdge.h:30
P4Helpers::deltaPhi
double deltaPhi(double phiA, double phiB)
delta Phi in range [-pi,pi[
Definition: P4Helpers.h:34
lumiFormat.i
int i
Definition: lumiFormat.py:85
FPGATrackSimGNNGraphConstructionTool::doMinMaxMask
bool doMinMaxMask(float val, float min, float max)
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:174
FPGATrackSimGNNEdge::setEdgeDZ
void setEdgeDZ(float v)
Definition: FPGATrackSimGNNEdge.h:31
FPGATrackSimGNNHit::getEta
float getEta() const
Definition: FPGATrackSimGNNHit.h:46
FPGATrackSimGNNEdge::setEdgeDR
void setEdgeDR(float v)
Definition: FPGATrackSimGNNEdge.h:29
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:59
FPGATrackSimGNNGraphConstructionTool::m_detamin_12
std::vector< float > m_detamin_12
Definition: FPGATrackSimGNNGraphConstructionTool.h:60
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:98
FPGATrackSimGNNHit::getZ
float getZ() const
Definition: FPGATrackSimGNNHit.h:43
FPGATrackSimGNNGraphConstructionTool::m_phislopemax_12
std::vector< float > m_phislopemax_12
Definition: FPGATrackSimGNNGraphConstructionTool.h:63
FPGATrackSimGNNGraphConstructionTool::m_moduleMapType
Gaudi::Property< std::string > m_moduleMapType
Definition: FPGATrackSimGNNGraphConstructionTool.h:46
FPGATrackSimGNNGraphConstructionTool::m_z0max_12
std::vector< float > m_z0max_12
Definition: FPGATrackSimGNNGraphConstructionTool.h:61
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:56
P4Helpers.h
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:228
FPGATrackSimGNNGraphConstructionTool::doModuleMap
void doModuleMap(const std::vector< std::shared_ptr< FPGATrackSimGNNHit >> &hits, std::vector< std::shared_ptr< FPGATrackSimGNNEdge >> &edges)
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:85
FPGATrackSimGNNGraphConstructionTool::m_detamax_12
std::vector< float > m_detamax_12
Definition: FPGATrackSimGNNGraphConstructionTool.h:64
FPGATrackSimGNNGraphConstructionTool::getEdges
virtual StatusCode getEdges(const std::vector< std::shared_ptr< FPGATrackSimGNNHit >> &hits, std::vector< std::shared_ptr< FPGATrackSimGNNEdge >> &edges)
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:33
FPGATrackSimGNNEdge::setEdgeRPhiSlope
void setEdgeRPhiSlope(float v)
Definition: FPGATrackSimGNNEdge.h:34
FPGATrackSimGNNGraphConstructionTool::m_mid1
std::vector< unsigned int > m_mid1
Definition: FPGATrackSimGNNGraphConstructionTool.h:55
FPGATrackSimGNNGraphConstructionTool::doMask
bool doMask(float val, float min, float max)
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:164
Pythia8_RapidityOrderMPI.val
val
Definition: Pythia8_RapidityOrderMPI.py:14
FPGATrackSimGNNGraphConstructionTool::m_moduleMapPath
Gaudi::Property< std::string > m_moduleMapPath
Definition: FPGATrackSimGNNGraphConstructionTool.h:49
FPGATrackSimGNNHit::getPhi
float getPhi() const
Definition: FPGATrackSimGNNHit.h:45
FPGATrackSimGNNGraphConstructionTool::featureSign
float featureSign(float feature)
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:186
AthAlgTool
Definition: AthAlgTool.h:26
dqBeamSpot.nEntries
int nEntries
Definition: dqBeamSpot.py:73
FPGATrackSimGNNEdge::setEdgeDEta
void setEdgeDEta(float v)
Definition: FPGATrackSimGNNEdge.h:32
FPGATrackSimGNNGraphConstructionTool::m_graphTool
Gaudi::Property< std::string > m_graphTool
Definition: FPGATrackSimGNNGraphConstructionTool.h:45
FPGATrackSimGNNGraphConstructionTool::loadDoubletModuleMap
void loadDoubletModuleMap()
Definition: FPGATrackSimGNNGraphConstructionTool.cxx:42
FPGATrackSimGNNGraphConstructionTool::m_dphimax_12
std::vector< float > m_dphimax_12
Definition: FPGATrackSimGNNGraphConstructionTool.h:62