ATLAS Offline Software
Loading...
Searching...
No Matches
FPGATrackSimGenScanTool.cxx
Go to the documentation of this file.
1// Copyright (C) 2002-2025 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());
83 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 ATH_CHECK(m_monitoring->registerHistograms(m_binnedhits.get()));
124
125 // write out the firmware LUTs
126 m_binnedhits->getBinTool().writeLUTs();
127
128 return StatusCode::SUCCESS;
129}
130
132// Main Algorithm
133
134StatusCode FPGATrackSimGenScanTool::getRoads(const std::vector<std::shared_ptr<const FPGATrackSimHit>> &hits,
135 std::vector<std::shared_ptr<const FPGATrackSimRoad>> &roads)
136{
137 ATH_MSG_DEBUG("In getRoads, Processing Event# " << ++m_evtsProcessed << " hit size = " << hits.size());
138
139
140 roads.clear();
141 m_roads.clear();
142 m_monitoring->resetDataFlowCounters();
143
144 // Currently assume that if less than 100 hits its a single track MC
145 m_monitoring->parseTruthInfo(getTruthTracks(),(hits.size() < 100));
146
147 // do the binning...
148 ATH_CHECK(m_binnedhits->fill(hits));
149 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 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 if (m_keepHitsStrategy > 0) {
197 for (std::unique_ptr<FPGATrackSimRoad>& r : m_roads) {
198 const std::vector<std::vector<std::shared_ptr<const FPGATrackSimHit>>>& theseHits = r->getAllHits();
199 layer_bitmask_t hitmask = r->getHitLayers();
200 std::vector<unsigned> toUse = PickHitsToUse(hitmask);
201
202 std::vector<std::vector<std::shared_ptr<const FPGATrackSimHit>>> vec(5); // even if not all layers have hits, they need to be in the vector as empty vectors
203 for (size_t ihit = 0; ihit < toUse.size(); ++ihit) {
204 unsigned int layer = toUse[ihit];
205 if (layer >= theseHits.size() || theseHits[layer].empty()) {
206 ATH_MSG_ERROR("Hit index out of range in keepHitsStrategy: layer=" << layer << ", hits.size()=" << theseHits.size());
207 return StatusCode::FAILURE;
208 }
209 vec[ihit].push_back(theseHits[layer][0]);
210 }
211 r->setHits(std::move(vec));
212 }
213 }
214
215
216 for (auto & r : m_roads) roads.push_back(std::move(r));
217 ATH_MSG_DEBUG("Roads = " << roads.size());
218
219 // clear previous event
220 // (reusing saves having to reallocate memory for each event)
221 m_binnedhits->resetBins();
222
224 return StatusCode::SUCCESS;
225}
226
227// Filter the bins above threshold into pairsets which output roads
229 std::vector<HitPairSet> &output_pairsets)
230{
231 ATH_MSG_VERBOSE("In pairThenGroupFilter");
232
233 // Organize Hits by Layer
234 std::vector<std::vector<const StoredHit *>> hitsByLayer(m_binnedhits->getNLayers());
235 ATH_CHECK(sortHitsByLayer(bindata, hitsByLayer));
236
237 // This is monitoring for each bin over threshold
238 // It's here so it can get the hitsByLayer
239 m_monitoring->fillHitsByLayer(hitsByLayer);
240
241 // Make Pairs
242 HitPairSet pairs;
243 ATH_CHECK(makePairs(hitsByLayer, pairs));
244
245 // Filter Pairs
246 HitPairSet filteredpairs;
247 ATH_CHECK(filterPairs(pairs, filteredpairs));
248
249 // Require road is still over threshold
250 bool passedPairFilter = (filteredpairs.lyrCnt() >= m_threshold);
251 m_monitoring->pairFilterCheck(pairs, filteredpairs, passedPairFilter);
252
253 // if passed Pair Filter proceed to group the filtered pairs into pairsets
254 if (passedPairFilter)
255 {
256 // Pair Set Grouping
257 std::vector<HitPairSet> pairsets;
258 ATH_CHECK(groupPairs(filteredpairs, pairsets, false));
259
260 // loop over pairsets and find those that are over thresold
261 for (const FPGATrackSimGenScanTool::HitPairSet& pairset : pairsets)
262 {
263 // if over threshold add it to the output
264 if (pairset.lyrCnt() >= m_threshold)
265 {
267 output_pairsets.push_back(pairset);
268 }
269 }
270 }
271 }
272
273 // set outputs if not all filters applied
274 if (!m_applyPairFilter) {
275 // output is just the filtered pairs
276 output_pairsets.push_back(std::move(pairs));
277 }
278 else if (passedPairFilter && !m_applyPairSetFilter)
279 {
280 // output is just the filtered pairs
281 output_pairsets.push_back(std::move(filteredpairs));
282 }
283
284 return StatusCode::SUCCESS;
285}
286
287
289 IntermediateState &outputstate,
290 unsigned lyridx,
291 const std::vector<const StoredHit *>& newhits)
292{
293 unsigned int allowed_missed_hits = m_binnedhits->getNLayers() - m_threshold;
294
295 std::vector<bool> pairset_used(inputstate.pairsets.size(),false);
296
297 for (auto &newhit : newhits) {
298
299 // don't make new pairs with hits that the new hit is already paired with in a group
300 std::set<const StoredHit *> vetoList;
301
302 // try adding hit to existing pair sets
303 for (unsigned ps_idx = 0; ps_idx < inputstate.pairsets.size(); ++ps_idx) {
304 auto &pairset = inputstate.pairsets[ps_idx];
305 HitPair nextpair(pairset.lastpair().second, newhit, m_reversePairDir);
306 if (pairMatchesPairSet(pairset, nextpair, false) || (m_applyPairSetFilter==false)) {
307 HitPairSet newset(pairset);
308 newset.addPair(nextpair);
309 pairset_used[ps_idx]=true;
310 outputstate.pairsets.push_back(std::move(newset));
311 // put inpair hits in list of hits not to pair again with the new hits
312 for (auto vetohit : pairset.hitlist) {
313 vetoList.insert(vetohit);
314 }
315 }
316 }
317
318 // make new pairsets with unpaired hits
319 for (auto prevhit : inputstate.unpairedHits) {
320 if (vetoList.count(prevhit) == 0) {
321 HitPair newpair(prevhit, newhit, m_reversePairDir);
322 if (pairPassesFilter(newpair) || (m_applyPairFilter == false)) {
323 HitPairSet newset;
324 newset.addPair(newpair);
325 outputstate.pairsets.push_back(std::move(newset));
326 }
327 }
328 }
329
330 // if this can be the start of a the start of a track and still
331 // have enough hits to make a track, add it to the unpaired hits list
332 if (lyridx <= allowed_missed_hits) {
333 outputstate.unpairedHits.push_back(newhit);
334 }
335 }
336
337 // Add groups to output without new hit if they have enough hits to skip
338 // this layer. Note expected hits at this point is lyridx+1, since we start
339 // counting lyridx from zero. Logic is then keep the pairset if
340 // expected hits <= hits in set + allows misses
341 for (unsigned ps_idx = 0; ps_idx < inputstate.pairsets.size(); ++ps_idx) {
342 auto &pairset = inputstate.pairsets[ps_idx];
343 if ((!pairset_used[ps_idx])&&(lyridx < (pairset.hitlist.size() + allowed_missed_hits))) {
344 outputstate.pairsets.push_back(pairset);
345 }
346 }
347
348 // Add hits to unpaired list if hits are still allowed to start a track
349 // ---------------------------------------------------------------------
350 // If you start a new track with a pair whose second element is
351 // layer N (layer numbering starting from zero), then you'll have missed N-1
352 // layers Minus 1 because the first hit was one of the N layers already
353 // passed.
354 //
355 // The next pass will be layer N=lyridx+1 where lyridx is the current value
356 // at this point in the code. You will have then missed lyridx layers, so...
357 // E.g. if you allow one missed layer then this stops putting hits in the
358 // unpairedHits list if lyridx>1, so tracks must start with layer 0 or 1,
359 // which makes sense
360 if (lyridx > allowed_missed_hits) {
361 // make no new track starts
362 outputstate.unpairedHits.clear();
363 } else {
364 // copy in any previous unpairedHits as well
365 for (auto prevhit : inputstate.unpairedHits) {
366 outputstate.unpairedHits.push_back(prevhit);
367 }
368 }
369}
370
371// Filter the bins above threshold into pairsets which output roads
373 std::vector<HitPairSet> &output_pairsets)
374{
375 ATH_MSG_VERBOSE("In buildGroupsWithPairs");
376
377 // Organize Hits by Layer
378 std::vector<std::vector<const StoredHit *>> hitsByLayer(m_binnedhits->getNLayers());
379 ATH_CHECK(sortHitsByLayer(bindata, hitsByLayer));
380
381 // This is monitoring for each bin over threshold
382 // It's here so it can get the hitsByLayer
383 m_monitoring->fillHitsByLayer(hitsByLayer);
384
385 std::vector<IntermediateState> states{m_binnedhits->getNLayers()+1};
386 for (unsigned lyridx = 0; lyridx < m_binnedhits->getNLayers(); lyridx++) {
387 unsigned lyr = m_pairingLayers[lyridx];
388 updateState(states[lyridx] , states[lyridx+1], lyridx, hitsByLayer[lyr]);
389 }
390
391 // this is a little ugly because it requires copying the output pairsets right now
392 output_pairsets=states[m_binnedhits->getNLayers()].pairsets;
393
394 m_monitoring->fillBuildGroupsWithPairs(states,m_binnedhits->getNLayers()-m_threshold);
395
396 return StatusCode::SUCCESS;
397}
398
399
400// 1st step of filter: sort hits by layer
402 std::vector<std::vector<const StoredHit *>>& hitsByLayer)
403{
404 ATH_MSG_DEBUG("In fillHitsByLayer");
405
406 for (const FPGATrackSimBinUtil::StoredHit &hit : bindata.hits) {
407 hitsByLayer.at(hit.layer).push_back(&hit);
408 }
409
410 return StatusCode::SUCCESS;
411}
412
413
414
415
416
417// 2nd step of filter: make pairs of hits from adjacent and next-to-adjacent layers
418StatusCode FPGATrackSimGenScanTool::makePairs(const std::vector<std::vector<const StoredHit *>>& hitsByLayer,
419 HitPairSet &pairs)
420{
421 ATH_MSG_VERBOSE("In makePairs");
422
423 std::vector<const StoredHit *> const * lastlyr = 0;
424 std::vector<const StoredHit *> const * lastlastlyr = 0;
425
426 // order here is designed so lower radius hits come first
427 for (unsigned lyr : m_pairingLayers) {
428 for (const StoredHit *const &ptr1 :
429 hitsByLayer[lyr]) {
430 if (lastlyr) {
431 for (const StoredHit *const &ptr2 : *lastlyr) {
432 pairs.addPair(HitPair(ptr2, ptr1,m_reversePairDir));
433 }
434 // Add Pairs that skip one layer
435 if (lastlastlyr) {
436 for (const StoredHit *const &ptr2 : *lastlastlyr) {
437 pairs.addPair(HitPair(ptr2, ptr1,m_reversePairDir));
438 }
439 }
440 }
441 }
442 lastlastlyr = lastlyr;
443 lastlyr = &hitsByLayer[lyr];
444 m_monitoring->fillPairingHits(lastlyr,lastlastlyr);
445 }
446
447 return StatusCode::SUCCESS;
448}
449
450// 3rd step of filter: make cuts on the pairs to ensure that are consisten with
451// the bin they are in
452
454 m_monitoring->fillPairFilterCuts(pair,m_rin,m_rout);
455 int lyr = std::min(pair.first->layer,pair.second->layer);
456 return (std::abs(pair.dPhi()) < m_pairFilterDeltaPhiCut[lyr]) &&
457 (std::abs(pair.dEta()) < m_pairFilterDeltaEtaCut[lyr]) &&
458 (std::abs(pair.PhiInExtrap(m_rin)) < m_pairFilterPhiExtrapCut[0]) &&
459 (std::abs(pair.PhiOutExtrap(m_rout)) < m_pairFilterPhiExtrapCut[1]) &&
460 (std::abs(pair.EtaInExtrap(m_rin)) < m_pairFilterEtaExtrapCut[0]) &&
461 (std::abs(pair.EtaOutExtrap(m_rout)) < m_pairFilterEtaExtrapCut[1]);
462}
463
465{
466 ATH_MSG_VERBOSE("In filterPairs");
467
468 for (const FPGATrackSimGenScanTool::HitPair &pair : pairs.pairList) {
469 if (pairPassesFilter(pair)) {
470 filteredpairs.addPair(pair);
471 }
472 }
473 return StatusCode::SUCCESS;
474}
475
476
477// 4th step of filter: group pairs into sets where they are all consistent with being from the same track
479 std::vector<HitPairSet> &pairsets,
480 bool verbose)
481{
482 ATH_MSG_VERBOSE("In groupPairs");
483 for (const FPGATrackSimGenScanTool::HitPair & pair : filteredpairs.pairList)
484 {
485 bool added = false;
486 for (FPGATrackSimGenScanTool::HitPairSet &pairset : pairsets)
487 {
488 // Only add skip pairs if skipped layer is not already hit
489 if ((std::abs(int(pair.second->layer) - int(pair.first->layer)) > 1) // gives if is a skip pair
490 && (pairset.hasLayer(std::min(pair.first->layer,pair.second->layer) + 1))) // gives true if skipped layer already in set
491 {
492 // if it matches mark as added so it doesn't start a new pairset
493 // false here is so it doesn't plot these either
494 if (pairMatchesPairSet(pairset, pair, verbose))
495 added = true;
496 if (!added) {
497 ATH_MSG_VERBOSE("Skip pair does not match non-skip pair");
498 }
499 }
500 else
501 {
502 if (pairMatchesPairSet(pairset, pair, verbose))
503 {
504 int size = pairset.addPair(pair);
505 if (verbose)
506 ATH_MSG_VERBOSE("addPair " << pairsets.size() << " " << pairset.pairList.size() << " " << size);
507 added = true;
508 }
509 }
510
511 }
512
513 if (!added)
514 {
515 HitPairSet newpairset;
516 newpairset.addPair(pair);
517 pairsets.push_back(std::move(newpairset));
518 }
519 }
520
521 return StatusCode::SUCCESS;
522
523}
524
525// used to determine if a pair is consistend with a pairset
527 const HitPair &pair,
528 bool verbose) {
529 // In order to make it easy to have a long list of possible cuts,
530 // a vector of cutvar structs is used to represent each cut
531 // then apply the AND of all the cuts is done with a std::count_if function
532
533 // define the struct (effectively mapping because variable, configured cut value,
534 // and histograms for plotting
535 struct cutvar {
536 cutvar(std::string name, double val, double cut, std::vector<TH1D *>& histset) :
537 m_name(std::move(name)), m_val(val), m_cut(cut), m_histset(histset) {}
538 bool passed() { return std::abs(m_val) < m_cut; }
539 void fill(unsigned cat) { m_histset[cat]->Fill(m_val); }
540 std::string m_name;
541 double m_val;
542 double m_cut;
543 std::vector<TH1D *> &m_histset;
544 };
545
546 // add the cuts to the list of all cuts
547 std::vector<cutvar> allcuts;
548 allcuts.push_back(cutvar("MatchPhi", pairset.MatchPhi(pair),
550 m_monitoring->m_pairSetMatchPhi));
551 allcuts.push_back(cutvar("MatchEta", pairset.MatchEta(pair),
553 m_monitoring->m_pairSetMatchEta));
554 allcuts.push_back(cutvar("DeltaDeltaPhi", pairset.DeltaDeltaPhi(pair),
556 m_monitoring->m_deltaDeltaPhi));
557 allcuts.push_back(cutvar("DeltaDeltaEta", pairset.DeltaDeltaEta(pair),
559 m_monitoring->m_deltaDeltaEta));
560 allcuts.push_back(cutvar("PhiCurvature", pairset.PhiCurvature(pair),
562 m_monitoring->m_phiCurvature));
563 allcuts.push_back(cutvar("EtaCurvature", pairset.EtaCurvature(pair),
565 m_monitoring->m_etaCurvature));
566 if (pairset.pairList.size() > 1) {
567 allcuts.push_back(cutvar(
568 "DeltaPhiCurvature", pairset.DeltaPhiCurvature(pair),
569 m_pairSetDeltaPhiCurvatureCut, m_monitoring->m_deltaPhiCurvature));
570 allcuts.push_back(cutvar(
571 "DeltaEtaCurvature", pairset.DeltaEtaCurvature(pair),
572 m_pairSetDeltaEtaCurvatureCut, m_monitoring->m_deltaEtaCurvature));
573 }
574 allcuts.push_back(cutvar(
575 "PhiInExtrapCurved", pairset.PhiInExtrapCurved(pair, m_rin),
576 m_pairSetPhiExtrapCurvedCut[0], m_monitoring->m_phiInExtrapCurved));
577 allcuts.push_back(cutvar(
578 "PhiOutExtrapCurved", pairset.PhiOutExtrapCurved(pair, m_rout),
579 m_pairSetPhiExtrapCurvedCut[1], m_monitoring->m_phiOutExtrapCurved));
580
581 // count number of cuts passed
582 unsigned passedCuts = std::count_if(allcuts.begin(), allcuts.end(),
583 [](cutvar& cut) { return cut.passed(); });
584 bool passedAll = (passedCuts == allcuts.size());
585
586 // monitoring
587 bool passedAllButOne = (passedCuts == allcuts.size() - 1);
588 for (cutvar& cut: allcuts) {
589 // the last value computes if an n-1 histogram should be filled
590 m_monitoring->fillPairSetFilterCut(cut.m_histset, cut.m_val, pair,
591 pairset.lastpair(),
592 (passedAll || (passedAllButOne && !cut.passed())));
593 }
594
595 if (verbose)
596 {
597 std::string s = "";
598 for (cutvar &cut : allcuts)
599 {
600 s += cut.m_name + " : (" + cut.passed() + ", " + cut.m_val + "), ";
601 }
602 ATH_MSG_DEBUG("PairSet test " << passedAll << " " << s);
603 ATH_MSG_DEBUG("Hits: \n " << *pairset.lastpair().first << "\n "
604 << *pairset.lastpair().second << "\n "
605 << *pair.first << "\n "
606 << *pair.second);
607 }
608
609 return passedAll;
610}
611
612// format final pairsets into expected output of getRoads
613void FPGATrackSimGenScanTool::addRoad(std::vector<const StoredHit *> const &hits, const FPGATrackSimBinUtil::IdxSet &idx)
614{
615 layer_bitmask_t hitLayers = 0;
616 std::vector<std::vector<std::shared_ptr<const FPGATrackSimHit>>>
617 sorted_hits(m_binnedhits->getNLayers(),std::vector<std::shared_ptr<const FPGATrackSimHit>>());
618 for (const FPGATrackSimBinUtil::StoredHit* hit : hits)
619 {
620 hitLayers |= 1 << hit->layer;
621 sorted_hits[hit->layer].push_back(hit->hitptr);
622 }
623
624 // "Fit" the track.
625 FPGATrackSimTrackPars fittedpars;
626 double chi2,chi2_phi,chi2_eta;
627 bool inBin = fitRoad(hits, idx, fittedpars, chi2, chi2_phi,chi2_eta);
628 if (!inBin && m_inBinFiltering) return;
629
630 m_roads.emplace_back(std::make_unique<FPGATrackSimRoad>());
631 FPGATrackSimRoad *r = m_roads.back().get();
632
633 r->setRoadID(m_roads.size() - 1);
634 // r.setPID(y * m_imageSize_y + x);
635 r->setHits(std::move(sorted_hits));
636
637 r->setBinIdx(idx);
638
639 FPGATrackSimBinUtil::ParSet binCenterPars = m_binnedhits->getBinTool().lastStep()->binCenter(idx);
640 FPGATrackSimTrackPars trackpars = m_binnedhits->getBinTool().binDesc()->parSetToTrackPars(binCenterPars);
641 r->setX(trackpars[FPGATrackSimTrackPars::IPHI]);
642 r->setY(trackpars[FPGATrackSimTrackPars::IHIP]);
643 r->setXBin(idx[3]);
644 r->setYBin(idx[4]);
645 r->setHitLayers(hitLayers);
646 r->setSubRegion(0);
647
648 // Store the fitted information on the track.
649 r->setFitParams(fittedpars);
650 r->setFitChi2(chi2);
651 r->setFitChi2_2d(chi2_phi,chi2_eta);
652}
653
654
656// HitPair and HitPairSet Storage
657
659{for (const FPGATrackSimBinUtil::StoredHit* hit : hitlist)
660 {
661 if (hit == newhit)
662 return true;
663 }
664 return false;
665}
666
668{
669 if (pairList.size() > 0)
670 {
671 // these need to be done before the pair is added
674 }
675
676 pairList.push_back(pair);
677
678 hitLayers |= (0x1 << pair.first->layer);
679 hitLayers |= (0x1 << pair.second->layer);
680 if (!hasHit(pair.first))
681 hitlist.push_back(pair.first);
682 if (!hasHit(pair.second))
683 hitlist.push_back(pair.second);
684 return this->pairList.size();
685}
686
687
689{
690 if ((lastpair().first->hitptr==pair.first->hitptr)||
691 (lastpair().first->hitptr==pair.second->hitptr)||
692 (lastpair().second->hitptr==pair.first->hitptr)||
693 (lastpair().second->hitptr==pair.second->hitptr))
694 return 0.0;
695
696 double dr = (lastpair().second->hitptr->getR() - pair.first->hitptr->getR()) / 2.0;
697 double lastpairextrap = lastpair().second->phiShift - lastpair().dPhi() / lastpair().dR() * dr;
698 double newpairextrap = pair.first->phiShift + pair.dPhi() / pair.dR() * dr;
699 return lastpairextrap - newpairextrap;
700}
701
703{
704 if ((lastpair().first->hitptr==pair.first->hitptr)||
705 (lastpair().first->hitptr==pair.second->hitptr)||
706 (lastpair().second->hitptr==pair.first->hitptr)||
707 (lastpair().second->hitptr==pair.second->hitptr))
708 return 0.0;
709
710 double dr = (lastpair().second->hitptr->getR() - pair.first->hitptr->getR()) / 2.0;
711 double lastpairextrap = lastpair().second->etaShift - lastpair().dEta() / lastpair().dR() * dr;
712 double newpairextrap = pair.first->etaShift + pair.dEta() / pair.dR() * dr;
713 return lastpairextrap - newpairextrap;
714}
715
716double
718 return (pair.dPhi() * lastpair().dR() - lastpair().dPhi() * pair.dR()) / (lastpair().dR() * pair.dR());
719}
720
722{
723 return (pair.dEta() * lastpair().dR() - lastpair().dEta() * pair.dR()) / (lastpair().dR() * pair.dR());
724}
725
727{
728 return 2 * DeltaDeltaPhi(pair) / (lastpair().dR() + pair.dR());
729}
731{
732 return 2 * DeltaDeltaEta(pair) / (lastpair().dR() + pair.dR());
733}
743{
744 double r = std::min(lastpair().first->hitptr->getR(),lastpair().second->hitptr->getR());
745 return lastpair().PhiInExtrap(r_in) + 0.5 * PhiCurvature(pair) * (r_in - r) * (r_in - r);
746}
748{
749 double r = std::max(lastpair().first->hitptr->getR(),lastpair().second->hitptr->getR());
750 return pair.PhiOutExtrap(r_out) + 0.5 * PhiCurvature(pair) * (r_out - r) * (r_out - r);
751}
752
753
754// format final pairsets into expected output of getRoads
755bool FPGATrackSimGenScanTool::fitRoad(std::vector<const StoredHit *> const &hits, const FPGATrackSimBinUtil::IdxSet &idx, FPGATrackSimTrackPars& trackpars, double& chi2, double& phi_chi2, double& eta_chi2) const
756{
757
758 double N =hits.size();
759 double sum_Phi = 0;
760 double sum_Phi2 = 0;
761 double sum_PhiR = 0;
762 double sum_PhiR2 = 0;
763 double sum_Eta = 0;
764 double sum_Eta2 = 0;
765 double sum_EtaR = 0;
766 double sum_R = 0;
767 double sum_R2 = 0;
768 double sum_R3 = 0;
769 double sum_R4 = 0;
770
771 for (const FPGATrackSimBinUtil::StoredHit* hit : hits)
772 {
773 // these are just relevant sums of variables (moments) needed for the chi2 calculation
774 // Calculate r^2, r^3, and r^4, we'll sum these up later below
775 double r = hit->hitptr->getR();
776 double r2 = r*r;
777 double r3 = r2*r;
778 double r4 = r3*r;
779 double dphi = hit->phiShift;
780 double dphi2 = dphi*dphi;
781 double deta = hit->etaShift;
782 double deta2 = deta*deta;
783
784 sum_Phi += dphi;
785 sum_Phi2 += dphi2;
786 sum_PhiR += dphi*r;
787 sum_PhiR2 += dphi*r2;
788 sum_Eta += deta;
789 sum_Eta2 += deta2;
790 sum_EtaR += deta*r;
791 sum_R += r;
792 sum_R2 += r2;
793 sum_R3 += r3;
794 sum_R4 += r4;
795 }
796
797 // phi var calculation
798 // phi(r) = phivars[0] + phivars[1]*r + phivars[2]*r^2
799 // math below is calculated by analytically minimizing the chi2
800 // and solving for the phivars.
801 // the terms below which recur in the phivar expression are just organized
802 // by the power of r (i.e. the dimension), but otherwise have no deep meaning
803 double r6_t0 = (-sum_R2 * sum_R4 + sum_R3*sum_R3);
804 double r5_t0 = (sum_R*sum_R4 - sum_R2*sum_R3);
805 double r4_t0 = (-sum_R * sum_R3 + sum_R2*sum_R2);
806 double r4_t1 = (-N*sum_R4 + sum_R2*sum_R2);
807 double r3_t0 = (N*sum_R3 - sum_R*sum_R2);
808 double r2_t0 = (-N*sum_R2 + sum_R*sum_R);
809
810 // all three phi var expresions use the same demoninator
811 const double denom_phi = N * r6_t0 + sum_R * r5_t0 + sum_R2 * r4_t0;
812 if (nearZero(denom_phi)){
813 ATH_MSG_ERROR("Divide by zero (phi) trapped in FPGATrackSimGenScanTool::fitRoad");
814 return false;
815 }
816
817 // phivar expresions from analytic chi2 minimization
818 std::vector<double> phivars(3);
819 phivars[0] = (sum_Phi*r6_t0 + sum_PhiR*r5_t0 + sum_PhiR2*r4_t0)/denom_phi;
820 phivars[1] = (sum_Phi*r5_t0 + sum_PhiR*r4_t1 + sum_PhiR2*r3_t0)/denom_phi;
821 phivars[2] = (sum_Phi*r4_t0 + sum_PhiR*r3_t0 + sum_PhiR2*r2_t0)/denom_phi;
822
823 // eta vars
824 // same as phi but with not curvature (r^2) term
825 const double denom_eta = N*sum_R2 - sum_R*sum_R;
826 if (nearZero(denom_eta)){
827 ATH_MSG_ERROR("Divide by zero (eta) trapped in FPGATrackSimGenScanTool::fitRoad");
828 return false;
829 }
830
831 std::vector<double> etavars(2);
832 etavars[0] = (-sum_R*sum_EtaR + sum_R2*sum_Eta)/denom_eta;
833 etavars[1] = (N*sum_EtaR - sum_R*sum_Eta)/denom_eta;
834
835 // bin center
836 FPGATrackSimBinUtil::ParSet parset = m_binnedhits->getBinTool().lastStep()->binCenter(idx);
837 double parshift[FPGATrackSimTrackPars::NPARS];
838 parshift[0] = etavars[0] + etavars[1]*m_rin; // z at r in
839 parshift[1]= etavars[0] + etavars[1]*m_rout; // z at r out
840 parshift[2]= -(phivars[0]/m_rin + phivars[1] + phivars[2]*m_rin) ; // phi at r in
841 parshift[3]= -(phivars[0]/m_rout + phivars[1] + phivars[2]*m_rout) ; // phi at r out
842 double y = (m_rout-m_rin); // midpoint between rin and rout from rin
843 parshift[4]= -phivars[2]/4.0*y*y ; // xm
844
845 bool inBin = true;
846 const auto& lastStep = m_binnedhits->getBinTool().lastStep();
847 for (int par = 0; par < FPGATrackSimTrackPars::NPARS; par++ ){
848 parset[par]+=parshift[par];
849 inBin = inBin && (std::abs(parshift[par]) < lastStep->binWidth(par)/2.0);
850 }
851
852 double ec0 = etavars[0];
853 double ec1 = etavars[1];
854 eta_chi2 = sum_Eta2 - 2*ec0*sum_Eta - 2*ec1*sum_EtaR + N*ec0*ec0 + 2*ec0*ec1*sum_R + ec1*ec1*sum_R2;
855
856 double pc0 = phivars[0];
857 double pc1 = phivars[1];
858 double pc2 = phivars[2];
859
860 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;
861
862 for (const FPGATrackSimBinUtil::StoredHit* hit : hits)
863 {
864 double r = hit->hitptr->getR();
865 ATH_MSG_VERBOSE("Fitted r= " << r << " phishift " << hit->phiShift << " =?= " << pc0+pc1*r+pc2*r*r << " etashift " << hit->etaShift << " =?= " << ec0+ec1*r);
866 }
867
868 ATH_MSG_DEBUG("Bin Info parset " << idx << " " << m_binnedhits->getBinTool().lastStep()->binCenter(idx));
869 ATH_MSG_DEBUG("Fitted parset inBin="<< inBin << " nhits=" << hits.size() << " " << parset << " chi2 " << eta_chi2 << "," << phi_chi2);
870
871 // Return by reference the "fitted" track parameters. shift q/pt dimension (GeV/MeV)
872 trackpars = m_binnedhits->getBinTool().binDesc()->parSetToTrackPars(parset);
873 trackpars[FPGATrackSimTrackPars::IHIP] = trackpars[FPGATrackSimTrackPars::IHIP] / 1000;
874 ATH_MSG_VERBOSE("Fitted track pars" << trackpars);
875
876 // and the summed chi2, which is assuming (right now) an even weighting between the two components.
877 // assume only options are 4 or 5 hits for now
878 chi2 = (hits.size() == 5) ?
879 m_etaWeight_5hits * eta_chi2 * eta_chi2 + m_phiWeight_5hits * phi_chi2 * phi_chi2 :
880 m_etaWeight_4hits * eta_chi2 * eta_chi2 + m_phiWeight_4hits * phi_chi2 * phi_chi2;
881
882
883 return inBin;
884}
885
886std::vector<unsigned> FPGATrackSimGenScanTool::PickHitsToUse(layer_bitmask_t hitmask) const
887{
888 std::vector<unsigned> toUse;
889 switch (m_keepHitsStrategy) {
890 case 1: // try and pick hits furthest apart, use only 3
891 {
892 if (hitmask == 0x1f) { // miss no hits
893 toUse = {0,2,4};
894 }
895 else if (hitmask == 0x1e) { // miss inner layer, ie layer 0
896 toUse = {1,3,4};
897 }
898 else if (hitmask == 0x1d) { // miss layer 1
899 toUse = {0,2,4};
900 }
901 else if (hitmask == 0x1b) { // miss layer 2
902 toUse = {0,3,4};
903 }
904 else if (hitmask == 0x17) { // miss layer 3
905 toUse = {0,2,4};
906 }
907 else if (hitmask == 0x0f) { // miss layer 4
908 toUse = {0,2,3};
909 }
910 }
911 break;
912 case 2: // pick inner hits, use only 3
913 {
914 if (hitmask == 0x1f) { // miss no hits
915 toUse = {0,1,2};
916 }
917 else if (hitmask == 0x1e) { // miss inner layer, ie layer 0
918 toUse = {1,2,3};
919 }
920 else if (hitmask == 0x1d) { // miss layer 1
921 toUse = {0,2,3};
922 }
923 else if (hitmask == 0x1b) { // miss layer 2
924 toUse = {0,1,3};
925 }
926 else if (hitmask == 0x17) { // miss layer 3
927 toUse = {0,1,2};
928 }
929 else if (hitmask == 0x0f) { // miss layer 4
930 toUse = {0,1,2};
931 }
932 }
933 break;
934 case 3: // pick outer hits, use only 3
935 {
936 if (hitmask == 0x1f) { // miss no hits
937 toUse = {2,3,4};
938 }
939 else if (hitmask == 0x1e) { // miss inner layer, ie layer 0
940 toUse = {2,3,4};
941 }
942 else if (hitmask == 0x1d) { // miss layer 1
943 toUse = {2,3,4};
944 }
945 else if (hitmask == 0x1b) { // miss layer 2
946 toUse = {1,3,4};
947 }
948 else if (hitmask == 0x17) { // miss layer 3
949 toUse = {1,2,4};
950 }
951 else if (hitmask == 0x0f) { // miss layer 4
952 toUse = {1,2,3};
953 }
954 }
955 break;
956 case 4: // keep 4 hits, choose middle one to drop if necessary
957 {
958 if (hitmask == 0x1f) { // miss no hits
959 toUse = {0,1,2,3};
960 }
961 else if (hitmask == 0x1e) { // miss inner layer, ie layer 0
962 toUse = {1,2,3,4};
963 }
964 else if (hitmask == 0x1d) { // miss layer 1
965 toUse = {0,2,3,4};
966 }
967 else if (hitmask == 0x1b) { // miss layer 2
968 toUse = {0,1,3,4};
969 }
970 else if (hitmask == 0x17) { // miss layer 3
971 toUse = {0,1,2,4};
972 }
973 else if (hitmask == 0x0f) { // miss layer 4
974 toUse = {0,1,2,3};
975 }
976 }
977 break;
978 }
979 return toUse;
980}
#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)
std::vector< size_t > vec
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)
virtual int getSubRegion() const override
Gaudi::Property< double > m_pairSetDeltaDeltaPhiCut
Gaudi::Property< double > m_pairSetDeltaEtaCurvatureCut
std::vector< unsigned > PickHitsToUse(layer_bitmask_t) const
ServiceHandle< IFPGATrackSimMappingSvc > m_FPGATrackSimMapping
Gaudi::Property< int > m_keepHitsStrategy
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< 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)
virtual StatusCode getRoads(const std::vector< std::shared_ptr< const FPGATrackSimHit > > &hits, std::vector< std::shared_ptr< const FPGATrackSimRoad > > &road) override
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
Gaudi::Property< double > m_rout
std::vector< std::unique_ptr< FPGATrackSimRoad > > m_roads
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)