ATLAS Offline Software
Loading...
Searching...
No Matches
FPGATrackSimNNPathfinderExtensionTool Class Reference

#include <FPGATrackSimNNPathfinderExtensionTool.h>

Inheritance diagram for FPGATrackSimNNPathfinderExtensionTool:
Collaboration diagram for FPGATrackSimNNPathfinderExtensionTool:

Public Member Functions

virtual StatusCode initialize () override
virtual StatusCode extendTracks (const std::vector< std::shared_ptr< const FPGATrackSimHit > > &hits, const FPGATrackSimTrackCollection &tracks, std::vector< FPGATrackSimRoad > &roads) override
virtual StatusCode setupSlices (FPGATrackSimLogicalEventInputHeader *slicedHitHeader) override

Private Member Functions

StatusCode fillInputTensorForNN (miniRoad &thisRoad, std::vector< float > &inputTensorValues)
StatusCode getPredictedHit (std::vector< float > &inputTensorValues, std::vector< float > &outputTensorValues, long &fineID)
StatusCode getPredictedHitBatched (const std::vector< std::vector< float > > &batchInputTensors, std::vector< std::vector< float > > &batchOutputTensors, std::vector< long > &batchFineIDs)
StatusCode addHitToRoad (miniRoad &newroad, miniRoad &currentRoad, const 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)
StatusCode getLastLayer (miniRoad &currentRoad, unsigned &lastHitLayer, std::shared_ptr< const FPGATrackSimHit > &lastHit)
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)
void printRoad (miniRoad &currentRoad)

Static Private Member Functions

static float getXScale ()
static float getYScale ()
static float getZScale ()
static float getRScale ()
static float getPhiScale ()

Private Attributes

ServiceHandle< IFPGATrackSimMappingSvcm_FPGATrackSimMapping {this, "FPGATrackSimMappingSvc", "FPGATrackSimMappingSvc"}
ServiceHandle< ITHistSvc > m_tHistSvc {this, "THistSvc", "THistSvc"}
ServiceHandle< IChronoStatSvc > m_chronoSvc {this, "ChronoStatSvc", "ChronoStatSvc"}
Gaudi::Property< unsigned > m_maxMiss { this, "threshold", 2, "Maximum number of missing hits to reject a road"}
Gaudi::Property< std::string > m_region { this, "OutputRegion", "", "region ID"}
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"}
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"}
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"}
Gaudi::Property< std::vector< int > > m_windowFineID { this, "windowFineID", {0}, "Fine ID indexing for windows"}
Gaudi::Property< float > m_lowPtValueForWindowRScaling { this, "lowPtValueWindowR", -1, "Value in MeV below which we scale the r window size"}
Gaudi::Property< float > m_lowPtWindowRScaling {this, "lowPtRScaling", 1.0, "Scaling factor for low pt in R"}
Gaudi::Property< float > m_lowPtValueForWindowZScaling { this, "lowPtValueWindowZ", -1, "Value in MeV below which we scale the r window size"}
Gaudi::Property< float > m_lowPtWindowZScaling {this, "lowPtZScaling", 1.0, "Scaling factor for low pt in Z"}
Gaudi::Property< float > m_lowPtValueForWindowPhiScaling { this, "lowPtValueWindowPhi", -1, "Value in MeV below which we scale the phi window size"}
Gaudi::Property< float > m_lowPtWindowPhiScaling {this, "lowPtPhiScaling", 1.0, "Scaling factor for low pt in Phi"}
Gaudi::Property< float > m_missedHitRScaling {this, "missedHitRScaling", -1, "Amount to scale R window if previous hit was missed. Negative means this is disabled"}
Gaudi::Property< float > m_missedHitZScaling {this, "missedHitZScaling", -1, "Amount to scale Z window if previous hit was missed. Negative means this is disabled"}
Gaudi::Property< float > m_missedHitPhiScaling {this, "missedHitPhiScaling", -1, "Amount to scale Phi window if previous hit was missed. Negative means this is disabled"}
Gaudi::Property< int > m_maxBranches { this, "maxBranches", -1, "Max number of branches before we stop, if negative this is disabled"}
Gaudi::Property< bool > m_doOutsideIn { this, "doOutsideIn", true, "Setup the tool so it's doing outside in extrap"}
Gaudi::Property< int > m_predictionWindowLength { this, "predictionWindowLength", 3, "Length of hits needed for prediction"}
Gaudi::Property< bool > m_useCartesian { this, "useCartesian", true, "If true, NNs use Cartestian coordinates. If false,they use cylindrical coordiantes"}
Gaudi::Property< unsigned int > m_batchSize { this, "batchSize", 1, "Batch size for NN inference (1 = sequential, >1 for true batching on GPU)"}
std::vector< FPGATrackSimRoadm_roads
unsigned m_nLayers_1stStage = 0
unsigned m_nLayers_2ndStage = 0
bool m_debugEvent = false
FPGATrackSimLogicalEventInputHeaderm_slicedHitHeader = nullptr
OnnxRuntimeBase m_extensionVolNN
OnnxRuntimeBase m_extensionHitNN

Detailed Description

Definition at line 128 of file FPGATrackSimNNPathfinderExtensionTool.h.

Member Function Documentation

◆ addHitToRoad()

StatusCode FPGATrackSimNNPathfinderExtensionTool::addHitToRoad ( miniRoad & newroad,
miniRoad & currentRoad,
const std::vector< std::shared_ptr< const FPGATrackSimHit > > & hits )
private

Definition at line 754 of file FPGATrackSimNNPathfinderExtensionTool.cxx.

755{
756 newroad.setHits(currentRoad.getHits());
757 newroad.addHits(hits);
758 return StatusCode::SUCCESS;
759}
void addHits(const std::vector< std::shared_ptr< const FPGATrackSimHit > > &hits)
std::vector< std::shared_ptr< const FPGATrackSimHit > > & getHits()
void setHits(std::vector< std::shared_ptr< const FPGATrackSimHit > > &hits)

◆ extendTracks()

StatusCode FPGATrackSimNNPathfinderExtensionTool::extendTracks ( const std::vector< std::shared_ptr< const FPGATrackSimHit > > & hits,
const FPGATrackSimTrackCollection & tracks,
std::vector< FPGATrackSimRoad > & roads )
overridevirtual

Definition at line 160 of file FPGATrackSimNNPathfinderExtensionTool.cxx.

162 {
163
164 // Reset the internal second stage roads storage.
165 roads.clear();
166 m_roads.clear();
167 const FPGATrackSimRegionMap* rmap_2nd = m_FPGATrackSimMapping->SubRegionMap_2nd();
168
169 // Create one "tower" per slice for this event.
170 // Note that there now might be only one "slice", at least for the time being.
171 if (m_slicedHitHeader) {
172 for (int ireg = 0; ireg < rmap_2nd->getNRegions(); ireg++) {
173 FPGATrackSimTowerInputHeader tower = FPGATrackSimTowerInputHeader(ireg);
174 m_slicedHitHeader->addTower(tower);
175 }
176 }
177
178 // Build spatial index (maps) once per event
179 HitSpatialIndex hitIndex;
180 {
181 Athena::Chrono chronoBuildIndex("NNPathfinder:BuildSpatialIndex", m_chronoSvc.get());
182 hitIndex.build(hits);
183 }
184
185 if(m_debugEvent) ATH_MSG_DEBUG("Got: "<<tracks.size()<<" tracks to extrapolate");
186 // Now, loop over the tracks.
187 for (const FPGATrackSimTrack& track : tracks) {
188 Athena::Chrono chronoTrackLoop("NNPathfinder:TrackLoop", m_chronoSvc.get());
189 if(m_debugEvent) ATH_MSG_DEBUG("\033[1;31m-------------------------- extraploating Track ------------------ \033[0m");
190 if (track.passedOR() == 0) {
191 continue;
192 }
193 const auto& hitsOnTrack = track.getFPGATrackSimHitPtrs();
194 miniRoad road;
195 float pt = track.getPt();
196
197 for (const auto& hit_ptr : hitsOnTrack) {
198 if (!hit_ptr) {
199 ATH_MSG_ERROR("Null hit pointer in track");
200 return StatusCode::FAILURE;
201 }
202 road.addHit(hit_ptr); // shared_ptr already points to active hit
203 }
204
205 if (m_debugEvent) {
206 ATH_MSG_DEBUG("-----------------Hits in event");
207 for (const auto& hit : hits) {
208 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());
209 }
210 }
211 // Using deque instead of vector for O(1) pop_front since we need operate only the first road each iteration
212 std::deque<miniRoad> roadsToExtrapolate;
213 roadsToExtrapolate.push_back(std::move(road));
214
215 std::vector<miniRoad> completedRoads;
216
217 int count = 0;
218 // FIXED: Add maximum iteration limit to prevent infinite loops
219 const int MAX_ROADS = 10000;
220
221 // Batching structures
222 std::vector<std::pair<miniRoad, std::vector<float>>> roadsAwaitingInference; // pairs of (road, input_tensor)
223
224 while(!roadsToExtrapolate.empty() && count < MAX_ROADS ) {
225 miniRoad currentRoad = std::move(roadsToExtrapolate.front());
226
227 // Erase this road from the deque
228 roadsToExtrapolate.pop_front();
229 count ++;
230 if(m_debugEvent) {
231 ATH_MSG_DEBUG("\033[1;31m-------------------------- extraploating road "<< count << "------------------ \033[0m");
232 printRoad(currentRoad);
233 }
234 // Check exit condition
235 if (currentRoad.getNHits() >= (m_nLayers_1stStage+m_nLayers_2ndStage))
236 {
237 completedRoads.push_back(std::move(currentRoad));
238 continue; // this one is done
239 }
240
241 // Prepare input tensor for NN
242 std::vector<float> inputTensorValues;
243 if (!fillInputTensorForNN(currentRoad, inputTensorValues)) {
244 ATH_MSG_WARNING("Failed to create input tensor for this road");
245 continue;
246 }
247
248 // Collect road for batching
249 roadsAwaitingInference.push_back(std::make_pair(std::move(currentRoad), inputTensorValues));
250
251 // When batch is full OR queue is empty, process batch
252 bool processBatch = (roadsAwaitingInference.size() >= (size_t)m_batchSize) || roadsToExtrapolate.empty();
253
254 if (processBatch && !roadsAwaitingInference.empty()) {
255 // Prepare batched inputs
256 std::vector<std::vector<float>> batchInputTensors;
257 std::vector<std::pair<miniRoad, std::vector<float>>> roadBatch = std::move(roadsAwaitingInference);
258 roadsAwaitingInference.clear();
259
260 for (const auto& [road, tensor] : roadBatch) {
261 batchInputTensors.push_back(tensor);
262 }
263
264 // Run batched inference
265 std::vector<std::vector<float>> batchOutputTensors;
266 std::vector<long> batchFineIDs;
267 {
268 Athena::Chrono chronoNN("NNPathfinder:NNInference", m_chronoSvc.get());
269 if (getPredictedHitBatched(batchInputTensors, batchOutputTensors, batchFineIDs) != StatusCode::SUCCESS) {
270 ATH_MSG_WARNING("Batched NN inference failed");
271 continue;
272 }
273 }
274
275 // Process each result in the batch
276 for (size_t batchIdx = 0; batchIdx < roadBatch.size(); ++batchIdx) {
277 miniRoad processRoad = std::move(roadBatch[batchIdx].first);
278 std::vector<float> predhit = batchOutputTensors[batchIdx];
279 long fineID = batchFineIDs[batchIdx];
280
281 // Check if exit conditions are met
282 if (m_doOutsideIn) {
283 // Make sure we are not predicting inside the inner most layer (x and y < 25, or r < 25)
284 // If we are, road is done
285 if ((m_useCartesian && (abs(predhit[0]) < 25 && abs(predhit[1]) < 25)) ||
286 (!m_useCartesian && abs(predhit[0]) < 25))
287 {
288 completedRoads.push_back(std::move(processRoad));
289 continue;
290 }
291 }
292 else {
293 // Make sure we are not predicting outside the outer most layer
294 // if we are, road is done
295 double radius = std::sqrt(predhit[0] * predhit[0] + predhit[1] * predhit[1]);
296 if ((m_useCartesian && (abs(predhit[0]) > 1024 || abs(predhit[1]) > 1024 || radius > 1024 || abs(predhit[2]) > 3000)) ||
297 (!m_useCartesian && (abs(predhit[0]) > 1024 || abs(predhit[2]) > 3000))) {
298 completedRoads.push_back(std::move(processRoad));
299 continue;
300 }
301 }
302 if(m_debugEvent) {
303 ATH_MSG_DEBUG("Predicted hit at: " << predhit[0] << " " << predhit[1] << " " << predhit[2]);
304 }
305
306 // Now search for the hits
307 bool foundhitForRoad = false;
308 if(fineID == 215){
309 ATH_MSG_DEBUG("Stopping condition reached");
310 completedRoads.push_back(std::move(processRoad));
311 continue;
312 }
313
314 // Get the last layer and hit in the road
315 unsigned lastLayerInRoad = 0;
316 std::shared_ptr<const FPGATrackSimHit> lastHit;
317 if(!getLastLayer(processRoad, lastLayerInRoad, lastHit) or !lastHit) {
318 ATH_MSG_WARNING("Failed to find last layer this road");
319 continue;
320 }
321 unsigned layer = lastLayerInRoad+1; // the layer we're looking to find
322 bool lastHitWasReal = lastHit->isReal();
323 float lastHitR = lastHit->getR();
324 if(layer >= (m_nLayers_1stStage + m_nLayers_2ndStage)) {
325 completedRoads.push_back(std::move(processRoad));
326 continue;
327 }
328 unsigned int hitsInWindow = 0;
329
330 // Cache predicted values and window parameters once per road, not per hit
331 const double predr = (m_useCartesian ? sqrt(predhit[0] * predhit[0] + predhit[1] * predhit[1]) : predhit[0]);
332 const double predphi = predhit[1];
333 const double predz = predhit[2];
334
335 // Calculate window parameters once
336 double windowR = m_windowR[0];
337 double windowPhi = m_windowPhi[0];
338 double windowZ = m_windowZ[0];
339 int fineID_index = 0;
340
341 if (m_windowZ.size() > 1 || m_windowR.size() > 1 || m_windowPhi.size() > 1) {
342 auto fineID_it = std::find(m_windowFineID.begin(), m_windowFineID.end(), fineID);
343 if (fineID_it == m_windowFineID.end()) {
344 ATH_MSG_DEBUG("No windows for predicted fineID " << fineID << ", using maximum in provided list instead!");
345 fineID_index = -1;
346 }
347 else {
348 fineID_index = fineID_it - m_windowFineID.begin();
349 }
350 }
351 if (m_windowR.size() > 1) {
352 windowR = (fineID_index == -1) ?
353 *std::max_element(m_windowR.begin(), m_windowR.end()) :
354 m_windowR[fineID_index];
355 }
356 if (m_windowZ.size() > 1) {
357 windowZ = (fineID_index == -1) ?
358 *std::max_element(m_windowZ.begin(), m_windowZ.end()) :
359 m_windowZ[fineID_index];
360 }
361 if (m_windowPhi.size() > 1) {
362 windowPhi = (fineID_index == -1) ?
363 *std::max_element(m_windowPhi.begin(), m_windowPhi.end()) :
364 m_windowPhi[fineID_index];
365 }
366
367 // Apply scaling factors
368 if (m_missedHitRScaling > 0 && !lastHitWasReal) windowR *= m_missedHitRScaling;
369 if (m_missedHitZScaling > 0 && !lastHitWasReal) windowZ *= m_missedHitZScaling;
370 if (m_missedHitPhiScaling > 0 && !lastHitWasReal) windowPhi *= m_missedHitPhiScaling;
371 if (pt < m_lowPtValueForWindowRScaling.value()) windowR *= m_lowPtWindowRScaling.value();
372 if (pt < m_lowPtValueForWindowZScaling.value()) windowZ *= m_lowPtWindowZScaling.value();
373 if (pt < m_lowPtValueForWindowPhiScaling.value()) windowPhi *= m_lowPtWindowPhiScaling.value();
374
375 // List of all the hits, with their distances to the predicted point
376 std::vector<std::vector<std::shared_ptr<const FPGATrackSimHit>>> listofHitsFound;
377
378 {
379 Athena::Chrono chronoHitSearch("NNPathfinder:HitSearchLoop", m_chronoSvc.get());
380
381 // Use spatial index (HitSpatialIndex) to get only candidate hits with matching fineID
382 const auto* candidateHits = hitIndex.getHits(fineID);
383 if (candidateHits) {
384 listofHitsFound.reserve(candidateHits->size());
385
386 for (const auto& hit_ptr : *candidateHits) {
387 const auto& hit = *hit_ptr;
388 // Apply direction filter
389 const float hitr = hit.getR();
390 if (m_doOutsideIn && hitr > lastHitR) continue;
391 if (!m_doOutsideIn && hitr < lastHitR) continue;
392
393 if(m_debugEvent) {
394 ATH_MSG_DEBUG("In the hit loop hit at x: " << hit.getX() << " y " << hit.getY() << " z " << hit.getZ() << " phi " << hit.getGPhi());
395 }
396
397 // Check if hit is within window
398 const double hitz = hit.getZ();
399 const double hitphi = hit.getGPhi();
400
401 const double dr = abs(hitr - predr);
402 const double dz = abs(hitz - predz);
403 double dphi = abs(hitphi - predphi);
404 while (dphi > pi) dphi -= pi;
405
406 const bool inWindow = (m_useCartesian && dr < windowR && dz < windowZ) ||
407 (!m_useCartesian && dphi < windowPhi && dz < windowZ && dr < windowR);
408
409 if (!inWindow) continue;
410
411 // Only create shared_ptr when hit passes all filters
412 std::vector<std::shared_ptr<const FPGATrackSimHit>> theseHits{ hit_ptr };
413 hitsInWindow = hitsInWindow + 1;
414
415 // Handle strip space points
416 if(hit.isStrip()) {
417 if (hit.getHitType() == HitType::spacepoint) {
418 // Use spatial index for strip matching
419 const float EPSILON = 0.00001f;
420 const float searchX = hit.getX();
421 const float searchY = hit.getY();
422 const float searchZ = hit.getZ();
423 const auto searchHash = hit.getIdentifierHash();
424 bool found = false;
425
426 const auto* coordCandidates = hitIndex.getHitsByCoord(searchX, searchY, searchZ);
427 if (coordCandidates) {
428 for (const auto& candidateHitPtr : *coordCandidates) {
429 const auto& candidateHit = *candidateHitPtr;
430 if (candidateHit.getIdentifierHash() != searchHash &&
431 abs(candidateHit.getX() - searchX) < EPSILON &&
432 abs(candidateHit.getY() - searchY) < EPSILON &&
433 abs(candidateHit.getZ() - searchZ) < EPSILON) {
434 theseHits.push_back(candidateHitPtr);
435 found = true;
436 break;
437 }
438 }
439 }
440 if (!found) {
441 ATH_MSG_WARNING("For a SP in layer " << layer << " Couldn't find a matching strip SP");
442 }
443 }
444 else {
445 std::shared_ptr<FPGATrackSimHit> guessedSecondHitPtr = std::make_shared<FPGATrackSimHit>();
446 guessedSecondHitPtr->setX(0);
447 guessedSecondHitPtr->setY(0);
448 guessedSecondHitPtr->setZ(0);
449 guessedSecondHitPtr->setPhysLayer(lastHit->getPhysLayer(true)+1);
450 guessedSecondHitPtr->setHitType(HitType::undefined);
451 if(isFineIDInStrip(fineID)) guessedSecondHitPtr->setDetType(SiliconTech::strip);
452 else guessedSecondHitPtr->setDetType(SiliconTech::pixel);
453
454 theseHits.push_back(guessedSecondHitPtr);
455 }
456 }
457 // Store the hits for now
458 listofHitsFound.push_back(std::move(theseHits));
459 }
460 }
461 } // end HitSearchLoop timing
462
463 {
464 Athena::Chrono chronoSort("NNPathfinder:HitSorting", m_chronoSvc.get());
465
466 // Use partial_sort instead of sort since we only need the top N hits
467 const size_t nToSort = (m_maxBranches.value() >= 0) ?
468 std::min(static_cast<size_t>(m_maxBranches.value()), listofHitsFound.size()) :
469 listofHitsFound.size();
470
471 // Sort the hit by the distance
472 if (m_useCartesian) {
473 auto comparatorFunc = [&](const auto& a, const auto& b) {
474 return cartesianComparator(a, b, predr, predz);
475 };
476 if (nToSort < listofHitsFound.size()) {
477 std::partial_sort(listofHitsFound.begin(),
478 listofHitsFound.begin() + nToSort,
479 listofHitsFound.end(),
480 comparatorFunc);
481 } else {
482 std::sort(listofHitsFound.begin(), listofHitsFound.end(), comparatorFunc);
483 }
484 }
485 else {
486 const double zScale = getZScale();
487 const double phiScale = getPhiScale();
488 const double rScale = getRScale();
489 const double zScale2 = zScale * zScale;
490 const double phiScale2 = phiScale * phiScale;
491 const double rScale2 = rScale * rScale;
492
493 auto comparatorFunc = [&](const auto& a, const auto& b) {
494 return polarComparator(a, b, predr, predphi, predz, zScale2, phiScale2, rScale2);
495 };
496 if (nToSort < listofHitsFound.size()) {
497 std::partial_sort(listofHitsFound.begin(),
498 listofHitsFound.begin() + nToSort,
499 listofHitsFound.end(),
500 comparatorFunc);
501 } else {
502 std::sort(listofHitsFound.begin(), listofHitsFound.end(), comparatorFunc);
503 }
504 }
505 } // end HitSorting timing
506
507 // Select the top N hits
508 std::vector<std::vector<std::shared_ptr<const FPGATrackSimHit>>> cleanHitsToGrow;
509
510 // If max branches are limited, pick only the top N from the list of hits found at each level
511 if (m_maxBranches.value() >= 0) {
512 int nHitsToChoose = std::min(int(m_maxBranches.value()), int(listofHitsFound.size()));
513 cleanHitsToGrow.reserve(nHitsToChoose);
514 std::copy(listofHitsFound.begin(), listofHitsFound.begin() + nHitsToChoose, std::back_inserter(cleanHitsToGrow));
515 }
516 else {
517 cleanHitsToGrow = std::move(listofHitsFound);
518 }
519
520 {
521 Athena::Chrono chronoRoadBuilding("NNPathfinder:RoadBuilding", m_chronoSvc.get());
522 for (auto& hitsFound: cleanHitsToGrow) {
523
524 // We got a hit, lets make a road
525 miniRoad newRoad;
526 if(!addHitToRoad(newRoad, processRoad, std::move(hitsFound))) {
527 ATH_MSG_WARNING("Failed to make a new road");
528 continue;
529 }
530 roadsToExtrapolate.push_back(newRoad);
531 foundhitForRoad = true;
532 if(m_debugEvent) {
533 ATH_MSG_DEBUG("------ road grown with hit from layer "<<layer<<" to");
534 printRoad(newRoad);
535 }
536 }
537 // If the hit wasn't found, push a fake hit
538 if (!foundhitForRoad) {
539 // 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
540 if (processRoad.getNWCLayers() >= m_maxMiss) {
541 // we don't want this road, so we continue
542 continue;
543 }
544 else {
545 std::vector<std::shared_ptr<const FPGATrackSimHit>> theseHits;
546 // first make the fake hit that we will add
547 if (!getFakeHit(processRoad, predhit, fineID, theseHits)) {
548 ATH_MSG_WARNING("Failed adding a guessed hit in extrapolation");
549 continue;
550 }
551
552 if((isFineIDInPixel(fineID) && theseHits.size() != 1) || (isFineIDInStrip(fineID) && theseHits.size() != 2)) {
553 continue;
554 }
555 // add the hit to the road
556 miniRoad newroad;
557 if (!addHitToRoad(newroad, processRoad, std::move(theseHits))) {
558 ATH_MSG_WARNING("Failed making a new road with fake hit");
559 continue;
560 }
561 roadsToExtrapolate.push_back(std::move(newroad));
562 }
563 }
564 } // end RoadBuilding timing
565 } // end batch processing for individual road
566 } // end batch inference processing
567 } // end while loop for all roads
568 {
569 Athena::Chrono chronoRoadConversion("NNPathfinder:RoadConversion", m_chronoSvc.get());
570 for (const auto &miniroad : completedRoads) {
571 FPGATrackSimRoad road;
572 road.setWCLayers(miniroad.getWCLayers());
573 road.setHitLayers(miniroad.getHitLayers());
574 road.setRoadID(m_roads.size() - 1);
575 // Set the "Hough x" and "Hough y" using the track parameters.
576 road.setX(track.getPhi());
577 road.setY(track.getQOverPt());
578 road.setXBin(track.getHoughXBin());
579 road.setYBin(track.getHoughYBin());
580 road.setSubRegion(track.getSubRegion());
581
582 // just force the right number of layers now, in case we find fewer than expected (needed downstream)
583 std::vector<std::vector<std::shared_ptr<const FPGATrackSimHit>>> roadhits = miniroad.getVecHits();
584 unsigned nexpected = m_nLayers_1stStage+m_nLayers_2ndStage;
585 if (roadhits.size() > nexpected) { // cut off the last ones
586 roadhits.resize(nexpected);
587 }
588 else if (roadhits.size() < nexpected) { // fill with missing hits
589 for (unsigned layer = roadhits.size(); layer < nexpected; layer++) {
590 std::shared_ptr<FPGATrackSimHit> emptyHitPtr = std::make_shared<FPGATrackSimHit>();
591 emptyHitPtr->setX(0);
592 emptyHitPtr->setY(0);
593 emptyHitPtr->setZ(0);
594 emptyHitPtr->setLayer(layer);
595 emptyHitPtr->setHitType(HitType::wildcard);
596
597 roadhits.emplace_back(1,emptyHitPtr);
598 layer_bitmask_t wclayers = road.getWCLayers();
599 wclayers |= (1 << layer);
600 road.setWCLayers(wclayers);
601 }
602 }
603 road.setHits(std::move(roadhits));
604
605 m_roads.push_back(std::move(road));
606 }
607 } // end RoadConversion timing
608 }
609 // Copy the roads we found into the output argument and return success.
610 roads.reserve(m_roads.size());
611 for (FPGATrackSimRoad & r : m_roads)
612 {
613 if (r.getNWCLayers() >= m_maxMiss) continue; // extra check on this
614 roads.emplace_back(r);
615 }
616 ATH_MSG_DEBUG("Found " << roads.size() << " new roads in second stage.");
617
618 return StatusCode::SUCCESS;
619}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
constexpr float EPSILON
bool isFineIDInStrip(long ID)
bool isFineIDInPixel(long ID)
uint32_t layer_bitmask_t
static Double_t a
size_t size() const
Number of registered mappings.
#define pi
ServiceHandle< IFPGATrackSimMappingSvc > m_FPGATrackSimMapping
StatusCode getPredictedHitBatched(const std::vector< std::vector< float > > &batchInputTensors, std::vector< std::vector< float > > &batchOutputTensors, std::vector< long > &batchFineIDs)
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:148
@ layer
Definition HitInfo.h:79
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
void partial_sort(DataModel_detail::iterator< DVL > beg, DataModel_detail::iterator< DVL > mid, DataModel_detail::iterator< DVL > end)
Specialization of partial_sort for DataVector/List.
void addHit(const std::shared_ptr< const FPGATrackSimHit > &hit)

◆ fillInputTensorForNN()

StatusCode FPGATrackSimNNPathfinderExtensionTool::fillInputTensorForNN ( miniRoad & thisRoad,
std::vector< float > & inputTensorValues )
private

Definition at line 621 of file FPGATrackSimNNPathfinderExtensionTool.cxx.

622{
623 std::vector<std::shared_ptr<const FPGATrackSimHit>> hitsR;
624
625 std::vector<std::shared_ptr<const FPGATrackSimHit>> hits = thisroad.getHits();
626
627
628 for (unsigned ihit = 0; ihit < hits.size(); ihit++) {
629 if (ihit < m_nLayers_1stStage && !(hits[ihit]->isReal())) continue; // skip guessed hits for the 1st stage
630 hitsR.push_back(hits[ihit]);
631 }
632 // Sort in increasing 3D distance.
633 std::sort(hitsR.begin(), hitsR.end(), [](auto& a, auto& b){
634 double dist_a = std::hypot(a->getX(), a->getY(), a->getZ());
635 double dist_b = std::hypot(b->getX(), b->getY(), b->getZ());
636 return dist_a < dist_b;
637 });
638
639 if(m_debugEvent) ATH_MSG_DEBUG("hitsR");
640 for (const auto& thit : hitsR)
641 {
642 if (m_debugEvent) ATH_MSG_DEBUG(thit->getX() << " " << thit->getY() << " " << thit->getZ() << " and phi = " << thit->getGPhi());
643 }
644
645 // Remove all the duplicate space points
646 std::vector<std::shared_ptr<const FPGATrackSimHit>> cleanHits;
647 bool skipHit = false;
648 for (auto thit : hitsR)
649 {
650 if(skipHit)
651 {
652 skipHit = false;
653 continue;
654 }
655 if (thit->isPixel())
656 {
657 cleanHits.push_back(std::move(thit));
658 }
659 else if (thit->isStrip() && (thit->getHitType() == HitType::spacepoint))
660 {
661 // This is a proper strips SP, push the first hit back and skip the next one since its a duplicate
662 cleanHits.push_back(std::move(thit));
663 skipHit = true;
664 }
665 else if (thit->isStrip() && (thit->getHitType() == HitType::guessed))
666 {
667 // this is a guessed strip SP, push the first hit back and skip the next one since its a duplicate
668 cleanHits.push_back(std::move(thit));
669 skipHit = true;
670 }
671 else if (thit->isStrip() && (thit->getHitType() == HitType::undefined))
672 {
673 // this is a fake hit, inserted for a unpaired SP, continue
674 continue;
675 }
676 else if (thit->isStrip() && thit->isReal())
677 {
678 // What is left here is a unpaired hit, push it back
679 cleanHits.push_back(std::move(thit));
680 }
681 else
682 {
683 ATH_MSG_WARNING("No clue how to deal with this hit in the NN predicition ");
684 continue;
685 }
686 }
687
688 // Reverse the hits as we changed the ordering before
689 if(!m_doOutsideIn)
690 {
691 std::reverse(cleanHits.begin(), cleanHits.end());
692 }
693
694 // Select the top N hits
695 std::vector<std::shared_ptr<const FPGATrackSimHit>> hitsToEncode;
696 std::copy(cleanHits.begin(), cleanHits.begin() + m_predictionWindowLength, std::back_inserter(hitsToEncode));
697
698 if(m_debugEvent) ATH_MSG_DEBUG("Clean hits");
699 for (const auto& thit : cleanHits)
700 {
701 if (m_debugEvent) ATH_MSG_DEBUG(thit->getX() << " " << thit->getY() << " " << thit->getZ() << " " << thit->isStrip() << " and gphi = " << thit->getGPhi());
702 }
703
704 // Reverse this vector so we can encode it the format as expected from the NN
705 std::reverse(hitsToEncode.begin(), hitsToEncode.end());
706
707 if (m_debugEvent) ATH_MSG_DEBUG("Input for NN prediction");
708 for (const auto& thit : hitsToEncode)
709 {
710 if (m_useCartesian) {
711 inputTensorValues.push_back(thit->getX() / getXScale());
712 inputTensorValues.push_back(thit->getY() / getYScale());
713 inputTensorValues.push_back(thit->getZ() / getZScale());
714 }
715 else {
716 inputTensorValues.push_back(thit->getR() / getRScale());
717 inputTensorValues.push_back(thit->getGPhi() / getPhiScale());
718 inputTensorValues.push_back(thit->getZ() / getZScale());
719 }
720
721 if (m_debugEvent) ATH_MSG_DEBUG(thit->getX() << " " << thit->getY() << " " << thit->getZ() << " and gphi = " << thit->getGPhi());
722 }
723
724 return StatusCode::SUCCESS;
725
726}
void reverse(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of reverse for DataVector/List.

◆ findHitinNextStripLayer()

StatusCode FPGATrackSimNNPathfinderExtensionTool::findHitinNextStripLayer ( std::shared_ptr< const FPGATrackSimHit > hit,
const std::vector< std::shared_ptr< const FPGATrackSimHit > > & hitList,
std::vector< std::shared_ptr< const FPGATrackSimHit > > & hits )
private

◆ getFakeHit()

StatusCode FPGATrackSimNNPathfinderExtensionTool::getFakeHit ( miniRoad & currentRoad,
std::vector< float > & predhit,
const long & fineID,
std::vector< std::shared_ptr< const FPGATrackSimHit > > & hits )
private

Definition at line 761 of file FPGATrackSimNNPathfinderExtensionTool.cxx.

761 {
762
763 unsigned guessedLayer(0);
764
765 std::shared_ptr<FPGATrackSimHit> guessedHitPtr = std::make_shared<FPGATrackSimHit>();
766 if (m_useCartesian) {
767 guessedHitPtr->setX(predhit[0]);
768 guessedHitPtr->setY(predhit[1]);
769 guessedHitPtr->setZ(predhit[2]);
770 }
771 else {
772 double r = predhit[0];
773 double phi = predhit[1];
774 guessedHitPtr->setX(r*cos(phi));
775 guessedHitPtr->setY(r*sin(phi));
776 guessedHitPtr->setZ(predhit[2]);
777 }
778
779
780 if (m_doOutsideIn) { // outside in
781 // TODO
782 }
783 else { // nope, inside out
784 guessedLayer = currentRoad.getNHits();
785 }
786
787 guessedHitPtr->setLayer(guessedLayer);
788 guessedHitPtr->setHitType(HitType::guessed);
789
790 if(isFineIDInStrip(fineID)) {
791 guessedHitPtr->setDetType(SiliconTech::strip);
792 guessedHitPtr->setPhysLayer(0); // it is just to set the side, it's not actually 0
793 }
794 else guessedHitPtr->setDetType(SiliconTech::pixel);
795
796
797 // Make sure that the hit is inside the boundaries of the layers
798 if(guessedLayer < (m_nLayers_1stStage + m_nLayers_2ndStage ))
799 {
800 hits.push_back(guessedHitPtr);
801 }
802 // if in strips, add the second hit into the list
803 if(isFineIDInStrip(fineID))
804 {
805 std::shared_ptr<FPGATrackSimHit> guessedSecondHitPtr = std::make_shared<FPGATrackSimHit>();
806 guessedSecondHitPtr->setX(0);
807 guessedSecondHitPtr->setY(0);
808 guessedSecondHitPtr->setZ(0);
809 guessedSecondHitPtr->setLayer( guessedLayer + 1 );
810 guessedSecondHitPtr->setHitType(HitType::guessed);
811
812 guessedSecondHitPtr->setDetType(SiliconTech::strip);
813 guessedHitPtr->setPhysLayer(1); // it is just to set the side, it's not actually 1
814 // Make sure that the hit is inside the boundaries of the layers
815 if((guessedLayer + 1) < (m_nLayers_1stStage + m_nLayers_2ndStage ))
816 {
817 hits.push_back(guessedSecondHitPtr);
818 }
819 }
820
821
822
823 return StatusCode::SUCCESS;
824
825}
Scalar phi() const
phi method

◆ getLastLayer()

StatusCode FPGATrackSimNNPathfinderExtensionTool::getLastLayer ( miniRoad & currentRoad,
unsigned & lastHitLayer,
std::shared_ptr< const FPGATrackSimHit > & lastHit )
private

Definition at line 862 of file FPGATrackSimNNPathfinderExtensionTool.cxx.

863{
864 lastHitLayer = currentRoad.getNHits()-1;
865 lastHit = currentRoad.getHit(lastHitLayer);
866 return StatusCode::SUCCESS;
867
868}
std::shared_ptr< const FPGATrackSimHit > getHit(size_t layer) const

◆ getPhiScale()

float FPGATrackSimNNPathfinderExtensionTool::getPhiScale ( )
inlinestaticprivate

Definition at line 188 of file FPGATrackSimNNPathfinderExtensionTool.h.

188{return 3.15;};

◆ getPredictedHit()

StatusCode FPGATrackSimNNPathfinderExtensionTool::getPredictedHit ( std::vector< float > & inputTensorValues,
std::vector< float > & outputTensorValues,
long & fineID )
private

Definition at line 729 of file FPGATrackSimNNPathfinderExtensionTool.cxx.

730{
731 std::vector<float> NNVoloutput = m_extensionVolNN.runONNXInference(inputTensorValues);
732 fineID = std::distance(NNVoloutput.begin(),std::max_element(NNVoloutput.begin(), NNVoloutput.end()));
733 // Insert the output for the second stage NN
734 inputTensorValues.insert(inputTensorValues.end(), NNVoloutput.begin(), NNVoloutput.end()); //now we append the first NN to the list of coordinates
735
736 // use the above to predict the next hit position
737 outputTensorValues = m_extensionHitNN.runONNXInference(inputTensorValues);
738
739 // now scale back
740 if (m_useCartesian) {
741 outputTensorValues[0] *= getXScale();
742 outputTensorValues[1] *= getYScale();
743 outputTensorValues[2] *= getZScale();
744 }
745 else {
746 outputTensorValues[0] *= getRScale();
747 outputTensorValues[1] *= getPhiScale();
748 outputTensorValues[2] *= getZScale();
749 }
750
751 return StatusCode::SUCCESS;
752}

◆ getPredictedHitBatched()

StatusCode FPGATrackSimNNPathfinderExtensionTool::getPredictedHitBatched ( const std::vector< std::vector< float > > & batchInputTensors,
std::vector< std::vector< float > > & batchOutputTensors,
std::vector< long > & batchFineIDs )
private

Definition at line 870 of file FPGATrackSimNNPathfinderExtensionTool.cxx.

873{
874 if (batchInputTensors.empty()) {
875 return StatusCode::SUCCESS;
876 }
877
878 size_t batchSize = batchInputTensors.size();
879 size_t featureSize = batchInputTensors[0].size();
880
881 // Convert to Eigen matrix format for proper batch dimension handling
882 // Rows = batch_size, Cols = features
883 NetworkBatchInput volInputMatrix(batchSize, featureSize);
884 for (size_t i = 0; i < batchSize; ++i) {
885 for (size_t j = 0; j < featureSize; ++j) {
886 volInputMatrix(i, j) = batchInputTensors[i][j];
887 }
888 }
889
890 // Run volume NN in batch using proper tensor format
891 auto volNNBatchedOutput = m_extensionVolNN.runONNXInference(volInputMatrix);
892
893 // Extract fineIDs from volume NN output (one per sample in batch)
894 // Output shape: [batch_size, num_classes]
895 batchFineIDs.reserve(batchSize);
896
897 for (size_t i = 0; i < volNNBatchedOutput.size(); ++i) {
898 const auto& output = volNNBatchedOutput[i];
899 auto maxIdx = std::distance(output.begin(), std::max_element(output.begin(), output.end()));
900 batchFineIDs.push_back(maxIdx);
901 }
902
903 // Prepare input for hit NN: concatenate original inputs with volume NN outputs
904 // Hit NN takes [original_features, volume_nn_output] as input
905 size_t numClasses = volNNBatchedOutput[0].size();
906 NetworkBatchInput hitInputMatrix(batchSize, featureSize + numClasses);
907
908 for (size_t i = 0; i < batchSize; ++i) {
909 // First part: original input features
910 for (size_t j = 0; j < featureSize; ++j) {
911 hitInputMatrix(i, j) = batchInputTensors[i][j];
912 }
913 // Second part: volume NN outputs
914 for (size_t j = 0; j < numClasses; ++j) {
915 hitInputMatrix(i, featureSize + j) = volNNBatchedOutput[i][j];
916 }
917 }
918
919 // Run hit NN in batch using proper tensor format
920 auto hitNNBatchedOutput = m_extensionHitNN.runONNXInference(hitInputMatrix);
921
922 batchOutputTensors.reserve(batchSize);
923
924 // Output (already in correct format: vector<vector<float>>) with batch_size samples
925 for (size_t i = 0; i < hitNNBatchedOutput.size(); ++i) {
926 std::vector<float> output = hitNNBatchedOutput[i];
927
928 // Scale back to original units
929 if (m_useCartesian) {
930 output[0] *= getXScale();
931 output[1] *= getYScale();
932 output[2] *= getZScale();
933 } else {
934 output[0] *= getRScale();
935 output[1] *= getPhiScale();
936 output[2] *= getZScale();
937 }
938
939 batchOutputTensors.push_back(std::move(output));
940 }
941
942 return StatusCode::SUCCESS;
943}
Eigen::Matrix< float, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > NetworkBatchInput
float j(const xAOD::IParticle &, const xAOD::TrackMeasurementValidation &hit, const Eigen::Matrix3d &jab_inv)
output
Definition merge.py:16

◆ getRScale()

float FPGATrackSimNNPathfinderExtensionTool::getRScale ( )
inlinestaticprivate

Definition at line 187 of file FPGATrackSimNNPathfinderExtensionTool.h.

187{return 1015.;};

◆ getXScale()

float FPGATrackSimNNPathfinderExtensionTool::getXScale ( )
inlinestaticprivate

Definition at line 184 of file FPGATrackSimNNPathfinderExtensionTool.h.

184{ return 1015.;};

◆ getYScale()

float FPGATrackSimNNPathfinderExtensionTool::getYScale ( )
inlinestaticprivate

Definition at line 185 of file FPGATrackSimNNPathfinderExtensionTool.h.

185{ return 1015.;};

◆ getZScale()

float FPGATrackSimNNPathfinderExtensionTool::getZScale ( )
inlinestaticprivate

Definition at line 186 of file FPGATrackSimNNPathfinderExtensionTool.h.

186{ return 3000.;};

◆ initialize()

StatusCode FPGATrackSimNNPathfinderExtensionTool::initialize ( )
overridevirtual

Definition at line 118 of file FPGATrackSimNNPathfinderExtensionTool.cxx.

118 {
119
120 // Retrieve the mapping service.
122
123 // hard code this for now but we may chance in the future
126
127 if (m_windowR.size() != 1 && m_windowR.size() != m_windowFineID.size()) {
128 ATH_MSG_ERROR("Window r size = " << m_windowR << " is not equal to 1 (for all layers) and not equal to " << m_windowFineID.size());
129 return StatusCode::FAILURE;
130 }
131
132 if (m_windowPhi.size() != 1 && m_windowPhi.size() != m_windowFineID.size()) {
133 ATH_MSG_ERROR("Window phi size = " << m_windowPhi << " is not equal to 1 (for all layers) and not equal to " << m_windowFineID.size());
134 return StatusCode::FAILURE;
135 }
136
137 if (m_windowZ.size() != 1 && m_windowZ.size() != m_windowFineID.size()) {
138 ATH_MSG_ERROR("Window z size = " << m_windowZ << " is not equal to 1 (for all layers) and not equal to " << m_windowFineID.size());
139 return StatusCode::FAILURE;
140 }
141
142 if (m_FPGATrackSimMapping->getExtensionNNVolMapString() != "" && m_FPGATrackSimMapping->getExtensionNNHitMapString() != "") {
143 ATH_MSG_INFO("Initializing extension hit NN with string = " << m_FPGATrackSimMapping->getExtensionNNHitMapString());
144 m_extensionHitNN.initialize(m_FPGATrackSimMapping->getExtensionNNHitMapString());
145 ATH_MSG_INFO("Initializing volume NN with string = " << m_FPGATrackSimMapping->getExtensionNNVolMapString());
146 m_extensionVolNN.initialize(m_FPGATrackSimMapping->getExtensionNNVolMapString());
147 }
148 else {
149 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.");
150 return StatusCode::FAILURE;
151 }
152
153 ATH_CHECK(m_tHistSvc.retrieve());
154 ATH_CHECK(m_chronoSvc.retrieve());
155
156 return StatusCode::SUCCESS;
157}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_INFO(x)

◆ printRoad()

void FPGATrackSimNNPathfinderExtensionTool::printRoad ( miniRoad & currentRoad)
private

Definition at line 827 of file FPGATrackSimNNPathfinderExtensionTool.cxx.

828{
829 if (!m_debugEvent) return;
830
831 // print this road
832 if (m_debugEvent)
833 {
834 std::vector<std::shared_ptr<const FPGATrackSimHit>> hitsR;
835 for (const auto& hit : currentRoad.getHits()) {
836 hitsR.push_back(hit);
837 }
838
839 // If outside in, sort in increasing R, otherwise, decreasing R
840 if (m_doOutsideIn)
841 {
842 std::sort(hitsR.begin(), hitsR.end(), [](const auto& a, const auto& b) {
843 if (a->getR() == b->getR()) return a->getLayer() < b->getLayer();
844 return a->getR() < b->getR();
845 });
846 }
847 else
848 {
849 std::sort(hitsR.begin(), hitsR.end(), [](const auto& a, const auto& b) {
850 return a->getR() > b->getR();
851 });
852 }
853
854 for (unsigned long i = 0; i < hitsR.size(); i++)
855 {
856 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());
857 }
858 }
859
860}

◆ setupSlices()

virtual StatusCode FPGATrackSimNNPathfinderExtensionTool::setupSlices ( FPGATrackSimLogicalEventInputHeader * slicedHitHeader)
inlineoverridevirtual

Definition at line 141 of file FPGATrackSimNNPathfinderExtensionTool.h.

141 {
142 m_slicedHitHeader = slicedHitHeader;
143 return StatusCode::SUCCESS;
144 };

Member Data Documentation

◆ m_batchSize

Gaudi::Property<unsigned int> FPGATrackSimNNPathfinderExtensionTool::m_batchSize { this, "batchSize", 1, "Batch size for NN inference (1 = sequential, >1 for true batching on GPU)"}
private

Definition at line 177 of file FPGATrackSimNNPathfinderExtensionTool.h.

177{ this, "batchSize", 1, "Batch size for NN inference (1 = sequential, >1 for true batching on GPU)"};

◆ m_chronoSvc

ServiceHandle<IChronoStatSvc> FPGATrackSimNNPathfinderExtensionTool::m_chronoSvc {this, "ChronoStatSvc", "ChronoStatSvc"}
private

Definition at line 150 of file FPGATrackSimNNPathfinderExtensionTool.h.

150{this, "ChronoStatSvc", "ChronoStatSvc"};

◆ m_debugEvent

bool FPGATrackSimNNPathfinderExtensionTool::m_debugEvent = false
private

Definition at line 190 of file FPGATrackSimNNPathfinderExtensionTool.h.

◆ m_doOutsideIn

Gaudi::Property<bool> FPGATrackSimNNPathfinderExtensionTool::m_doOutsideIn { this, "doOutsideIn", true, "Setup the tool so it's doing outside in extrap"}
private

Definition at line 174 of file FPGATrackSimNNPathfinderExtensionTool.h.

174{ this, "doOutsideIn", true, "Setup the tool so it's doing outside in extrap"};

◆ m_extensionHitNN

OnnxRuntimeBase FPGATrackSimNNPathfinderExtensionTool::m_extensionHitNN
private

Definition at line 197 of file FPGATrackSimNNPathfinderExtensionTool.h.

◆ m_extensionVolNN

OnnxRuntimeBase FPGATrackSimNNPathfinderExtensionTool::m_extensionVolNN
private

Definition at line 196 of file FPGATrackSimNNPathfinderExtensionTool.h.

◆ m_FPGATrackSimMapping

ServiceHandle<IFPGATrackSimMappingSvc> FPGATrackSimNNPathfinderExtensionTool::m_FPGATrackSimMapping {this, "FPGATrackSimMappingSvc", "FPGATrackSimMappingSvc"}
private

Definition at line 148 of file FPGATrackSimNNPathfinderExtensionTool.h.

148{this, "FPGATrackSimMappingSvc", "FPGATrackSimMappingSvc"};

◆ m_lowPtValueForWindowPhiScaling

Gaudi::Property<float> FPGATrackSimNNPathfinderExtensionTool::m_lowPtValueForWindowPhiScaling { this, "lowPtValueWindowPhi", -1, "Value in MeV below which we scale the phi window size"}
private

Definition at line 167 of file FPGATrackSimNNPathfinderExtensionTool.h.

167{ this, "lowPtValueWindowPhi", -1, "Value in MeV below which we scale the phi window size"};

◆ m_lowPtValueForWindowRScaling

Gaudi::Property<float> FPGATrackSimNNPathfinderExtensionTool::m_lowPtValueForWindowRScaling { this, "lowPtValueWindowR", -1, "Value in MeV below which we scale the r window size"}
private

Definition at line 163 of file FPGATrackSimNNPathfinderExtensionTool.h.

163{ this, "lowPtValueWindowR", -1, "Value in MeV below which we scale the r window size"};

◆ m_lowPtValueForWindowZScaling

Gaudi::Property<float> FPGATrackSimNNPathfinderExtensionTool::m_lowPtValueForWindowZScaling { this, "lowPtValueWindowZ", -1, "Value in MeV below which we scale the r window size"}
private

Definition at line 165 of file FPGATrackSimNNPathfinderExtensionTool.h.

165{ this, "lowPtValueWindowZ", -1, "Value in MeV below which we scale the r window size"};

◆ m_lowPtWindowPhiScaling

Gaudi::Property<float> FPGATrackSimNNPathfinderExtensionTool::m_lowPtWindowPhiScaling {this, "lowPtPhiScaling", 1.0, "Scaling factor for low pt in Phi"}
private

Definition at line 168 of file FPGATrackSimNNPathfinderExtensionTool.h.

168{this, "lowPtPhiScaling", 1.0, "Scaling factor for low pt in Phi"};

◆ m_lowPtWindowRScaling

Gaudi::Property<float> FPGATrackSimNNPathfinderExtensionTool::m_lowPtWindowRScaling {this, "lowPtRScaling", 1.0, "Scaling factor for low pt in R"}
private

Definition at line 164 of file FPGATrackSimNNPathfinderExtensionTool.h.

164{this, "lowPtRScaling", 1.0, "Scaling factor for low pt in R"};

◆ m_lowPtWindowZScaling

Gaudi::Property<float> FPGATrackSimNNPathfinderExtensionTool::m_lowPtWindowZScaling {this, "lowPtZScaling", 1.0, "Scaling factor for low pt in Z"}
private

Definition at line 166 of file FPGATrackSimNNPathfinderExtensionTool.h.

166{this, "lowPtZScaling", 1.0, "Scaling factor for low pt in Z"};

◆ m_maxBranches

Gaudi::Property<int> FPGATrackSimNNPathfinderExtensionTool::m_maxBranches { this, "maxBranches", -1, "Max number of branches before we stop, if negative this is disabled"}
private

Definition at line 173 of file FPGATrackSimNNPathfinderExtensionTool.h.

173{ this, "maxBranches", -1, "Max number of branches before we stop, if negative this is disabled"};

◆ m_maxMiss

Gaudi::Property<unsigned> FPGATrackSimNNPathfinderExtensionTool::m_maxMiss { this, "threshold", 2, "Maximum number of missing hits to reject a road"}
private

Definition at line 152 of file FPGATrackSimNNPathfinderExtensionTool.h.

152{ this, "threshold", 2, "Maximum number of missing hits to reject a road"};

◆ m_missedHitPhiScaling

Gaudi::Property<float> FPGATrackSimNNPathfinderExtensionTool::m_missedHitPhiScaling {this, "missedHitPhiScaling", -1, "Amount to scale Phi window if previous hit was missed. Negative means this is disabled"}
private

Definition at line 172 of file FPGATrackSimNNPathfinderExtensionTool.h.

172{this, "missedHitPhiScaling", -1, "Amount to scale Phi window if previous hit was missed. Negative means this is disabled"};

◆ m_missedHitRScaling

Gaudi::Property<float> FPGATrackSimNNPathfinderExtensionTool::m_missedHitRScaling {this, "missedHitRScaling", -1, "Amount to scale R window if previous hit was missed. Negative means this is disabled"}
private

Definition at line 170 of file FPGATrackSimNNPathfinderExtensionTool.h.

170{this, "missedHitRScaling", -1, "Amount to scale R window if previous hit was missed. Negative means this is disabled"};

◆ m_missedHitZScaling

Gaudi::Property<float> FPGATrackSimNNPathfinderExtensionTool::m_missedHitZScaling {this, "missedHitZScaling", -1, "Amount to scale Z window if previous hit was missed. Negative means this is disabled"}
private

Definition at line 171 of file FPGATrackSimNNPathfinderExtensionTool.h.

171{this, "missedHitZScaling", -1, "Amount to scale Z window if previous hit was missed. Negative means this is disabled"};

◆ m_nLayers_1stStage

unsigned FPGATrackSimNNPathfinderExtensionTool::m_nLayers_1stStage = 0
private

Definition at line 181 of file FPGATrackSimNNPathfinderExtensionTool.h.

◆ m_nLayers_2ndStage

unsigned FPGATrackSimNNPathfinderExtensionTool::m_nLayers_2ndStage = 0
private

Definition at line 182 of file FPGATrackSimNNPathfinderExtensionTool.h.

◆ m_predictionWindowLength

Gaudi::Property<int> FPGATrackSimNNPathfinderExtensionTool::m_predictionWindowLength { this, "predictionWindowLength", 3, "Length of hits needed for prediction"}
private

Definition at line 175 of file FPGATrackSimNNPathfinderExtensionTool.h.

175{ this, "predictionWindowLength", 3, "Length of hits needed for prediction"};

◆ m_region

Gaudi::Property<std::string> FPGATrackSimNNPathfinderExtensionTool::m_region { this, "OutputRegion", "", "region ID"}
private

Definition at line 153 of file FPGATrackSimNNPathfinderExtensionTool.h.

153{ this, "OutputRegion", "", "region ID"};

◆ m_roads

std::vector<FPGATrackSimRoad> FPGATrackSimNNPathfinderExtensionTool::m_roads
private

Definition at line 180 of file FPGATrackSimNNPathfinderExtensionTool.h.

◆ m_slicedHitHeader

FPGATrackSimLogicalEventInputHeader* FPGATrackSimNNPathfinderExtensionTool::m_slicedHitHeader = nullptr
private

Definition at line 194 of file FPGATrackSimNNPathfinderExtensionTool.h.

◆ m_tHistSvc

ServiceHandle<ITHistSvc> FPGATrackSimNNPathfinderExtensionTool::m_tHistSvc {this, "THistSvc", "THistSvc"}
private

Definition at line 149 of file FPGATrackSimNNPathfinderExtensionTool.h.

149{this, "THistSvc", "THistSvc"};

◆ m_useCartesian

Gaudi::Property<bool> FPGATrackSimNNPathfinderExtensionTool::m_useCartesian { this, "useCartesian", true, "If true, NNs use Cartestian coordinates. If false,they use cylindrical coordiantes"}
private

Definition at line 176 of file FPGATrackSimNNPathfinderExtensionTool.h.

176{ this, "useCartesian", true, "If true, NNs use Cartestian coordinates. If false,they use cylindrical coordiantes"};

◆ m_windowFineID

Gaudi::Property<std::vector<int> > FPGATrackSimNNPathfinderExtensionTool::m_windowFineID { this, "windowFineID", {0}, "Fine ID indexing for windows"}
private

Definition at line 162 of file FPGATrackSimNNPathfinderExtensionTool.h.

162{ this, "windowFineID", {0}, "Fine ID indexing for windows"};

◆ m_windowPhi

Gaudi::Property<std::vector<float> > FPGATrackSimNNPathfinderExtensionTool::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"}
private

Definition at line 161 of file FPGATrackSimNNPathfinderExtensionTool.h.

161{ this, "windowPhi", {0.4}, "Window Size to search in for phi, either pass one value for all layers or use the number of layers"};

◆ m_windowR

Gaudi::Property<std::vector<float> > FPGATrackSimNNPathfinderExtensionTool::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"}
private

Definition at line 159 of file FPGATrackSimNNPathfinderExtensionTool.h.

159{ this, "windowR", {20.0}, "Window Size to search in for r, either pass one value for all layers or use the number of layers"};

◆ m_windowZ

Gaudi::Property<std::vector<float> > FPGATrackSimNNPathfinderExtensionTool::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"}
private

Definition at line 160 of file FPGATrackSimNNPathfinderExtensionTool.h.

160{ this, "windowZ", {20.0}, "Window Size to search in for z, either pass one value for all layers or use the number of layers"};

The documentation for this class was generated from the following files: