ATLAS Offline Software
Loading...
Searching...
No Matches
AmbiguityProcessorBase.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
8#include "GaudiKernel/RenounceToolInputsVisitor.h"
9#include "GaudiKernel/ToolVisitor.h"
10
11namespace Trk {
13 const std::string& n,
14 const IInterface* p)
15 : AthAlgTool(t, n, p),
16 m_etaBounds{0.8, 1.6, 2.5, 4.0},
17 m_stat(m_etaBounds) {}
18//
20 const Trk::Track& track) const {
21 return m_tryBremFit and
22 not(track.info().trackProperties(Trk::TrackInfo::BremFit)) and
23 (track.trackParameters()->front()->pT() > m_pTminBrem) and
24 ((not m_caloSeededBrem) or
25 track.info().patternRecoInfo(Trk::TrackInfo::TrackInCaloROI));
26}
27
29 // init ClusterSplitProb input and output container handles and renounce the
30 // output handle from its tools. the latter works because the ClusterSplitProb
31 // output container is created and recorded before the tools are called.
33 !m_clusterSplitProbContainerIn.key().empty()));
35 !m_clusterSplitProbContainerOut.key().empty()));
36 if (!m_clusterSplitProbContainerOut.key().empty()) {
37 auto visitor = [this](const IAlgTool* tool) {
38 const AlgTool* alg_tool = dynamic_cast<const AlgTool*>(tool);
39 for (Gaudi::DataHandle* handle : alg_tool->inputHandles()) {
40 this->msg(MSG::DEBUG) << " input Handle " << tool->name() << " . "
41 << handle->objKey() << endmsg;
42 }
43 for (const auto& elm : alg_tool->inputDataObjs()) {
44 this->msg(MSG::DEBUG)
45 << " input object " << tool->name() << " . " << elm.key() << endmsg;
46 }
47 };
48 if (msgLvl(MSG::DEBUG)) {
49 ToolVisitor::visit(tools(), visitor);
50 }
52 ATH_MSG_DEBUG(" Renounce " << m_clusterSplitProbContainerOut.objKey());
53 auto logger = RenounceToolInputsVisitor::createLogger(
54 [this](const std::string_view& tool_name, const std::string_view& key) {
55 this->msg(MSG::INFO)
56 << " Renounce " << tool_name << " . " << key << endmsg;
57 });
58 RenounceToolInputsVisitor renounce(
59 std::vector<DataObjID>{m_clusterSplitProbContainerOut.fullKey()},
60 logger);
61 ToolVisitor::visit(tools(), renounce);
62 ATH_MSG_DEBUG(" renounced " << m_clusterSplitProbContainerOut.objKey());
63 if (msgLvl(MSG::DEBUG)) {
64 ToolVisitor::visit(tools(), visitor);
65 }
66 }
67 return StatusCode::SUCCESS;
68}
69
72 const EventContext& ctx) const {
74 if (!m_clusterSplitProbContainerIn.key().empty()) {
75 splitProbContainerIn = SG::ReadHandle(m_clusterSplitProbContainerIn, ctx);
76 if (!splitProbContainerIn.isValid()) {
77 ATH_MSG_ERROR("Failed to get input cluster split probability container "
79 }
80 }
81 std::unique_ptr<Trk::ClusterSplitProbabilityContainer> newSplitProbContainer(
83 ? std::make_unique<Trk::ClusterSplitProbabilityContainer>(
84 *splitProbContainerIn)
85 : std::make_unique<Trk::ClusterSplitProbabilityContainer>());
87 splitProbContainerHandle;
88 if (!m_clusterSplitProbContainerOut.key().empty()) {
89 splitProbContainerHandle =
92 if (splitProbContainerHandle.record(std::move(newSplitProbContainer))
93 .isFailure()) {
95 "Failed to record output cluster split probability container "
97 }
98 // newSplitProbContainer owned by storegate -> no cleanup
99 return {splitProbContainerHandle.ptr(),
101 } else {
102 // when not recording the split prob container in storegate the container
103 // must be deleted once going out of scope
104 return {newSplitProbContainer.release(),
105 [](Trk::ClusterSplitProbabilityContainer* ptr) { delete ptr; }};
106 }
107}
108
110 const Trk::Track& track, const TrackParameters* pPar) const {
111 return m_tryBremFit and (pPar->pT() > m_pTminBrem) and
112 ((not m_caloSeededBrem) or
113 track.info().patternRecoInfo(Trk::TrackInfo::TrackInCaloROI));
114}
115//
117 const Trk::Track* track,
118 Trk::PRDtoTrackMap& prdToTrackMap,
119 Counter& stat, int trackId,
120 int subtrackId) const {
121 std::unique_ptr<Trk::Track> newTrack;
122 if (!m_suppressTrackFit) {
123 if (m_refitPrds) {
124 // simple case, fit PRD directly
125 ATH_MSG_VERBOSE("Refit track " << track << " from PRDs");
126 newTrack.reset(refitPrds(ctx, track, prdToTrackMap, stat));
127 } else {
128 // ok, we fit ROTs
129 ATH_MSG_VERBOSE("Refit track " << track << " from ROTs");
130 newTrack.reset(refitRots(ctx, track, stat));
131 }
132 } else {
134 }
135 if (newTrack) {
136 if (m_observerTool.isEnabled()) {
139 m_observerTool->addSubTrack(subtrackId, trackId, *newTrack);
140 }
141 ATH_MSG_DEBUG("New track " << newTrack.get() << " successfully fitted from "
142 << track);
143 } else {
144 if (m_observerTool.isEnabled()) {
147 }
148 ATH_MSG_DEBUG("Fit failed !");
149 }
150 return newTrack.release();
151}
152//
154 const EventContext& ctx,
155 Trk::Track* in_track, const bool fitted, TrackScoreMap& trackScoreTrackMap,
156 std::vector<std::unique_ptr<const Trk::Track> >& trackDustbin,
157 Counter& stat, int parentTrackId) const {
158
159 std::unique_ptr<Trk::Track> atrack(in_track);
160 // compute score
161 TrackScore score = Trk::TrackScore(0);;
162 bool suppressHoleSearch = fitted ? m_suppressHoleSearch : true;
163 bool passBasicSelections = m_scoringTool->passBasicSelections(*atrack);
164 if(passBasicSelections){
165 if (m_trackSummaryTool.isEnabled()) {
166 m_trackSummaryTool->computeAndReplaceTrackSummary(ctx, *atrack,
167 suppressHoleSearch);
168 }
169 bool recheckBasicSel = false;
170 score = m_scoringTool->score(*atrack, recheckBasicSel);
171 }
172
173 if (m_observerTool.isEnabled()) {
174 m_observerTool->updateScore(parentTrackId, static_cast<double>(score));
175 }
176
177 // do we accept the track ?
178 if (score != 0) {
179 ATH_MSG_DEBUG("Track (" << atrack.get() << ") has score " << score);
180 // statistic
181 stat.incrementCounterByRegion(CounterIndex::kNscoreOk, atrack.get());
182 // add track to map, map is sorted small to big !
183 if (m_observerTool.isEnabled()) {
184 trackScoreTrackMap.emplace(
185 -score, TrackPtr(atrack.release(), fitted, parentTrackId));
186 } else {
187 trackScoreTrackMap.emplace(-score, TrackPtr(atrack.release(), fitted));
188 }
189 return;
190 }
191 // do we try to recover the track ?
192 if (fitted and shouldTryBremRecovery(*atrack)) {
193 ATH_MSG_DEBUG("Track score is zero, try to recover it via brem fit");
194 // run track fit using electron hypothesis
195 auto bremTrack(doBremRefit(ctx, *atrack));
196 if (!bremTrack) {
197 ATH_MSG_DEBUG("Brem refit failed, drop track");
198 if (m_observerTool.isEnabled()) {
199 m_observerTool->rejectTrack(parentTrackId,
202 }
203 // statistic
204 stat.incrementCounterByRegion(CounterIndex::kNscoreZeroBremRefitFailed,
205 atrack.get());
206 stat.incrementCounterByRegion(CounterIndex::kNfailedFits, atrack.get());
207 // clean up
208 trackDustbin.push_back(std::move(atrack));
209 } else {
210 int newTrackId = AmbiguityProcessor::getUid();
211 if (m_observerTool.isEnabled()) {
212 m_observerTool->rejectTrack(
213 parentTrackId, xAOD::RejectionStep::addTrack,
215 m_observerTool->addSubTrack(newTrackId, parentTrackId, *bremTrack);
216 }
217 // statistic
218 stat.incrementCounterByRegion(CounterIndex::kNgoodFits, bremTrack.get());
219 // rerun score
220 score = Trk::TrackScore(0);
221 passBasicSelections = m_scoringTool->passBasicSelections(*bremTrack);
222 if(passBasicSelections){
223 if (m_trackSummaryTool.isEnabled()) {
224 m_trackSummaryTool->computeAndReplaceTrackSummary(ctx, *bremTrack,
225 suppressHoleSearch);
226 }
227 bool recheckBasicSel = false;
228 score = m_scoringTool->score(*bremTrack, recheckBasicSel);
229 }
230 if (m_observerTool.isEnabled()) {
231 m_observerTool->updateScore(newTrackId, static_cast<double>(score));
232 }
233 // put original track in the bin, ready to preserve a new Brem track
234 trackDustbin.push_back(std::move(atrack));
235 // do we accept the track ?
236 if (score != 0) {
237 ATH_MSG_DEBUG("Brem refit successful, recovered track ("
238 << bremTrack.get() << ") has score " << score);
239 // statistics
240 stat.incrementCounterByRegion(CounterIndex::kNscoreZeroBremRefit,
241 bremTrack.get());
242 // add track to map, map is sorted small to big !
243 if (m_observerTool.isEnabled()) {
244 m_observerTool->addSubTrack(newTrackId, parentTrackId, *bremTrack);
245 trackScoreTrackMap.emplace(
246 -score, TrackPtr(bremTrack.release(), fitted, newTrackId));
247 } else {
248 trackScoreTrackMap.emplace(-score,
249 TrackPtr(bremTrack.release(), fitted));
250 }
251 return;
252 } else {
253 ATH_MSG_DEBUG("Brem refit gave still track score zero, reject it");
254 if (m_observerTool.isEnabled()) {
255 m_observerTool->rejectTrack(
258 }
259 // statistic
260 stat.incrementCounterByRegion(
262 }
263 }
264 } else {
265 ATH_MSG_DEBUG("Track score is zero, reject it");
266 if (m_observerTool.isEnabled()) {
267 m_observerTool->rejectTrack(parentTrackId, xAOD::RejectionStep::addTrack,
269 }
270 // statistic
271 stat.incrementCounterByRegion(CounterIndex::kNscoreZero, atrack.get());
272 trackDustbin.push_back(std::move(atrack));
273 }
274}
275
277 const Trk::Track* track) const {
278 const TrackParameters* par = track->perigeeParameters();
279 if (not par) {
280 ATH_MSG_DEBUG("Track (" << track << ") has no perigee! Try any other ?");
281 par = track->trackParameters()->front();
282 }
283 if (not par)
284 ATH_MSG_DEBUG("Track (" << track
285 << ") has no Track Parameters ! No refit !");
286 return par;
287}
288
289//==================================================================================================
290
292 const Trk::Track* track,
293 Counter& stat) const {
294 ATH_MSG_VERBOSE("Refit track " << track);
295 // refit using first parameter, do outliers
296 std::unique_ptr<Trk::Track> newTrack{};
297 if (m_tryBremFit and track->info().trackProperties(Trk::TrackInfo::BremFit)) {
298 stat.incrementCounterByRegion(CounterIndex::kNbremFits, track);
299 ATH_MSG_VERBOSE("Brem track, refit with electron brem fit");
300 newTrack = doBremRefit(ctx, *track);
301 } else {
302 stat.incrementCounterByRegion(CounterIndex::kNfits, track);
303 ATH_MSG_VERBOSE("Normal track, refit");
304 newTrack = fit(ctx, *track, true, m_particleHypothesis);
305 if ((not newTrack) and shouldTryBremRecovery(*track)) {
306 stat.incrementCounterByRegion(CounterIndex::kNrecoveryBremFits, track);
307 ATH_MSG_VERBOSE("Normal fit failed, try brem recovery");
308 newTrack = doBremRefit(ctx, *track);
309 }
310 }
311
312 if (newTrack) {
313 stat.incrementCounterByRegion(CounterIndex::kNgoodFits, newTrack.get());
314 // keeping the track of previously accumulated TrackInfo
315 const Trk::TrackInfo& originalInfo = track->info();
316 newTrack->info().addPatternReco(originalInfo);
317 } else {
318 stat.incrementCounterByRegion(CounterIndex::kNfailedFits, track);
319 }
320 return newTrack.release();
321}
322} // namespace Trk
#define endmsg
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_DEBUG(x)
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)
bool msgLvl(const MSG::Level lvl) const
MsgStream & msg() const
virtual bool isValid() override final
Can the handle be successfully dereferenced?
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
pointer_type ptr()
Dereference the pointer.
AmbiguityProcessorBase(const std::string &, const std::string &, const IInterface *)
AmbiguityProcessorBase::UniqueClusterSplitProbabilityContainerPtr createAndRecordClusterSplitProbContainer(const EventContext &ctx) const
Create a new cluster splitting probability container and (optionally) record it in storegate The new ...
AmbiCounter< CounterIndex > Counter
void addTrack(const EventContext &ctx, Trk::Track *in_track, const bool fitted, TrackScoreMap &trackScoreTrackMap, std::vector< std::unique_ptr< const Trk::Track > > &trackDustbin, Counter &stat, int parentTrackId) const
ToolHandle< Trk::IExtendedTrackSummaryTool > m_trackSummaryTool
SG::WriteHandleKey< Trk::ClusterSplitProbabilityContainer > m_clusterSplitProbContainerOut
Track * refitTrack(const EventContext &ctx, const Trk::Track *track, Trk::PRDtoTrackMap &prdToTrackMap, Counter &stat, int trackId, int subtrackId) const
refit track
virtual std::unique_ptr< Trk::Track > doBremRefit(const EventContext &ctx, const Trk::Track &track) const =0
StatusCode initializeClusterSplitProbContainer()
Initialize read and write handles for ClusterSplitProbabilityContainers.
Trk::ParticleHypothesis m_particleHypothesis
virtual std::unique_ptr< Trk::Track > fit(const EventContext &ctx, const Track &track, bool flag, Trk::ParticleHypothesis hypo) const =0
virtual Trk::Track * refitRots(const EventContext &ctx, const Trk::Track *track, Counter &stat) const
ToolHandle< ITrackScoringTool > m_scoringTool
Scoring tool This tool is used to 'score' the tracks, i.e.
virtual Trk::Track * refitPrds(const EventContext &ctx, const Trk::Track *track, Trk::PRDtoTrackMap &prdToTrackMap, Counter &stat) const =0
bool m_tryBremFit
brem recovery mode with brem fit ?
std::vector< float > m_etaBounds
eta intervals for internal monitoring
std::unique_ptr< Trk::ClusterSplitProbabilityContainer, void(*)(Trk::ClusterSplitProbabilityContainer *)> UniqueClusterSplitProbabilityContainerPtr
const TrackParameters * getTrackParameters(const Trk::Track *track) const
std::multimap< TrackScore, TrackPtr > TrackScoreMap
SG::ReadHandleKey< Trk::ClusterSplitProbabilityContainer > m_clusterSplitProbContainerIn
bool shouldTryBremRecovery(const Trk::Track &track) const
PublicToolHandle< Trk::ITrkObserverTool > m_observerTool
Observer tool This tool is used to observe the tracks and their 'score'.
Container to associate Cluster with cluster splitting probabilities.
double pT() const
Access method for transverse momentum.
Contains information about the 'fitter' of this track.
@ BremFit
A brem fit was performed on this track.
static Root::TMsgLogger logger("iLumiCalc")
std::unique_ptr< Trk::Track > createNewFitQualityTrack(const Trk::Track &track)
Ensure that the ATLAS eigen extensions are properly loaded.
float TrackScore
Definition TrackScore.h:10
ParametersBase< TrackParametersDim, Charged > TrackParameters
@ bremRefitTrackScoreZero
@ bremRefitSubtrackCreated