ATLAS Offline Software
Loading...
Searching...
No Matches
FPGATrackSimGenScanTool.cxx
Go to the documentation of this file.
1// Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
2
9
21
22#include <sstream>
23#include <cmath>
24#include <algorithm>
25#include <limits>
26
27
28#include <nlohmann/json.hpp>
29
30#include "TH1.h"
31
32namespace {
33 bool
34 nearZero(const double & v){
35 return std::abs(v)<=std::numeric_limits<double>::min();
36 }
37}
38
39
40
42// Debug Print Tools
43template<class T>
44static inline std::string to_string(const std::vector<T>& v)
45{
46 std::ostringstream oss;
47 oss << "[";
48 if (!v.empty())
49 {
50 std::copy(v.begin(), v.end() - 1, std::ostream_iterator<T>(oss, ", "));
51 oss << v.back();
52 }
53 oss << "]";
54 return oss.str();
55}
56
57
59// AthAlgTool
60
61FPGATrackSimGenScanTool::FPGATrackSimGenScanTool(const std::string& algname, const std::string &name, const IInterface *ifc) :
62 base_class(algname, name, ifc)
63{
64 declareInterface<IFPGATrackSimRoadFinderTool>(this);
65}
66
67
69{
70 // Dump the configuration to make sure it propagated through right
71 const std::vector<Gaudi::Details::PropertyBase*> props = this->getProperties();
72 for( Gaudi::Details::PropertyBase* prop : props ) {
73 if (prop->ownerTypeName()==this->type()) {
74 ATH_MSG_DEBUG("Property:\t" << prop->name() << "\t : \t" << prop->toString());
75 }
76 }
77
78 // Retrieve info
80 ATH_MSG_INFO("Map specifies :" << m_binnedhits->getNLayers());
81 ATH_CHECK(m_binnedhits.retrieve());
82 ATH_CHECK(m_monitoring.retrieve(EnableTool{m_enableMonitoring}));
83 if (m_enableMonitoring) ATH_MSG_INFO("Monitoring Dir :" << m_monitoring->dir());
84
85 // Setup layer configuration if not already set from layerMap
86 if (m_binnedhits->getNLayers()==0){
87 auto nLogicalLayers = m_FPGATrackSimMapping->PlaneMap_1st(getSubRegion())->getNLogiLayers();
88 if (nLogicalLayers == 0){
89 ATH_MSG_ERROR("Number of logical layers is zero in FPGATrackSimGenScanTool::initialize");
90 return StatusCode::FAILURE;
91 }
92 m_binnedhits->setNLayers(nLogicalLayers);
93 }
94 // This is the layers they get paired with previous layers
95 for (unsigned lyr = 0; lyr < m_binnedhits->getNLayers(); ++lyr) m_pairingLayers.push_back(lyr);
96 if (m_reversePairDir) {
98 }
99 ATH_MSG_INFO("Pairing Layers: " << m_pairingLayers);
100
101 // Check inputs
102 bool ok = false;
103 const int signedSize = static_cast<int>(m_binnedhits->getNLayers()) - 1;
104 if (std::ssize(m_pairFilterDeltaPhiCut) != signedSize)
105 ATH_MSG_FATAL("initialize() pairFilterDeltaPhiCut must have size nLayers-1=" << signedSize << " found " << m_pairFilterDeltaPhiCut.size());
106 else if (std::ssize(m_pairFilterDeltaEtaCut) != signedSize)
107 ATH_MSG_FATAL("initialize() pairFilterDeltaEtaCut must have size nLayers-1=" << signedSize << " found " << m_pairFilterDeltaEtaCut.size());
108 else if (m_pairFilterPhiExtrapCut.size() != 2)
109 ATH_MSG_FATAL("initialize() pairFilterPhiExtrapCut must have size 2 found " << m_pairFilterPhiExtrapCut.size());
110 else if (m_pairFilterEtaExtrapCut.size() != 2)
111 ATH_MSG_FATAL("initialize() pairFilterEtaExtrapCut must have size 2 found " << m_pairFilterEtaExtrapCut.size());
112 else if (m_pairSetPhiExtrapCurvedCut.size() != 2)
113 ATH_MSG_FATAL("initialize() PairSetPhiExtrapCurvedCut must have size 2found " << m_pairSetPhiExtrapCurvedCut.size());
114 else if ((m_rin < 0.0) || (m_rout < 0.0))
115 ATH_MSG_FATAL("Radii not set");
116 else
117 ok = true;
118 if (!ok)
119 return StatusCode::FAILURE;
120
121
122 // register histograms
123 if (m_enableMonitoring) ATH_CHECK(m_monitoring->registerHistograms(m_binnedhits.get()));
124 // write out the firmware LUTs
125 m_binnedhits->getBinTool().writeLUTs();
126
127 return StatusCode::SUCCESS;
128}
129
131// Main Algorithm
132
133StatusCode FPGATrackSimGenScanTool::getRoads(const std::vector<std::shared_ptr<const FPGATrackSimHit>> &hits,
134 std::vector<FPGATrackSimRoad> &roads)
135{
136 ATH_MSG_DEBUG("In getRoads, Processing Event# " << ++m_evtsProcessed << " hit size = " << hits.size());
137
138
139 roads.clear();
140 m_roads.clear();
141 if (m_enableMonitoring) {
142 m_monitoring->resetDataFlowCounters();
143 // Currently assume that if less than 100 hits its a single track MC
144 m_monitoring->parseTruthInfo(getTruthTracks(),(hits.size() < 100));
145 }
146
147 // do the binning...
148 ATH_CHECK(m_binnedhits->fill(hits));
149 if (m_enableMonitoring) m_monitoring->fillBinningSummary(hits);
150
151 // scan over image building pairs for bins over threshold
153 {
154 // apply threshold
155 if (bin.data().lyrCnt() < m_threshold) continue;
156 ATH_MSG_DEBUG("Bin passes threshold " << bin.data().lyrCnt() << " "
157 << bin.idx());
158
159 // Monitor contents of bins passing threshold
160 if (m_enableMonitoring) m_monitoring->fillBinLevelOutput(bin.idx(), bin.data());
161 if (m_binningOnly) continue;
162
163 // pass hits for bin to filterRoad and get back pairs of hits grouped into pairsets
164 std::vector<HitPairSet> pairsets;
165 if (m_binFilter=="PairThenGroup") {
166 ATH_CHECK(pairThenGroupFilter(bin.data(), pairsets));
167 } else if (m_binFilter=="IncrementalBuild") {
168 ATH_CHECK(incrementalBuildFilter(bin.data(), pairsets));
169 } else {
170 ATH_MSG_FATAL("Unknown bin filter" << m_binFilter);
171 }
172 ATH_MSG_DEBUG("grouped PairSets " << pairsets.size());
173
174 // convert the group pairsets to FPGATrackSimRoads
175 for (const FPGATrackSimGenScanTool::HitPairSet& pairset: pairsets)
176 {
177 addRoad(pairset.hitlist, bin.idx());
178 ATH_MSG_DEBUG("Output road size=" <<pairset.hitlist.size());
179
180 // debug statement if more than one road found in bin
181 // not necessarily a problem
182 if (pairsets.size() >1) {
183 std::string s = "";
184 for (const FPGATrackSimBinUtil::StoredHit* const hit : pairset.hitlist)
185 {
186 s += "(" + std::to_string(hit->layer) + "," + std::to_string(hit->hitptr->getR()) + "), ";
187 }
188 ATH_MSG_DEBUG("Duplicate Group " << s);
189 }
190 }
191 }
192
193 // copy roads to output vector
194 roads.reserve(m_roads.size());
195
196 for (auto & r : m_roads) roads.push_back(std::move(r));
197 ATH_MSG_DEBUG("Roads = " << roads.size());
198
199 // clear previous event
200 // (reusing saves having to reallocate memory for each event)
201 m_binnedhits->resetBins();
202
204 return StatusCode::SUCCESS;
205}
206
207// Filter the bins above threshold into pairsets which output roads
209 std::vector<HitPairSet> &output_pairsets)
210{
211 ATH_MSG_VERBOSE("In pairThenGroupFilter");
212
213 // Organize Hits by Layer
214 std::vector<std::vector<const StoredHit *>> hitsByLayer(m_binnedhits->getNLayers());
215 ATH_CHECK(sortHitsByLayer(bindata, hitsByLayer));
216
217 // This is monitoring for each bin over threshold
218 // It's here so it can get the hitsByLayer
219 if (m_enableMonitoring) m_monitoring->fillHitsByLayer(hitsByLayer);
220
221 // Make Pairs
222 HitPairSet pairs;
223 ATH_CHECK(makePairs(hitsByLayer, pairs));
224
225 // Filter Pairs
226 HitPairSet filteredpairs;
227 ATH_CHECK(filterPairs(pairs, filteredpairs));
228
229 // Require road is still over threshold
230 bool passedPairFilter = (filteredpairs.lyrCnt() >= m_threshold);
231 if (m_enableMonitoring) m_monitoring->pairFilterCheck(pairs, filteredpairs, passedPairFilter);
232
233 // if passed Pair Filter proceed to group the filtered pairs into pairsets
234 if (passedPairFilter)
235 {
236 // Pair Set Grouping
237 std::vector<HitPairSet> pairsets;
238 ATH_CHECK(groupPairs(filteredpairs, pairsets, false));
239
240 // loop over pairsets and find those that are over thresold
241 for (const FPGATrackSimGenScanTool::HitPairSet& pairset : pairsets)
242 {
243 // if over threshold add it to the output
244 if (pairset.lyrCnt() >= m_threshold)
245 {
247 output_pairsets.push_back(pairset);
248 }
249 }
250 }
251 }
252
253 // set outputs if not all filters applied
254 if (!m_applyPairFilter) {
255 // output is just the filtered pairs
256 output_pairsets.push_back(std::move(pairs));
257 }
258 else if (passedPairFilter && !m_applyPairSetFilter)
259 {
260 // output is just the filtered pairs
261 output_pairsets.push_back(std::move(filteredpairs));
262 }
263
264 return StatusCode::SUCCESS;
265}
266
267
269 IntermediateState &outputstate,
270 unsigned lyridx,
271 const std::vector<const StoredHit *>& newhits)
272{
273 unsigned int allowed_missed_hits = m_binnedhits->getNLayers() - m_threshold;
274
275 std::vector<bool> pairset_used(inputstate.pairsets.size(),false);
276
277 for (auto &newhit : newhits) {
278
279 // don't make new pairs with hits that the new hit is already paired with in a group
280 std::set<const StoredHit *> vetoList;
281
282 // try adding hit to existing pair sets
283 for (unsigned ps_idx = 0; ps_idx < inputstate.pairsets.size(); ++ps_idx) {
284 auto &pairset = inputstate.pairsets[ps_idx];
285 HitPair nextpair(pairset.lastpair().second, newhit, m_reversePairDir);
286 if (pairMatchesPairSet(pairset, nextpair, false) || (m_applyPairSetFilter==false)) {
287 HitPairSet newset(pairset);
288 newset.addPair(nextpair);
289 pairset_used[ps_idx]=true;
290 outputstate.pairsets.push_back(std::move(newset));
291 // put inpair hits in list of hits not to pair again with the new hits
292 for (auto vetohit : pairset.hitlist) {
293 vetoList.insert(vetohit);
294 }
295 }
296 }
297
298 // make new pairsets with unpaired hits
299 for (auto prevhit : inputstate.unpairedHits) {
300 if (vetoList.count(prevhit) == 0) {
301 HitPair newpair(prevhit, newhit, m_reversePairDir);
302 if (pairPassesFilter(newpair) || (m_applyPairFilter == false)) {
303 HitPairSet newset;
304 newset.addPair(newpair);
305 outputstate.pairsets.push_back(std::move(newset));
306 }
307 }
308 }
309
310 // if this can be the start of a the start of a track and still
311 // have enough hits to make a track, add it to the unpaired hits list
312 if (lyridx <= allowed_missed_hits) {
313 outputstate.unpairedHits.push_back(newhit);
314 }
315 }
316
317 // Add groups to output without new hit if they have enough hits to skip
318 // this layer. Note expected hits at this point is lyridx+1, since we start
319 // counting lyridx from zero. Logic is then keep the pairset if
320 // expected hits <= hits in set + allows misses
321 for (unsigned ps_idx = 0; ps_idx < inputstate.pairsets.size(); ++ps_idx) {
322 auto &pairset = inputstate.pairsets[ps_idx];
323 if ((!pairset_used[ps_idx])&&(lyridx < (pairset.hitlist.size() + allowed_missed_hits))) {
324 outputstate.pairsets.push_back(pairset);
325 }
326 }
327
328 // Add hits to unpaired list if hits are still allowed to start a track
329 // ---------------------------------------------------------------------
330 // If you start a new track with a pair whose second element is
331 // layer N (layer numbering starting from zero), then you'll have missed N-1
332 // layers Minus 1 because the first hit was one of the N layers already
333 // passed.
334 //
335 // The next pass will be layer N=lyridx+1 where lyridx is the current value
336 // at this point in the code. You will have then missed lyridx layers, so...
337 // E.g. if you allow one missed layer then this stops putting hits in the
338 // unpairedHits list if lyridx>1, so tracks must start with layer 0 or 1,
339 // which makes sense
340 if (lyridx > allowed_missed_hits) {
341 // make no new track starts
342 outputstate.unpairedHits.clear();
343 } else {
344 // copy in any previous unpairedHits as well
345 for (auto prevhit : inputstate.unpairedHits) {
346 outputstate.unpairedHits.push_back(prevhit);
347 }
348 }
349}
350
351// Filter the bins above threshold into pairsets which output roads
353 std::vector<HitPairSet> &output_pairsets)
354{
355 ATH_MSG_VERBOSE("In buildGroupsWithPairs");
356
357 // Organize Hits by Layer
358 std::vector<std::vector<const StoredHit *>> hitsByLayer(m_binnedhits->getNLayers());
359 ATH_CHECK(sortHitsByLayer(bindata, hitsByLayer));
360
361 // This is monitoring for each bin over threshold
362 // It's here so it can get the hitsByLayer
363 if (m_enableMonitoring) m_monitoring->fillHitsByLayer(hitsByLayer);
364
365 std::vector<IntermediateState> states{m_binnedhits->getNLayers()+1};
366 for (unsigned lyridx = 0; lyridx < m_binnedhits->getNLayers(); lyridx++) {
367 unsigned lyr = m_pairingLayers[lyridx];
368 updateState(states[lyridx] , states[lyridx+1], lyridx, hitsByLayer[lyr]);
369 }
370
371 // this is a little ugly because it requires copying the output pairsets right now
372 output_pairsets=states[m_binnedhits->getNLayers()].pairsets;
373
374 if (m_enableMonitoring) m_monitoring->fillBuildGroupsWithPairs(states,m_binnedhits->getNLayers()-m_threshold);
375
376 return StatusCode::SUCCESS;
377}
378
379
380// 1st step of filter: sort hits by layer
382 std::vector<std::vector<const StoredHit *>>& hitsByLayer)
383{
384 ATH_MSG_DEBUG("In fillHitsByLayer");
385
386 for (const FPGATrackSimBinUtil::StoredHit &hit : bindata.hits) {
387 hitsByLayer.at(hit.layer).push_back(&hit);
388 }
389
390 return StatusCode::SUCCESS;
391}
392
393
394
395
396
397// 2nd step of filter: make pairs of hits from adjacent and next-to-adjacent layers
398StatusCode FPGATrackSimGenScanTool::makePairs(const std::vector<std::vector<const StoredHit *>>& hitsByLayer,
399 HitPairSet &pairs)
400{
401 ATH_MSG_VERBOSE("In makePairs");
402
403 std::vector<const StoredHit *> const * lastlyr = 0;
404 std::vector<const StoredHit *> const * lastlastlyr = 0;
405
406 // order here is designed so lower radius hits come first
407 for (unsigned lyr : m_pairingLayers) {
408 for (const StoredHit *const &ptr1 :
409 hitsByLayer[lyr]) {
410 if (lastlyr) {
411 for (const StoredHit *const &ptr2 : *lastlyr) {
412 pairs.addPair(HitPair(ptr2, ptr1,m_reversePairDir));
413 }
414 // Add Pairs that skip one layer
415 if (lastlastlyr) {
416 for (const StoredHit *const &ptr2 : *lastlastlyr) {
417 pairs.addPair(HitPair(ptr2, ptr1,m_reversePairDir));
418 }
419 }
420 }
421 }
422 lastlastlyr = lastlyr;
423 lastlyr = &hitsByLayer[lyr];
424 if (m_enableMonitoring) m_monitoring->fillPairingHits(lastlyr,lastlastlyr);
425 }
426
427 return StatusCode::SUCCESS;
428}
429
430// 3rd step of filter: make cuts on the pairs to ensure that are consisten with
431// the bin they are in
432
434 if (m_enableMonitoring) m_monitoring->fillPairFilterCuts(pair,m_rin,m_rout);
435 int lyr = std::min(pair.first->layer,pair.second->layer);
436 return (std::abs(pair.dPhi()) < m_pairFilterDeltaPhiCut[lyr]) &&
437 (std::abs(pair.dEta()) < m_pairFilterDeltaEtaCut[lyr]) &&
438 (std::abs(pair.PhiInExtrap(m_rin)) < m_pairFilterPhiExtrapCut[0]) &&
439 (std::abs(pair.PhiOutExtrap(m_rout)) < m_pairFilterPhiExtrapCut[1]) &&
440 (std::abs(pair.EtaInExtrap(m_rin)) < m_pairFilterEtaExtrapCut[0]) &&
441 (std::abs(pair.EtaOutExtrap(m_rout)) < m_pairFilterEtaExtrapCut[1]);
442}
443
445{
446 ATH_MSG_VERBOSE("In filterPairs");
447
448 for (const FPGATrackSimGenScanTool::HitPair &pair : pairs.pairList) {
449 if (pairPassesFilter(pair)) {
450 filteredpairs.addPair(pair);
451 }
452 }
453 return StatusCode::SUCCESS;
454}
455
456
457// 4th step of filter: group pairs into sets where they are all consistent with being from the same track
459 std::vector<HitPairSet> &pairsets,
460 bool verbose)
461{
462 ATH_MSG_VERBOSE("In groupPairs");
463 for (const FPGATrackSimGenScanTool::HitPair & pair : filteredpairs.pairList)
464 {
465 bool added = false;
466 for (FPGATrackSimGenScanTool::HitPairSet &pairset : pairsets)
467 {
468 // Only add skip pairs if skipped layer is not already hit
469 if ((std::abs(int(pair.second->layer) - int(pair.first->layer)) > 1) // gives if is a skip pair
470 && (pairset.hasLayer(std::min(pair.first->layer,pair.second->layer) + 1))) // gives true if skipped layer already in set
471 {
472 // if it matches mark as added so it doesn't start a new pairset
473 // false here is so it doesn't plot these either
474 if (pairMatchesPairSet(pairset, pair, verbose))
475 added = true;
476 if (!added) {
477 ATH_MSG_VERBOSE("Skip pair does not match non-skip pair");
478 }
479 }
480 else
481 {
482 if (pairMatchesPairSet(pairset, pair, verbose))
483 {
484 int size = pairset.addPair(pair);
485 if (verbose)
486 ATH_MSG_VERBOSE("addPair " << pairsets.size() << " " << pairset.pairList.size() << " " << size);
487 added = true;
488 }
489 }
490
491 }
492
493 if (!added)
494 {
495 HitPairSet newpairset;
496 newpairset.addPair(pair);
497 pairsets.push_back(std::move(newpairset));
498 }
499 }
500
501 return StatusCode::SUCCESS;
502
503}
504
505// used to determine if a pair is consistend with a pairset
507 const HitPair &pair,
508 bool verbose) {
509 // Compute all the cut values needed for filtering
510 bool matchPhiPassed = (std::abs(pairset.MatchPhi(pair)) < m_pairSetMatchPhiCut);
511 bool matchEtaPassed = (std::abs(pairset.MatchEta(pair)) < m_pairSetMatchEtaCut);
512 bool deltaDeltaPhiPassed = (std::abs(pairset.DeltaDeltaPhi(pair)) < m_pairSetDeltaDeltaPhiCut);
513 bool deltaDeltaEtaPassed = (std::abs(pairset.DeltaDeltaEta(pair)) < m_pairSetDeltaDeltaEtaCut);
514 bool phiCurvaturePassed = (std::abs(pairset.PhiCurvature(pair)) < m_pairSetPhiCurvatureCut);
515 bool etaCurvaturePassed = (std::abs(pairset.EtaCurvature(pair)) < m_pairSetEtaCurvatureCut);
516 bool phiInExtrapPassed = (std::abs(pairset.PhiInExtrapCurved(pair, m_rin)) < m_pairSetPhiExtrapCurvedCut[0]);
517 bool phiOutExtrapPassed = (std::abs(pairset.PhiOutExtrapCurved(pair, m_rout)) < m_pairSetPhiExtrapCurvedCut[1]);
518
519 std::vector<bool> allcutsPassed;
520 allcutsPassed.push_back(matchPhiPassed);
521 allcutsPassed.push_back(matchEtaPassed);
522 allcutsPassed.push_back(deltaDeltaPhiPassed);
523 allcutsPassed.push_back(deltaDeltaEtaPassed);
524 allcutsPassed.push_back(phiCurvaturePassed);
525 allcutsPassed.push_back(etaCurvaturePassed);
526
527 bool deltaPhiCurvaturePassed = true;
528 bool deltaEtaCurvaturePassed = true;
529 if (pairset.pairList.size() > 1) {
530 deltaPhiCurvaturePassed = (std::abs(pairset.DeltaPhiCurvature(pair)) < m_pairSetDeltaPhiCurvatureCut);
531 deltaEtaCurvaturePassed = (std::abs(pairset.DeltaEtaCurvature(pair)) < m_pairSetDeltaEtaCurvatureCut);
532 allcutsPassed.push_back(deltaPhiCurvaturePassed);
533 allcutsPassed.push_back(deltaEtaCurvaturePassed);
534 }
535 allcutsPassed.push_back(phiInExtrapPassed);
536 allcutsPassed.push_back(phiOutExtrapPassed);
537
538 // count number of cuts passed
539 unsigned passedCuts = std::count_if(allcutsPassed.begin(), allcutsPassed.end(),
540 [](bool cut) { return cut; });
541 bool passedAll = (passedCuts == allcutsPassed.size());
542
543 // Detailed monitoring with histograms (only if enabled)
544 if (m_enableMonitoring) {
545 // In order to make it easy to have a long list of possible cuts,
546 // a vector of cutvar structs is used to represent each cut
547 // then apply the AND of all the cuts is done with a std::count_if function
548
549 // define the struct (effectively mapping because variable, configured cut value,
550 // and histograms for plotting
551 struct cutvar {
552 cutvar(std::string name, double val, double cut, std::vector<TH1D *>& histset) :
553 m_name(std::move(name)), m_val(val), m_cut(cut), m_histset(histset) {}
554 bool passed() { return std::abs(m_val) < m_cut; }
555 void fill(unsigned cat) { m_histset[cat]->Fill(m_val); }
556 std::string m_name;
557 double m_val;
558 double m_cut;
559 std::vector<TH1D *> &m_histset;
560 };
561
562 // add the cuts to the list of all cuts
563 std::vector<cutvar> allcuts;
564 allcuts.push_back(cutvar("MatchPhi", pairset.MatchPhi(pair),
566 m_monitoring->m_pairSetMatchPhi));
567 allcuts.push_back(cutvar("MatchEta", pairset.MatchEta(pair),
569 m_monitoring->m_pairSetMatchEta));
570 allcuts.push_back(cutvar("DeltaDeltaPhi", pairset.DeltaDeltaPhi(pair),
572 m_monitoring->m_deltaDeltaPhi));
573 allcuts.push_back(cutvar("DeltaDeltaEta", pairset.DeltaDeltaEta(pair),
575 m_monitoring->m_deltaDeltaEta));
576 allcuts.push_back(cutvar("PhiCurvature", pairset.PhiCurvature(pair),
578 m_monitoring->m_phiCurvature));
579 allcuts.push_back(cutvar("EtaCurvature", pairset.EtaCurvature(pair),
581 m_monitoring->m_etaCurvature));
582 if (pairset.pairList.size() > 1) {
583 allcuts.push_back(cutvar(
584 "DeltaPhiCurvature", pairset.DeltaPhiCurvature(pair),
585 m_pairSetDeltaPhiCurvatureCut, m_monitoring->m_deltaPhiCurvature));
586 allcuts.push_back(cutvar(
587 "DeltaEtaCurvature", pairset.DeltaEtaCurvature(pair),
588 m_pairSetDeltaEtaCurvatureCut, m_monitoring->m_deltaEtaCurvature));
589 }
590 allcuts.push_back(cutvar(
591 "PhiInExtrapCurved", pairset.PhiInExtrapCurved(pair, m_rin),
592 m_pairSetPhiExtrapCurvedCut[0], m_monitoring->m_phiInExtrapCurved));
593 allcuts.push_back(cutvar(
594 "PhiOutExtrapCurved", pairset.PhiOutExtrapCurved(pair, m_rout),
595 m_pairSetPhiExtrapCurvedCut[1], m_monitoring->m_phiOutExtrapCurved));
596
597 // monitoring
598 unsigned monitoringPassedCuts = std::count_if(allcuts.begin(), allcuts.end(),
599 [](cutvar& cut) { return cut.passed(); });
600 bool monitoringPassedAll = (monitoringPassedCuts == allcuts.size());
601 bool monitoringPassedAllButOne = (monitoringPassedCuts == allcuts.size() - 1);
602 for (cutvar& cut: allcuts) {
603 // the last value computes if an n-1 histogram should be filled
604 m_monitoring->fillPairSetFilterCut(cut.m_histset, cut.m_val, pair,
605 pairset.lastpair(),
606 (monitoringPassedAll || (monitoringPassedAllButOne && !cut.passed())));
607 }
608
609 if (verbose)
610 {
611 std::string s = "";
612 for (cutvar &cut : allcuts)
613 {
614 s += cut.m_name + " : (" + cut.passed() + ", " + cut.m_val + "), ";
615 }
616 ATH_MSG_DEBUG("PairSet test " << monitoringPassedAll << " " << s);
617 ATH_MSG_DEBUG("Hits: \n " << *pairset.lastpair().first << "\n "
618 << *pairset.lastpair().second << "\n "
619 << *pair.first << "\n "
620 << *pair.second);
621 }
622 }
623
624 return passedAll;
625}
626
627// format final pairsets into expected output of getRoads
628void FPGATrackSimGenScanTool::addRoad(std::vector<const StoredHit *> const &hits, const FPGATrackSimBinUtil::IdxSet &idx)
629{
630 layer_bitmask_t hitLayers = 0;
631 std::vector<std::vector<std::shared_ptr<const FPGATrackSimHit>>>
632 sorted_hits(m_binnedhits->getNLayers(),std::vector<std::shared_ptr<const FPGATrackSimHit>>());
633 for (const FPGATrackSimBinUtil::StoredHit* hit : hits)
634 {
635 hitLayers |= 1 << hit->layer;
636 sorted_hits[hit->layer].push_back(hit->hitptr);
637 }
638
639 // "Fit" the track.
640 FPGATrackSimTrackPars fittedpars;
641 double chi2,chi2_phi,chi2_eta;
642 bool inBin = fitRoad(hits, idx, fittedpars, chi2, chi2_phi,chi2_eta);
643 if (!inBin && m_inBinFiltering) return;
644
645 m_roads.emplace_back();
646 FPGATrackSimRoad &r = m_roads.back();
647
648 r.setRoadID(static_cast<int>(m_roads.size()) - 1);
649 // r.setPID(y * m_imageSize_y + x);
650 r.setHits(std::move(sorted_hits));
651
652 r.setBinIdx(idx);
653
654 FPGATrackSimBinUtil::ParSet binCenterPars = m_binnedhits->getBinTool().lastStep()->binCenter(idx);
655 FPGATrackSimTrackPars trackpars = m_binnedhits->getBinTool().binDesc()->parSetToTrackPars(binCenterPars);
656 r.setX(trackpars[FPGATrackSimTrackPars::IPHI]);
657 r.setY(trackpars[FPGATrackSimTrackPars::IHIP]);
658 r.setXBin(idx[3]);
659 r.setYBin(idx[4]);
660 r.setHitLayers(hitLayers);
661 r.setSubRegion(0);
662
663 // Store the fitted information on the track.
664 r.setFitParams(fittedpars);
665 r.setFitChi2(chi2);
666 r.setFitChi2_2d(chi2_phi,chi2_eta);
667}
668
669
671// HitPair and HitPairSet Storage
672
674{for (const FPGATrackSimBinUtil::StoredHit* hit : hitlist)
675 {
676 if (hit == newhit)
677 return true;
678 }
679 return false;
680}
681
683{
684 if (pairList.size() > 0)
685 {
686 // these need to be done before the pair is added
689 }
690
691 pairList.push_back(pair);
692
693 hitLayers |= (0x1 << pair.first->layer);
694 hitLayers |= (0x1 << pair.second->layer);
695 if (!hasHit(pair.first))
696 hitlist.push_back(pair.first);
697 if (!hasHit(pair.second))
698 hitlist.push_back(pair.second);
699 return this->pairList.size();
700}
701
702
704{
705 if ((lastpair().first->hitptr==pair.first->hitptr)||
706 (lastpair().first->hitptr==pair.second->hitptr)||
707 (lastpair().second->hitptr==pair.first->hitptr)||
708 (lastpair().second->hitptr==pair.second->hitptr))
709 return 0.0;
710
711 double dr = (lastpair().second->hitptr->getR() - pair.first->hitptr->getR()) / 2.0;
712 double lastpairextrap = lastpair().second->phiShift - lastpair().dPhi() / lastpair().dR() * dr;
713 double newpairextrap = pair.first->phiShift + pair.dPhi() / pair.dR() * dr;
714 return lastpairextrap - newpairextrap;
715}
716
718{
719 if ((lastpair().first->hitptr==pair.first->hitptr)||
720 (lastpair().first->hitptr==pair.second->hitptr)||
721 (lastpair().second->hitptr==pair.first->hitptr)||
722 (lastpair().second->hitptr==pair.second->hitptr))
723 return 0.0;
724
725 double dr = (lastpair().second->hitptr->getR() - pair.first->hitptr->getR()) / 2.0;
726 double lastpairextrap = lastpair().second->etaShift - lastpair().dEta() / lastpair().dR() * dr;
727 double newpairextrap = pair.first->etaShift + pair.dEta() / pair.dR() * dr;
728 return lastpairextrap - newpairextrap;
729}
730
731double
733 return (pair.dPhi() * lastpair().dR() - lastpair().dPhi() * pair.dR()) / (lastpair().dR() * pair.dR());
734}
735
737{
738 return (pair.dEta() * lastpair().dR() - lastpair().dEta() * pair.dR()) / (lastpair().dR() * pair.dR());
739}
740
742{
743 return 2 * DeltaDeltaPhi(pair) / (lastpair().dR() + pair.dR());
744}
746{
747 return 2 * DeltaDeltaEta(pair) / (lastpair().dR() + pair.dR());
748}
758{
759 double r = std::min(lastpair().first->hitptr->getR(),lastpair().second->hitptr->getR());
760 return lastpair().PhiInExtrap(r_in) + 0.5 * PhiCurvature(pair) * (r_in - r) * (r_in - r);
761}
763{
764 double r = std::max(lastpair().first->hitptr->getR(),lastpair().second->hitptr->getR());
765 return pair.PhiOutExtrap(r_out) + 0.5 * PhiCurvature(pair) * (r_out - r) * (r_out - r);
766}
767
768
769// format final pairsets into expected output of getRoads
770bool FPGATrackSimGenScanTool::fitRoad(std::vector<const StoredHit *> const &hits, const FPGATrackSimBinUtil::IdxSet &idx, FPGATrackSimTrackPars& trackpars, double& chi2, double& phi_chi2, double& eta_chi2) const
771{
772
773 double N =hits.size();
774 double sum_Phi = 0;
775 double sum_Phi2 = 0;
776 double sum_PhiR = 0;
777 double sum_PhiR2 = 0;
778 double sum_Eta = 0;
779 double sum_Eta2 = 0;
780 double sum_EtaR = 0;
781 double sum_R = 0;
782 double sum_R2 = 0;
783 double sum_R3 = 0;
784 double sum_R4 = 0;
785
786 for (const FPGATrackSimBinUtil::StoredHit* hit : hits)
787 {
788 // these are just relevant sums of variables (moments) needed for the chi2 calculation
789 // Calculate r^2, r^3, and r^4, we'll sum these up later below
790 double r = hit->hitptr->getR();
791 double r2 = r*r;
792 double r3 = r2*r;
793 double r4 = r3*r;
794 double dphi = hit->phiShift;
795 double dphi2 = dphi*dphi;
796 double deta = hit->etaShift;
797 double deta2 = deta*deta;
798
799 sum_Phi += dphi;
800 sum_Phi2 += dphi2;
801 sum_PhiR += dphi*r;
802 sum_PhiR2 += dphi*r2;
803 sum_Eta += deta;
804 sum_Eta2 += deta2;
805 sum_EtaR += deta*r;
806 sum_R += r;
807 sum_R2 += r2;
808 sum_R3 += r3;
809 sum_R4 += r4;
810 }
811
812 // phi var calculation
813 // phi(r) = phivars[0] + phivars[1]*r + phivars[2]*r^2
814 // math below is calculated by analytically minimizing the chi2
815 // and solving for the phivars.
816 // the terms below which recur in the phivar expression are just organized
817 // by the power of r (i.e. the dimension), but otherwise have no deep meaning
818 double r6_t0 = (-sum_R2 * sum_R4 + sum_R3*sum_R3);
819 double r5_t0 = (sum_R*sum_R4 - sum_R2*sum_R3);
820 double r4_t0 = (-sum_R * sum_R3 + sum_R2*sum_R2);
821 double r4_t1 = (-N*sum_R4 + sum_R2*sum_R2);
822 double r3_t0 = (N*sum_R3 - sum_R*sum_R2);
823 double r2_t0 = (-N*sum_R2 + sum_R*sum_R);
824
825 // all three phi var expresions use the same demoninator
826 const double denom_phi = N * r6_t0 + sum_R * r5_t0 + sum_R2 * r4_t0;
827 if (nearZero(denom_phi)){
828 ATH_MSG_ERROR("Divide by zero (phi) trapped in FPGATrackSimGenScanTool::fitRoad");
829 return false;
830 }
831
832 // phivar expresions from analytic chi2 minimization
833 std::vector<double> phivars(3);
834 phivars[0] = (sum_Phi*r6_t0 + sum_PhiR*r5_t0 + sum_PhiR2*r4_t0)/denom_phi;
835 phivars[1] = (sum_Phi*r5_t0 + sum_PhiR*r4_t1 + sum_PhiR2*r3_t0)/denom_phi;
836 phivars[2] = (sum_Phi*r4_t0 + sum_PhiR*r3_t0 + sum_PhiR2*r2_t0)/denom_phi;
837
838 // eta vars
839 // same as phi but with not curvature (r^2) term
840 const double denom_eta = N*sum_R2 - sum_R*sum_R;
841 if (nearZero(denom_eta)){
842 ATH_MSG_ERROR("Divide by zero (eta) trapped in FPGATrackSimGenScanTool::fitRoad");
843 return false;
844 }
845
846 std::vector<double> etavars(2);
847 etavars[0] = (-sum_R*sum_EtaR + sum_R2*sum_Eta)/denom_eta;
848 etavars[1] = (N*sum_EtaR - sum_R*sum_Eta)/denom_eta;
849
850 // bin center
851 FPGATrackSimBinUtil::ParSet parset = m_binnedhits->getBinTool().lastStep()->binCenter(idx);
852 double parshift[FPGATrackSimTrackPars::NPARS];
853 parshift[0] = etavars[0] + etavars[1]*m_rin; // z at r in
854 parshift[1]= etavars[0] + etavars[1]*m_rout; // z at r out
855 parshift[2]= -(phivars[0]/m_rin + phivars[1] + phivars[2]*m_rin) ; // phi at r in
856 parshift[3]= -(phivars[0]/m_rout + phivars[1] + phivars[2]*m_rout) ; // phi at r out
857 double y = (m_rout-m_rin); // midpoint between rin and rout from rin
858 parshift[4]= -phivars[2]/4.0*y*y ; // xm
859
860 bool inBin = true;
861 const auto& lastStep = m_binnedhits->getBinTool().lastStep();
862 for (int par = 0; par < FPGATrackSimTrackPars::NPARS; par++ ){
863 parset[par]+=parshift[par];
864 inBin = inBin && (std::abs(parshift[par]) < lastStep->binWidth(par)/2.0);
865 }
866
867 double ec0 = etavars[0];
868 double ec1 = etavars[1];
869 eta_chi2 = sum_Eta2 - 2*ec0*sum_Eta - 2*ec1*sum_EtaR + N*ec0*ec0 + 2*ec0*ec1*sum_R + ec1*ec1*sum_R2;
870
871 double pc0 = phivars[0];
872 double pc1 = phivars[1];
873 double pc2 = phivars[2];
874
875 phi_chi2 = sum_Phi2 - 2*pc0*sum_Phi - 2*pc1*sum_PhiR - 2*pc2*sum_PhiR2 + N*pc0*pc0 + 2*pc0*pc1*sum_R + 2*pc0*pc2*sum_R2 + 2*pc1*pc2*sum_R3 + pc1*pc1*sum_R2 + pc2*pc2*sum_R4;
876
877 for (const FPGATrackSimBinUtil::StoredHit* hit : hits)
878 {
879 double r = hit->hitptr->getR();
880 ATH_MSG_VERBOSE("Fitted r= " << r << " phishift " << hit->phiShift << " =?= " << pc0+pc1*r+pc2*r*r << " etashift " << hit->etaShift << " =?= " << ec0+ec1*r);
881 }
882
883 ATH_MSG_DEBUG("Bin Info parset " << idx << " " << m_binnedhits->getBinTool().lastStep()->binCenter(idx));
884 ATH_MSG_DEBUG("Fitted parset inBin="<< inBin << " nhits=" << hits.size() << " " << parset << " chi2 " << eta_chi2 << "," << phi_chi2);
885
886 // Return by reference the "fitted" track parameters. shift q/pt dimension (GeV/MeV)
887 trackpars = m_binnedhits->getBinTool().binDesc()->parSetToTrackPars(parset);
888 trackpars[FPGATrackSimTrackPars::IHIP] = trackpars[FPGATrackSimTrackPars::IHIP] / 1000;
889 ATH_MSG_VERBOSE("Fitted track pars" << trackpars);
890
891 // and the summed chi2, which is assuming (right now) an even weighting between the two components.
892 // assume only options are 4 or 5 hits for now
893 chi2 = (hits.size() == 5) ?
894 m_etaWeight_5hits * eta_chi2 + m_phiWeight_5hits * phi_chi2 :
895 m_etaWeight_4hits * eta_chi2 + m_phiWeight_4hits * phi_chi2 ;
896
897 return inBin;
898}
899
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_DEBUG(x)
Binning Utilities for GenScanTool.
This is the monitoring for the FPGATrackSimGenScanTool.
static std::string to_string(const std::vector< T > &v)
Implements a generalized 5d scan which filters pair of hits in bins linearized about nominal trajecto...
: FPGATrackSim-specific class to represent an hit in the detector.
Maps physical layers to logical layers.
Maps ITK module indices to FPGATrackSim regions.
Structs that store the 5 track parameters.
uint32_t layer_bitmask_t
bool passed(DecisionID id, const DecisionIDContainer &)
checks if required decision ID is in the set of IDs in the container
#define y
Gaudi::Property< double > m_pairSetMatchPhiCut
bool pairMatchesPairSet(const HitPairSet &pairset, const HitPair &pair, bool verbose)
bool fitRoad(std::vector< const StoredHit * > const &hits, const FPGATrackSimBinUtil::IdxSet &idx, FPGATrackSimTrackPars &trackpars, double &chi2, double &chi2_phi, double &chi2_eta) const
StatusCode groupPairs(HitPairSet &filteredpairs, std::vector< HitPairSet > &clusters, bool verbose)
Gaudi::Property< double > m_etaWeight_4hits
FPGATrackSimGenScanTool(const std::string &algname, const std::string &name, const IInterface *ifc)
Gaudi::Property< std::vector< double > > m_pairFilterDeltaEtaCut
virtual StatusCode initialize() override
bool pairPassesFilter(const HitPair &pair)
std::vector< FPGATrackSimRoad > m_roads
virtual int getSubRegion() const override
Gaudi::Property< double > m_pairSetDeltaDeltaPhiCut
Gaudi::Property< double > m_pairSetDeltaEtaCurvatureCut
ServiceHandle< IFPGATrackSimMappingSvc > m_FPGATrackSimMapping
Gaudi::Property< double > m_phiWeight_5hits
Gaudi::Property< double > m_rin
Gaudi::Property< std::string > m_binFilter
Gaudi::Property< std::vector< double > > m_pairSetPhiExtrapCurvedCut
Gaudi::Property< double > m_pairSetDeltaDeltaEtaCut
Gaudi::Property< bool > m_enableMonitoring
Gaudi::Property< unsigned > m_threshold
Gaudi::Property< double > m_pairSetEtaCurvatureCut
Gaudi::Property< double > m_pairSetMatchEtaCut
Gaudi::Property< bool > m_applyPairFilter
void addRoad(std::vector< const StoredHit * > const &hits, const FPGATrackSimBinUtil::IdxSet &idx)
std::vector< unsigned int > m_pairingLayers
Gaudi::Property< bool > m_binningOnly
StatusCode sortHitsByLayer(const BinEntry &bindata, std::vector< std::vector< const StoredHit * > > &hitsByLayer)
StatusCode pairThenGroupFilter(const BinEntry &bindata, std::vector< HitPairSet > &output_pairset)
Gaudi::Property< double > m_etaWeight_5hits
Gaudi::Property< double > m_phiWeight_4hits
StatusCode filterPairs(HitPairSet &pairs, HitPairSet &filteredpairs)
Gaudi::Property< std::vector< double > > m_pairFilterPhiExtrapCut
Gaudi::Property< std::vector< double > > m_pairFilterDeltaPhiCut
StatusCode makePairs(const std::vector< std::vector< const StoredHit * > > &hitsByLayer, HitPairSet &pairs)
Gaudi::Property< bool > m_reversePairDir
FPGATrackSimBinnedHits::BinEntry BinEntry
virtual StatusCode getRoads(const std::vector< std::shared_ptr< const FPGATrackSimHit > > &hits, std::vector< FPGATrackSimRoad > &road) override
Gaudi::Property< double > m_rout
StatusCode incrementalBuildFilter(const BinEntry &bindata, std::vector< HitPairSet > &output_pairset)
FPGATrackSimBinUtil::StoredHit StoredHit
Gaudi::Property< std::vector< double > > m_pairFilterEtaExtrapCut
Gaudi::Property< bool > m_inBinFiltering
ToolHandle< FPGATrackSimBinnedHits > m_binnedhits
void updateState(const IntermediateState &inputstate, IntermediateState &outputstate, unsigned lyridx, const std::vector< const StoredHit * > &newhits)
Gaudi::Property< bool > m_applyPairSetFilter
Gaudi::Property< double > m_pairSetPhiCurvatureCut
Gaudi::Property< double > m_pairSetDeltaPhiCurvatureCut
ToolHandle< FPGATrackSimGenScanMonitoring > m_monitoring
float getR() const
STL class.
double chi2(TH1 *h0, TH1 *h1)
int r
Definition globals.cxx:22
bool verbose
Definition hcg.cxx:73
void reverse(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of reverse for DataVector/List.
std::shared_ptr< const FPGATrackSimHit > hitptr
std::vector< FPGATrackSimBinUtil::StoredHit > hits
double PhiCurvature(const HitPair &pair) const
double MatchPhi(const HitPair &pair) const
double DeltaDeltaPhi(const HitPair &pair) const
bool hasHit(const StoredHit *hit) const
double PhiInExtrapCurved(const HitPair &pair, double r_in) const
double DeltaDeltaEta(const HitPair &pair) const
double EtaCurvature(const HitPair &pair) const
double MatchEta(const HitPair &pair) const
double DeltaEtaCurvature(const HitPair &pair) const
double PhiOutExtrapCurved(const HitPair &pair, double r_out) const
double DeltaPhiCurvature(const HitPair &pair) const
std::vector< const StoredHit * > hitlist
void fill(H5::Group &out_file, size_t iterations)