ATLAS Offline Software
Loading...
Searching...
No Matches
FPGATrackSimNNPathfinderExtensionTool.cxx
Go to the documentation of this file.
1// Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
2
3
13
14
16
17#include <cmath>
18#include <algorithm>
19#include "CLHEP/Units/SystemOfUnits.h"
20
21using CLHEP::pi;
22
24
25 // Retrieve the mapping service.
27
28 // hard code this for now but we may chance in the future
31
32 if (m_windowR.size() != 1 && m_windowR.size() != m_windowFineID.size()) {
33 ATH_MSG_ERROR("Window r size = " << m_windowR << " is not equal to 1 (for all layers) and not equal to " << m_windowFineID.size());
34 return StatusCode::FAILURE;
35 }
36
37 if (m_windowPhi.size() != 1 && m_windowPhi.size() != m_windowFineID.size()) {
38 ATH_MSG_ERROR("Window phi size = " << m_windowPhi << " is not equal to 1 (for all layers) and not equal to " << m_windowFineID.size());
39 return StatusCode::FAILURE;
40 }
41
42 if (m_windowZ.size() != 1 && m_windowZ.size() != m_windowFineID.size()) {
43 ATH_MSG_ERROR("Window z size = " << m_windowZ << " is not equal to 1 (for all layers) and not equal to " << m_windowFineID.size());
44 return StatusCode::FAILURE;
45 }
46
47 if (m_FPGATrackSimMapping->getExtensionNNVolMapString() != "" && m_FPGATrackSimMapping->getExtensionNNHitMapString() != "") {
48 ATH_MSG_INFO("Initializing extension hit NN with string = " << m_FPGATrackSimMapping->getExtensionNNHitMapString());
49 m_extensionHitNN.initialize(m_FPGATrackSimMapping->getExtensionNNHitMapString());
50 ATH_MSG_INFO("Initializing volume NN with string = " << m_FPGATrackSimMapping->getExtensionNNVolMapString());
51 m_extensionVolNN.initialize(m_FPGATrackSimMapping->getExtensionNNVolMapString());
52 }
53 else {
54 ATH_MSG_ERROR("Path to NN-based track extension ONNX file is empty! If you want to run this pipeline, you need to provide an input file.");
55 return StatusCode::FAILURE;
56 }
57
58 ATH_CHECK(m_tHistSvc.retrieve());
59
60 return StatusCode::SUCCESS;
61}
62
63
64StatusCode FPGATrackSimNNPathfinderExtensionTool::extendTracks(const std::vector<std::shared_ptr<const FPGATrackSimHit>> & hits,
65 const std::vector<std::shared_ptr<const FPGATrackSimTrack>> & tracks,
66 std::vector<std::shared_ptr<const FPGATrackSimRoad>> & roads) {
67
68 // Reset the internal second stage roads storage.
69 roads.clear();
70 m_roads.clear();
71 const FPGATrackSimRegionMap* rmap_2nd = m_FPGATrackSimMapping->SubRegionMap_2nd();
72
73 // Create one "tower" per slice for this event.
74 // Note that there now might be only one "slice", at least for the time being.
76 for (int ireg = 0; ireg < rmap_2nd->getNRegions(); ireg++) {
78 m_slicedHitHeader->addTower(tower);
79 }
80 }
81 if(m_debugEvent) ATH_MSG_DEBUG("Got: "<<tracks.size()<<" tracks to extrapolate");
82 // Now, loop over the tracks.
83 for (std::shared_ptr<const FPGATrackSimTrack> track : tracks) {
84 if(m_debugEvent) ATH_MSG_DEBUG("\033[1;31m-------------------------- extraploating Track ------------------ \033[0m");
85 if (track->passedOR() == 0) {
86 continue;
87 }
88 const std::vector<FPGATrackSimHit> hitsOnTrack = track->getFPGATrackSimHits();
89 miniRoad road;
90 float pt = track->getPt();
91
92 for (const auto &thit : hitsOnTrack) {
93 road.addHit(std::make_shared<const FPGATrackSimHit>(thit)); // add all hits, we check if WC later
94 }
95
96 if(m_debugEvent) {
97 ATH_MSG_DEBUG("-----------------Hits in event");
98 for (const std::shared_ptr<const FPGATrackSimHit>& hit: hits) {
99 ATH_MSG_DEBUG("Hit " << " X: " << hit->getX() << " Y: " << hit->getY() << " Z: " << hit->getZ() << " R: " << hit->getR() << "phi = " << hit->getGPhi() << " hitType: " << hit->getHitType() << " getDetType: " << hit->getDetType());
100 }
101 }
102 std::vector<miniRoad> roadsToExtrapolate;
103 roadsToExtrapolate.push_back(std::move(road));
104
105 std::vector<miniRoad> completedRoads;
106
107 int count = 0;
108 // FIXED: Add maximum iteration limit to prevent infinite loops
109 const int MAX_ROADS = 10000;
110 while(roadsToExtrapolate.size() > 0 && count < MAX_ROADS ) {
111 miniRoad currentRoad = *roadsToExtrapolate.begin();
112
113 // Erase this road from the vector
114 roadsToExtrapolate.erase(roadsToExtrapolate.begin());
115 count ++;
116 if(m_debugEvent) {
117 ATH_MSG_DEBUG("\033[1;31m-------------------------- extraploating road "<< count << "------------------ \033[0m");
118 printRoad(currentRoad);
119 }
120 // Check exit condition
121 if (currentRoad.getNHits() >= (m_nLayers_1stStage+m_nLayers_2ndStage))
122 {
123 completedRoads.push_back(std::move(currentRoad));
124 continue; // this one is done
125 }
126 // Other try to find the next hit in this road
127 std::vector<float> inputTensorValues;
128 if (!fillInputTensorForNN(currentRoad, inputTensorValues)) {
129 ATH_MSG_WARNING("Failed to create input tensor for this road");
130 continue;
131 }
132 std::vector<float> predhit;
133 long fineID;
134 if (!getPredictedHit(inputTensorValues, predhit, fineID)) {
135 ATH_MSG_WARNING("Failed to predict hit for this road");
136 continue;
137 }
138 // Check if exist conditions are there
139 if (m_doOutsideIn) {
140 // Make sure we are not predicting inside the inner most layer (x and y < 25, or r < 25)
141 // If we are, road is done
142 if ((m_useCartesian && (abs(predhit[0]) < 25 && abs(predhit[1]) < 25)) ||
143 (!m_useCartesian && abs(predhit[0]) < 25))
144 {
145 completedRoads.push_back(std::move(currentRoad));
146 continue;
147 }
148 }
149 else {
150 // Make sure we are not predicting outside the outer most layer
151 // if we are, road is done
152 double radius = std::sqrt(predhit[0] * predhit[0] + predhit[1] * predhit[1]);
153 if ((m_useCartesian && (abs(predhit[0]) > 1024 || abs(predhit[1]) > 1024 || radius > 1024 || abs(predhit[2]) > 3000)) ||
154 (!m_useCartesian && (abs(predhit[0]) > 1024 || abs(predhit[2]) > 3000))) {
155 completedRoads.push_back(std::move(currentRoad));
156 continue;
157 }
158 }
159 if(m_debugEvent) {
160 ATH_MSG_DEBUG("Predicted hit at: " << predhit[0] << " " << predhit[1] << " " << predhit[2]);
161 }
162
163 // Now search for the hits
164 bool foundhitForRoad = false;
165 if(fineID == 215){
166 ATH_MSG_DEBUG("Stopping condition reached");
167 completedRoads.push_back(std::move(currentRoad));
168 continue;
169 }
170 // Get the last layer and hit in the road
171 unsigned lastLayerInRoad = 0;
172 std::shared_ptr<const FPGATrackSimHit> lastHit;
173 if(!getLastLayer(currentRoad, lastLayerInRoad, lastHit) or !lastHit) {
174 ATH_MSG_WARNING("Failed to find last layer this road");
175 continue;
176 }
177 unsigned layer = lastLayerInRoad+1; // the layer we're looking to find
178 bool lastHitWasReal = lastHit->isReal();
179 float lastHitR = lastHit->getR();
180 if(layer >= (m_nLayers_1stStage + m_nLayers_2ndStage)) {
181 completedRoads.push_back(std::move(currentRoad));
182 continue;
183 }
184 unsigned int hitsInWindow = 0;
185
186 // List of all the hits, with their distances to the predicted point
187 std::vector<std::vector<std::shared_ptr<const FPGATrackSimHit>>> listofHitsFound;
188
189 for (const std::shared_ptr<const FPGATrackSimHit>& hit: hits) {
190 if (m_doOutsideIn && (hit->getR() > lastHitR)) continue;
191 if (!m_doOutsideIn && (hit->getR() < lastHitR)) continue;
192 if ((hit->getHitType() == HitType::spacepoint) && ((hit->getPhysLayer(true)) %2 == 1)) continue; // ignore outer parts of SP, they get added separately
193 if(m_debugEvent) {
194 ATH_MSG_DEBUG("In the hit loop hit at x: " << hit->getX() << " y " << hit->getY() << " z " << hit->getZ() << " phi " << hit->getGPhi());
195 }
196 // loop over hits in that layer
197 if (getFineID(*hit) == fineID && hit->isReal()) {
198 // a hit is in the right fine ID == layer
199 double hitz = hit->getZ();
200 double hitr = hit->getR();
201 double hitphi = hit->getGPhi();
202 double predr = (m_useCartesian ? sqrt(predhit[0] * predhit[0] + predhit[1] * predhit[1]) : predhit[0]);
203 double predphi = predhit[1]; // only for non cartesian, ie cylinndrical
204 double predz = predhit[2];
205 double windowR = m_windowR[0]; // default for all layers
206 double windowPhi = m_windowPhi[0]; // default for all layers
207 double windowZ = m_windowZ[0]; // default for all layers
208 int fineID_index = 0;
209 // But if available pick up per-window values
210 if (m_windowZ.size() > 1 || m_windowR.size() > 1 || m_windowPhi.size() > 1) {
211
212 auto fineID_it = std::find(m_windowFineID.begin(), m_windowFineID.end(), fineID);
213 if (fineID_it == m_windowFineID.end()){
214 ATH_MSG_DEBUG("No windows for predicted fineID " << fineID << ", using maximum in provided list instead!");
215 fineID_index = -1;
216 }
217 else {
218 fineID_index = fineID_it - m_windowFineID.begin();
219 }
220 }
221 if (m_windowR.size() > 1) {
222 if (fineID_index == -1) {
223 windowR = *std::max_element(m_windowR.begin(), m_windowR.end());
224 }
225 else {
226 windowR = m_windowR[fineID_index];
227 }
228 }
229 if (m_windowZ.size() > 1) {
230 if (fineID_index == -1) {
231 windowZ = *std::max_element(m_windowZ.begin(), m_windowZ.end());
232 }
233 else {
234 windowZ = m_windowZ[fineID_index];
235 }
236 }
237 if (m_windowPhi.size() > 1) {
238 if (fineID_index == -1) {
239 windowPhi = *std::max_element(m_windowPhi.begin(), m_windowPhi.end());
240 }
241 else {
242 windowPhi = m_windowPhi[fineID_index];
243 }
244 }
245
246 // If last hit was not real and we want to, scale the window
247 if (m_missedHitRScaling > 0 && !lastHitWasReal) windowR *= m_missedHitRScaling;
248 if (m_missedHitZScaling > 0 && !lastHitWasReal) windowZ *= m_missedHitZScaling;
249 if (m_missedHitPhiScaling > 0 && !lastHitWasReal) windowPhi *= m_missedHitPhiScaling;
250 // now scale windows for low pt, if desired
251 if (pt < m_lowPtValueForWindowRScaling.value()) windowR *= m_lowPtWindowRScaling.value();
252 if (pt < m_lowPtValueForWindowZScaling.value()) windowZ *= m_lowPtWindowZScaling.value();
253 if (pt < m_lowPtValueForWindowPhiScaling.value()) windowPhi *= m_lowPtWindowPhiScaling.value();
254
255 double dr = abs(hitr - predr);
256 double dz = abs(hitz - predz);
257 double dphi = abs(hitphi - predphi);
258 while (dphi > pi) dphi -= pi;
259 if ((m_useCartesian && dr < windowR && dz < windowZ) ||
260 (!m_useCartesian && dphi < windowPhi && dz < windowZ && dr < windowR))
261 {
262 std::vector<std::shared_ptr<const FPGATrackSimHit>> theseHits {hit};
263 hitsInWindow = hitsInWindow + 1;
264 // If the hit is a space point, skip the next layer, as it will be a duplicated space point and we have already taken care of that in the adding of the hits
265 if(hit->isStrip()) {
266 if (hit->getHitType() == HitType::spacepoint) {
267 // find the hit in the next layer
268 if(!findHitinNextStripLayer(hit, hits, theseHits)) {
269 ATH_MSG_WARNING("For a SP in layer " << layer << " Couldn't find a matching strip SP");
270 }
271 }
272 else {
273 std::shared_ptr<FPGATrackSimHit> guessedSecondHitPtr = std::make_shared<FPGATrackSimHit>();
274 guessedSecondHitPtr->setX(0);
275 guessedSecondHitPtr->setY(0);
276 guessedSecondHitPtr->setZ(0);
277 guessedSecondHitPtr->setPhysLayer(lastHit->getPhysLayer(true)+1);
278 guessedSecondHitPtr->setHitType(HitType::undefined);
279 if(isFineIDInStrip(fineID)) guessedSecondHitPtr->setDetType(SiliconTech::strip);
280 else guessedSecondHitPtr->setDetType(SiliconTech::pixel);
281
282 theseHits.push_back(guessedSecondHitPtr);
283 }
284 }
285 // Store the hits for now
286 listofHitsFound.push_back(std::move(theseHits));
287 }
288 }
289 }
290 // Sort the hit by the distance
291 if (m_useCartesian) {
292 std::sort(listofHitsFound.begin(), listofHitsFound.end(), [&predhit](auto& a, auto& b){
293 double predr = sqrt(predhit[0] * predhit[0] + predhit[1] * predhit[1]);
294 double predz = predhit[2];
295
296 // HitA
297 double hitz = a[0]->getZ();
298 double hitr = a[0]->getR();
299 float distance_a = (hitr - predr)*(hitr - predr) + (hitz - predz)*(hitz - predz);
300
301 // HitB
302 hitz = b[0]->getZ();
303 hitr = b[0]->getR();
304 float distance_b = (hitr - predr)*(hitr - predr) + (hitz - predz)*(hitz - predz);
305
306 return distance_a < distance_b;
307 });
308 }
309 else {
310 std::sort(listofHitsFound.begin(), listofHitsFound.end(), [&predhit](auto& a, auto& b){
311
312 double predr = predhit[0];
313 double predphi = predhit[1];
314 double predz = predhit[2];
315 // HitA
316 double hitr = a[0]->getR();
317 double hitphi = a[0]->getGPhi();
318 double hitz = a[0]->getZ();
319 double dz = abs(hitz - predz);
320 double dr = abs(hitr - predr);
321 double dphi = abs(hitphi - predphi);
322 while (dphi > pi) dphi -= pi;
323
324 // scaled distance because z and phi and r are not in same units
325 float distance_a = dz*dz/(getZScale()*getZScale()) + dphi*dphi/(getPhiScale()*getPhiScale()) + dr*dr/(getRScale()*getRScale());
326
327 // HitB
328 hitr = b[0]->getR();
329 hitphi = b[0]->getGPhi();
330 hitz = b[0]->getZ();
331 dz = abs(hitz - predz);
332 dr = abs(hitr - predr);
333 dphi = abs(hitphi - predphi);
334 while (dphi > pi) dphi -= pi;
335
336 // scaled distance because z and phi and r are not in same units
337 float distance_b = dz*dz/(getZScale()*getZScale()) + dphi*dphi/(getPhiScale()*getPhiScale()) + dr*dr/(getRScale()*getRScale());
338
339 return distance_a < distance_b;
340 });
341 }
342 // Select the top N hits
343 std::vector<std::vector<std::shared_ptr<const FPGATrackSimHit>>> cleanHitsToGrow;
344
345 // If max branches are limited, pick only the top N from the list of hits found at each level
346 if (m_maxBranches.value() >= 0) {
347 int nHitsToChoose = std::min(int(m_maxBranches.value()), int(listofHitsFound.size()));
348 cleanHitsToGrow.reserve(nHitsToChoose);
349 std::copy(listofHitsFound.begin(), listofHitsFound.begin() + nHitsToChoose, std::back_inserter(cleanHitsToGrow));
350 }
351 else {
352 cleanHitsToGrow = std::move(listofHitsFound);
353 }
354
355 for (auto& hitsFound: cleanHitsToGrow) {
356
357 // We got a hit, lets make a road
358 miniRoad newRoad;
359 if(!addHitToRoad(newRoad, currentRoad, std::move(hitsFound))) {
360 ATH_MSG_WARNING("Failed to make a new road");
361 continue;
362 }
363 roadsToExtrapolate.push_back(newRoad);
364 foundhitForRoad = true;
365 if(m_debugEvent) {
366 ATH_MSG_DEBUG("------ road grown with hit from layer "<<layer<<" to");
367 printRoad(newRoad);
368 }
369 }
370 // If the hit wasn't found, push a fake hit
371 if (!foundhitForRoad) {
372 // did not find a hit to extrapolate to, check if we need to delete this road. if not, add a guessed hit if still useful
373 if (currentRoad.getNWCLayers() >= m_maxMiss) {
374 // we don't want this road, so we continue
375 continue;
376 }
377 else {
378 std::vector<std::shared_ptr<const FPGATrackSimHit>> theseHits;
379 // first make the fake hit that we will add
380 if (!getFakeHit(currentRoad, predhit, fineID, theseHits)) {
381 ATH_MSG_WARNING("Failed adding a guessed hit in extrapolation");
382 continue;
383 }
384
385 if((isFineIDInPixel(fineID) && theseHits.size() != 1) || (isFineIDInStrip(fineID) && theseHits.size() != 2)) {
386 continue;
387 }
388 // add the hit to the road
389 miniRoad newroad;
390 if (!addHitToRoad(newroad, currentRoad, std::move(theseHits))) {
391 ATH_MSG_WARNING("Failed making a new road with fake hit");
392 continue;
393 }
394 roadsToExtrapolate.push_back(std::move(newroad));
395 }
396 }
397 }
398 // This track has been extrapolated, copy the completed tracks to the full list with full road objects
399 for (const auto &miniroad : completedRoads) {
400 FPGATrackSimRoad road;
401 road.setWCLayers(miniroad.getWCLayers());
402 road.setHitLayers(miniroad.getHitLayers());
403 road.setRoadID(m_roads.size() - 1);
404 // Set the "Hough x" and "Hough y" using the track parameters.
405 road.setX(track->getPhi());
406 road.setY(track->getQOverPt());
407 road.setXBin(track->getHoughXBin());
408 road.setYBin(track->getHoughYBin());
409 road.setSubRegion(track->getSubRegion());
410
411 // just force the right number of layers now, in case we find fewer than expected (needed downstream)
412 std::vector<std::vector<std::shared_ptr<const FPGATrackSimHit>>> roadhits = miniroad.getVecHits();
413 unsigned nexpected = m_nLayers_1stStage+m_nLayers_2ndStage;
414 if (roadhits.size() > nexpected) { // cut off the last ones
415 roadhits.resize(nexpected);
416 }
417 else if (roadhits.size() < nexpected) { // fill with missing hits
418 for (unsigned layer = roadhits.size(); layer < nexpected; layer++) {
419 std::shared_ptr<FPGATrackSimHit> emptyHitPtr = std::make_shared<FPGATrackSimHit>();
420 emptyHitPtr->setX(0);
421 emptyHitPtr->setY(0);
422 emptyHitPtr->setZ(0);
423 emptyHitPtr->setLayer(layer);
424 emptyHitPtr->setHitType(HitType::wildcard);
425
426 roadhits.emplace_back(1,emptyHitPtr);
427 layer_bitmask_t wclayers = road.getWCLayers();
428 wclayers |= (1 << layer);
429 road.setWCLayers(wclayers);
430 }
431 }
432 road.setHits(std::move(roadhits));
433
434 m_roads.push_back(std::move(road));
435 }
436 }
437 // Copy the roads we found into the output argument and return success.
438 roads.reserve(m_roads.size());
439 for (FPGATrackSimRoad & r : m_roads)
440 {
441 if (r.getNWCLayers() >= m_maxMiss) continue; // extra check on this
442 roads.emplace_back(std::make_shared<const FPGATrackSimRoad>(r));
443 }
444 ATH_MSG_DEBUG("Found " << roads.size() << " new roads in second stage.");
445
446 return StatusCode::SUCCESS;
447}
448
449StatusCode FPGATrackSimNNPathfinderExtensionTool::fillInputTensorForNN(miniRoad& thisroad, std::vector<float>& inputTensorValues)
450{
451 std::vector<std::shared_ptr<const FPGATrackSimHit>> hitsR;
452
453 std::vector<std::shared_ptr<const FPGATrackSimHit>> hits = thisroad.getHits();
454
455
456 for (unsigned ihit = 0; ihit < hits.size(); ihit++) {
457 if (ihit < m_nLayers_1stStage && !(hits[ihit]->isReal())) continue; // skip guessed hits for the 1st stage
458 hitsR.push_back(hits[ihit]);
459 }
460 // Sort in increasing 3D distance.
461 std::sort(hitsR.begin(), hitsR.end(), [](auto& a, auto& b){
462 double dist_a = std::hypot(a->getX(), a->getY(), a->getZ());
463 double dist_b = std::hypot(b->getX(), b->getY(), b->getZ());
464 return dist_a < dist_b;
465 });
466
467 if(m_debugEvent) ATH_MSG_DEBUG("hitsR");
468 for (auto thit : hitsR)
469 {
470 if(m_debugEvent) ATH_MSG_DEBUG(thit->getX()<<" "<<thit->getY()<<" "<<thit->getZ() << " and phi = " << thit->getGPhi());
471 }
472
473 // Remove all the duplicate space points
474 std::vector<std::shared_ptr<const FPGATrackSimHit>> cleanHits;
475 bool skipHit = false;
476 for (auto thit : hitsR)
477 {
478 if(skipHit)
479 {
480 skipHit = false;
481 continue;
482 }
483 if (thit->isPixel())
484 {
485 cleanHits.push_back(std::move(thit));
486 }
487 else if (thit->isStrip() && (thit->getHitType() == HitType::spacepoint))
488 {
489 // This is a proper strips SP, push the first hit back and skip the next one since its a duplicate
490 cleanHits.push_back(std::move(thit));
491 skipHit = true;
492 }
493 else if (thit->isStrip() && (thit->getHitType() == HitType::guessed))
494 {
495 // this is a guessed strip SP, push the first hit back and skip the next one since its a duplicate
496 cleanHits.push_back(std::move(thit));
497 skipHit = true;
498 }
499 else if (thit->isStrip() && (thit->getHitType() == HitType::undefined))
500 {
501 // this is a fake hit, inserted for a unpaired SP, continue
502 continue;
503 }
504 else if (thit->isStrip() && thit->isReal())
505 {
506 // What is left here is a unpaired hit, push it back
507 cleanHits.push_back(std::move(thit));
508 }
509 else
510 {
511 ATH_MSG_WARNING("No clue how to deal with this hit in the NN predicition ");
512 continue;
513 }
514 }
515
516 // Reverse the hits as we changed the ordering before
517 if(!m_doOutsideIn)
518 {
519 std::reverse(cleanHits.begin(), cleanHits.end());
520 }
521
522 // Select the top N hits
523 std::vector<std::shared_ptr<const FPGATrackSimHit>> hitsToEncode;
524 std::copy(cleanHits.begin(), cleanHits.begin() + m_predictionWindowLength, std::back_inserter(hitsToEncode));
525
526 if(m_debugEvent) ATH_MSG_DEBUG("Clean hits");
527 for (auto thit : cleanHits)
528 {
529 if(m_debugEvent) ATH_MSG_DEBUG(thit->getX()<<" "<<thit->getY()<<" "<<thit->getZ()<<" "<<thit->isStrip() << " and gphi = " << thit->getGPhi());
530 }
531
532 // Reverse this vector so we can encode it the format as expected from the NN
533 std::reverse(hitsToEncode.begin(), hitsToEncode.end());
534
535 if(m_debugEvent) ATH_MSG_DEBUG("Input for NN prediction");
536 for (auto thit : hitsToEncode)
537 {
538 if (m_useCartesian) {
539 inputTensorValues.push_back(thit->getX()/ getXScale());
540 inputTensorValues.push_back(thit->getY()/ getYScale());
541 inputTensorValues.push_back(thit->getZ()/ getZScale());
542 }
543 else {
544 inputTensorValues.push_back(thit->getR()/ getRScale());
545 inputTensorValues.push_back(thit->getGPhi()/ getPhiScale());
546 inputTensorValues.push_back(thit->getZ()/ getZScale());
547 }
548
549 if(m_debugEvent) ATH_MSG_DEBUG(thit->getX()<<" "<<thit->getY()<<" "<<thit->getZ() << " and gphi = " << thit->getGPhi());
550 }
551
552 return StatusCode::SUCCESS;
553
554}
555
556
557StatusCode FPGATrackSimNNPathfinderExtensionTool::getPredictedHit(std::vector<float>& inputTensorValues, std::vector<float>& outputTensorValues, long& fineID)
558{
559 std::vector<float> NNVoloutput = m_extensionVolNN.runONNXInference(inputTensorValues);
560 fineID = std::distance(NNVoloutput.begin(),std::max_element(NNVoloutput.begin(), NNVoloutput.end()));
561 // Insert the output for the second stage NN
562 inputTensorValues.insert(inputTensorValues.end(), NNVoloutput.begin(), NNVoloutput.end()); //now we append the first NN to the list of coordinates
563
564 // use the above to predict the next hit position
565 outputTensorValues = m_extensionHitNN.runONNXInference(inputTensorValues);
566
567 // now scale back
568 if (m_useCartesian) {
569 outputTensorValues[0] *= getXScale();
570 outputTensorValues[1] *= getYScale();
571 outputTensorValues[2] *= getZScale();
572 }
573 else {
574 outputTensorValues[0] *= getRScale();
575 outputTensorValues[1] *= getPhiScale();
576 outputTensorValues[2] *= getZScale();
577 }
578
579 return StatusCode::SUCCESS;
580}
581
582StatusCode FPGATrackSimNNPathfinderExtensionTool::addHitToRoad(miniRoad& newroad, miniRoad& currentRoad, const std::vector<std::shared_ptr<const FPGATrackSimHit>>& hits)
583{
584 newroad.setHits(currentRoad.getHits());
585 newroad.addHits(hits);
586 return StatusCode::SUCCESS;
587}
588
589StatusCode FPGATrackSimNNPathfinderExtensionTool::getFakeHit(miniRoad& currentRoad, std::vector<float>& predhit, const long& fineID, std::vector<std::shared_ptr<const FPGATrackSimHit>>& hits) {
590
591 unsigned guessedLayer(0);
592
593 std::shared_ptr<FPGATrackSimHit> guessedHitPtr = std::make_shared<FPGATrackSimHit>();
594 if (m_useCartesian) {
595 guessedHitPtr->setX(predhit[0]);
596 guessedHitPtr->setY(predhit[1]);
597 guessedHitPtr->setZ(predhit[2]);
598 }
599 else {
600 double r = predhit[0];
601 double phi = predhit[1];
602 guessedHitPtr->setX(r*cos(phi));
603 guessedHitPtr->setY(r*sin(phi));
604 guessedHitPtr->setZ(predhit[2]);
605 }
606
607
608 if (m_doOutsideIn) { // outside in
609 // TODO
610 }
611 else { // nope, inside out
612 guessedLayer = currentRoad.getNHits();
613 }
614
615 guessedHitPtr->setLayer(guessedLayer);
616 guessedHitPtr->setHitType(HitType::guessed);
617
618 if(isFineIDInStrip(fineID)) {
619 guessedHitPtr->setDetType(SiliconTech::strip);
620 guessedHitPtr->setPhysLayer(0); // it is just to set the side, it's not actually 0
621 }
622 else guessedHitPtr->setDetType(SiliconTech::pixel);
623
624
625 // Make sure that the hit is inside the boundaries of the layers
626 if(guessedLayer < (m_nLayers_1stStage + m_nLayers_2ndStage ))
627 {
628 hits.push_back(guessedHitPtr);
629 }
630 // if in strips, add the second hit into the list
631 if(isFineIDInStrip(fineID))
632 {
633 std::shared_ptr<FPGATrackSimHit> guessedSecondHitPtr = std::make_shared<FPGATrackSimHit>();
634 guessedSecondHitPtr->setX(0);
635 guessedSecondHitPtr->setY(0);
636 guessedSecondHitPtr->setZ(0);
637 guessedSecondHitPtr->setLayer( guessedLayer + 1 );
638 guessedSecondHitPtr->setHitType(HitType::guessed);
639
640 guessedSecondHitPtr->setDetType(SiliconTech::strip);
641 guessedHitPtr->setPhysLayer(1); // it is just to set the side, it's not actually 1
642 // Make sure that the hit is inside the boundaries of the layers
643 if((guessedLayer + 1) < (m_nLayers_1stStage + m_nLayers_2ndStage ))
644 {
645 hits.push_back(guessedSecondHitPtr);
646 }
647 }
648
649
650
651 return StatusCode::SUCCESS;
652
653}
654
655StatusCode FPGATrackSimNNPathfinderExtensionTool::findHitinNextStripLayer(std::shared_ptr<const FPGATrackSimHit> hitToSearch, const std::vector<std::shared_ptr<const FPGATrackSimHit>>& hitList, std::vector<std::shared_ptr<const FPGATrackSimHit>>& hits)
656{
657 float EPSILON = 0.00001;
658 for (const std::shared_ptr<const FPGATrackSimHit>& hit: hitList)
659 {
660 if (abs(hit->getX() - hitToSearch->getX()) < EPSILON && abs(hit->getY() - hitToSearch->getY()) < EPSILON && abs(hit->getZ() - hitToSearch->getZ()) < EPSILON && hit->getIdentifierHash() != hitToSearch->getIdentifierHash())
661 {
662 hits.push_back(hit);
663 return StatusCode::SUCCESS;
664 }
665
666 }
667
668 ATH_MSG_WARNING("Didn't find a matching space point");
669
670 return StatusCode::FAILURE;
671}
672
674{
675 if(!m_debugEvent) return;
676
677 // print this road
678 if(m_debugEvent)
679 {
680 std::vector<std::shared_ptr<const FPGATrackSimHit>> hitsR;
681 for (auto &hit : currentRoad.getHits()) {
682 hitsR.push_back(hit);
683 }
684
685 // If outside in, sort in increasing R, otherwise, decreasing R
686 if(m_doOutsideIn)
687 {
688 std::sort(hitsR.begin(), hitsR.end(), [](auto& a, auto& b){
689 if(a->getR() == b->getR()) return a->getLayer() < b->getLayer();
690 return a->getR() < b->getR();
691 });
692 }
693 else
694 {
695 std::sort(hitsR.begin(), hitsR.end(), [](auto& a, auto& b){
696 return a->getR() > b->getR();
697 });
698 }
699
700 for (unsigned long i = 0; i < hitsR.size(); i++)
701 {
702 ATH_MSG_DEBUG("Hit i "<<i<<" X: "<<hitsR[i]->getX()<<" Y: "<<hitsR[i]->getY()<<" Z: "<<hitsR[i]->getZ()<<" R: "<<hitsR[i]->getR()<< " hitType: "<<hitsR[i]->getHitType()<<" getDetType: "<<hitsR[i]->getDetType() << "phi = " << hitsR[i]->getGPhi());
703 }
704 }
705
706}
707
708StatusCode FPGATrackSimNNPathfinderExtensionTool::getLastLayer(miniRoad& currentRoad, unsigned& lastHitLayer, std::shared_ptr<const FPGATrackSimHit>& lastHit)
709{
710 lastHitLayer = currentRoad.getNHits()-1;
711 lastHit = currentRoad.getHit(lastHitLayer);
712 return StatusCode::SUCCESS;
713
714}
Scalar phi() const
phi method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
constexpr float EPSILON
bool isFineIDInStrip(long ID)
long getFineID(const FPGATrackSimHit &hit)
bool isFineIDInPixel(long ID)
Default track extension algorithm to produce "second stage" roads.
This file declares a class that stores the module IDs of the sectors.
uint32_t layer_bitmask_t
static Double_t a
#define pi
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)
void setHitLayers(layer_bitmask_t hit_layers)
void setY(float v)
void setHits(std::vector< std::vector< std::shared_ptr< const FPGATrackSimHit > > > &&hits)
void setXBin(unsigned v)
void setYBin(unsigned v)
void setRoadID(int roadID)
void setSubRegion(int v)
layer_bitmask_t getWCLayers() const
void setX(float v)
void setWCLayers(layer_bitmask_t wc_layers)
int r
Definition globals.cxx:22
int count(std::string s, const std::string &regx)
count how many occurances of a regx are in a string
Definition hcg.cxx:146
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
void reverse(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of reverse for DataVector/List.
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::shared_ptr< const FPGATrackSimHit > getHit(size_t layer) const
void setHits(std::vector< std::shared_ptr< const FPGATrackSimHit > > &hits)