ATLAS Offline Software
Loading...
Searching...
No Matches
SiSPSeededTrackFinderRoI.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
12
13#include "xAODTracking/Vertex.h"
16
17#include <set>
18#include <vector>
19#include <sstream>
20#include <fstream>
21#include <string>
22
23namespace {
25// Track quality calculation
27
28double trackQuality(const Trk::Track* Tr) {
29
30 double quality = 0.;
31 double baseScorePerHit = 17.;
33 for (const Trk::TrackStateOnSurface* m : *(Tr->trackStateOnSurfaces())) {
36 continue;
38 const Trk::FitQualityOnSurface fq = m->fitQualityOnSurface();
39 if (!fq)
40 continue;
41
42 double x2 = fq.chiSquared();
43 double hitQualityScore;
46 if (fq.numberDoF() == 2)
47 hitQualityScore = (1.2 * (baseScorePerHit - x2 * .5)); // pix
48 else
49 hitQualityScore = (baseScorePerHit - x2); // sct
50 if (hitQualityScore < 0.)
51 hitQualityScore =
52 0.; // do not allow a bad hit to decrement the overall score
53 quality += hitQualityScore;
54 }
57 quality *= 0.7;
58
59 return quality;
60}
61} // namespace
62 //
64// Constructor
66
68(const std::string& name,ISvcLocator* pSvcLocator) : AthReentrantAlgorithm(name, pSvcLocator)
69{
70}
71
73// Initialisation
75
77{
78
79 // Initialize read/write handles
80 ATH_CHECK( m_outputTracksKey.initialize() );
81 ATH_CHECK( m_vxOutputKey.initialize() );
82 ATH_CHECK( m_prdToTrackMap.initialize( !m_prdToTrackMap.key().empty() ) );
83 ATH_CHECK( m_beamSpotKey.initialize() );
84 ATH_CHECK(m_evtKey.initialize());
85
86 // Retrieve Tools
87 ATH_CHECK( m_ZWindowRoISeedTool.retrieve() );
88 ATH_CHECK( m_RandomRoISeedTool.retrieve( DisableTool{(not m_doRandomSpot) or m_RandomRoISeedTool.name().empty()} ) );
89 ATH_CHECK( m_seedsmaker.retrieve() );
90 ATH_CHECK( m_trackmaker.retrieve() );
91 ATH_CHECK( m_trackSummaryTool.retrieve( DisableTool{ m_trackSummaryTool.name().empty()} ));
92 if (not m_etaDependentCutsSvc.name().empty()) ATH_CHECK(m_etaDependentCutsSvc.retrieve());
93 magneticFieldInit(); //init magnetic field for simplified propagation
94
96 if (msgLvl(MSG::DEBUG)) {
97 dump(MSG::DEBUG, nullptr);
98 }
100 m_problemsTotal = 0;
101
102 return StatusCode::SUCCESS;
103}
104
105
107// Extended local EDM
109
110namespace InDet {
112 {
113 public:
120 private:
122 };
123}
124
126// Execute
128
129StatusCode InDet::SiSPSeededTrackFinderRoI::execute(const EventContext& ctx) const
130{
132 ATH_CHECK(outputTracks.record(std::make_unique<TrackCollection>()));
133
135 ATH_CHECK( beamSpotHandle.isValid() );
136 Trk::PerigeeSurface beamPosPerigee(beamSpotHandle->beamPos());
137
138 bool PIX = true ;
139 bool SCT = true ;
140 bool ERR = false;
141
142 //Create container for RoI information
143 auto theVertexContainer = std::make_unique<xAOD::VertexContainer>();
144 auto theVertexAuxContainer = std::make_unique<xAOD::VertexAuxContainer>();
145 theVertexContainer->setStore( theVertexAuxContainer.get() );
146
147 // Find reference point of the event and create z boundary region
148 //
149 std::vector<InDet::IZWindowRoISeedTool::ZWindow> listRoIs;
150 listRoIs = m_ZWindowRoISeedTool->getRoIs(ctx);
151 double ZBoundary[2] = {0.0, 0.0};
152 //if no RoI found; no need to go further
153 if ( listRoIs.empty() ) {
154 ATH_MSG_DEBUG("no selectedRoIs" );
155 } else {
156 //listRoIs[0].zReference is the midpoint
157 ZBoundary[0] = listRoIs[0].zWindow[0];
158 ZBoundary[1] = listRoIs[0].zWindow[1];
159 ATH_MSG_DEBUG("selectedRoIs " << ZBoundary[0] <<" " << ZBoundary[1]);
160 }
161
162 //Store RoI information in a xAOD::Vertex object
163 static const SG::AuxElement::Accessor<float> vtxDecor_boundaryLow("boundaryLow");
164 static const SG::AuxElement::Accessor<float> vtxDecor_boundaryHigh("boundaryHigh");
165 static const SG::AuxElement::Accessor<float> vtxDecor_perigeeZ0Lead("perigeeZ0Lead");
166 static const SG::AuxElement::Accessor<float> vtxDecor_perigeeZ0Sublead("perigeeZ0Sublead");
167 static const SG::AuxElement::Accessor<float> vtxDecor_isHS("isHS");
168
169 for( size_t r = 0; r < listRoIs.size(); r++ ){
170
171 theVertexContainer->push_back( new xAOD::Vertex() );
172
173 theVertexContainer->back()->setZ( listRoIs[r].zReference );
174 vtxDecor_boundaryLow(*theVertexContainer->back()) = listRoIs[r].zWindow[0];;
175 vtxDecor_boundaryHigh(*theVertexContainer->back()) = listRoIs[r].zWindow[1];
176 vtxDecor_perigeeZ0Lead(*theVertexContainer->back()) = listRoIs[r].zPerigeePos[0];
177 vtxDecor_perigeeZ0Sublead(*theVertexContainer->back()) = listRoIs[r].zPerigeePos[1];
178 vtxDecor_isHS(*theVertexContainer->back()) = 1;
179 }
180
181 // Analyses that want to run low-pt tracking with a region of interest care
182 // about the beam conditions near a collision of interest. Validation of the
183 // beam conditions elsewhere in the beamspot (regarding low-pt tracks) will be
184 // needed to establish meaningful uncertainties. Choosing a random position
185 // allows for this check. Run with RAWtoESD section of postexec:
186 // ToolSvc.InDetSiSpTrackFinder_LowPtRoI.doRandomSpot = True
187 double RandZBoundary[2];
188 std::vector<InDet::IZWindowRoISeedTool::ZWindow> listRandRoIs;
189 if(m_doRandomSpot and not listRoIs.empty()){
190 //Finding Random Spot in beamspot
191 listRandRoIs = m_RandomRoISeedTool->getRoIs(ctx);
192
193 while( std::abs( listRoIs[0].zReference - listRandRoIs[0].zReference ) < 5. || std::abs(listRandRoIs[0].zReference) > 250.0 ){
194 listRandRoIs.clear();
195 listRandRoIs = m_RandomRoISeedTool->getRoIs(ctx);
196 }
197
198 RandZBoundary[0] = listRandRoIs[0].zWindow[0];
199 RandZBoundary[1] = listRandRoIs[0].zWindow[1];
200 for( size_t r = 0; r < listRandRoIs.size(); r++ ){
201
202 theVertexContainer->push_back( new xAOD::Vertex() );
203
204 theVertexContainer->back()->setZ( listRandRoIs[r].zReference );
205 vtxDecor_boundaryLow(*theVertexContainer->back()) = listRoIs[r].zWindow[0];;
206 vtxDecor_boundaryHigh(*theVertexContainer->back()) = listRoIs[r].zWindow[1];
207 vtxDecor_perigeeZ0Lead(*theVertexContainer->back()) = listRoIs[r].zPerigeePos[0];
208 vtxDecor_perigeeZ0Sublead(*theVertexContainer->back()) = listRoIs[r].zPerigeePos[1];
209 vtxDecor_isHS(*theVertexContainer->back()) = 0;
210 }
211 }
212
213 // Record the RoI information
215 CHECK( vxOut_h.record ( std::move(theVertexContainer), std::move(theVertexAuxContainer) ) );
216
217 // Find seeds that point within the RoI region in z
218 //
220 m_seedsmaker->newEvent(ctx, seedEventData, -1);
221 if (not listRoIs.empty()) {
222 std::list<Trk::Vertex> VZ;
223 if(m_useRoIWidth) m_seedsmaker->find3Sp(ctx, seedEventData, VZ, ZBoundary);
224 //If you want to disable the RoI but still have a separate container for low-pt tracks
225 // The RoI "vertex" container will still be there in case you want to use
226 // that information for whatever reason (ie where the RoI would have been centered).
227 else m_seedsmaker->find3Sp(ctx, seedEventData, VZ);
228 if(m_doRandomSpot) m_seedsmaker->find3Sp(ctx, seedEventData, VZ, RandZBoundary);
229 }
230
232 m_trackmaker->newEvent(ctx, trackEventData, PIX, SCT);
233
234
235 // Get the value of the seed maker validation ntuple writing switch
236 bool doWriteNtuple = m_seedsmaker->getWriteNtupleBoolProperty();
237 unsigned long EvNumber = 0; //Event number variable to be used for the validation ntuple
238
239 if (doWriteNtuple) {
241 EvNumber = !eventInfo.isValid() ? 0 : eventInfo->eventNumber();
242 }
243
244 // Loop through all seed and create track candidates
245 //
246 ERR = false;
247 Counter_t counter{};
248 std::multimap<double,Trk::Track*> qualitySortedTrackCandidates;
249 const InDet::SiSpacePointsSeed* seed = nullptr;
250
251 while((seed = m_seedsmaker->next(ctx, seedEventData))) {
252 ++counter[kNSeeds];
253 const std::list<Trk::Track*> trackList = m_trackmaker->getTracks(ctx, trackEventData, seed->spacePoints());
254 for(Trk::Track* t: trackList) {
255 qualitySortedTrackCandidates.insert(std::make_pair( -trackQuality(t), t ));
256 }
257 if(doWriteNtuple) {
258 // Note: determining if "pixel" or "strips" or "mixed" seed depending on innermost/outermost radii - hardcoded boundaries for Run 1-3 detector
260 if (seed->r3() > 200.) {
261 if (seed->r1() < 200.) {
263 } else {
265 }
266 }
267 m_seedsmaker->writeNtuple(seed, !trackList.empty() ? trackList.front() : nullptr, seedType, EvNumber) ;
268 }
269 if( counter[kNSeeds] >= m_maxNumberSeeds) {
270 ERR = true;
272 break;
273 }
274 }
275
276 m_trackmaker->endEvent(trackEventData);
277
278 // Remove shared tracks with worse quality
279 //
280 filterSharedTracks(qualitySortedTrackCandidates);
281
282 // Save good tracks in track collection
283 //
284 for (const std::pair<const double, Trk::Track*> & qualityAndTrack: qualitySortedTrackCandidates) {
285 ++counter[kNTracks];
286 if (m_trackSummaryTool.isEnabled()) {
287 m_trackSummaryTool->computeAndReplaceTrackSummary(*(qualityAndTrack.second),
288 false /* DO NOT suppress hole search*/);
289 }
290 outputTracks->push_back(qualityAndTrack.second);
291 }
292 m_counterTotal[kNSeeds] += counter[kNSeeds];
294
295 // In case of errors, clear output tracks
296 //
297 if (ERR) {
298 outputTracks->clear();
299 } else {
300 m_counterTotal[kNTracks] += counter[kNTracks];
301 }
302
303 // Print common event information
304 //
305 if (msgLvl(MSG::DEBUG)) {
306 dump(MSG::DEBUG, &counter);
307 }
308
309 return StatusCode::SUCCESS;
310
311}
312
314// Finalize
316
318{
319 dump(MSG::INFO, &m_counterTotal);
320
321 return StatusCode::SUCCESS;
322}
323
325// Dumps relevant information into the MsgStream
327
328MsgStream& InDet::SiSPSeededTrackFinderRoI::dump(MSG::Level assign_level, const InDet::SiSPSeededTrackFinderRoI::Counter_t* counter) const
329{
330 msg(assign_level) <<std::endl;
331 MsgStream& out_msg=msg();
332 if (counter) dumpevent(out_msg ,*counter);
333 else dumptools(out_msg);
334 out_msg << endmsg;
335 return out_msg;
336}
337
339// Dumps conditions information into the MsgStream
341
342MsgStream& InDet::SiSPSeededTrackFinderRoI::dumptools( MsgStream& out ) const
343{
344 int n;
345 n = 65-m_seedsmaker.type().size();
346 std::string s2; for(int i=0; i<n; ++i) s2.append(" "); s2.append("|");
347 n = 65-m_trackmaker.type().size();
348 std::string s3; for(int i=0; i<n; ++i) s3.append(" "); s3.append("|");
349 n = 65-m_outputTracksKey.key().size();
350 std::string s4; for(int i=0; i<n; ++i) s4.append(" "); s4.append("|");
351
352 out<<"|----------------------------------------------------------------"
353 <<"----------------------------------------------------|"
354 <<std::endl;
355 out<<"| Tool for space points seeds finding | "<<m_seedsmaker.type()<<s2
356 <<std::endl;
357 out<<"| Tool for space points seeded track finding | "<<m_trackmaker.type()<<s3
358 <<std::endl;
359 out<<"| Location of output tracks | "<<m_outputTracksKey.key()<<s4
360 <<std::endl;
361 out<<"|----------------------------------------------------------------"
362 <<"----------------------------------------------------|"
363 <<std::endl;
364 return out;
365}
366
368// Dumps event information into the ostream
370
372{
373 int ns = counter[kNSeeds];
374 int nt = counter[kNTracks];
375
376 out<<"|-------------------------------------------------------------------|" <<std::endl;
377 out<<"| Investigated "
378 <<std::setw(9)<<ns<<" space points seeds and found ";
379 out<<std::setw(9)<<nt<<" tracks using RoI-z strategy |"<<std::endl;
380
381 out<<"|-------------------------------------------------------------------|" <<std::endl;
382 if(m_problemsTotal > 0) {
383 out<<"| Events "
384 <<std::setw(7)<<m_neventsTotal <<" |"
385 <<std::endl;
386 out<<"| Problems "
387 <<std::setw(7)<<m_problemsTotal <<" |"
388 <<std::endl;
389 out<<"|-------------------------------------------------------------------|" <<std::endl;
390 }
391 return out;
392}
393
395// Filer shared tracks
397
399(std::multimap<double,Trk::Track*>& qualitySortedTracks) const
400{
401 std::set<const Trk::PrepRawData*> clusters;
402
403 std::vector<const Trk::PrepRawData*> freeClusters;
404 freeClusters.reserve(15);
405
406 std::multimap<double, Trk::Track*>::iterator it_qualityAndTrack = qualitySortedTracks.begin();
407
409 while (it_qualityAndTrack!=qualitySortedTracks.end()) {
410 freeClusters.clear();
411
412 std::set<const Trk::PrepRawData*>::iterator it_clustersEnd = clusters.end();
413
414 int nClusters = 0;
415 int nPixels = 0;
416
418 for (const Trk::TrackStateOnSurface* tsos: *((*it_qualityAndTrack).second->trackStateOnSurfaces())) {
420 if(!tsos->type(Trk::TrackStateOnSurface::Measurement)) continue;
421 const Trk::FitQualityOnSurface fq = tsos->fitQualityOnSurface();
422 if(!fq) continue;
423 if(fq.numberDoF() == 2) ++nPixels;
424
425 const Trk::MeasurementBase* mb = tsos->measurementOnTrack();
426 const Trk::RIO_OnTrack* ri = dynamic_cast<const Trk::RIO_OnTrack*>(mb);
427 if(!ri) continue;
428 const Trk::PrepRawData* pr = ri->prepRawData();
429 if (not pr) continue;
431 ++nClusters;
433 if (clusters.find(pr)==it_clustersEnd) {
435 freeClusters.push_back(pr);
436 }
437 }
438
440 int nFreeClusters = static_cast<int>(freeClusters.size());
441 if (nFreeClusters >= m_nfreeCut || nFreeClusters==nClusters) {
444 clusters.insert(freeClusters.begin(), freeClusters.end());
445
446 if (m_etaDependentCutsSvc.name().empty()) {
447 ++it_qualityAndTrack;
448 } else {
449 //eta-dependent cuts
450 int nFreeClusters = static_cast<int>(freeClusters.size());
451 if( passEtaDepCuts( (*it_qualityAndTrack).second, nClusters, nFreeClusters, nPixels) ){
453 ++it_qualityAndTrack;
454 } else {
456 delete (*it_qualityAndTrack).second;
457 qualitySortedTracks.erase(it_qualityAndTrack++);
458 }
459 } //eta-dependent cuts
460 } else {
462 delete (*it_qualityAndTrack).second;
463 qualitySortedTracks.erase(it_qualityAndTrack++);
464 }
465 }
466}
467
469// Callback function - get the magnetic field /
471
473{
474 // Build MagneticFieldProperties
475 //
476 if(m_fieldmode == "NoField") {
478 } else {
480 }
481}
482
484 int nClusters,
485 int nFreeClusters,
486 int nPixels) const
487{
488 Trk::TrackStates::const_iterator m = track->trackStateOnSurfaces()->begin();
489 const Trk::TrackParameters* par = (*m)->trackParameters();
490 if(!par) return false;
491
492 double eta = std::abs(par->eta());
493 if(nClusters < m_etaDependentCutsSvc->getMinSiHitsAtEta(eta)) return false;
494 if(nFreeClusters < m_etaDependentCutsSvc->getMinSiNotSharedAtEta(eta)) return false;
495 if(nClusters-nFreeClusters > m_etaDependentCutsSvc->getMaxSharedAtEta(eta)) return false;
496 if(nPixels < m_etaDependentCutsSvc->getMinPixelHitsAtEta(eta)) return false;
497
498 if(par->pT() < m_etaDependentCutsSvc->getMinPtAtEta(eta)) return false;
499 if(!(*m)->type(Trk::TrackStateOnSurface::Perigee)) return true ;
500 if(std::abs(par->localPosition()[0]) > m_etaDependentCutsSvc->getMaxPrimaryImpactAtEta(eta)) return false;
501 return true;
502}
503
504
505
Scalar eta() const
pseudorapidity method
#define endmsg
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_DEBUG(x)
#define CHECK(...)
Evaluate an expression and check for errors.
bool msgLvl(const MSG::Level lvl) const
An algorithm that can be simultaneously executed in multiple threads.
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
ExtendedSiTrackMakerEventData_xk(const SG::ReadHandleKey< Trk::PRDtoTrackMap > &key)
SG::ReadHandle< Trk::PRDtoTrackMap > m_prdToTrackMap
std::atomic_int m_neventsTotal
Number events.
ToolHandle< Trk::IExtendedTrackSummaryTool > m_trackSummaryTool
SG::WriteHandleKey< xAOD::VertexContainer > m_vxOutputKey
virtual StatusCode initialize() override
void filterSharedTracks(std::multimap< double, Trk::Track * > &) const
cleans up the collection of quality filtered tracks.
std::atomic_int m_problemsTotal
Number events with number seeds > maxNumber.
MsgStream & dumpevent(MsgStream &out, const SiSPSeededTrackFinderRoI::Counter_t &counter) const
SG::ReadHandleKey< xAOD::EventInfo > m_evtKey
bool passEtaDepCuts(const Trk::Track *track, int nClusters, int nFreeClusters, int nPixels) const
apply eta-dependent selections
ToolHandle< ISiSpacePointsSeedMaker > m_seedsmaker
SiSPSeededTrackFinderRoI(const std::string &name, ISvcLocator *pSvcLocator)
MsgStream & dumptools(MsgStream &out) const
SG::WriteHandleKey< TrackCollection > m_outputTracksKey
ServiceHandle< IInDetEtaDependentCutsSvc > m_etaDependentCutsSvc
Trk::MagneticFieldProperties m_fieldprop
SG::ReadHandleKey< Trk::PRDtoTrackMap > m_prdToTrackMap
ToolHandle< ISiTrackMaker > m_trackmaker
virtual StatusCode execute(const EventContext &ctx) const override
MsgStream & dump(MSG::Level lvl, const SiSPSeededTrackFinderRoI::Counter_t *) const
ToolHandle< IZWindowRoISeedTool > m_RandomRoISeedTool
SG::ReadCondHandleKey< InDet::BeamSpotData > m_beamSpotKey
ToolHandle< IZWindowRoISeedTool > m_ZWindowRoISeedTool
virtual StatusCode finalize() override
InDet::SiSpacePointsSeedMakerEventData holds event dependent data used by ISiSpacePointsSeedMaker.
InDet::SiTrackMakerEventData_xk holds event dependent data used by ISiTrackMaker.
void setPRDtoTrackMap(const Trk::PRDtoTrackMap *prd_to_track_map)
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
Property holding a SG store/key/clid from which a ReadHandle is made.
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.
int numberDoF() const
returns the number of degrees of freedom of the overall track or vertex fit as integer
Definition FitQuality.h:60
double chiSquared() const
returns the of the overall track fit
Definition FitQuality.h:56
magnetic field properties to steer the behavior of the extrapolation
This class is the pure abstract base class for all fittable tracking measurements.
Class describing the Line to which the Perigee refers to.
Class to handle RIO On Tracks ROT) for InDet and Muons, it inherits from the common MeasurementBase.
Definition RIO_OnTrack.h:70
virtual const Trk::PrepRawData * prepRawData() const =0
returns the PrepRawData (also known as RIO) object to which this RIO_OnTrack is associated.
bool trackProperties(const TrackProperties &property) const
Access methods for track properties.
@ BremFit
A brem fit was performed on this track.
represents the track state (measurement, material, fit parameters and quality) at a surface.
@ Measurement
This is a measurement, and will at least contain a Trk::MeasurementBase.
@ Perigee
This represents a perigee, and so will contain a Perigee object only.
const Trk::TrackStates * trackStateOnSurfaces() const
return a pointer to a const DataVector of const TrackStateOnSurfaces.
const TrackInfo & info() const
Returns a const ref to info of a const tracks.
int r
Definition globals.cxx:22
Primary Vertex Finder.
@ FastField
call the fast field access method of the FieldSvc
@ NoField
Field is set to 0., 0., 0.,.
ParametersBase< TrackParametersDim, Charged > TrackParameters
-event-from-file
Vertex_v1 Vertex
Define the latest version of the vertex class.