ATLAS Offline Software
Loading...
Searching...
No Matches
InDetConversionFinderTools.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5/***************************************************************************
6 InDetConversionFinderTools.cxx - Description
7 -------------------
8 begin : 28-08-2006
9 authors : Tatjana Lenz
10 email : tatjana.lenz@cern.ch
11 changes :
12***************************************************************************/
14
15#include "AthLinks/ElementLink.h"
18
19#include <vector>
20#include <utility>
21
22namespace InDet
23{
24
26 const std::string& n,
27 const IInterface* p)
28 : AthAlgTool(t, n, p)
29{
30 declareInterface<IVertexFinder>(this);
31}
32
34
36 {
37 StatusCode sc = AthAlgTool::initialize();
38 if (sc.isFailure()) {
39 ATH_MSG_FATAL("Unable to initialize InDetConversionFinderTools");
40 return StatusCode::FAILURE;
41 }
42 /* Get the right vertex fitting tool from ToolSvc */
43 if (m_iVertexFitter.retrieve().isFailure()) {
44 ATH_MSG_FATAL("Failed to retrieve tool " << m_iVertexFitter);
45 return StatusCode::FAILURE;
46 }
47 ATH_MSG_DEBUG("Retrieved tool " << m_iVertexFitter);
48
49 /* Get the track pairs selector tool from ToolSvc */
50 if (m_trackPairsSelector.retrieve().isFailure()) {
51 ATH_MSG_FATAL("Failed to retrieve tool " << m_trackPairsSelector);
52 return StatusCode::FAILURE;
53 }
54 ATH_MSG_DEBUG("Retrieved tool " << m_trackPairsSelector);
55
56 /* Get the vertex point estimator tool from ToolSvc */
57 if (m_vertexEstimator.retrieve().isFailure()) {
58 ATH_MSG_FATAL("Failed to retrieve tool " << m_vertexEstimator);
59 return StatusCode::FAILURE;
60 }
61 ATH_MSG_DEBUG("Retrieved tool " << m_vertexEstimator);
62
63 /* Get the postselector tool from ToolSvc */
64 if (m_postSelector.retrieve().isFailure()) {
65 ATH_MSG_FATAL("Failed to retrieve tool " << m_postSelector);
66 return StatusCode::FAILURE;
67 }
68 ATH_MSG_DEBUG("Retrieved tool " << m_postSelector);
69
70 /* Get the single track conversion tool from ToolSvc */
71 if (m_singleTrkConvTool.retrieve().isFailure()) {
72 ATH_MSG_FATAL("Failed to retrieve tool " << m_singleTrkConvTool);
73 return StatusCode::FAILURE;
74 }
75 ATH_MSG_DEBUG("Retrieved tool " << m_singleTrkConvTool);
76
77 /* Get the track selector tool from ToolSvc */
78 if (m_trkSelector.retrieve().isFailure()) {
79 ATH_MSG_FATAL("Failed to retrieve tool " << m_trkSelector);
80 return StatusCode::FAILURE;
81 }
82 ATH_MSG_DEBUG("Retrieved tool " << m_trkSelector);
83
84 ATH_MSG_DEBUG("Initialization successful");
85 return StatusCode::SUCCESS;
86 }
87
89 {
90 return StatusCode::SUCCESS;
91 }
92
93 //TrackCollection
94 std::pair<xAOD::VertexContainer*, xAOD::VertexAuxContainer*>
95 InDetConversionFinderTools::findVertex(const EventContext& /*ctx*/,
96 const TrackCollection* /*trk_coll*/) const
97 {
98
99 ATH_MSG_ERROR("Using Track Container not currently supported returning an empty conatiner");
100
101 // Make collection for conversions.
102 xAOD::VertexContainer* InDetConversionContainer = new xAOD::VertexContainer();
103 xAOD::VertexAuxContainer* InDetConversionContainerAux = new xAOD::VertexAuxContainer();
104 InDetConversionContainer->setStore( InDetConversionContainerAux );
105
106 return std::make_pair(InDetConversionContainer,InDetConversionContainerAux);
107 }
108
109 // TrackParticle Collection
110 std::pair<xAOD::VertexContainer*, xAOD::VertexAuxContainer*>
112 const xAOD::TrackParticleContainer* trk_coll) const
113 {
114 // Make collection for conversions.
115 xAOD::VertexContainer* InDetConversionContainer =
117 xAOD::VertexAuxContainer* InDetConversionContainerAux =
119 InDetConversionContainer->setStore(InDetConversionContainerAux);
120 // Count for number of successful conversions
121 int numConversions = 0;
122
123 // have to be used for the vertex fit
124 Amg::Vector3D pos(0., 0., 0.);
125 Amg::Vector3D initPos(0., 0., 0.);
126
127 // Make seperate lists of positive and negative tracks (after presection
128 // cuts)
129 std::vector<const xAOD::TrackParticle*> negSelectedTracks;
130 negSelectedTracks.clear();
131 std::vector<const xAOD::TrackParticle*> posSelectedTracks;
132 posSelectedTracks.clear();
133 std::vector<int> negIndx;
134 std::vector<int> posIndx;
135
137 // track preselection -->pt-cut
139 for (iter = (*trk_coll).begin(); iter != (*trk_coll).end(); ++iter) {
140 if (m_trkSelector->decision(
141 **iter,
142 nullptr)) { // Only save if track passes the pt, d0, z0 and TRT PID cuts
143 ATH_MSG_DEBUG("Track passed preselection");
144 if ((*iter)->charge() < 0) {
145 negSelectedTracks.push_back(*iter);
146 negIndx.push_back(0);
147 } else {
148 posSelectedTracks.push_back(*iter);
149 posIndx.push_back(0);
150 }
151 } else
152 ATH_MSG_DEBUG("Track failed preselection");
153 } // end pt,d0.z0-cuts
154
155 // Make track pairs. To be used for double leg conversions
156 std::vector<const xAOD::TrackParticle*>::const_iterator iter_pos;
157 std::vector<const xAOD::TrackParticle*>::const_iterator iter_neg;
158 std::vector<Amg::Vector3D> positionList;
159 positionList.clear();
160 std::vector<const xAOD::TrackParticle*> trackParticleList;
161 trackParticleList.clear();
162 std::vector<const xAOD::TrackParticle*> singleTrackConvList;
163 singleTrackConvList.clear();
164 int ipos = -1;
165 int ineg = -1;
166 // Outer loop: Loop over positive tracks
167 for (iter_pos = posSelectedTracks.begin();
168 iter_pos != posSelectedTracks.end();
169 ++iter_pos) {
170 ineg = -1;
171 ipos++;
172 // Inner loop: Loop over negative tracks
173 for (iter_neg = negSelectedTracks.begin();
174 iter_neg != negSelectedTracks.end();
175 ++iter_neg) {
176 ineg++;
177 int flag = 0;
178
179 std::map<std::string, float> intersectionDecors;
180 if (!passPreSelection(cache,
181 *iter_pos,
182 *iter_neg,
183 positionList,
184 initPos,
185 flag,
186 intersectionDecors)) {
187 positionList.clear();
188 continue;
189 }
190
191 // Do the fit
192 if (positionList.size() < 2) {
193 ATH_MSG_DEBUG("No tracks to fit ");
194 positionList.clear();
195 continue;
196 }
197
198 trackParticleList.push_back(*iter_pos);
199 trackParticleList.push_back(*iter_neg);
200
201 std::unique_ptr<xAOD::Vertex> myVertex =
202 m_iVertexFitter->fit(ctx, trackParticleList, initPos);
203 trackParticleList.clear();
204
205 // We have a new vertex
206 if (myVertex) {
207 ATH_MSG_DEBUG("VertexFit successful!");
208 int type = -1;
209 if ((m_isConversion && m_postSelector->selectConversionCandidate(
210 myVertex.get(), flag, positionList)) ||
211 (!m_isConversion &&
212 m_postSelector->selectSecVtxCandidate(
213 myVertex.get(), flag, positionList, type))) {
214
215 ATH_MSG_DEBUG(" Conversion passed postselection cuts");
216 // Remove old element links
217 myVertex->clearTracks();
218
219 // If we do not have a valid type just reset
220 if (!m_isConversion && !(type == 101) && !(type == 110) &&
221 !(type == 11)) {
222 myVertex.reset();
223 }
224
225 // If we have the right type (not reset above) lets fill information
226 // and then push to the containers
227 if (myVertex) {
228 if (m_decorateVertices) {
230 "Decorating vertex with values used in track pair selector");
231 for (const auto& kv :
233 SG::Accessor<float> acc (kv.first);
234 acc (*myVertex) = kv.second;
235 }
236 ATH_MSG_DEBUG("Decorating vertex with values used in vertex "
237 "point estimator");
238 for (const auto& kv : intersectionDecors) {
239 SG::Accessor<float> acc (kv.first);
240 acc (*myVertex) = kv.second;
241 }
242 }
243
245 *trk_coll);
247 *trk_coll);
248 myVertex->addTrackAtVertex(newLinkPos);
249 myVertex->addTrackAtVertex(newLinkNeg);
250
251 // Now fill in the containers depending on the 2 possible
252 // cases
253 if (m_isConversion) {
254 myVertex->setVertexType(xAOD::VxType::ConvVtx);
255 InDetConversionContainer->push_back(std::move(myVertex));
256 }
257 else if (type == 101 || type == 110 || type == 11) { // V0
258 myVertex->setVertexType(xAOD::VxType::V0Vtx);
259 InDetConversionContainer->push_back(std::move(myVertex));
260 }
261
262 } // End if on right type
263
264 negIndx[ineg] = 1;
265 posIndx[ipos] = 1;
266 ++numConversions;
267 } else {
268 ATH_MSG_DEBUG("VxCandidate failed the post selection cuts!");
269 myVertex.reset();
270 }
271 } else {
272 ATH_MSG_DEBUG("VertexFit was NOT successful!");
273 }
274 positionList.clear();
275 } // neg loop
276 } // pos loop
277 ATH_MSG_DEBUG("Number of conversions found passing post selection cuts: "
278 << numConversions);
279
280 // single track conversions
281 if (m_isConversion) {
282 for (int ip = 0; ip < int(posIndx.size()); ++ip) {
283 if (posIndx[ip] == 0)
284 singleTrackConvList.push_back(posSelectedTracks[ip]);
285 }
286 for (int in = 0; in < int(negIndx.size()); ++in) {
287 if (negIndx[in] == 0)
288 singleTrackConvList.push_back(negSelectedTracks[in]);
289 }
290
291 std::vector<const xAOD::TrackParticle*>::iterator itk;
292 std::vector<const xAOD::TrackParticle*>::iterator itke = singleTrackConvList.end();
293 int numSingle = 0;
294 for (itk = singleTrackConvList.begin(); itk != itke; ++itk) {
295 if (!m_singleTrkConvTool->selectSingleTrackParticleConversion((*itk)))
296 ATH_MSG_DEBUG("Track failed single track conversion selection");
297 else {
298 xAOD::Vertex* sConver(nullptr);
299 sConver = m_singleTrkConvTool->buildSingleTrackParticleConversion(
300 (*itk), InDetConversionContainer);
301 if (sConver) {
302 sConver->clearTracks();
303
305 newLink.toContainedElement(*trk_coll, *itk);
306 sConver->addTrackAtVertex(newLink);
308 numSingle++;
309
310 if (m_decorateVertices) {
311 ATH_MSG_DEBUG("Decorating single track vertex with dummy values "
312 "used in track pair selector");
313 for (const auto& kv : InDet::TrackPairsSelector::getLastValues(cache)) {
314 SG::Accessor<float> acc (kv.first);
315 acc (*sConver) = 0.;
316 }
317
318 ATH_MSG_DEBUG("Decorating single track vertex with dummy values "
319 "used in vertex point estimator");
320 for (const std::string& k : InDet::VertexPointEstimator::decorKeys()) {
321 SG::Accessor<float> acc (k);
322 acc (*sConver) = 0.;
323 }
324
325 ATH_MSG_DEBUG("Decorating single track vertex with dummy values "
326 "used in post selector");
327 InDet::ConversionPostSelector::decorateVertex(*sConver, 0., 0., 0., 0., 0.);
328 }
329 }
330 }
331 }
332 ATH_MSG_DEBUG("Number successful reconstructed single track conversion: "
333 << numSingle);
334 }
335
336 return std::make_pair(InDetConversionContainer,
337 InDetConversionContainerAux);
338 }
339
340 bool
343 const xAOD::TrackParticle* track_pos,
344 const xAOD::TrackParticle* track_neg,
345 std::vector<Amg::Vector3D>& trackList,
346 Amg::Vector3D& initPos,
347 int& flag,
348 std::map<std::string, float>& intersectionDecors) const
349 {
350 //Track summary information
351
352 uint8_t nclusPos(0);
353 uint8_t dummy(0);
354 if(track_pos->summaryValue(dummy,xAOD::numberOfSCTHits)){
355 nclusPos += dummy;
356 }
357 if(track_pos->summaryValue(dummy,xAOD::numberOfPixelHits)){
358 nclusPos += dummy;
359 }
360
361 uint8_t nclusNeg(0);
362 if(track_neg->summaryValue(dummy,xAOD::numberOfSCTHits)){
363 nclusNeg += dummy;
364 }
365 if(track_neg->summaryValue(dummy,xAOD::numberOfPixelHits)){
366 nclusNeg += dummy;
367 }
368
369
370 if(nclusNeg>0 && nclusPos>0) flag= 0;
371 if((nclusNeg>0 && nclusPos==0) || (nclusNeg==0 && nclusPos>0)) flag = 1;
372 if(nclusNeg==0 && nclusPos==0) flag = 2;
373 if(m_removeTrt && (flag==1 || flag==2)) return false;
374
375 if (m_trackPairsSelector->selectTrackParticlePair( track_pos,track_neg,cache)){
376
377 const Trk::Perigee& perPos = track_pos->perigeeParameters();
378 const Trk::Perigee& perNeg = track_neg->perigeeParameters();
379 int errorcode = 0;
380 Amg::Vector3D startingPoint(
381 m_vertexEstimator->getCirclesIntersectionPoint(
382 &perPos, &perNeg, flag, errorcode, intersectionDecors));
383 if(m_isConversion && errorcode != 0) return false;
384 if(!m_isConversion){
385 Amg::Vector3D v_direction = perPos.momentum() + perNeg.momentum();
386 double d_alpha = (startingPoint.adjoint()*v_direction)[0]/(startingPoint.mag()*v_direction.mag());
387 if(d_alpha<m_MinFlightAngle) return false;
388 }
389 initPos = startingPoint;
390 double startingPointR = startingPoint.perp();
391 if(startingPointR<800.) {
392 //Get measured track parameters at first track measurement.
393
394 //Position of first hit in track particle
395 Trk::CurvilinearParameters trkPar_pos;
396 Trk::CurvilinearParameters trkPar_neg;
397
398 int index(-1);
399 for(unsigned int i(0); i< track_pos->numberOfParameters() ; ++i ){
400 if( xAOD::FirstMeasurement == track_pos->parameterPosition(i) ){
401 index = i;
402 break;
403 }
404 }
405 if(index!=-1){
406 trkPar_pos = track_pos->curvilinearParameters(index);
407 } else {
408 ATH_MSG_WARNING("Track Particle does not contain first Measurement track parameters");
409 return false;
410 }
411
412 index = -1;
413 for(unsigned int i(0); i< track_neg->numberOfParameters() ; ++i ){
414 if( xAOD::FirstMeasurement == track_neg->parameterPosition(i) ){
415 index = i;
416 break;
417 }
418 }
419 if(index!=-1){
420 trkPar_neg = track_neg->curvilinearParameters(index);
421 } else {
422 ATH_MSG_WARNING("Track Particle does not contain first Measurement track parameters");
423 return false;
424 }
425
426
427 if(!m_isConversion){
428 double posR = trkPar_pos.position().perp();
429 double negR = trkPar_neg.position().perp();
430 double diffR = 1000.;
431 if((startingPointR-posR)<(startingPointR-negR)) diffR = startingPointR-posR;
432 else diffR = startingPointR-negR;
433 if(startingPointR<m_MinInitVtxR) return false;
434 if(diffR<m_mindR || diffR>m_maxdR) return false;
435 }
436
437 //Alternative way to estimate the guess vertex
438 initPos = trkPar_pos.position() + trkPar_neg.position();
439 initPos *= 0.5;
440 // Not correct but will do for now
441 trackList.push_back(trkPar_pos.position());
442 trackList.push_back(trkPar_neg.position());
443
444 }
445 }
446 return true;
447 }
448}// end namespace InDet
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
Helper class to provide type-safe access to aux data.
static Double_t sc
DataVector< Trk::Track > TrackCollection
This typedef represents a collection of Trk::Track objects.
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
value_type push_back(value_type pElem)
Add an element to the end of the collection.
static void decorateVertex(xAOD::Vertex &vertex, float inv_mass, float pt1, float pt2, float fR, float deltaPhiVtxTrk)
Decorate vertices with values used in post selector.
ToolHandle< Trk::ITrackSelectorTool > m_trkSelector
Track Selector Tool.
virtual std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > findVertex(const EventContext &ctx, const TrackCollection *trk_coll) const override
Find vertex from Trk::TrackCollection.
BooleanProperty m_decorateVertices
Conversion candidate reconstruction for Trk::Tracks.
bool passPreSelection(TrackPairsSelector::Cache &cache, const xAOD::TrackParticle *track_pos, const xAOD::TrackParticle *track_neg, std::vector< Amg::Vector3D > &trackList, Amg::Vector3D &initPos, int &flag, std::map< std::string, float > &intersectionDecors) const
ToolHandle< InDet::VertexPointEstimator > m_vertexEstimator
Initial conversion vertex estimator tool.
ToolHandle< InDet::TrackPairsSelector > m_trackPairsSelector
Initial conversion vertex estimator tool.
ToolHandle< Trk::IVertexFitter > m_iVertexFitter
Vertex fitter interface.
InDetConversionFinderTools(const std::string &t, const std::string &n, const IInterface *p)
ToolHandle< InDet::ConversionPostSelector > m_postSelector
Conversion post-fit selector tool.
ToolHandle< InDet::SingleTrackConversionTool > m_singleTrkConvTool
Single track conversion tool.
static std::map< std::string, float > getLastValues(const Cache &cache)
Return a map with the values calculated for the last pair to decorate the vertex once it is created.
static std::vector< std::string > decorKeys()
Return list of keys used for decorations.
Helper class to provide type-safe access to aux data.
const Amg::Vector3D & momentum() const
Access method for the momentum.
const Amg::Vector3D & position() const
Access method for the position.
const Trk::Perigee & perigeeParameters() const
Returns the Trk::MeasuredPerigee track parameters.
const Trk::CurvilinearParameters curvilinearParameters(unsigned int index) const
Returns a curvilinear representation of the parameters at 'index'.
xAOD::ParameterPosition parameterPosition(unsigned int index) const
Return the ParameterPosition of the parameters at 'index'.
bool summaryValue(uint8_t &value, const SummaryType &information) const
Accessor for TrackSummary values.
size_t numberOfParameters() const
Returns the number of additional parameters stored in the TrackParticle.
void clearTracks()
Remove all tracks from the vertex.
void addTrackAtVertex(const ElementLink< TrackParticleContainer > &tr, float weight=1.0)
Add a new track to the vertex.
void setVertexType(VxType::VertexType vType)
Set the type of the vertex.
Eigen::Matrix< double, 3, 1 > Vector3D
Primary Vertex Finder.
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
CurvilinearParametersT< TrackParametersDim, Charged, PlaneSurface > CurvilinearParameters
Definition index.py:1
@ V0Vtx
Vertex from V0 decay.
@ ConvVtx
Conversion vertex.
VertexAuxContainer_v1 VertexAuxContainer
Definition of the current jet auxiliary container.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".
@ numberOfSCTHits
number of hits in SCT [unit8_t].
@ numberOfPixelHits
these are the pixel hits, including the b-layer [unit8_t].
@ FirstMeasurement
Parameter defined at the position of the 1st measurement.