ATLAS Offline Software
Loading...
Searching...
No Matches
InDet::InDetAdaptiveMultiSecVtxFinderTool Class Reference

#include <InDetAdaptiveMultiSecVtxFinderTool.h>

Inheritance diagram for InDet::InDetAdaptiveMultiSecVtxFinderTool:
Collaboration diagram for InDet::InDetAdaptiveMultiSecVtxFinderTool:

Public Member Functions

 InDetAdaptiveMultiSecVtxFinderTool (const std::string &t, const std::string &n, const IInterface *p)
virtual ~InDetAdaptiveMultiSecVtxFinderTool ()=default
 Destructor.
StatusCode initialize () override
StatusCode finalize () override
std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > findVertex (const xAOD::TrackParticleContainer *trackParticles) override
void setPrimaryVertexPosition (double, double, double) override

Private Member Functions

std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > doVertexing (const std::vector< Trk::ITrackLink * > &trackVector)
float findCompatibleTracks (const EventContext &ctx, Amg::Vector3D &seedVertex, Trk::ITrackLink *trkIter) const
void countTracksAndNdf (xAOD::Vertex *myxAODVertex, float &ndf, int &ntracks) const
bool checkFit (xAOD::Vertex *actualCandidate) const
int removeTracksFromSeeds (xAOD::Vertex *actualCandidate, std::vector< Trk::ITrackLink * > &seedTracks) const
void removeClosestTrack (const EventContext &ctx, Amg::Vector3D &seedVertex, std::vector< Trk::ITrackLink * > &seedTracks, int &nFound) const
bool V0check (const std::vector< Amg::Vector3D > &momenta, const Amg::Vector3D &posi) const
const std::vector< Amg::Vector3DgetVertexMomenta (xAOD::Vertex *myxAODVertex) const

Private Attributes

ToolHandle< Trk::AdaptiveMultiVertexFitterm_VertexFitter {this, "VertexFitterTool", "Trk::AdaptiveMultiVertexFitter", "Multi Vertex Fitter"}
ToolHandle< InDet::IInDetTrackSelectionToolm_trkFilter {this, "BaseTrackSelector", "InDet::DetailedTrackSelectToolRelax", "Base track selection tool"}
ToolHandle< InDet::IInDetTrackSelectionToolm_SVtrkFilter {this, "SecVtxTrackSelector", "InDet::SecVtxTrackSelector", "SV track selection tool"}
ToolHandle< Trk::IVertexSeedFinderm_SeedFinder {this, "SeedFinder", "Trk::IndexedCrossDistancesSeedFinder", "Seed finder"}
ToolHandle< Trk::IImpactPoint3dEstimatorm_ImpactPoint3dEstimator {this, "ImpactPoint3dEstimator", "Trk::ImpactPoint3dEstimator", "Impact point estimator"}
FloatProperty m_privtxRef {this, "MomentumProjectionOnDirection", -999999.9, "pri vtx ref"}
DoubleProperty m_significanceCutSeeding {this, "significanceCutSeeding", 10, "significanceCutSeeding"}
DoubleProperty m_minWghtAtVtx {this, "minTrackWeightAtVtx", 0., "minTrackWeightAtVtx"}
DoubleProperty m_maxIterations {this, "maxVertices", 25, "max iterations"}
Amg::Vector3D m_privtx

Detailed Description

Definition at line 53 of file InDetAdaptiveMultiSecVtxFinderTool.h.

Constructor & Destructor Documentation

◆ InDetAdaptiveMultiSecVtxFinderTool()

InDet::InDetAdaptiveMultiSecVtxFinderTool::InDetAdaptiveMultiSecVtxFinderTool ( const std::string & t,
const std::string & n,
const IInterface * p )

Definition at line 56 of file InDetAdaptiveMultiSecVtxFinderTool.cxx.

57 :
58 base_class(t, n, p) {}

◆ ~InDetAdaptiveMultiSecVtxFinderTool()

virtual InDet::InDetAdaptiveMultiSecVtxFinderTool::~InDetAdaptiveMultiSecVtxFinderTool ( )
virtualdefault

Destructor.

Member Function Documentation

◆ checkFit()

bool InDet::InDetAdaptiveMultiSecVtxFinderTool::checkFit ( xAOD::Vertex * actualCandidate) const
private

Definition at line 484 of file InDetAdaptiveMultiSecVtxFinderTool.cxx.

484 {
485 int ntracks = 0;
486 if (not actualCandidate) return false;
487 float ndf = actualCandidate->numberDoF();
488
489 static const xAOD::Vertex::Decorator<std::vector<Trk::VxTrackAtVertex*>> VTAV("VTAV");
490
491 for (Trk::VxTrackAtVertex* trkAtVtxIter : VTAV(*actualCandidate)) {
492 if ((trkAtVtxIter)->weight() > m_minWghtAtVtx) { ntracks += 1; }
493 }
494
495 ATH_MSG_DEBUG(" xAOD::Vertex : " << (actualCandidate != nullptr ? 1 : 0) << ", #dof = " << ndf
496 << ", #tracks (weight>0.01) = " << ntracks);
497
498 return ( ndf > 0 && ntracks >= 2);
499 }
#define ATH_MSG_DEBUG(x)

◆ countTracksAndNdf()

void InDet::InDetAdaptiveMultiSecVtxFinderTool::countTracksAndNdf ( xAOD::Vertex * myxAODVertex,
float & ndf,
int & ntracks ) const
private

Definition at line 442 of file InDetAdaptiveMultiSecVtxFinderTool.cxx.

442 {
443 ndf = -3.0;
444 ntrk = 0;
445
446 if (myxAODVertex) {
447 ndf = myxAODVertex->numberDoF();
448
449 static const xAOD::Vertex::Decorator<std::vector<Trk::VxTrackAtVertex*>> VTAV("VTAV");
450
451 for (Trk::VxTrackAtVertex* trkAtVtxIter : VTAV(*myxAODVertex)) {
452 if ((trkAtVtxIter)->weight() > m_minWghtAtVtx) { ntrk += 1; }
453 }
454 }
455 }
float numberDoF() const
Returns the number of degrees of freedom of the vertex fit as float.

◆ doVertexing()

std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > InDet::InDetAdaptiveMultiSecVtxFinderTool::doVertexing ( const std::vector< Trk::ITrackLink * > & trackVector)
private

prepare decorators to hold multi vertex fit information

place all tracks in the origTrack vector and initially all in the seedTrack vector

Definition at line 122 of file InDetAdaptiveMultiSecVtxFinderTool.cxx.

122 {
123 xAOD::VertexContainer* theVertexContainer = new xAOD::VertexContainer;
124 xAOD::VertexAuxContainer* theVertexAuxContainer = new xAOD::VertexAuxContainer;
125 theVertexContainer->setStore(theVertexAuxContainer);
126
128 static const xAOD::Vertex::Decorator<Trk::MvfFitInfo*> MvfFitInfo("MvfFitInfo");
129 static const xAOD::Vertex::Decorator<bool> isInitialized("isInitialized");
130 static const xAOD::Vertex::Decorator<std::vector<Trk::VxTrackAtVertex*>> VTAV("VTAV");
131
132 const EventContext& ctx = Gaudi::Hive::currentContext();
133
134 std::vector<xAODVertex_pair> myxAODVertices;
135
136 std::vector<Trk::TrackToVtxLink*> myTrackToVtxLinks;
137
139 std::vector<Trk::ITrackLink*> origTracks = trackVector;
140 std::vector<Trk::ITrackLink*> seedTracks = trackVector;
141
142 // create a map between ITrackLink* and TrackToVtxLink*
143 std::map<Trk::ITrackLink*, Trk::TrackToVtxLink*> TrackLinkOf;
144
145 // fill vector of track parameters used for fitting
146 // create map
147
148 for (Trk::ITrackLink* trkIter : origTracks) {
149 Trk::TrackToVtxLink* newTrkToVtxLink(new Trk::TrackToVtxLink(new std::vector<xAOD::Vertex*>));
150
151 TrackLinkOf[trkIter] = newTrkToVtxLink;
152 myTrackToVtxLinks.push_back(newTrkToVtxLink);
153 }
154
155 int iteration = 0;
156 unsigned int seedtracknumber = seedTracks.size();
157
158 do {
159 if (seedtracknumber == 0) { ATH_MSG_DEBUG("New iteration. No tracks available after track selection for seeding."); }
160
161 iteration += 1;
162 ATH_MSG_DEBUG("Iteration number " << iteration << " and tracks left for seeding " << seedtracknumber);
163
164 std::vector<const Trk::TrackParameters*> perigeeList;
165
166 perigeeList.reserve(seedTracks.size());
167for (const Trk::ITrackLink* seedtrkAtVtxIter : seedTracks) { perigeeList.push_back((seedtrkAtVtxIter)->parameters()); }
168
169 ATH_MSG_DEBUG("Going to seed finder");
170
171 std::unique_ptr<Trk::IMode3dInfo> info;
172 Amg::Vector3D seedVertex;
173 seedVertex = m_SeedFinder->findSeed(m_privtx.x(), m_privtx.y(), info, perigeeList);
174
175 ATH_MSG_DEBUG("Found seed at x: " << seedVertex.x() << " at y: " << seedVertex.y() << " at z: " << seedVertex.z());
176
177 xAOD::Vertex* seededxAODVertex = new xAOD::Vertex;
178 theVertexContainer->push_back(seededxAODVertex);
179 seededxAODVertex->setPosition(seedVertex);
180 Amg::MatrixX looseConstraintCovariance(3, 3);
181 looseConstraintCovariance.setIdentity();
182 looseConstraintCovariance = looseConstraintCovariance * 1e+8;
183 seededxAODVertex->setCovariancePosition(looseConstraintCovariance);
184 seededxAODVertex->vxTrackAtVertex() = std::vector<Trk::VxTrackAtVertex>();
185 seededxAODVertex->setVertexType(xAOD::VxType::NotSpecified);
186
187 if (seedVertex.z() == 0) {
188 ATH_MSG_DEBUG("No good seed found: no further vertices in event");
189 ATH_MSG_DEBUG("Number of input tracks: " << perigeeList.size() << ", but no good seed returned");
190 break;
191 }
192
193 xAOD::Vertex* constraintVertex = nullptr;
194 looseConstraintCovariance.setIdentity();
195 looseConstraintCovariance = looseConstraintCovariance * 1e+8;
196 constraintVertex = new xAOD::Vertex();
197 constraintVertex->makePrivateStore();
198 constraintVertex->setPosition(seedVertex);
199 constraintVertex->setCovariancePosition(looseConstraintCovariance);
200 constraintVertex->setFitQuality(0., -3.);
201 constraintVertex->setVertexType(xAOD::VxType::NotSpecified);
202
203 xAOD::Vertex* actualCandidate = new xAOD::Vertex;
204 actualCandidate->makePrivateStore();
206
207 MvfFitInfo(*actualCandidate) =
208 new Trk::MvfFitInfo(constraintVertex, new Amg::Vector3D(seedVertex), new Amg::Vector3D(seedVertex));
209 isInitialized(*actualCandidate) = false;
210 std::vector<Trk::VxTrackAtVertex*> vectorOfTracks(0);
211 VTAV(*actualCandidate) = vectorOfTracks;
212
213 for (Trk::ITrackLink* trkIter : origTracks) {
214 // now fill perigeesToFit list of track parameters from origTracks
215 float doe = findCompatibleTracks(ctx, seedVertex, trkIter);
216 if (doe < m_significanceCutSeeding) {
217 Trk::TrackToVtxLink* actualLink = TrackLinkOf[trkIter];
218 std::vector<xAOD::Vertex*>* actualvtxlink = actualLink->vertices();
219 // adding vertex to candidates of track
220 actualvtxlink->push_back(actualCandidate);
221 VTAV(*actualCandidate).push_back(new Trk::MVFVxTrackAtVertex((trkIter)->clone(), actualLink));
222 }
223 }
224
225 ATH_MSG_DEBUG(" Considering n. " << VTAV(*actualCandidate).size() << " tracks for the fit. ");
226
227 if (VTAV(*actualCandidate).size() < 2) {
228 ATH_MSG_DEBUG("No tracks found near seed, while at least two tracks were expected.");
229
230 if (VTAV.isAvailable(*actualCandidate)) {
231 for (auto *tav : VTAV(*actualCandidate)) {
232 if (tav == nullptr) continue;
233
234 (static_cast<Trk::MVFVxTrackAtVertex*>(tav))->setLinkToVertices(nullptr);
235 delete tav;
236 tav = nullptr;
237 }
238 VTAV(*actualCandidate).clear();
239 }
240 if (MvfFitInfo.isAvailable(*actualCandidate) && MvfFitInfo(*actualCandidate) != nullptr) {
241 delete MvfFitInfo(*actualCandidate);
242 MvfFitInfo(*actualCandidate) = nullptr;
243 }
244 delete actualCandidate;
245 actualCandidate = nullptr;
246
247 break;
248 }
249
250 ATH_MSG_DEBUG("Going to fitter.");
251
252 m_VertexFitter->addVtxTofit(ctx, actualCandidate);
253
254 ATH_MSG_DEBUG("Deleting tracks with really good fit to vertex from seeding tracks.");
255 int nFound = removeTracksFromSeeds(actualCandidate, seedTracks);
256
257 ATH_MSG_DEBUG("Found and deleted " << nFound << " tracks from seeding tracks.");
258 if (nFound == 0) {
259 ATH_MSG_DEBUG("All tracks used for fitting came from fiting tracks, removing closest from seeding.");
260 // all the tracks used for the fit came from fitting track list
261 //-> so remove the closest track to seed from seeding, otherwise you'll keep finding the same seed position
262
263 removeClosestTrack(ctx, seedVertex, seedTracks, nFound);
264 }
265
266 if (nFound == 0) {
267 ATH_MSG_DEBUG("You still have not removed any tracks from seeds! Aborting.");
268 break;
269 }
270
271 ATH_MSG_DEBUG("Checking goodness of fit.");
272 bool goodVertex = checkFit(actualCandidate);
273
274 if (!goodVertex) {
275 ATH_MSG_DEBUG("Bad vertex, deleting the vertex and clearing all pointers");
276
277 seededxAODVertex->setVertexType(xAOD::VxType::KinkVtx);
278
279 if (actualCandidate) {
280 if (VTAV.isAvailable(*actualCandidate)) {
281 for (auto *tav : VTAV(*actualCandidate)) {
282 if (tav == nullptr) continue;
283
284 (static_cast<Trk::MVFVxTrackAtVertex*>(tav))->setLinkToVertices(nullptr);
285 delete tav;
286 tav = nullptr;
287 }
288 VTAV(*actualCandidate).clear();
289 }
290 if (MvfFitInfo.isAvailable(*actualCandidate) && MvfFitInfo(*actualCandidate) != nullptr) {
291 delete MvfFitInfo(*actualCandidate);
292 MvfFitInfo(*actualCandidate) = nullptr;
293 }
294
295 delete actualCandidate;
296 actualCandidate = nullptr;
297 }
298
299 } else {
300 ATH_MSG_DEBUG("I have found a good vertex!");
301
302 seededxAODVertex->setVertexType(xAOD::VxType::NoVtx);
303 actualCandidate->setVertexType(xAOD::VxType::SecVtx);
304 myxAODVertices.emplace_back(0, actualCandidate);
305 }
306 seedtracknumber = seedTracks.size();
307 } while (seedTracks.size() > 1 && iteration < m_maxIterations);
308
309 if (iteration >= m_maxIterations) {
310 ATH_MSG_DEBUG("Maximum number of iterations ("
311 << m_maxIterations << ") reached; to reconstruct more vertices, set maxIterations to a higher value.");
312 }
313
314 ATH_MSG_DEBUG("Secondary vertex finding complete with " << iteration << " iterations and " << myxAODVertices.size()
315 << " vertices found.");
316
317 for (const xAODVertex_pair& vtxIter : myxAODVertices) {
318 xAOD::Vertex* fittedVert = vtxIter.second;
319
320 xAOD::Vertex* cand = new xAOD::Vertex;
321 theVertexContainer->push_back(cand);
322 cand->setPosition(fittedVert->position());
323 cand->setCovariancePosition(fittedVert->covariancePosition());
324 cand->setFitQuality(fittedVert->chiSquared(), fittedVert->numberDoF());
326
327 std::vector<Trk::VxTrackAtVertex>* tracksOfVertex = &(cand->vxTrackAtVertex());
328 tracksOfVertex->clear();
329
330 for (Trk::VxTrackAtVertex* MVFtrkIter : VTAV(*fittedVert)) {
331 if ((*MVFtrkIter).initialPerigee()) { (*MVFtrkIter).setPerigeeAtVertex(((*MVFtrkIter).initialPerigee())->clone()); }
332 tracksOfVertex->push_back(*MVFtrkIter);
333 }
334 }
335
336 for (const xAODVertex_pair& vtxIter : myxAODVertices) {
337 xAOD::Vertex* cand = vtxIter.second;
338
339 for (Trk::VxTrackAtVertex* MVFtrkIter : VTAV(*cand)) {
340 (static_cast<Trk::MVFVxTrackAtVertex*>(MVFtrkIter))->setLinkToVertices(nullptr);
341 delete MVFtrkIter;
342 MVFtrkIter = nullptr;
343 }
344
345 delete MvfFitInfo(*cand);
346 }
347
348 ATH_MSG_DEBUG("Looping over vertex container");
349
350 for (xAOD::Vertex* vxIter : *theVertexContainer) {
351 std::vector<Trk::VxTrackAtVertex>* myVxTracksAtVtx = &((vxIter)->vxTrackAtVertex());
352 if (!myVxTracksAtVtx) continue;
353
354 for (Trk::VxTrackAtVertex& tracksIter : *myVxTracksAtVtx) {
355 Trk::LinkToXAODTrackParticle* linkToXAODTP = nullptr;
356 Trk::ITrackLink* tmpLink = (tracksIter).trackOrParticleLink();
357 if (tmpLink->type() == Trk::ITrackLink::ToxAODTrackParticle) {
358 linkToXAODTP = static_cast<Trk::LinkToXAODTrackParticle*>(tmpLink);
359 }
360
361 if (linkToXAODTP) { (vxIter)->addTrackAtVertex(*linkToXAODTP, (tracksIter).weight()); }
362 }
363
364 int ntrk = myVxTracksAtVtx->size();
365 if (ntrk == 2) {
366 ATH_MSG_DEBUG("Could do a V0 search");
367
368 bool isV0 = V0check(getVertexMomenta(vxIter), (&(*vxIter))->position());
369 if (isV0) {
370 ATH_MSG_DEBUG("Labeling as V0");
371 (vxIter)->setVertexType(xAOD::VxType::V0Vtx);
372 }
373 }
374 }
375
376 // delete all TrackToVtxLink objects
377 for (Trk::TrackToVtxLink* iterator : myTrackToVtxLinks) { delete iterator; }
378
379 if (!theVertexContainer->empty()) {
380 xAOD::Vertex* secVtx = theVertexContainer->front();
381 if (!secVtx->vxTrackAtVertex().empty()) {
383 xAOD::Vertex* dummyxAODVertex = new xAOD::Vertex;
384 theVertexContainer->push_back(dummyxAODVertex); // have to add vertex to container here first so it can use its aux store
385 dummyxAODVertex->setPosition(secVtx->position());
386 dummyxAODVertex->setCovariancePosition(secVtx->covariancePosition());
387 dummyxAODVertex->vxTrackAtVertex() = std::vector<Trk::VxTrackAtVertex>();
388 dummyxAODVertex->setVertexType(xAOD::VxType::NoVtx);
389 } else {
391 }
392 }
393
394 else if (theVertexContainer->empty()) {
395 xAOD::Vertex* dummyxAODVertex = new xAOD::Vertex;
396 theVertexContainer->push_back(dummyxAODVertex); // have to add vertex to container here first so it can use its aux store
397 dummyxAODVertex->setPosition(Amg::Vector3D(0, 0, 0));
398 Amg::MatrixX looseConstraintCovariance(3, 3);
399 looseConstraintCovariance.setIdentity();
400 looseConstraintCovariance = looseConstraintCovariance * 1e+8;
401 dummyxAODVertex->setCovariancePosition(looseConstraintCovariance);
402 dummyxAODVertex->vxTrackAtVertex() = std::vector<Trk::VxTrackAtVertex>();
403 dummyxAODVertex->setVertexType(xAOD::VxType::NoVtx);
404 }
405
406 int noVtx = 0;
407 int kinkVtx = 0;
408 int notSpec = 0;
409 int secVtx = 0;
410 int V0vtx = 0;
411 for (unsigned int i = 0; i < theVertexContainer->size(); i++) {
413 vtxType = static_cast<xAOD::VxType::VertexType>((*theVertexContainer)[i]->vertexType());
414 switch (vtxType) {
415 case xAOD::VxType::NoVtx: noVtx++; break;
416 case xAOD::VxType::KinkVtx: kinkVtx++; break;
417 case xAOD::VxType::NotSpecified: notSpec++; break;
418 case xAOD::VxType::V0Vtx: V0vtx++; break;
419 case xAOD::VxType::SecVtx: secVtx++; break;
420 default: ATH_MSG_DEBUG("Unfamiliar vertex type");
421 }
422
423 ATH_MSG_DEBUG(" Vtx: " << i << " x= " << (*theVertexContainer)[i]->position().x() << " y= "
424 << (*theVertexContainer)[i]->position().y() << " z= " << (*theVertexContainer)[i]->position().z()
425 << " ntracks= " << (*theVertexContainer)[i]->vxTrackAtVertex().size() << " chi2= "
426 << (*theVertexContainer)[i]->chiSquared() << " #dof = " << (*theVertexContainer)[i]->numberDoF());
427 }
428
429 ATH_MSG_DEBUG("Done finding " << theVertexContainer->size() << " vertices and cleaning the container.");
430 ATH_MSG_DEBUG("Seeds good/bad/all : " << noVtx << "/" << kinkVtx << "/" << notSpec);
431 ATH_MSG_DEBUG("'Good' secondaries : " << secVtx << " and V0: " << V0vtx);
432
433 return std::make_pair(theVertexContainer, theVertexAuxContainer);
434 }
size_t size() const
Number of registered mappings.
#define y
#define x
#define z
value_type push_back(value_type pElem)
Add an element to the end of the collection.
ToolHandle< Trk::AdaptiveMultiVertexFitter > m_VertexFitter
const std::vector< Amg::Vector3D > getVertexMomenta(xAOD::Vertex *myxAODVertex) const
int removeTracksFromSeeds(xAOD::Vertex *actualCandidate, std::vector< Trk::ITrackLink * > &seedTracks) const
bool V0check(const std::vector< Amg::Vector3D > &momenta, const Amg::Vector3D &posi) const
float findCompatibleTracks(const EventContext &ctx, Amg::Vector3D &seedVertex, Trk::ITrackLink *trkIter) const
void removeClosestTrack(const EventContext &ctx, Amg::Vector3D &seedVertex, std::vector< Trk::ITrackLink * > &seedTracks, int &nFound) const
void setCovariancePosition(const AmgSymMatrix(3)&covariancePosition)
Sets the vertex covariance matrix.
void setVertexType(VxType::VertexType vType)
Set the type of the vertex.
void setPosition(const Amg::Vector3D &position)
Sets the 3-position.
float chiSquared() const
Returns the of the vertex fit as float.
std::vector< Trk::VxTrackAtVertex > & vxTrackAtVertex()
Non-const access to the VxTrackAtVertex vector.
void setFitQuality(float chiSquared, float numberDoF)
Set the 'Fit Quality' information.
const Amg::Vector3D & position() const
Returns the 3-pos.
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Eigen::Matrix< double, 3, 1 > Vector3D
float chiSquared(const U &p)
void iterator
const Amg::Vector3D & position() const
Method to retrieve the position of the Intersection.
VertexType
Vertex types.
@ KinkVtx
Kink vertex.
@ V0Vtx
Vertex from V0 decay.
@ NotSpecified
Default value, no explicit type set.
@ SecVtx
Secondary vertex.
@ NoVtx
Dummy vertex. TrackParticle was not used in vertex fit.
VertexAuxContainer_v1 VertexAuxContainer
Definition of the current jet auxiliary container.
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.

◆ finalize()

StatusCode InDet::InDetAdaptiveMultiSecVtxFinderTool::finalize ( )
override

Definition at line 440 of file InDetAdaptiveMultiSecVtxFinderTool.cxx.

440{ return StatusCode::SUCCESS; }

◆ findCompatibleTracks()

float InDet::InDetAdaptiveMultiSecVtxFinderTool::findCompatibleTracks ( const EventContext & ctx,
Amg::Vector3D & seedVertex,
Trk::ITrackLink * trkIter ) const
private

Definition at line 540 of file InDetAdaptiveMultiSecVtxFinderTool.cxx.

540 {
541 double distance = 0.;
542
543 try {
544 std::unique_ptr<Trk::PlaneSurface> mySurface =
545 m_ImpactPoint3dEstimator->Estimate3dIP(ctx, (trkIter)->parameters(), &seedVertex, distance);
546 ATH_MSG_VERBOSE(" ImpactPoint3dEstimator done ");
547 } catch (error::ImpactPoint3dEstimatorProblem err) {
548 ATH_MSG_DEBUG(" ImpactPoint3dEstimator failed to find minimum distance between track and vertex seed: " << err.p);
549 }
550
551 if (distance < 0) { ATH_MSG_DEBUG(" Distance between track and seed vtx is negative: " << distance); }
552
553 const Trk::TrackParameters* myPerigee = ((trkIter)->parameters());
554 double doe = 99999999.9;
555 double error = 0.;
556
557 if (myPerigee && myPerigee->covariance()) {
558 error = std::sqrt((*myPerigee->covariance())(Trk::d0, Trk::d0) + (*myPerigee->covariance())(Trk::z0, Trk::z0));
559 } // end of the security check
560
561 if (error == 0.) {
562 ATH_MSG_ERROR(" Error is zero! " << distance);
563 error = 1.;
564 }
565
566 doe = distance / error;
567
568 ATH_MSG_VERBOSE("Distance between track and seed vtx: " << distance << " d/s(d) = " << distance / error << " err " << error);
569
570 return doe;
571 }
#define ATH_MSG_ERROR(x)
#define ATH_MSG_VERBOSE(x)
ToolHandle< Trk::IImpactPoint3dEstimator > m_ImpactPoint3dEstimator
float distance(const Amg::Vector3D &p1, const Amg::Vector3D &p2)
calculates the distance between two point in 3D space
@ d0
Definition ParamDefs.h:63
@ z0
Definition ParamDefs.h:64
ParametersBase< TrackParametersDim, Charged > TrackParameters

◆ findVertex()

std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > InDet::InDetAdaptiveMultiSecVtxFinderTool::findVertex ( const xAOD::TrackParticleContainer * trackParticles)
override

Definition at line 81 of file InDetAdaptiveMultiSecVtxFinderTool.cxx.

82 {
83 ATH_MSG_DEBUG(" Number of input tracks before track selection: " << trackParticles->size());
84
85 std::vector<Trk::ITrackLink*> selectedTracks;
86 bool selectionPassed;
87 xAOD::TrackParticle::Decorator<bool> trackPass("TrackPassedSelection");
88
89 for (const xAOD::TrackParticle* itr : *trackParticles) {
90 xAOD::Vertex null;
91 null.makePrivateStore();
92 null.setPosition(Amg::Vector3D(0, 0, 0));
93 AmgSymMatrix(3) vertexError;
94 vertexError.setZero();
95 null.setCovariancePosition(vertexError);
96 selectionPassed = static_cast<bool>(m_trkFilter->accept(*itr, &null));
97 if (selectionPassed) selectionPassed = static_cast<bool>(m_SVtrkFilter->accept(*itr, &null));
98
99 if (selectionPassed) {
100 trackPass(*itr) = true;
101
102 Amg::VectorX par = (itr)->definingParameters();
103 par[0] = (itr)->hitPattern();
104
105 ElementLink<xAOD::TrackParticleContainer> link;
106 link.setElement(itr);
107 Trk::LinkToXAODTrackParticle* linkTT = new Trk::LinkToXAODTrackParticle(link);
108
109 linkTT->setStorableObject(*trackParticles);
110 selectedTracks.push_back(linkTT);
111 }
112 }
113
114 ATH_MSG_DEBUG("Of " << trackParticles->size() << " tracks " << selectedTracks.size() << " survived the preselection.");
115
116 std::pair<xAOD::VertexContainer*, xAOD::VertexAuxContainer*> returnContainers = doVertexing(selectedTracks);
117
118 return returnContainers;
119 }
#define AmgSymMatrix(dim)
if(pathvar)
size_type size() const noexcept
Returns the number of elements in the collection.
ToolHandle< InDet::IInDetTrackSelectionTool > m_trkFilter
std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > doVertexing(const std::vector< Trk::ITrackLink * > &trackVector)
ToolHandle< InDet::IInDetTrackSelectionTool > m_SVtrkFilter
Eigen::Matrix< double, Eigen::Dynamic, 1 > VectorX
Dynamic Vector - dynamic allocation.
pointer & link(pointer p) const
Return a reference to the link for an element.
TrackParticle_v1 TrackParticle
Reference the current persistent version:

◆ getVertexMomenta()

const std::vector< Amg::Vector3D > InDet::InDetAdaptiveMultiSecVtxFinderTool::getVertexMomenta ( xAOD::Vertex * myxAODVertex) const
private

Definition at line 457 of file InDetAdaptiveMultiSecVtxFinderTool.cxx.

457 {
458 std::vector<Amg::Vector3D> TrkAtVtxMomenta;
459
460 std::vector<Trk::VxTrackAtVertex>* tracksAtVertex = &(myxAODVertex->vxTrackAtVertex());
461
462 ATH_MSG_DEBUG(" getVertexMomenta ... #Tracks associated at vertex : " << tracksAtVertex->size());
463
464 for (const Trk::VxTrackAtVertex& tracksAtVertexIter : *tracksAtVertex) {
465 if ((tracksAtVertexIter).weight() <= m_minWghtAtVtx) continue;
466 {
467 const Trk::TrackParameters* sv_perigee = (tracksAtVertexIter).perigeeAtVertex();
468 if (!sv_perigee) {
469 ATH_MSG_DEBUG("perigeeAtVertex not available!!");
470 continue;
471 }
472
473 double qp = 1. / (std::fabs(sv_perigee->parameters()[Trk::qOverP]));
474 double theta = sv_perigee->parameters()[Trk::theta];
475 double phi = sv_perigee->parameters()[Trk::phi];
476
477 TrkAtVtxMomenta.emplace_back(qp * sin(theta) * cos(phi), qp * sin(theta) * sin(phi), qp * cos(theta));
478 }
479 }
480
481 return TrkAtVtxMomenta;
482 }
Scalar phi() const
phi method
Scalar theta() const
theta method
@ theta
Definition ParamDefs.h:66
@ qOverP
perigee
Definition ParamDefs.h:67
@ phi
Definition ParamDefs.h:75

◆ initialize()

StatusCode InDet::InDetAdaptiveMultiSecVtxFinderTool::initialize ( )
override

Definition at line 60 of file InDetAdaptiveMultiSecVtxFinderTool.cxx.

60 {
61 /* Get the right vertex fitting tool from ToolSvc */
62 ATH_CHECK(m_VertexFitter.retrieve());
63 ATH_CHECK(m_SeedFinder.retrieve());
65 ATH_CHECK(m_trkFilter.retrieve());
66
67 ATH_MSG_DEBUG("Initialization successful");
68 return StatusCode::SUCCESS;
69 }
#define ATH_CHECK
Evaluate an expression and check for errors.

◆ removeClosestTrack()

void InDet::InDetAdaptiveMultiSecVtxFinderTool::removeClosestTrack ( const EventContext & ctx,
Amg::Vector3D & seedVertex,
std::vector< Trk::ITrackLink * > & seedTracks,
int & nFound ) const
private

Definition at line 573 of file InDetAdaptiveMultiSecVtxFinderTool.cxx.

575 {
576 const Trk::ITrackLink* nearestTrack = nullptr;
577 double dist = 1e8;
578
579 for (Trk::ITrackLink* trkIter : seedTracks) {
580 double distance = 0.;
581 try {
582 std::unique_ptr<Trk::PlaneSurface> mySurface =
583 m_ImpactPoint3dEstimator->Estimate3dIP(ctx, (trkIter)->parameters(), &seedVertex, distance);
584 } catch (error::ImpactPoint3dEstimatorProblem err) {
585 ATH_MSG_DEBUG(" ImpactPoint3dEstimator failed to find minimum distance between this track and vertex seed: " << err.p);
586 }
587 ATH_MSG_DEBUG("Seed to track dist: " << distance);
588 if (distance < 0) { ATH_MSG_DEBUG("Distance was negative!"); }
589
590 if (distance > 0 && !nearestTrack) {
591 dist = distance;
592 nearestTrack = trkIter;
593 }
594 if (distance > 0 && distance < dist) {
595 dist = distance;
596 nearestTrack = trkIter;
597 }
598 }
599 if (nearestTrack) {
600 ATH_MSG_DEBUG("Found closest track to seed and deleting.");
601 std::vector<Trk::ITrackLink*>::iterator seedBegin = seedTracks.begin();
602 std::vector<Trk::ITrackLink*>::iterator seedEnd = seedTracks.end();
603 nFound += 1;
604 std::vector<Trk::ITrackLink*>::iterator foundTrack = std::find(seedBegin, seedEnd, nearestTrack);
605 if (foundTrack != seedEnd) {
606 seedTracks.erase(foundTrack);
607 seedBegin = seedTracks.begin();
608 seedEnd = seedTracks.end();
609 } else {
610 ATH_MSG_DEBUG("The nearest track was not found!");
611 }
612 } else {
613 ATH_MSG_DEBUG("What else can I try?");
614 }
615 }

◆ removeTracksFromSeeds()

int InDet::InDetAdaptiveMultiSecVtxFinderTool::removeTracksFromSeeds ( xAOD::Vertex * actualCandidate,
std::vector< Trk::ITrackLink * > & seedTracks ) const
private

Definition at line 501 of file InDetAdaptiveMultiSecVtxFinderTool.cxx.

502 {
503 if (not actualCandidate) return 0;
504 static const xAOD::Vertex::Decorator<std::vector<Trk::VxTrackAtVertex*>> VTAV("VTAV");
505
506 std::vector<Trk::ITrackLink*>::iterator seedBegin = seedTracks.begin();
507 std::vector<Trk::ITrackLink*>::iterator seedEnd = seedTracks.end();
508
509 bool goodVertex = checkFit(actualCandidate);
510
511 int nFound = 0;
512
513 for (Trk::VxTrackAtVertex* trkAtVtxIter : VTAV(*actualCandidate)) {
514 // delete the pointer to this vertex if the vertex was bad
515 if (!goodVertex) {
516 (static_cast<Trk::MVFVxTrackAtVertex*>(trkAtVtxIter))->linkToVertices()->vertices()->pop_back();
517 }
518
519 std::vector<Trk::ITrackLink*>::iterator foundTrack = seedEnd;
520 for (std::vector<Trk::ITrackLink*>::iterator seedtrkiter = seedBegin; seedtrkiter != seedEnd; ++seedtrkiter) {
521 if ((*seedtrkiter)->parameters() == (trkAtVtxIter)->trackOrParticleLink()->parameters() &&
522 (trkAtVtxIter)->weight() > m_minWghtAtVtx) {
523 foundTrack = seedtrkiter;
524 }
525 }
526
527 if (foundTrack != seedEnd) {
528 seedTracks.erase(foundTrack);
529
530 nFound += 1;
531
532 seedBegin = seedTracks.begin();
533 seedEnd = seedTracks.end();
534 }
535 }
536
537 return nFound;
538 }

◆ setPrimaryVertexPosition()

void InDet::InDetAdaptiveMultiSecVtxFinderTool::setPrimaryVertexPosition ( double vx,
double vy,
double vz )
override

Definition at line 436 of file InDetAdaptiveMultiSecVtxFinderTool.cxx.

436 {
437 m_privtx = Amg::Vector3D(vx, vy, vz);
438 }

◆ V0check()

bool InDet::InDetAdaptiveMultiSecVtxFinderTool::V0check ( const std::vector< Amg::Vector3D > & momenta,
const Amg::Vector3D & posi ) const
private

Definition at line 617 of file InDetAdaptiveMultiSecVtxFinderTool.cxx.

617 {
618 int ntrk = momenta.size();
619
620 if (ntrk < 2) {
621 ATH_MSG_DEBUG(" ntrk < 2 , Meaningless to test mass ");
622 return false;
623 }
624
625 std::vector<double> Pv(ntrk);
626 double vx = 0., vy = 0., vz = 0., eK0 = 0.;
628
629 for (int t = 0; t < ntrk; t++) {
630 Amg::Vector3D trk = momenta[t];
631
632 vz += trk.z();
633 vx += trk.x();
634 vy += trk.y();
635 Pv[t] = trk.x() * trk.x() + trk.y() * trk.y() + trk.z() * trk.z();
636 eK0 += std::sqrt(Pv[t] + pi2);
637 }
638
639 double mnt2 = vx * vx + vy * vy + vz * vz;
640 double mass = eK0 * eK0 - mnt2;
641 mass = 0.001 * (std::sqrt(std::abs(mass)));
642
643 Amg::Vector3D vdif = posi - m_privtx;
644 Amg::Vector3D vmoment = Amg::Vector3D(vx, vy, vz);
645
646 double modir = vmoment.dot(vdif) / std::sqrt(mnt2);
647
648 // borrowed from InnerDetector/InDetRecAlgs/InDetV0Finder/InDetV0FinderTool
649 double a0z = (vdif + vmoment * vmoment.dot(vdif) / (mnt2 + 0.00001)).z();
650 double Rxy = vdif.perp();
651
652 ATH_MSG_DEBUG(" V0kine : a0z = " << a0z << " Rxy = " << Rxy << " direction " << modir);
653
654 if (ntrk != 2) {
655 ATH_MSG_DEBUG(" ntrk != 2 , Meaningless to test V0 ");
656 return false;
657 }
658
659 if (a0z > 15. || Rxy > 500.) { return false; }
660
661 // 1 eV^(-1) of time = hbar / eV = 6.582173*10^(-16) second, for energy-time in natural unit
662 // double planck = 6.582173 ;
663
665 double mGam = eGam * eGam - mnt2;
666
668 double eLam = Pv[0] > Pv[1] ? std::sqrt(Pv[0] + prtn2) + std::sqrt(Pv[1] + pi2) : std::sqrt(Pv[0] + pi2) + std::sqrt(Pv[1] + prtn2);
669 double mLam = eLam * eLam - mnt2;
670
671 ATH_MSG_DEBUG(" V0 masses : " << mass << " " << std::sqrt(std::abs(mGam)) << " " << std::sqrt(std::abs(mLam)));
672
673 return ((fabs(mass - ParticleConstants::KZeroMassInMeV) < 100.) // K short
674 || (mGam > 0 && std::sqrt(mGam) < 40.) // gamma conversion ;
675 || (mLam > 0 && std::abs(std::sqrt(mLam) - ParticleConstants::lambdaMassInMeV) < 200.) // Lambda
676 );
677 }
constexpr double protonMassInMeV
the mass of the proton (in MeV)
constexpr double KZeroMassInMeV
the mass of the neutral kaon (K0) (in MeV)
constexpr double chargedPionMassInMeV
the mass of the charged pion (in MeV)
constexpr double electronMassInMeV
the mass of the electron (in MeV)
constexpr double lambdaMassInMeV
the mass of the lambda baryon (in MeV)

Member Data Documentation

◆ m_ImpactPoint3dEstimator

ToolHandle<Trk::IImpactPoint3dEstimator> InDet::InDetAdaptiveMultiSecVtxFinderTool::m_ImpactPoint3dEstimator {this, "ImpactPoint3dEstimator", "Trk::ImpactPoint3dEstimator", "Impact point estimator"}
private

Definition at line 95 of file InDetAdaptiveMultiSecVtxFinderTool.h.

95{this, "ImpactPoint3dEstimator", "Trk::ImpactPoint3dEstimator", "Impact point estimator"};

◆ m_maxIterations

DoubleProperty InDet::InDetAdaptiveMultiSecVtxFinderTool::m_maxIterations {this, "maxVertices", 25, "max iterations"}
private

Definition at line 101 of file InDetAdaptiveMultiSecVtxFinderTool.h.

101{this, "maxVertices", 25, "max iterations"};

◆ m_minWghtAtVtx

DoubleProperty InDet::InDetAdaptiveMultiSecVtxFinderTool::m_minWghtAtVtx {this, "minTrackWeightAtVtx", 0., "minTrackWeightAtVtx"}
private

Definition at line 100 of file InDetAdaptiveMultiSecVtxFinderTool.h.

100{this, "minTrackWeightAtVtx", 0., "minTrackWeightAtVtx"};

◆ m_privtx

Amg::Vector3D InDet::InDetAdaptiveMultiSecVtxFinderTool::m_privtx
private

Definition at line 102 of file InDetAdaptiveMultiSecVtxFinderTool.h.

◆ m_privtxRef

FloatProperty InDet::InDetAdaptiveMultiSecVtxFinderTool::m_privtxRef {this, "MomentumProjectionOnDirection", -999999.9, "pri vtx ref"}
private

Definition at line 98 of file InDetAdaptiveMultiSecVtxFinderTool.h.

98{this, "MomentumProjectionOnDirection", -999999.9, "pri vtx ref"};

◆ m_SeedFinder

ToolHandle<Trk::IVertexSeedFinder> InDet::InDetAdaptiveMultiSecVtxFinderTool::m_SeedFinder {this, "SeedFinder", "Trk::IndexedCrossDistancesSeedFinder", "Seed finder"}
private

Definition at line 94 of file InDetAdaptiveMultiSecVtxFinderTool.h.

94{this, "SeedFinder", "Trk::IndexedCrossDistancesSeedFinder", "Seed finder"};

◆ m_significanceCutSeeding

DoubleProperty InDet::InDetAdaptiveMultiSecVtxFinderTool::m_significanceCutSeeding {this, "significanceCutSeeding", 10, "significanceCutSeeding"}
private

Definition at line 99 of file InDetAdaptiveMultiSecVtxFinderTool.h.

99{this, "significanceCutSeeding", 10, "significanceCutSeeding"};

◆ m_SVtrkFilter

ToolHandle<InDet::IInDetTrackSelectionTool> InDet::InDetAdaptiveMultiSecVtxFinderTool::m_SVtrkFilter {this, "SecVtxTrackSelector", "InDet::SecVtxTrackSelector", "SV track selection tool"}
private

Definition at line 92 of file InDetAdaptiveMultiSecVtxFinderTool.h.

92{this, "SecVtxTrackSelector", "InDet::SecVtxTrackSelector", "SV track selection tool"};

◆ m_trkFilter

ToolHandle<InDet::IInDetTrackSelectionTool> InDet::InDetAdaptiveMultiSecVtxFinderTool::m_trkFilter {this, "BaseTrackSelector", "InDet::DetailedTrackSelectToolRelax", "Base track selection tool"}
private

Definition at line 91 of file InDetAdaptiveMultiSecVtxFinderTool.h.

91{this, "BaseTrackSelector", "InDet::DetailedTrackSelectToolRelax", "Base track selection tool"};

◆ m_VertexFitter

ToolHandle<Trk::AdaptiveMultiVertexFitter> InDet::InDetAdaptiveMultiSecVtxFinderTool::m_VertexFitter {this, "VertexFitterTool", "Trk::AdaptiveMultiVertexFitter", "Multi Vertex Fitter"}
private

Definition at line 90 of file InDetAdaptiveMultiSecVtxFinderTool.h.

90{this, "VertexFitterTool", "Trk::AdaptiveMultiVertexFitter", "Multi Vertex Fitter"};

The documentation for this class was generated from the following files: