ATLAS Offline Software
Loading...
Searching...
No Matches
TrackFindingAlg.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7// Athena
12
13// ACTS
14#include "Acts/Geometry/TrackingGeometry.hpp"
15#include "Acts/Geometry/GeometryIdentifier.hpp"
16#include "Acts/Utilities/TrackHelpers.hpp"
17#include "Acts/TrackFitting/MbfSmoother.hpp"
18#include "Acts/Utilities/Logger.hpp"
19#include "ActsInterop/Logger.h"
20
21// ActsTrk
29
30// STL
31#include <Acts/Propagator/StandardAborters.hpp>
32#include <sstream>
33#include <functional>
34#include <stdexcept>
35#include <utility>
36#include <algorithm>
37#include <variant>
38
39namespace {
40 static std::size_t sourceLinkHash(const Acts::SourceLink& slink) {
41 const ActsTrk::ATLASUncalibSourceLink &atlasSourceLink = slink.get<ActsTrk::ATLASUncalibSourceLink>();
42 const xAOD::UncalibratedMeasurement &uncalibMeas = ActsTrk::getUncalibratedMeasurement(atlasSourceLink);
43 return uncalibMeas.identifier();
44 }
45
46 static bool sourceLinkEquality(const Acts::SourceLink& a, const Acts::SourceLink& b) {
49 return uncalibMeas_a.identifier() == uncalibMeas_b.identifier();
50 }
51
52 static std::optional<ActsTrk::detail::RecoTrackStateContainerProxy> getFirstMeasurementFromTrack(typename ActsTrk::detail::RecoTrackContainer::TrackProxy trackProxy) {
53 std::optional<ActsTrk::detail::RecoTrackStateContainerProxy> firstMeasurement {std::nullopt};
54 for (auto st : trackProxy.trackStatesReversed()) {
55 // We are excluding non measurement states and outlier here. Those can
56 // decrease resolution because only the smoothing corrected the very
57 // first prediction as filtering is not possible.
58 if (not st.typeFlags().hasMeasurement()) continue;
59 if (st.typeFlags().isOutlier()) continue;
60 firstMeasurement = st;
61 }
62 return firstMeasurement;
63 }
64
65}
66
67
68
69namespace ActsTrk
70{
72
73 TrackFindingAlg::TrackFindingAlg(const std::string &name, ISvcLocator *pSvcLocator)
74 : TrackFindingBaseAlg(name, pSvcLocator) {}
75
76 // === initialize ==========================================================
77
79 {
80 ATH_MSG_INFO("Initializing " << name() << " ... ");
81
94
96 ATH_CHECK(m_seedContainerKeys.initialize());
97 ATH_CHECK(m_detEleCollKeys.initialize());
100 ATH_CHECK(m_detElStatus.initialize());
101 ATH_CHECK(m_beamSpotKey.initialize());
102
103 m_storeDestinies = not m_seedDestiny.empty();
105
106 if (m_seedContainerKeys.size() != m_detEleCollKeys.size())
107 {
108 ATH_MSG_FATAL("There are " << m_detEleCollKeys.size() << " DetectorElementsKeys, but " << m_seedContainerKeys.size() << " SeedContainerKeys");
109 return StatusCode::FAILURE;
110 }
111
112 if (m_detEleCollKeys.size() != m_seedLabels.size())
113 {
114 ATH_MSG_FATAL("There are " << m_seedLabels.size() << " SeedLabels, but " << m_detEleCollKeys.size() << " DetectorElementsKeys");
115 return StatusCode::FAILURE;
116 }
117
118 if (m_useTopSpRZboundary.size() != 2)
119 {
120 ATH_MSG_FATAL("useTopSpRZboundary must have 2 elements, but has " << m_useTopSpRZboundary.size());
121 return StatusCode::FAILURE;
122 }
123
124 if (m_ambiStrategy != 0u /* OUTSIDE_TF */) m_showResolvedStats = true;
125 if (m_ambiStrategy == 1u /* END_OF_TF */) {
126 Acts::GreedyAmbiguityResolution::Config cfg;
127 cfg.maximumSharedHits = m_maximumSharedHits;
128 cfg.maximumIterations = m_maximumIterations;
129 cfg.nMeasurementsMin = m_nMeasurementsMin;
130
131 m_ambi.emplace(std::move(cfg), makeActsAthenaLogger(this, "Acts"));
132 }
133
134 if (m_storeDestinies) {
135 if (m_seedDestiny.size() != m_seedContainerKeys.size()) {
136 ATH_MSG_ERROR("There are " << m_seedDestiny.size() << " seed destiny collections, but " << m_seedContainerKeys.size() << " seed collections");
137 return StatusCode::FAILURE;
138 }
139 }
140
141 return StatusCode::SUCCESS;
142 }
143
144 // === finalize ============================================================
145
148 return StatusCode::SUCCESS;
149 }
150
151 // === execute =============================================================
152
153 StatusCode TrackFindingAlg::execute(const EventContext &ctx) const
154 {
155 ATH_MSG_DEBUG("Executing " << name() << " ... ");
156
157 auto timer = Monitored::Timer<std::chrono::milliseconds>("TIME_execute");
158 auto mon_nTracks = Monitored::Scalar<int>("nTracks");
159 auto mon = Monitored::Group(m_monTool, timer, mon_nTracks);
160
161 // ================================================== //
162 // ===================== INPUTS ===================== //
163 // ================================================== //
164
165 // SEED TRIPLETS
166 std::vector<const ActsTrk::SeedContainer *> seedContainers;
167 std::size_t total_seeds = 0;
168 ATH_CHECK(getContainersFromKeys(ctx, m_seedContainerKeys, seedContainers, total_seeds));
169
170 // DESTINIES
171 std::vector< std::unique_ptr< std::vector<int> > > destinies {};
172 if (m_storeDestinies) {
173 destinies.reserve( seedContainers.size() );
174 for (std::size_t i(0); i<seedContainers.size(); ++i) {
175 destinies.push_back( std::make_unique< std::vector<int> >( seedContainers.at(i)->size(), DestinyType::UNKNOWN) );
176 }
177 }
178
179 // MEASUREMENTS
180 std::vector<const xAOD::UncalibratedMeasurementContainer *> uncalibratedMeasurementContainers;
181 std::size_t total_measurements = 0;
182 ATH_CHECK(getContainersFromKeys(ctx, m_uncalibratedMeasurementContainerKeys, uncalibratedMeasurementContainers, total_measurements));
183
184 // map detector element status to volume ids
186 volumeIdToDetectorElementCollMap(m_volumeIdToDetectorElementCollMapKey,ctx);
187 ATH_CHECK(volumeIdToDetectorElementCollMap.isValid());
188 std::vector< const InDet::SiDetectorElementStatus *> det_el_status_arr;
189 const std::vector<const InDetDD::SiDetectorElementCollection*> &det_el_collections =volumeIdToDetectorElementCollMap->collections();
190 det_el_status_arr.resize( det_el_collections.size(), nullptr);
192 SG::ReadHandle<InDet::SiDetectorElementStatus> det_el_status(det_el_status_key,ctx);
193 ATH_CHECK( det_el_status.isValid());
194 const std::vector<const InDetDD::SiDetectorElementCollection*>::const_iterator
195 det_el_col_iter = std::find(det_el_collections.begin(),
196 det_el_collections.end(),
197 &det_el_status->getDetectorElements());
198 det_el_status_arr.at(det_el_col_iter - det_el_collections.begin()) = det_el_status.cptr();
199 }
200
201 detail::MeasurementIndex measurementIndex(uncalibratedMeasurementContainers.size());
202 for (std::size_t icontainer = 0; icontainer < uncalibratedMeasurementContainers.size(); ++icontainer) {
203 measurementIndex.addMeasurements(*uncalibratedMeasurementContainers[icontainer]);
204 }
205
206 detail::TrackFindingMeasurements measurements(uncalibratedMeasurementContainers.size());
207 for (std::size_t icontainer = 0; icontainer < uncalibratedMeasurementContainers.size(); ++icontainer) {
208 ATH_MSG_DEBUG("Create " << uncalibratedMeasurementContainers[icontainer]->size() << " source links from measurements in " << m_uncalibratedMeasurementContainerKeys[icontainer].key());
209 measurements.addMeasurements(icontainer,
210 *uncalibratedMeasurementContainers[icontainer],
211 *m_trackingGeometryTool->surfaceIdMap(),
212 m_forceTrackOnSeed ? &measurementIndex : nullptr);
213 }
214
215 ATH_MSG_DEBUG("measurement index size = " << measurementIndex.size());
216
217 ATH_CHECK( propagateDetectorElementStatusToMeasurements(*(volumeIdToDetectorElementCollMap.cptr()), det_el_status_arr, measurements) );
218
219 if (m_trackStatePrinter.isSet()) {
220 m_trackStatePrinter->printMeasurements(ctx, uncalibratedMeasurementContainers, measurements.measurementOffsets());
221 }
222
223 detail::DuplicateSeedDetector duplicateSeedDetector(total_seeds,
224 m_seedMeasOffset.value(),
226 for (std::size_t icontainer = 0; icontainer < seedContainers.size(); ++icontainer)
227 {
228 duplicateSeedDetector.addSeeds(icontainer, *seedContainers[icontainer], measurementIndex,
229 m_paramEstimationTool->spacePointIndicesFun(),
230 [this,icontainer](const ActsTrk::Seed& seed) -> bool {
231 const bool reverseSearch = m_autoReverseSearch && shouldReverseSearch(seed);
232 const bool refitSeeds = icontainer < m_refitSeeds.size() && m_refitSeeds[icontainer];
233 const bool useTopSp = reverseSearch && !refitSeeds;
234 return useTopSp;
235 });
236 }
237
238 // Get Beam pos and make pSurface
240 ATH_CHECK( beamSpotHandle.isValid() );
241 const InDet::BeamSpotData* beamSpotData = beamSpotHandle.cptr();
242
243 // Beam Spot Position
244 Acts::Vector3 beamPos( beamSpotData->beamPos().x() * Acts::UnitConstants::mm,
245 beamSpotData->beamPos().y() * Acts::UnitConstants::mm,
246 0 );
247
248 // Construct a perigee surface as the target surface
249 std::shared_ptr<Acts::PerigeeSurface> pSurface = Acts::Surface::makeShared<Acts::PerigeeSurface>(beamPos);
250
251 // ================================================== //
252 // ===================== CONDS ====================== //
253 // ================================================== //
254
255 std::vector<const InDetDD::SiDetectorElementCollection*> detElementsCollections;
256 std::size_t total_detElems = 0;
257 ATH_CHECK(getContainersFromKeys(ctx, m_detEleCollKeys, detElementsCollections, total_detElems));
258
259 // ================================================== //
260 // ===================== COMPUTATION ================ //
261 // ================================================== //
262 Acts::VectorTrackContainer actsTrackBackend;
263 Acts::VectorMultiTrajectory actsTrackStateBackend;
264 {
265 std::lock_guard<std::mutex> lock( m_mutex );
266 actsTrackBackend.reserve(m_nTrackReserve);
267 actsTrackStateBackend.reserve(m_nTrackStateReserve);
268 }
269 detail::RecoTrackContainer actsTracksContainer(actsTrackBackend,
270 actsTrackStateBackend);
271
272 if (m_addCounts) {
273 addCounts(actsTracksContainer);
274 }
275
276 detail::ExpectedLayerPatternHelper::add(actsTracksContainer);
277
278 EventStats event_stat;
279 event_stat.resize(m_stat.size());
280
281 DetectorContextHolder detContext {
282 .geometry = m_trackingGeometryTool->getGeometryContext(ctx).context(),
283 .magField = m_extrapolationTool->getMagneticFieldContext(ctx),
284 // CalibrationContext converter not implemented yet.
285 .calib = getCalibrationContext(ctx)
286 };
287
288 detail::SharedHitCounter sharedHits;
289 std::optional<std::vector<unsigned int>> trackCategories;
290 if (m_ambi) trackCategories.emplace(); // only needed if m_ambiStrategy == END_OF_TF
291
292 // Perform the track finding for all initial parameters.
293 for (std::size_t icontainer = 0; icontainer < seedContainers.size(); ++icontainer)
294 {
296 detContext,
297 measurements,
298 measurementIndex,
299 sharedHits,
300 duplicateSeedDetector,
301 *seedContainers.at(icontainer),
302 *detElementsCollections.at(icontainer),
303 actsTracksContainer,
304 icontainer,
305 icontainer < m_seedLabels.size() ? m_seedLabels[icontainer].c_str() : m_seedContainerKeys[icontainer].key().c_str(),
306 event_stat,
307 m_storeDestinies ? destinies.at(icontainer).get() : nullptr,
308 *pSurface.get(),
309 trackCategories));
310 }
311
312 ATH_MSG_DEBUG(" \\__ Created " << actsTracksContainer.size() << " tracks");
313
314 mon_nTracks = actsTracksContainer.size();
315
316
317 // ================================================== //
318 // ===================== OUTPUTS ==================== //
319 // ================================================== //
320
321 // Save the seed destinies
322 if (m_storeDestinies) {
323 for (std::size_t i(0); i<destinies.size(); ++i) {
324 const SG::WriteHandleKey< std::vector<int> >& writeKey = m_seedDestiny.at(i);
325 // make the handle and record
326 SG::WriteHandle< std::vector<int> > destinyHandle = SG::makeHandle( writeKey, ctx );
327 ATH_CHECK( destinyHandle.record( std::move( destinies.at(i) ) ) );
328 }
329 }
330
331 {
332 std::lock_guard<std::mutex> lock( m_mutex );
333 // update the reserve space
334 if (actsTrackBackend.size() > m_nTrackReserve) {
335 m_nTrackReserve = static_cast<std::size_t>( std::ceil(m_memorySafetyMargin * actsTrackBackend.size()) );
336 }
337 if (actsTrackStateBackend.size() > m_nTrackStateReserve) {
338 m_nTrackStateReserve = static_cast<std::size_t>( std::ceil(m_memorySafetyMargin * actsTrackStateBackend.size()) );
339 }
340 }
341
342 // handle the ambiguity
343 // we potentially need to short list the track candidates and make some copies
344 if (not m_ambi) {
345 copyStats(event_stat);
346 // no need to shortlist anything. just use the actsTracksContainer
347 ATH_MSG_DEBUG(" \\__ Created " << actsTracksContainer.size() << " resolved tracks");
349 std::move(actsTrackBackend),
350 std::move(actsTrackStateBackend) ) );
351 return StatusCode::SUCCESS;
352 }
353
354 // we have asked for the ambi
355 // we start by shortlisting the container
356 Acts::VectorTrackContainer resolvedTrackBackend;
357 Acts::VectorMultiTrajectory resolvedTrackStateBackend;
358 resolvedTrackBackend.reserve( actsTrackBackend.size() );
359 resolvedTrackStateBackend.reserve( actsTrackStateBackend.size() );
360 detail::RecoTrackContainer resolvedTracksContainer(resolvedTrackBackend, resolvedTrackStateBackend);
361 detail::ExpectedLayerPatternHelper::add(resolvedTracksContainer);
362
363 if (m_addCounts) {
364 addCounts(resolvedTracksContainer);
365 }
366
367 // Start ambiguity resolution
368 Acts::GreedyAmbiguityResolution::State state;
369 m_ambi->computeInitialState(actsTracksContainer, state, &sourceLinkHash,
370 &sourceLinkEquality);
371 m_ambi->resolve(state);
372
373 // Copy the resolved tracks into the output container
374 // We need a different sharedHits counter here because it saves the track index
375 // and since I ran the resolving the track indices changed.
376 detail::SharedHitCounter sharedHits_forFinalAmbi;
377
378 // shotlist
379 for (auto iTrack : state.selectedTracks) {
380 int actsTrackIndex = state.trackTips.at(iTrack);
381 auto destProxy = resolvedTracksContainer.makeTrack();
382 destProxy.copyFrom(actsTracksContainer.getTrack(actsTrackIndex));
383
384 unsigned int category_i = trackCategories->at(actsTrackIndex);
385 ++event_stat[category_i][kNResolvedTracks];
386
387 if (m_countSharedHits) {
388 auto [nShared, nBadTrackMeasurements] = sharedHits_forFinalAmbi.computeSharedHits(destProxy, resolvedTracksContainer, measurementIndex);
389 if (nBadTrackMeasurements > 0)
390 ATH_MSG_ERROR("computeSharedHits: " << nBadTrackMeasurements << " track measurements not found in input track");
391 }
392 } // loop on tracks
393
394 ATH_MSG_DEBUG(" \\__ Created " << resolvedTracksContainer.size() << " resolved tracks");
395 copyStats(event_stat);
396
398 std::move(resolvedTrackBackend),
399 std::move(resolvedTrackStateBackend)) );
400
401 return StatusCode::SUCCESS;
402 }
403
404 StatusCode TrackFindingAlg::storeTrackCollectionToStoreGate(const EventContext& ctx,
405 Acts::VectorTrackContainer&& originalTrackBackend,
406 Acts::VectorMultiTrajectory&& originalTrackStateBackend) const
407 {
408 // convert to const
409 Acts::ConstVectorTrackContainer constTrackBackend( std::move(originalTrackBackend) );
410 Acts::ConstVectorMultiTrajectory constTrackStateBackend( std::move(originalTrackStateBackend) );
411 std::unique_ptr< ActsTrk::TrackContainer> constTracksContainer = std::make_unique< ActsTrk::TrackContainer >( std::move(constTrackBackend),
412 std::move(constTrackStateBackend) );
413
415 ATH_MSG_DEBUG(" \\__ Tracks Container `" << m_trackContainerKey.key() << "` created ...");
416 ATH_CHECK(trackContainerHandle.record(std::move(constTracksContainer)));
417 return StatusCode::SUCCESS;
418 }
419
421 const xAOD::SpacePoint* bottom_sp = seed.sp().front();
422
423 const double r = bottom_sp->radius();
424 const double z = std::abs(bottom_sp->z());
425
426 const double rBoundary = m_useTopSpRZboundary.value()[0];
427 const double zBoundary = m_useTopSpRZboundary.value()[1];
428
429 return r > rBoundary || z > zBoundary;
430 }
431
432 // === findTracks ==========================================================
433
434 StatusCode
435 TrackFindingAlg::findTracks(const EventContext &ctx,
436 const DetectorContextHolder& detContext,
437 const detail::TrackFindingMeasurements &measurements,
438 const detail::MeasurementIndex &measurementIndex,
439 detail::SharedHitCounter &sharedHits,
440 detail::DuplicateSeedDetector &duplicateSeedDetector,
441 const ActsTrk::SeedContainer &seeds,
442 const InDetDD::SiDetectorElementCollection& detElements,
443 detail::RecoTrackContainer &actsTracksContainer,
444 std::size_t typeIndex,
445 const char *seedType,
446 EventStats &event_stat,
447 std::vector<int>* destiny,
448 const Acts::PerigeeSurface& pSurface,
449 std::optional<std::vector<unsigned int>>& trackCategories) const
450 {
451 ATH_MSG_DEBUG(name() << "::" << __FUNCTION__);
452
453 auto [options, secondOptions, measurementSelector] = getDefaultOptions(ctx, detContext, measurements, &pSurface);
454
455 // ActsTrk::MutableTrackContainer tracksContainerTemp;
456 Acts::VectorTrackContainer trackBackend;
457 Acts::VectorMultiTrajectory trackStateBackend;
458 detail::RecoTrackContainer tracksContainerTemp(trackBackend, trackStateBackend);
459
460 if (m_addCounts) {
461 addCounts(tracksContainerTemp);
462 }
463
464 detail::ExpectedLayerPatternHelper::add(tracksContainerTemp);
465
466 std::size_t category_i = 0;
467 const auto &trackSelectorCfg = trackFinder().trackSelector.config();
468 auto stopBranchProxy = [&](const detail::RecoTrackContainer::TrackProxy &track,
469 const detail::RecoTrackContainer::TrackStateProxy &trackState) -> BranchStopperResult {
470 return stopBranch(track, trackState, trackSelectorCfg, detContext.geometry, measurementIndex, typeIndex, event_stat[category_i]);
471 };
472 options.extensions.branchStopper.connect(stopBranchProxy);
473
474 Acts::PropagatorOptions<detail::Stepper::Options, detail::Navigator::Options,
475 Acts::ActorList<Acts::MaterialInteractor>>
476 extrapolationOptions(detContext.geometry, detContext.magField);
477
478 Acts::TrackExtrapolationStrategy extrapolationStrategy =
479 Acts::TrackExtrapolationStrategy::first;
480
481 // Perform the track finding for all initial parameters
482 ATH_MSG_DEBUG("Invoke track finding with " << seeds.size() << ' ' << seedType << " seeds.");
483
484
485 std::size_t nPrinted = 0;
486
487 // Function for Estimate Track Parameters
488 auto retrieveSurfaceFunction =
489 [this, &detElements] (const ActsTrk::Seed& seed, bool useTopSp) -> const Acts::Surface& {
490 const xAOD::SpacePoint* sp = useTopSp ? seed.sp().back() : seed.sp().front();
491 const InDetDD::SiDetectorElement* element = detElements.getDetectorElement(useTopSp ? sp->elementIdList().back()
492 : sp->elementIdList().front());
493 const Trk::Surface& atlas_surface = element->surface();
494 return *m_ATLASConverterTool->trkSurfaceToActsSurface(atlas_surface);
495 };
496
497
498
499 // Loop over the track finding results for all initial parameters
500 for (unsigned int iseed = 0; iseed < seeds.size(); ++iseed)
501 {
502 // Get the seed
503 const ActsTrk::Seed seed = seeds[iseed];
504
505 category_i = typeIndex * (m_statEtaBins.size() + 1);
506 tracksContainerTemp.clear();
507
508 const bool reverseSearch = m_autoReverseSearch && shouldReverseSearch(seed);
509 const bool refitSeeds = typeIndex < m_refitSeeds.size() && m_refitSeeds[typeIndex];
510 const bool useTopSp = reverseSearch && !refitSeeds;
511
512 // Check if the seed is a duplicate seed
513 const bool isDupSeed = duplicateSeedDetector.isDuplicate(typeIndex, iseed);
514
515 if (isDupSeed) {
516 ATH_MSG_DEBUG("skip " << seedType << " seed " << iseed << " - already found");
517 category_i = getSeedCategory(typeIndex, seed, useTopSp);
518 ++event_stat[category_i][kNTotalSeeds];
519 ++event_stat[category_i][kNDuplicateSeeds];
520 if (m_storeDestinies) destiny->at(iseed) = DestinyType::DUPLICATE;
521 if (!m_trackStatePrinter.isSet()) continue; // delay continue to estimate track parms for TrackStatePrinter?
522 }
523
524 // Set the option accordingly - we change the direction and the target surface accordingly
525 options.propagatorPlainOptions.direction = reverseSearch ? Acts::Direction::Backward() : Acts::Direction::Forward();
526 secondOptions.propagatorPlainOptions.direction = options.propagatorPlainOptions.direction.invert();
527 options.targetSurface = reverseSearch ? &pSurface : nullptr;
528 secondOptions.targetSurface = reverseSearch ? nullptr : &pSurface;
529 // TODO since the second pass is strictly an extension we should have a separate branch stopper which never drops and always extrapolates to the target surface
530
531 // Get first estimate of parameters from the seed
532 std::optional<Acts::BoundTrackParameters> optTrackParams =
533 m_paramEstimationTool->estimateTrackParameters(seed,
534 useTopSp,
535 detContext.geometry,
536 detContext.magField,
537 retrieveSurfaceFunction);
538
539 if (!optTrackParams) {
540 ATH_MSG_DEBUG("Failed to estimate track parameters for seed " << iseed);
541 if (!isDupSeed) {
542 category_i = getSeedCategory(typeIndex, seed, useTopSp);
543 ++event_stat[category_i][kNTotalSeeds];
544 ++event_stat[category_i][kNNoEstimatedParams];
545 if (m_storeDestinies) destiny->at(iseed) = DestinyType::FAILURE;
546 }
547 continue;
548 }
549
550 Acts::BoundTrackParameters *initialParameters = &(*optTrackParams);
551 printSeed(iseed, detContext, seeds, *initialParameters, measurementIndex, nPrinted, seedType);
552 if (isDupSeed) continue; // skip now if not done before
553
554 double etaInitial = -std::log(std::tan(0.5 * initialParameters->theta()));
555 category_i = getStatCategory(typeIndex, etaInitial);
556 ++event_stat[category_i][kNTotalSeeds]; // also updated for duplicate seeds
557 ++event_stat[category_i][kNUsedSeeds];
558
559 // Optional refit track parameters to get more refined value
560 std::unique_ptr<Acts::BoundTrackParameters> refitSeedParameters;
561 if (refitSeeds) {
562 refitSeedParameters = doRefit(seed, *initialParameters, detContext, reverseSearch);
563 if (refitSeedParameters.get() == nullptr) {
564 ++event_stat[category_i][kNRejectedRefinedSeeds];
565 if (m_storeDestinies) destiny->at(iseed) = DestinyType::FAILURE;
566 continue;
567 }
568 if (refitSeedParameters.get() != initialParameters) {
569 initialParameters = refitSeedParameters.get();
570 printSeed(iseed, detContext, seeds, *initialParameters, measurementIndex, nPrinted, seedType, true);
571 }
572 }
573
574 auto measurementRangesForced =
575 m_forceTrackOnSeed ? measurements.createMeasurementRangesForced(seed, measurementIndex)
576 : std::unique_ptr<ActsTrk::detail::MeasurementRangeListFlat>();
577 measurementSelector->setMeasurementRangesForced(measurementRangesForced.get());
578 if (measurementRangesForced)
579 event_stat[category_i][kNForcedSeedMeasurements] += measurementRangesForced->size();
580
581 // Get the Acts tracks, given this seed
582 Acts::Result<std::vector<TrkProxy> > result =
583 trackFinder().ckf.findTracks(*initialParameters, options, tracksContainerTemp);
584
585 // The result for this seed
586 if (not result.ok()) {
587 ATH_MSG_WARNING("Track finding failed for " << seedType << " seed " << iseed << " with error" << result.error());
588 if (m_storeDestinies) destiny->at(iseed) = DestinyType::FAILURE;
589 continue;
590 }
591 auto &tracksForSeed = result.value();
592
593
594
595 std::size_t ntracks = 0ul;
596
597 // loop on the tracks we have just found from the seed
598 std::size_t nfirst = 0;
599 for (TrkProxy &firstTrack : tracksForSeed) {
600 // smoothing
601 auto smoothingResult = Acts::smoothTrack(detContext.geometry, firstTrack, logger(), Acts::MbfSmoother());
602 if (!smoothingResult.ok()) {
603 ATH_MSG_DEBUG("Smoothing for seed "
604 << iseed << " and first track " << firstTrack.index()
605 << " failed with error " << smoothingResult.error());
606 continue;
607 }
608
609 // if no two way, just add the track and move on
610 if (not m_doTwoWay) {
611 // add the track to the collection
612 ATH_CHECK( addTrack(detContext,
613 firstTrack,
614 pSurface,
615 extrapolationStrategy,
616 sharedHits,
617 actsTracksContainer,
618 measurementIndex,
619 tracksContainerTemp,
620 duplicateSeedDetector,
621 destiny,
622 event_stat,
623 ntracks,
624 iseed,
625 category_i,
626 seedType,
627 trackCategories) );
628 ++nfirst;
629 continue;
630 }
631
632 // TWO WAY STARTS HERE
633 // We need the first measurement of the track
634 std::optional<detail::RecoTrackStateContainerProxy> firstMeas = getFirstMeasurementFromTrack(firstTrack);
635 // we are supposed to find a measurement
636 if (not firstMeas.has_value()) {
637 ATH_MSG_ERROR("Could not retrieve first measurement from track proxy. Is it ill-formed?");
638 return StatusCode::FAILURE;
639 }
640 detail::RecoTrackStateContainerProxy& firstMeasurement = firstMeas.value();
641
642 // Get the tracks from the second track finding
643 std::vector<typename detail::RecoTrackContainer::TrackProxy> secondTracksForSeed =
644 doTwoWayTrackFinding(firstMeasurement,
645 firstTrack,
646 tracksContainerTemp,
647 secondOptions);
648
649 if ( secondTracksForSeed.empty() ) {
650 ATH_MSG_DEBUG("No viable result from second track finding for " << seedType << " seed " << iseed << " track " << nfirst);
651 ++event_stat[category_i][kNoSecond];
652 ATH_CHECK( addTrack(detContext,
653 firstTrack,
654 pSurface,
655 extrapolationStrategy,
656 sharedHits,
657 actsTracksContainer,
658 measurementIndex,
659 tracksContainerTemp,
660 duplicateSeedDetector,
661 destiny,
662 event_stat,
663 ntracks,
664 iseed,
665 category_i,
666 seedType,
667 trackCategories) );
668 }
669
670 // need to add tracks here
671 // do the stiching
672 // store the original previous state to restore it later
673 auto originalFirstMeasurementPrevious = firstMeasurement.previous();
674 for (auto &secondTrack : secondTracksForSeed) {
675 secondTrack.reverseTrackStates(true);
676
677 firstMeasurement.previous() = secondTrack.outermostTrackState().index();
678 secondTrack.tipIndex() = firstTrack.tipIndex();
679
680 if (reverseSearch) {
681 // smooth the full track
682 auto secondSmoothingResult = Acts::smoothTrack(detContext.geometry,
683 secondTrack,
684 logger());
685 if ( not secondSmoothingResult.ok() ) {
686 continue;
687 }
688 secondTrack.reverseTrackStates(true);
689 }
690
691 // Add track to collection
692 ATH_CHECK( addTrack(detContext,
693 secondTrack,
694 pSurface,
695 extrapolationStrategy,
696 sharedHits,
697 actsTracksContainer,
698 measurementIndex,
699 tracksContainerTemp,
700 duplicateSeedDetector,
701 destiny,
702 event_stat,
703 ntracks,
704 iseed,
705 category_i,
706 seedType,
707 trackCategories) );
708 } // loop on tracks
709
710 // finish the stiching
711 // restore the original previous state for the first track
712 firstMeasurement.previous() = originalFirstMeasurementPrevious;
713
714 nfirst++;
715 } // loop on tracks from seed
716
717 if (m_storeDestinies) {
718 if (ntracks == 0) {
719 destiny->at(iseed) = DestinyType::FAILURE;
720 } else {
721 destiny->at(iseed) = DestinyType::SUCCEED;
722 }
723 }
724
725 if (ntracks == 0) {
726 ATH_MSG_DEBUG("Track finding found no track candidates for " << seedType << " seed " << iseed);
727 ++event_stat[category_i][kNoTrack];
728 } else if (ntracks >= 2) {
729 ++event_stat[category_i][kMultipleBranches];
730 }
731
732 if (m_trackStatePrinter.isSet())
733 std::cout << std::flush;
734 } // loop on seeds
735
736 ATH_MSG_DEBUG("Completed " << seedType << " track finding with " << computeStatSum(typeIndex, kNOutputTracks, event_stat) << " track candidates.");
737
738 return StatusCode::SUCCESS;
739 }
740
741 void
744 detail::DuplicateSeedDetector &duplicateSeedDetector,
745 const detail::MeasurementIndex &measurementIndex) const {
746
747 const auto lastMeasurementIndex = track.tipIndex();
748 duplicateSeedDetector.newTrajectory();
749
750 tracksContainer.trackStateContainer().visitBackwards(
751 lastMeasurementIndex,
752 [&duplicateSeedDetector,&measurementIndex](const detail::RecoTrackStateContainer::ConstTrackStateProxy &state) -> void
753 {
754 // Check there is a source link
755 if (not state.hasUncalibratedSourceLink())
756 return;
757
758 // Fill the duplicate selector
759 auto sl = state.getUncalibratedSourceLink().template get<ATLASUncalibSourceLink>();
760 duplicateSeedDetector.addMeasurement(sl, measurementIndex);
761 }); // end visitBackwards
762 }
763
765 const std::vector< const InDet::SiDetectorElementStatus *> &det_el_status_arr,
766 detail::TrackFindingMeasurements &measurements) const {
767 const Acts::TrackingGeometry *
768 acts_tracking_geometry = m_trackingGeometryTool->trackingGeometry().get();
769 ATH_CHECK(acts_tracking_geometry != nullptr);
770
771 using Counter = struct { unsigned int n_volumes, n_volumes_with_status, n_missing_detector_elements, n_detector_elements, n_disabled_detector_elements;};
772 Counter counter {0u,0u,0u,0u,0u};
773 acts_tracking_geometry->visitVolumes([&counter,
774 &volume_id_to_det_el_coll,
775 &det_el_status_arr,
776 &measurements,
777 this](const Acts::TrackingVolume *volume_ptr) {
778 ++counter.n_volumes;
779 if (!volume_ptr) return;
780
782 det_el_status = det_el_status_arr.at(volume_id_to_det_el_coll.collecionMap().at(volume_ptr->geometryId().volume()));
783 if (det_el_status) {
784 ++counter.n_volumes_with_status;
785 volume_ptr->visitSurfaces([&counter, det_el_status, &measurements,this](const Acts::Surface *surface_ptr) {
786 if (!surface_ptr) return;
787 const Acts::Surface &surface = *surface_ptr;
788 const Acts::SurfacePlacementBase* detector_element = surface.surfacePlacement();
789 if (detector_element) {
790 ++counter.n_detector_elements;
791 const ActsDetectorElement *acts_detector_element = static_cast<const ActsDetectorElement*>(detector_element);
792 if (!det_el_status->isGood( acts_detector_element->identifyHash() )) {
793 ActsTrk::detail::MeasurementRange old_range = measurements.markSurfaceInsensitive(surface_ptr->geometryId());
794 if (!old_range.empty()) {
795 auto geoid_to_string = [](const Acts::GeometryIdentifier &id) -> std::string {
796 std::stringstream amsg;
797 amsg << id;
798 return amsg.str();
799 };
800 std::string a_msg ( geoid_to_string(surface_ptr->geometryId()));
801 ATH_MSG_WARNING("Reject " << (old_range.elementEndIndex() - old_range.elementBeginIndex())
802 << " measurements because surface " << a_msg);
803 }
804 ++counter.n_disabled_detector_elements;
805 }
806 }
807 }, true /*only sensitive surfaces*/);
808 }
809 else {
810 ++counter.n_missing_detector_elements;
811 }
812 });
813 ATH_MSG_DEBUG("Volumes with detector element status " << counter.n_volumes_with_status << " / " << counter.n_volumes
814 << " disabled detector elements " << counter.n_disabled_detector_elements
815 << " / " << counter.n_detector_elements
816 << " missing detector elements "
817 << counter.n_missing_detector_elements);
818 return StatusCode::SUCCESS;
819 }
820
821 std::size_t TrackFindingAlg::getSeedCategory(std::size_t typeIndex,
822 const ActsTrk::Seed& seed,
823 bool useTopSp) const
824 {
825 const xAOD::SpacePoint* sp = useTopSp ? seed.sp().back() : seed.sp().front();
826 const xAOD::SpacePoint::ConstVectorMap pos = sp->globalPosition();
827 double etaSeed = std::atanh(pos[2] / pos.norm());
828 return getStatCategory(typeIndex, etaSeed);
829 }
830
831 void TrackFindingAlg::printSeed(unsigned int iseed,
832 const DetectorContextHolder& detContext,
833 const ActsTrk::SeedContainer& seeds,
834 const Acts::BoundTrackParameters &seedParameters,
835 const detail::MeasurementIndex &measurementIndex,
836 std::size_t& nPrinted,
837 const char *seedType,
838 bool isKF) const
839 {
840 if (not m_trackStatePrinter.isSet()) return;
841
842 if (nPrinted == 0) {
843 ATH_MSG_INFO("CKF results for " << seeds.size() << ' ' << seedType << " seeds:");
844 }
845 ++nPrinted;
846 m_trackStatePrinter->printSeed(detContext.geometry, seeds[iseed], seedParameters, measurementIndex, iseed, isKF);
847 }
848
849namespace {
850struct Collector {
851 using result_type = TrackFindingAlg::ExpectedLayerPattern*;
852
853 template <typename propagator_state_t, typename stepper_t,
854 typename navigator_t>
855 Acts::Result<void> act(propagator_state_t& state, const stepper_t& /*stepper*/,
856 const navigator_t& navigator, result_type& result,
857 const Acts::Logger& /*logger*/) const {
858 const Acts::Surface* currentSurface = navigator.currentSurface(state.navigation);
859 if (currentSurface == nullptr) {
860 return Acts::Result<void>::success();
861 }
862
863 assert(result != nullptr && "Result type is nullptr");
864
865 if (currentSurface->surfacePlacement() != nullptr) {
866 const auto* detElem = dynamic_cast<const ActsDetectorElement*>(currentSurface->surfacePlacement());
867 if(detElem != nullptr) {
868 detail::addToExpectedLayerPattern(*result, *detElem);
869 }
870 }
871
872 return Acts::Result<void>::success();
873 }
874};
875}
876
878 const DetectorContextHolder& detContext,
880 const Acts::Surface &referenceSurface,
881 const detail::Extrapolator &propagator,
882 Acts::TrackExtrapolationStrategy strategy,
883 ExpectedLayerPattern& expectedLayerPattern) const {
884
885 Acts::PropagatorOptions<detail::Stepper::Options, detail::Navigator::Options,
886 Acts::ActorList<Acts::MaterialInteractor, Collector>>
887 options(detContext.geometry, detContext.magField);
888
889 auto findResult = findTrackStateForExtrapolation(
890 options.geoContext, track, referenceSurface, strategy, logger());
891
892 if (!findResult.ok()) {
893 ATH_MSG_WARNING("Failed to find track state for extrapolation");
894 return findResult.error();
895 }
896
897 auto &[trackState, distance] = *findResult;
898
899 options.direction = Acts::Direction::fromScalarZeroAsPositive(distance);
900
901 Acts::BoundTrackParameters parameters = track.createParametersFromState(trackState);
902 ATH_MSG_VERBOSE("Extrapolating track to reference surface at distance "
903 << distance << " with direction " << options.direction
904 << " with starting parameters " << parameters);
905
906 auto state = propagator.makeState<decltype(options), Acts::ForcedSurfaceReached>(referenceSurface, options);
907 ExpectedLayerPattern*& collectorResult = state.get<TrackFindingAlg::ExpectedLayerPattern*>();
908 collectorResult = &expectedLayerPattern;
909
910 auto initRes = propagator.initialize(state, parameters);
911 if(!initRes.ok()) {
912 ATH_MSG_WARNING("Failed to initialize propagation state: " << initRes.error().message());
913 return initRes.error();
914 }
915
916
917 auto propagateOnlyResult =
918 propagator.propagate(state);
919
920 if (!propagateOnlyResult.ok()) {
921 ATH_MSG_WARNING("Failed to extrapolate track: " << propagateOnlyResult.error().message());
922 return propagateOnlyResult.error();
923 }
924
925 auto propagateResult = propagator.makeResult(
926 std::move(state), propagateOnlyResult, options, true, &referenceSurface);
927
928 if (!propagateResult.ok()) {
929 ATH_MSG_WARNING("Failed to extrapolate track: " << propagateResult.error().message());
930 return propagateResult.error();
931 }
932
933 track.setReferenceSurface(referenceSurface.getSharedPtr());
934 track.parameters() = propagateResult->endParameters.value().parameters();
935 track.covariance() =
936 propagateResult->endParameters.value().covariance().value();
937
938 return Acts::Result<void>::success();
939 }
940
943 const Acts::Surface& pSurface,
944 const Acts::TrackExtrapolationStrategy& extrapolationStrategy,
945 detail::SharedHitCounter &sharedHits,
946 detail::RecoTrackContainer &actsTracksContainer,
947 const detail::MeasurementIndex& measurementIndex,
948 const detail::RecoTrackContainer& tracksContainerTemp,
949 detail::DuplicateSeedDetector& duplicateSeedDetector,
950 std::vector<int>* destiny,
951 EventStats& event_stat,
952 std::size_t& ntracks,
953 std::size_t iseed,
954 std::size_t category_i,
955 const char *seedType,
956 std::optional<std::vector<unsigned int>>& trackCategories) const
957 {
958
959 std::array<unsigned int, 4> expectedLayerPattern{};
960
961 // if the the perigeeSurface was not hit (in particular the case for the inside-out pass,
962 // the track has no reference surface and the extrapolation to the perigee has not been done
963 // yet.
964 if (not track.hasReferenceSurface()) {
965 auto extrapolationResult =
966 extrapolateTrackToReferenceSurface(detContext, track,
967 pSurface,
968 trackFinder().extrapolator,
969 extrapolationStrategy,
970 expectedLayerPattern);
971
972 if (not extrapolationResult.ok()) {
973 ATH_MSG_WARNING("Extrapolation for seed "
974 << iseed << " and " << track.index()
975 << " failed with error " << extrapolationResult.error()
976 << " dropping track candidate.");
977 if (m_storeDestinies) destiny->at(iseed) = DestinyType::FAILURE;
978 return StatusCode::SUCCESS;
979 }
980 }
981
982 // Before trimming, inspect encountered surfaces from all track states
983 for(const auto ts : track.trackStatesReversed()) {
984 const auto& surface = ts.referenceSurface();
985 if(surface.surfacePlacement() != nullptr) {
986 const auto* detElem = dynamic_cast<const ActsDetectorElement*>(surface.surfacePlacement());
987 if(detElem != nullptr) {
988 detail::addToExpectedLayerPattern(expectedLayerPattern, *detElem);
989 }
990 }
991 }
992
993 // Trim tracks
994 // - trimHoles
995 // - trimOutliers
996 // - trimMaterial
997 // - trimOtherNoneMeasurement
998 Acts::trimTrack(track, true, true, true, true);
999 Acts::calculateTrackQuantities(track);
1000 if (m_addCounts) {
1001 initCounts(track);
1002 for (const auto trackState : track.trackStatesReversed()) {
1003 updateCounts(track, trackState.typeFlags(), measurementType(trackState));
1004 }
1005 if (m_checkCounts) {
1006 checkCounts(track);
1007 }
1008 }
1009
1010 ++ntracks;
1011 ++event_stat[category_i][kNOutputTracks];
1012
1013 if ( not trackFinder().trackSelector.isValidTrack(track) or
1014 not selectCountsFinal(track)) {
1015 ATH_MSG_DEBUG("Track " << ntracks << " from " << seedType << " seed " << iseed << " failed track selection");
1016 if ( m_trackStatePrinter.isSet() ) {
1017 m_trackStatePrinter->printTrack(detContext.geometry, tracksContainerTemp, track, measurementIndex, true);
1018 }
1019 return StatusCode::SUCCESS;
1020 }
1021
1022 ++event_stat[category_i][kNSelectedTracks];
1023
1024 // Fill the track infos into the duplicate seed detector
1026 storeSeedInfo(tracksContainerTemp, track, duplicateSeedDetector, measurementIndex);
1027 }
1028
1029 auto actsDestProxy = actsTracksContainer.makeTrack();
1030 actsDestProxy.copyFrom(track); // make sure we copy track states!
1031
1032 detail::ExpectedLayerPatternHelper::set(actsDestProxy, expectedLayerPattern);
1033
1034 auto setTrackCategory = [&]() {
1035 if (!trackCategories) return;
1036 if (!(actsDestProxy.index() < trackCategories->size())) trackCategories->resize(actsDestProxy.index()+1);
1037 trackCategories->at(actsDestProxy.index()) = category_i;
1038 };
1039
1040 if (not m_countSharedHits) {
1041 return StatusCode::SUCCESS;
1042 }
1043
1044 auto [nShared, nBadTrackMeasurements] = sharedHits.computeSharedHits(actsDestProxy, actsTracksContainer, measurementIndex);
1045
1046 if (nBadTrackMeasurements > 0) {
1047 ATH_MSG_ERROR("computeSharedHits: " << nBadTrackMeasurements << " track measurements not found in input for " << seedType << " seed " << iseed << " track");
1048 }
1049
1050 ATH_MSG_DEBUG("found " << actsDestProxy.nSharedHits() << " shared hits in " << seedType << " seed " << iseed << " track");
1051
1052 event_stat[category_i][kNTotalSharedHits] += nShared;
1053
1054 if (m_ambiStrategy == 2u) { // run the ambiguity during track selection
1055
1056 if (actsDestProxy.nSharedHits() <= m_maximumSharedHits) {
1057 setTrackCategory();
1058 ++event_stat[category_i][kNResolvedTracks];
1059 }
1060 else { // track fails the shared hit selection
1061
1062 ATH_MSG_DEBUG("found " << actsDestProxy.nSharedHits() << " shared hits in " << seedType << " seed " << iseed << " track");
1063 // Reset the original track shared hits by running coumputeSharedHits
1064 // with removeSharedHits flag to true
1065 // nSharedRemoved contains the total shared hits that will be removed
1066 auto [nSharedRemoved, nRemoveBadTrackMeasurements] = sharedHits.computeSharedHits(actsDestProxy, actsTracksContainer, measurementIndex, true);
1067
1068 ATH_MSG_DEBUG("Removed " << nSharedRemoved << " shared hits in " << seedType << " seed " << iseed << " track and the matching track");
1069
1070 if (nRemoveBadTrackMeasurements > 0) {
1071 ATH_MSG_ERROR("computeSharedHits with remove flag ON: " << nRemoveBadTrackMeasurements <<
1072 " track measurements not found in input for " << seedType << " seed " << iseed << " track");
1073 }
1074
1075 if (actsDestProxy.nSharedHits() != 0) {
1076 ATH_MSG_ERROR("computeSharedHits with remove flag ON returned " <<
1077 actsDestProxy.nSharedHits()<< " while expecting 0 for" <<
1078 seedType << " seed " << iseed << " track");
1079 }
1080
1081 // Remove the track from the container
1082 actsTracksContainer.removeTrack(actsDestProxy.index());
1083 ATH_MSG_DEBUG("Track " << ntracks << " from " << seedType << " seed " << iseed << " failed shared hit selection");
1084 }
1085 }
1086 else {
1087 // run ambi later
1088 setTrackCategory();
1089 if (m_trackStatePrinter.isSet()) {
1090 m_trackStatePrinter->printTrack(detContext.geometry, actsTracksContainer, actsDestProxy, measurementIndex);
1091 }
1092 }
1093
1094 return StatusCode::SUCCESS;
1095 }
1096
1097} // namespace
#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_WARNING(x)
#define ATH_MSG_DEBUG(x)
virtual void lock()=0
Interface to allow an object to lock itself when made const in SG.
static Double_t sp
static Double_t a
Header file to be included by clients of the Monitored infrastructure.
size_t size() const
Number of registered mappings.
std::unique_ptr< const Acts::Logger > makeActsAthenaLogger(IMessageSvc *svc, const std::string &name, int level, std::optional< std::string > parent_name)
#define z
IdentifierHash identifyHash() const
Identifier hash.
SG::WriteHandleKeyArray< std::vector< int > > m_seedDestiny
Gaudi::Property< bool > m_countSharedHits
Gaudi::Property< unsigned int > m_maximumSharedHits
Gaudi::Property< float > m_memorySafetyMargin
SG::ReadHandleKeyArray< xAOD::UncalibratedMeasurementContainer > m_uncalibratedMeasurementContainerKeys
SG::ReadHandleKeyArray< ActsTrk::SeedContainer > m_seedContainerKeys
std::size_t getSeedCategory(std::size_t typeIndex, const ActsTrk::Seed &seed, bool useTopSp) const
SG::ReadCondHandleKey< ActsTrk::ActsVolumeIdToDetectorElementCollectionMap > m_volumeIdToDetectorElementCollMapKey
SG::ReadCondHandleKeyArray< InDetDD::SiDetectorElementCollection > m_detEleCollKeys
TrackFindingAlg(const std::string &name, ISvcLocator *pSvcLocator)
void printSeed(unsigned int iseed, const DetectorContextHolder &detContext, const ActsTrk::SeedContainer &seeds, const Acts::BoundTrackParameters &seedParameters, const detail::MeasurementIndex &measurementIndex, std::size_t &nPrinted, const char *seedType, bool isKF=false) const
Gaudi::Property< std::vector< bool > > m_refitSeeds
void storeSeedInfo(const detail::RecoTrackContainer &tracksContainer, const detail::RecoTrackContainerProxy &track, detail::DuplicateSeedDetector &duplicateSeedDetector, const detail::MeasurementIndex &measurementIndex) const
std::optional< Acts::GreedyAmbiguityResolution > m_ambi
StatusCode propagateDetectorElementStatusToMeasurements(const ActsTrk::ActsVolumeIdToDetectorElementCollectionMap &volume_id_to_det_el_coll, const std::vector< const InDet::SiDetectorElementStatus * > &det_el_status_arr, detail::TrackFindingMeasurements &measurements) const
Gaudi::Property< unsigned int > m_nMeasurementsMin
ToolHandle< ActsTrk::ITrackParamsEstimationTool > m_paramEstimationTool
Gaudi::Property< unsigned int > m_seedMeasOffset
virtual StatusCode initialize() override
StatusCode findTracks(const EventContext &ctx, const DetectorContextHolder &detContext, const detail::TrackFindingMeasurements &measurements, const detail::MeasurementIndex &measurementIndex, detail::SharedHitCounter &sharedHits, detail::DuplicateSeedDetector &duplicateSeedDetector, const ActsTrk::SeedContainer &seeds, const InDetDD::SiDetectorElementCollection &detElements, detail::RecoTrackContainer &actsTracksContainer, std::size_t seedCollectionIndex, const char *seedType, EventStats &event_stat, std::vector< int > *destiny, const Acts::PerigeeSurface &pSurface, std::optional< std::vector< unsigned int > > &trackCategories) const
invoke track finding procedure
Gaudi::Property< std::vector< double > > m_useTopSpRZboundary
Acts::Result< void > extrapolateTrackToReferenceSurface(const DetectorContextHolder &detContext, detail::RecoTrackContainerProxy &track, const Acts::Surface &referenceSurface, const detail::Extrapolator &propagator, Acts::TrackExtrapolationStrategy strategy, ExpectedLayerPattern &expectedLayerPattern) const
Gaudi::Property< bool > m_autoReverseSearch
StatusCode storeTrackCollectionToStoreGate(const EventContext &ctx, Acts::VectorTrackContainer &&originalTrackBackend, Acts::VectorMultiTrajectory &&originalTrackStateBackend) const
Gaudi::Property< unsigned int > m_maximumIterations
std::array< unsigned int, 4 > ExpectedLayerPattern
Gaudi::Property< bool > m_forceTrackOnSeed
Gaudi::Property< bool > m_skipDuplicateSeeds
SG::ReadCondHandleKey< InDet::BeamSpotData > m_beamSpotKey
virtual StatusCode execute(const EventContext &ctx) const override
virtual StatusCode finalize() override
StatusCode addTrack(const DetectorContextHolder &detContext, detail::RecoTrackContainerProxy &track, const Acts::Surface &pSurface, const Acts::TrackExtrapolationStrategy &extrapolationStrategy, detail::SharedHitCounter &sharedHits, detail::RecoTrackContainer &actsTracksContainer, const detail::MeasurementIndex &measurementIndex, const detail::RecoTrackContainer &tracksContainerTemp, detail::DuplicateSeedDetector &duplicateSeedDetector, std::vector< int > *destiny, EventStats &event_stat, std::size_t &ntracks, std::size_t iseed, std::size_t category_i, const char *seedType, std::optional< std::vector< unsigned int > > &trackCategories) const
Gaudi::Property< std::size_t > m_ambiStrategy
bool shouldReverseSearch(const ActsTrk::Seed &seed) const
SG::ReadHandleKeyArray< InDet::SiDetectorElementStatus > m_detElStatus
Gaudi::Property< std::vector< float > > m_statEtaBins
Gaudi::Property< bool > m_doTwoWay
Gaudi::Property< bool > m_dumpAllStatEtaBins
static xAOD::UncalibMeasType measurementType(const detail::RecoTrackContainer::TrackStateProxy &trackState)
ToolHandle< ActsTrk::TrackStatePrinterTool > m_trackStatePrinter
SG::WriteHandleKey< ActsTrk::TrackContainer > m_trackContainerKey
detail::RecoTrackContainer::TrackProxy TrkProxy
void copyStats(const EventStats &event_stat) const
std::size_t getStatCategory(std::size_t seed_collection, float eta) const
ToolHandle< GenericMonitoringTool > m_monTool
static void addCounts(detail::RecoTrackContainer &tracksContainer)
std::vector< typename detail::RecoTrackContainer::TrackProxy > doTwoWayTrackFinding(const detail::RecoTrackStateContainerProxy &firstMeasurement, const TrkProxy &trackProxy, detail::RecoTrackContainer &tracksContainerTemp, const TrackFinderOptions &options) const
Perform two-way track finding.
TrackFindingBaseAlg(const std::string &name, ISvcLocator *pSvcLocator)
PublicToolHandle< ActsTrk::ITrackingGeometryTool > m_trackingGeometryTool
ToolHandle< ActsTrk::IActsToTrkConverterTool > m_ATLASConverterTool
ToolHandle< ActsTrk::IExtrapolationTool > m_extrapolationTool
std::size_t computeStatSum(std::size_t seed_collection, EStat counter_i, const EventStats &stat) const
virtual StatusCode initialize() override
StatusCode getContainersFromKeys(const EventContext &ctx, HandleArrayKeyType &handleKeyArray, std::vector< const ContainerType * > &outputContainers, std::size_t &sum) const
Take the array of handle keys and for each key retrieve containers, then append them to the output ve...
void checkCounts(const detail::RecoTrackContainer::TrackProxy &track) const
const Acts::Logger & logger() const
Private access to the logger.
BranchStopperResult stopBranch(const detail::RecoTrackContainer::TrackProxy &track, const detail::RecoTrackContainer::TrackStateProxy &trackState, const Acts::TrackSelector::EtaBinnedConfig &trackSelectorCfg, const Acts::GeometryContext &tgContext, const detail::MeasurementIndex &measurementIndex, const std::size_t typeIndex, EventStats::value_type &event_stat_category_i) const
Branch stopper.
static void initCounts(const detail::RecoTrackContainer::TrackProxy &track)
Gaudi::Property< bool > m_checkCounts
TrackFindingDefaultOptions getDefaultOptions(const EventContext &ctx, const DetectorContextHolder &detContext, const detail::TrackFindingMeasurements &measurements, const Acts::PerigeeSurface *pSurface) const
Get CKF options for first and second pass + pointer to MeasurementSelector.
std::unique_ptr< Acts::BoundTrackParameters > doRefit(const MeasurementSource &measurement, const Acts::BoundTrackParameters &initialParameters, const DetectorContextHolder &detContext, const bool paramsAtOutermostSurface) const
Perform Kalman Filter fit and update given initialParameters.
virtual StatusCode finalize() override
std::vector< std::array< unsigned int, kNStat > > EventStats
static void updateCounts(const detail::RecoTrackContainer::TrackProxy &track, Acts::ConstTrackStateTypeMap typeFlags, xAOD::UncalibMeasType detType)
Gaudi::Property< bool > m_addCounts
Acts::CombinatorialKalmanFilterBranchStopperResult BranchStopperResult
bool selectCountsFinal(const detail::RecoTrackContainer::TrackProxy &track) const
Gaudi::Property< std::vector< std::string > > m_seedLabels
void addMeasurement(const ActsTrk::ATLASUncalibSourceLink &sl, const MeasurementIndex &measurementIndex)
bool isDuplicate(std::size_t typeIndex, index_t iseed)
void addSeeds(std::size_t typeIndex, const ActsTrk::SeedContainer &seeds, const MeasurementIndex &measurementIndex)
void addMeasurements(const xAOD::UncalibratedMeasurementContainer &clusterContainer)
auto computeSharedHits(typename track_container_t::TrackProxy &track, track_container_t &tracks, const MeasurementIndex &measurementIndex, bool removeSharedHits=false) -> ReturnSharedAndBad
const std::vector< std::size_t > & measurementOffsets() const
std::unique_ptr< MeasurementRangeListFlat > createMeasurementRangesForced(const ActsTrk::Seed &seed, const MeasurementIndex &measurementIndex) const
void addMeasurements(std::size_t typeIndex, const xAOD::UncalibratedMeasurementContainer &clusterContainer, const DetectorElementToActsGeometryIdMap &detectorElementToGeoid, const MeasurementIndex *measurementIndex=nullptr)
MeasurementRange markSurfaceInsensitive(const Acts::GeometryIdentifier &identifier)
Class to hold the SiDetectorElement objects to be put in the detector store.
const SiDetectorElement * getDetectorElement(const IdentifierHash &hash) const
Class to hold geometrical description of a silicon detector element.
Trk::Surface & surface()
Element Surface.
const Amg::Vector3D & beamPos() const noexcept
bool isGood(IdentifierHash hash) const
Group of local monitoring quantities and retain correlation when filling histograms
Declare a monitored scalar variable.
A monitored timer.
const_pointer_type cptr()
Property holding a SG store/key/clid from which a ReadHandle is made.
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type cptr()
Dereference the pointer.
Property holding a SG store/key/clid from which a WriteHandle is made.
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
Abstract Base Class for tracking surfaces.
Definition Surface.h:79
float z() const
Eigen::Map< const Eigen::Matrix< float, 3, 1 > > ConstVectorMap
float radius() const
DetectorIdentType identifier() const
Returns the full Identifier of the measurement.
int ts
Definition globals.cxx:24
int r
Definition globals.cxx:22
T * get(TKey *tobj)
get a TObject* from a TKey* (why can't a TObject be a TKey?)
Definition hcg.cxx:132
RecoTrackStateContainer::TrackStateProxy RecoTrackStateContainerProxy
Acts::TrackContainer< Acts::VectorTrackContainer, Acts::VectorMultiTrajectory > RecoTrackContainer
void addToExpectedLayerPattern(std::array< unsigned int, 4 > &pattern, const ActsDetectorElement &detElement)
The AlignStoreProviderAlg loads the rigid alignment corrections and pipes them through the readout ge...
const xAOD::UncalibratedMeasurement & getUncalibratedMeasurement(const ATLASUncalibSourceLink &source_link)
Acts::CalibrationContext getCalibrationContext(const EventContext &ctx)
The Acts::Calibration context is piped through the Acts fitters to (re)calibrate the Acts::SourceLink...
const xAOD::UncalibratedMeasurement * ATLASUncalibSourceLink
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
UncalibratedMeasurement_v1 UncalibratedMeasurement
Define the version of the uncalibrated measurement class.
std::size_t size() const noexcept
static void add(track_container_t &trackContainer)
static void set(track_proxy_t &track, std::array< unsigned int, 4 > values)