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

#include <InDetIterativeSecVtxFinderTool.h>

Inheritance diagram for InDet::InDetIterativeSecVtxFinderTool:
Collaboration diagram for InDet::InDetIterativeSecVtxFinderTool:

Public Member Functions

 InDetIterativeSecVtxFinderTool (const std::string &t, const std::string &n, const IInterface *p)
virtual ~InDetIterativeSecVtxFinderTool ()=default
 Destructor.
StatusCode initialize () override
StatusCode finalize () override
std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > findVertex (const TrackCollection *trackTES) override
std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > findVertex (const Trk::TrackParticleBaseCollection *trackTES) override
std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > findVertex (const xAOD::TrackParticleContainer *trackParticles) override
 Finding method.
std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > findVertex (const std::vector< const xAOD::IParticle * > &inputTracks) override
void setPriVtxPosition (double, double, double) override
int getModes1d (std::vector< int > *, std::vector< int > *, std::vector< int > *) const

Private Member Functions

std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > findVertex (const std::vector< Trk::ITrackLink * > &trackVector)
bool passHitsFilter (const Trk::TrackParameters *, float vtxR, float absvz) const
bool V0kine (const std::vector< Amg::Vector3D > &, const Amg::Vector3D &position, float &, float &) const
const std::vector< Amg::Vector3DgetVertexMomenta (xAOD::Vertex *myxAODVertex) const
float removeTracksInBadSeed (xAOD::Vertex *myxAODVertex, std::vector< const Trk::TrackParameters * > &) const
void FillXcheckdefauls ()
void removeCompatibleTracks (xAOD::Vertex *myxAODVertex, std::vector< const Trk::TrackParameters * > &perigeesToFit, std::vector< Trk::ITrackLink * > &seedTracks) const
void removeAllFrom (std::vector< const Trk::TrackParameters * > &perigeesToFit, std::vector< Trk::ITrackLink * > &seedTracks) const
double compatibility (const Trk::TrackParameters &measPerigee, const AmgSymMatrix(3) &covariancePosition, const Amg::Vector3D &position) const
void countTracksAndNdf (xAOD::Vertex *myxAODVertex, float &ndf, int &ntracks) const
void SGError (const std::string &errService)
virtual void printParameterSettings ()

Static Private Member Functions

static double VrtVrtDist (xAOD::Vertex *v1, xAOD::Vertex *v2)

Private Attributes

Amg::Vector3D m_privtx
std::vector< Amg::Vector3Dm_TrkAtVtxMomenta
ToolHandle< Trk::AdaptiveVertexFitterm_iVertexFitter {this, "VertexFitterTool", "Trk::AdaptiveVertexFitter","Vertex Fitter"}
ToolHandle< InDet::IInDetTrackSelectionToolm_trkFilter {this, "BaseTrackSelector", "InDet::InDetTrackSelection", "base track selector"}
ToolHandle< InDet::IInDetTrackSelectionToolm_SVtrkFilter {this, "SecVtxTrackSelector", "InDet::InDetSecVtxTrackSelection", "SV track selector"}
ToolHandle< Trk::IVertexSeedFinderm_SeedFinder {this, "SeedFinder", "Trk::IndexedCrossDistancesSeedFinder", "seed finder"}
ToolHandle< Trk::IImpactPoint3dEstimatorm_ImpactPoint3dEstimator {this, "ImpactPoint3dEstimator", "Trk::ImpactPoint3dEstimator", "impact point estimator"}
ToolHandle< Trk::IVertexLinearizedTrackFactorym_LinearizedTrackFactory {this, "LinearizedTrackFactory", "Trk::FullLinearizedTrackFactory", "linearized track factory"}
Gaudi::Property< int > m_filterLevel {this,"VertexFilterLevel",0,""}
Gaudi::Property< double > m_significanceCutSeeding {this, "significanceCutSeeding", 9., ""}
Gaudi::Property< double > m_maximumChi2cutForSeeding {this, "maxCompatibilityCutSeeding",18., ""}
Gaudi::Property< double > m_minWghtAtVtx {this, "minTrackWeightAtVtx",0.02, ""}
Gaudi::Property< double > m_maxVertices {this, "maxVertices",20, ""}
Gaudi::Property< double > m_CutHitsFilter {this, "TrackInnerOuterFraction",0.95, ""}
Gaudi::Property< float > m_privtxRef {this, "MomentumProjectionOnDirection",-999.9,""}
Gaudi::Property< float > m_minVtxDist {this, "SeedsMinimumDistance",0.1,""}
Gaudi::Property< bool > m_createSplitVertices {this, "createSplitVertices", false, ""}
Gaudi::Property< int > m_splitVerticesTrkInvFraction {this, "splitVerticesTrkInvFraction",2,""}
 Integer: 1./fraction of tracks to be assigned to the tag split vertex.
Gaudi::Property< bool > m_reassignTracksAfterFirstFit {this, "reassignTracksAfterFirstFit", true, ""}
Gaudi::Property< bool > m_doMaxTracksCut {this, "doMaxTracksCut", false, ""}
Gaudi::Property< unsigned int > m_maxTracks {this, "MaxTracks", 5000,""}
TTree * m_OTree {}
long int m_evtNum {}
int m_iterations {}
std::vector< int > * m_leastmodes {}
std::vector< std::vector< float > > * m_sdFsmwX {}
std::vector< std::vector< float > > * m_sdFsmwY {}
std::vector< std::vector< float > > * m_sdFsmwZ {}
std::vector< std::vector< float > > * m_sdcrsWght {}
std::vector< int > * m_nperiseed {}
std::vector< float > * m_seedX {}
std::vector< float > * m_seedY {}
std::vector< float > * m_seedZ {}
std::vector< float > * m_seedXYdist {}
std::vector< float > * m_seedZdist {}
std::vector< int > * m_seedac {}
float m_v0mass {}
float m_v0ee {}
float m_dir {}
float m_ndf {}
float m_hif {}
int m_ntracks {}
std::vector< Amg::VectorXm_trkdefiPars

Detailed Description

Author
Lianyou SHAN ( IHEP )

developed upon the InDetPriVxFinderTool by :

Author
Giacinto Piacquadio (Freiburg University)

This class provides an implementation for a secondary vertex finding tool, which uses the Adaptive Vertex Fitter to reject outliers not belonging to the vertex interaction.

Definition at line 66 of file InDetIterativeSecVtxFinderTool.h.

Constructor & Destructor Documentation

◆ InDetIterativeSecVtxFinderTool()

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

Definition at line 67 of file InDetIterativeSecVtxFinderTool.cxx.

68 : base_class(t,n,p){}

◆ ~InDetIterativeSecVtxFinderTool()

virtual InDet::InDetIterativeSecVtxFinderTool::~InDetIterativeSecVtxFinderTool ( )
virtualdefault

Destructor.

Member Function Documentation

◆ compatibility()

double InDet::InDetIterativeSecVtxFinderTool::compatibility ( const Trk::TrackParameters & measPerigee,
const AmgSymMatrix(3) & covariancePosition,
const Amg::Vector3D & position ) const
private

Definition at line 1696 of file InDetIterativeSecVtxFinderTool.cxx.

1698{
1699 double returnValue = 0. ;
1700
1701 Trk::LinearizedTrack* myLinearizedTrack=m_LinearizedTrackFactory->linearizedTrack(
1702 &measPerigee, position );
1703
1704// Amg::Vector3D mapca = myLinearizedTrack->expectedMomentumAtPCA() ;
1705
1706 AmgMatrix(2,2) weightReduced=myLinearizedTrack->expectedCovarianceAtPCA().block<2,2>(0,0);
1707
1708 AmgMatrix(2,2) errorVertexReduced=(myLinearizedTrack->positionJacobian()*
1709 ( covariancePosition*myLinearizedTrack->positionJacobian().transpose())).block<2,2>(0,0);
1710
1711 weightReduced+=errorVertexReduced;
1712
1713 weightReduced = weightReduced.inverse().eval();
1714 Amg::Vector2D trackParameters2D = myLinearizedTrack->expectedParametersAtPCA().block<2,1>(0,0);
1715 returnValue += trackParameters2D.dot(weightReduced*trackParameters2D);
1716
1717 delete myLinearizedTrack;
1718 myLinearizedTrack=nullptr;
1719
1720 return returnValue;
1721}
#define AmgMatrix(rows, cols)
ToolHandle< Trk::IVertexLinearizedTrackFactory > m_LinearizedTrackFactory

◆ countTracksAndNdf()

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

Definition at line 1768 of file InDetIterativeSecVtxFinderTool.cxx.

1770{
1771 ndf = -3.0 ;
1772 ntrk = 0 ;
1773
1774 if ( myxAODVertex )
1775 {
1776 ndf = myxAODVertex->numberDoF();
1777 std::vector<Trk::VxTrackAtVertex> myVxTracksAtVtx = myxAODVertex->vxTrackAtVertex();
1778 std::vector<Trk::VxTrackAtVertex>::iterator tracksBegin=myVxTracksAtVtx.begin();
1779 std::vector<Trk::VxTrackAtVertex>::iterator tracksEnd=myVxTracksAtVtx.end();
1780
1781 for (std::vector<Trk::VxTrackAtVertex>::iterator tracksIter=tracksBegin;
1782 tracksIter!=tracksEnd;++tracksIter)
1783 {
1784
1785 if ( (*tracksIter).weight() > m_minWghtAtVtx )
1786 {
1787 ntrk+=1;
1788 }
1789 }
1790 }
1791}
float numberDoF() const
Returns the number of degrees of freedom of the vertex fit as float.
std::vector< Trk::VxTrackAtVertex > & vxTrackAtVertex()
Non-const access to the VxTrackAtVertex vector.

◆ FillXcheckdefauls()

void InDet::InDetIterativeSecVtxFinderTool::FillXcheckdefauls ( )
private

Definition at line 1502 of file InDetIterativeSecVtxFinderTool.cxx.

1503{
1504 m_seedac->push_back( 0 ) ;
1505 m_nperiseed->push_back( 0 ) ;
1506}

◆ finalize()

StatusCode InDet::InDetIterativeSecVtxFinderTool::finalize ( )
override

Definition at line 1648 of file InDetIterativeSecVtxFinderTool.cxx.

1649{
1650#ifdef MONITORTUNES
1651 m_OTree->Write();
1652
1653 delete m_leastmodes ;
1654 delete m_sdFsmwX ;
1655 delete m_sdFsmwY ;
1656 delete m_sdFsmwZ ;
1657 delete m_sdcrsWght ;
1658 delete m_nperiseed ;
1659 delete m_seedX ;
1660 delete m_seedY ;
1661 delete m_seedZ ;
1662 delete m_seedXYdist ;
1663 delete m_seedZdist ;
1664
1665 delete m_seedac ;
1666
1667 m_leastmodes = nullptr ;
1668 m_sdFsmwX = nullptr ;
1669 m_sdFsmwY = nullptr ;
1670 m_sdFsmwZ = nullptr ;
1671 m_sdcrsWght = nullptr ;
1672 m_nperiseed = nullptr ;
1673 m_seedX = nullptr ;
1674 m_seedY = nullptr ;
1675 m_seedZ = nullptr ;
1676 m_seedXYdist = nullptr ;
1677 m_seedZdist = nullptr ;
1678 m_seedac = nullptr ;
1679
1680#endif
1681
1682 return StatusCode::SUCCESS;
1683}
std::vector< std::vector< float > > * m_sdcrsWght
std::vector< std::vector< float > > * m_sdFsmwY
std::vector< std::vector< float > > * m_sdFsmwZ
std::vector< std::vector< float > > * m_sdFsmwX

◆ findVertex() [1/5]

std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > InDet::InDetIterativeSecVtxFinderTool::findVertex ( const std::vector< const xAOD::IParticle * > & inputTracks)
override

Definition at line 125 of file InDetIterativeSecVtxFinderTool.cxx.

126{
127
128 ATH_MSG_DEBUG(" Number of input tracks before track selection: " << inputTracks.size());
129
130 m_evtNum ++ ;
131
132 std::vector<Trk::ITrackLink*> selectedTracks;
133
134
135 bool selectionPassed;
136 m_trkdefiPars.clear() ;
137 xAOD::Vertex null;
138 null.makePrivateStore();
139 null.setPosition(Amg::Vector3D(0,0,0));
140 AmgSymMatrix(3) vertexError;
141 vertexError.setZero();
142 null.setCovariancePosition(vertexError);
143 std::vector<const xAOD::IParticle*>::const_iterator trk_iter;
144 for (trk_iter= inputTracks.begin(); trk_iter != inputTracks.end(); ++trk_iter)
145 {
146 const xAOD::TrackParticle * tmp=dynamic_cast<const xAOD::TrackParticle *> ((*trk_iter));
147
148 selectionPassed=static_cast<bool>(m_trkFilter->accept( *tmp, &null));
149 if ( selectionPassed ) selectionPassed =static_cast<bool>(m_SVtrkFilter->accept(*tmp,&null));
150
151 if (selectionPassed)
152 {
153 Amg::VectorX par = tmp->definingParameters();
154 par[0] = 1.0*tmp->hitPattern() ;
155 m_trkdefiPars.push_back( par ) ;
156
157 ElementLink<xAOD::TrackParticleContainer> linkTP;
158 linkTP.setElement(tmp);
159
160 Trk::LinkToXAODTrackParticle* link= new Trk::LinkToXAODTrackParticle(linkTP);
161
162 selectedTracks.push_back(link);
163 }
164
165 }
166
167 std::pair<xAOD::VertexContainer*, xAOD::VertexAuxContainer*> returnContainers
168 =findVertex( selectedTracks );
169
170 return returnContainers;
171
172}
#define ATH_MSG_DEBUG(x)
#define AmgSymMatrix(dim)
if(febId1==febId2)
ToolHandle< InDet::IInDetTrackSelectionTool > m_SVtrkFilter
std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > findVertex(const TrackCollection *trackTES) override
ToolHandle< InDet::IInDetTrackSelectionTool > m_trkFilter
void makePrivateStore()
Create a new (empty) private store for this object.
void setPosition(const Amg::Vector3D &position)
Sets the 3-position.
Eigen::Matrix< double, 3, 1 > Vector3D
Eigen::Matrix< double, Eigen::Dynamic, 1 > VectorX
Dynamic Vector - dynamic allocation.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
Vertex_v1 Vertex
Define the latest version of the vertex class.

◆ findVertex() [2/5]

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

Definition at line 176 of file InDetIterativeSecVtxFinderTool.cxx.

177{
178
179 //two things need to be added
180 //1) the seeding mechanism
181 //2) the iterative removal of tracks
182 std::vector<Trk::ITrackLink*> origTracks=trackVector;
183 std::vector<Trk::ITrackLink*> seedTracks=trackVector;
184// in the iteration from the below do { ... } while () loop,
185// the container of seedTracks is updated/dynamic : some tracks are moved out
186// once a vertex is successfully found
187//
188
189 xAOD::VertexContainer* theVertexContainer = new xAOD::VertexContainer;
190 xAOD::VertexAuxContainer* theVertexAuxContainer = new xAOD::VertexAuxContainer;
191 theVertexContainer->setStore( theVertexAuxContainer );
192 const auto invalidResponse = std::make_pair(theVertexContainer, theVertexAuxContainer);
193
194 //bail out early with only Dummy vertex if multiplicity cut is applied and exceeded
195 if (m_doMaxTracksCut && (trackVector.size() > m_maxTracks)){
196 ATH_MSG_WARNING( trackVector.size() << " tracks - exceeds maximum (" << m_maxTracks << "), skipping vertexing and returning only dummy..." );
197 return std::make_pair(theVertexContainer, theVertexAuxContainer);
198 }
199
200 m_iterations = -1 ;
201 unsigned int seedtracknumber=seedTracks.size();
202
203 //used to store seed info
204 Amg::Vector3D seedVertex;
205
206 //prepare iterators for tracks only necessary for seeding
207 std::vector<Trk::ITrackLink*>::iterator seedBegin;
208 std::vector<Trk::ITrackLink*>::iterator seedEnd;
209
210// fortunately the tool InDetIterativeSecVtxFinderTool::findVertex is called only once per event
211
212#ifdef MONITORTUNES
213
214 m_leastmodes->clear() ;
215 m_sdFsmwX->clear() ;
216 m_sdFsmwY->clear() ;
217 m_sdFsmwZ->clear() ;
218 m_sdcrsWght->clear() ;
219
220 m_nperiseed->clear() ;
221 m_seedX->clear() ;
222 m_seedY->clear() ;
223 m_seedZ->clear() ;
224 m_seedXYdist->clear() ;
225 m_seedZdist->clear() ;
226 m_seedac->clear() ;
227
228 SG::AuxElement::Decorator<std::vector<float> > mDecor_trkWght( "trkWeight" ) ;
229 SG::AuxElement::Decorator<std::vector<float> > mDecor_trkDOE( "trkDistOverError" ) ;
230 SG::AuxElement::Decorator<float> mDecor_direction( "MomentaDirection" ) ;
231 SG::AuxElement::Decorator< float > mDecor_HitsFilter( "radiiPattern" );
232 SG::AuxElement::Decorator< int > mDecor_NumTrk( "NumTrkAtVtx" );
233#endif
234
235 SG::AuxElement::Decorator<float> mDecor_sumPt2( "sumPt2" );
236 SG::AuxElement::Decorator<float> mDecor_mass( "mass" );
237 SG::AuxElement::Decorator<float> mDecor_energy( "ee" );
238 SG::AuxElement::Decorator<int> mDecor_nrobbed( "nrobbed" );
239 SG::AuxElement::Decorator<int> mDecor_intrk( "NumInputTrk" );
240
241 do
242 {
243 seedBegin=seedTracks.begin();
244 seedEnd=seedTracks.end();
245
246 if (seedtracknumber ==0)
247 {
248 ATH_MSG_DEBUG( " New iteration. No tracks available after track selection for seeding. No finding done." );
249 break;
250 }
251
252 m_iterations += 1;
253 ATH_MSG_DEBUG( "ITERATION NUMBER " << m_iterations );
254
255 //now find a new SEED
256
257 std::vector<const Trk::TrackParameters*> perigeeList;
258 for (std::vector<Trk::ITrackLink*>::iterator seedtrkAtVtxIter=seedBegin;
259 seedtrkAtVtxIter!=seedEnd;++seedtrkAtVtxIter) {
260 perigeeList.push_back( (*seedtrkAtVtxIter)->parameters() );
261 }
262
263 int ncandi = 0 ;
264 xAOD::Vertex theconstraint;
265
266
267 ATH_MSG_DEBUG( " goto seed finder " );
268
269 std::unique_ptr<Trk::IMode3dInfo> info;
270 seedVertex = m_SeedFinder->findSeed(m_privtx.x(), m_privtx.y(),
271 info, perigeeList);
272
273 ATH_MSG_DEBUG( " seedFinder finished " );
274
275#ifdef MONITORTUNES
276 std::vector<float> FsmwX, FsmwY, FsmwZ, wght ;
277
278
279 double cXY = -9.9, cZ = -9.9 ;
280 info->getCorrelationDistance( cXY, cZ ) ;
281
282 m_sdFsmwX->push_back( FsmwX ) ;
283 m_sdFsmwY->push_back( FsmwY ) ;
284 m_sdFsmwZ->push_back( FsmwZ ) ;
285 m_sdcrsWght->push_back( wght ) ;
286
287 m_seedX->push_back( seedVertex.x() ) ;
288 m_seedY->push_back( seedVertex.y() ) ;
289 m_seedZ->push_back( seedVertex.z() ) ;
290 m_seedXYdist->push_back( cXY ) ;
291 m_seedZdist->push_back( cZ ) ;
292#endif
293
294
295
296 Amg::MatrixX looseConstraintCovariance(3,3);
297 looseConstraintCovariance.setIdentity();
298 looseConstraintCovariance = looseConstraintCovariance * 1e+8;
299 theconstraint = xAOD::Vertex();
300 theconstraint.setPosition( seedVertex );
301 theconstraint.setCovariancePosition( looseConstraintCovariance );
302 theconstraint.setFitQuality( -99.9, 1.0 );
303
304// on the other side VKalFitter : Tracking/TrkVertexFitter/TrkVKalVrtFitter/src/VKalVrtFitSvc
305
306 ATH_MSG_DEBUG( " seed at x: " << seedVertex.x() <<
307 " at y: " << seedVertex.y() <<
308 " at z: " << seedVertex.z() );
309
310
311 if ( seedVertex.z()==0. ) {
312
313 ATH_MSG_DEBUG( "No seed found: no further vertices in event" );
314 ATH_MSG_DEBUG( "Number of input tracks: " << perigeeList.size() << " but no seed returned." );
315
316
317#ifdef MONITORTUNES
319#endif
320
321 break;
322 }
323
324#ifdef MONITORTUNES
325 m_nperiseed->push_back( ncandi ) ;
326 m_seedac->push_back( 1 ) ;
327#endif
328
329 //now question (remove tracks which are too far away??? I think so...)
330 std::vector<const Trk::TrackParameters*> perigeesToFit;
331 std::vector<const Trk::TrackParameters*> perigeesToFitSplitVertex;
332
333 int numberOfTracks( perigeeList.size() );
334
335 std::vector<const Trk::TrackParameters*>::const_iterator perigeeListBegin=perigeeList.begin();
336 std::vector<const Trk::TrackParameters*>::const_iterator perigeeListEnd=perigeeList.end();
337
338 int counter=0;
339 for (std::vector<const Trk::TrackParameters*>::const_iterator perigeeListIter=perigeeListBegin;
340 perigeeListIter!=perigeeListEnd;++perigeeListIter)
341 {
342
343 if (numberOfTracks<=2)
344 {
345 perigeesToFit.push_back(*perigeeListIter);
346 counter+=1;
347 }
348 else if ( numberOfTracks <= 3 && !m_createSplitVertices )
349 {
350 perigeesToFit.push_back(*perigeeListIter);
351 counter+=1;
352 }
353 else if (numberOfTracks<=4*m_splitVerticesTrkInvFraction && m_createSplitVertices)
354 {
355 // few tracks are left, put them into the fit regardless their position!
356 if (counter % m_splitVerticesTrkInvFraction == 0)
357 {
358 perigeesToFit.push_back(*perigeeListIter);
359 counter+=1;
360 }
361 else
362 {
363 perigeesToFitSplitVertex.push_back(*perigeeListIter);
364 counter+=1;
365 }
366 }
367 else
368 { //check first whether it is not too far away!
369
370 double distance=0.;
371 try
372 {
373 std::unique_ptr<Trk::PlaneSurface> mySurface=m_ImpactPoint3dEstimator->Estimate3dIP(*perigeeListIter,&seedVertex,distance);
374 ATH_MSG_VERBOSE( " ImpactPoint3dEstimator done " );
375 }
376 catch (error::ImpactPoint3dEstimatorProblem err)
377 {
378 ATH_MSG_WARNING( " ImpactPoint3dEstimator failed to find minimum distance between track and vertex seed: " << err.p );
379 }
380
381 if (distance<0)
382 {
383 ATH_MSG_WARNING( " Distance between track and seed vtx is negative: " << distance );
384 }
385
386 const Trk::TrackParameters* myPerigee=(*perigeeListIter);
387
388 //very approximate error
389 double doe = 99999999.9 ;
390 double error= 0.;
391
392 if( myPerigee && myPerigee->covariance() )
393 {
394 error = std::sqrt((*myPerigee->covariance())(Trk::d0,Trk::d0)+
395 (*myPerigee->covariance())(Trk::z0,Trk::z0));
396 }//end of the security check
397
398 if (error==0.)
399 {
400 ATH_MSG_ERROR( " Error is zero! " << distance );
401 error=1.;
402 }
403 doe = distance/error ;
404
405
406 ATH_MSG_VERBOSE( " Distance between track and seed vtx: " << distance << " d/s(d) = " <<
407 distance/error << " err " << error );
408
409
410 if ( doe < m_significanceCutSeeding )
411 {
413 {
414 perigeesToFit.push_back(*perigeeListIter);
415 counter+=1;
416 }
417 else
418 {
419 perigeesToFitSplitVertex.push_back(*perigeeListIter);
420 counter+=1;
421 }
422 }
423
424 }
425 } // end of loop for filling perigeeListIter into perigeesToFit
426
427// here is while-okay
428
429
430 ATH_MSG_VERBOSE( " Considering n. " << perigeesToFit.size() << " tracks for the fit. " );
432 {
433 ATH_MSG_VERBOSE( " and n. " << perigeesToFitSplitVertex.size() << " tracks for split vertex fit. " );
434 }
435
436
437 m_ndf = -3. ;
438 m_ntracks = 0;
439 m_v0mass = -299.9 ;
440 m_v0ee = -299.9 ;
441 m_dir = -99999.9 ;
442 m_hif = -1.0 ;
443
444 if (perigeesToFit.empty())
445 {
446
447 ATH_MSG_DEBUG( " No good seed found. Exiting search for vertices..." );
448
449#ifdef MONITORTUNES
450 xAOD::Vertex * seededxAODVertex = new xAOD::Vertex;
451 theVertexContainer->push_back( seededxAODVertex ); //
452
453 seededxAODVertex->setPosition( seedVertex );
454 Amg::MatrixX dummyCovariance(3,3);
455 dummyCovariance.setIdentity();
456 seededxAODVertex->setCovariancePosition( dummyCovariance );
457 seededxAODVertex->vxTrackAtVertex() = std::vector<Trk::VxTrackAtVertex>();
458 seededxAODVertex->setFitQuality( 99.9, m_ndf ) ;
459 seededxAODVertex->setVertexType( xAOD::VxType::NotSpecified );
460
461 mDecor_mass( *seededxAODVertex ) = m_v0mass ;
462 mDecor_energy( *seededxAODVertex ) = m_v0ee ;
463 mDecor_NumTrk( *seededxAODVertex ) = m_ntracks ;
464 mDecor_HitsFilter( *seededxAODVertex ) = m_hif ;
465 mDecor_direction( *seededxAODVertex ) = m_dir ;
466 mDecor_intrk( *seededxAODVertex ) = numberOfTracks ;
467#endif
468
469 break;
470 }
471 // one has to break off the loop if no seed is found
472
473 //now you have perigeeToFit and perigeeToFitSplitVertices
474 //AND HERE YOU DO THE FIT
475 //to reassign vertices you look ino what is already in myVxCandidate
476 //you do it only ONCE!
477
478 xAOD::Vertex * myxAODVertex = nullptr;
479
480
481 if ( perigeesToFit.size()>1) // our main use case
482 {
483 myxAODVertex=m_iVertexFitter->fit( perigeesToFit, seedVertex );
484 }
485
486 //
487 // The fitter neither tell whether the vertex is good or not, but the chi2/dof
488 //
489
490 countTracksAndNdf( myxAODVertex, m_ndf, m_ntracks);
491
492 m_v0mass = -199.9 ;
493
494 bool goodVertex = myxAODVertex != nullptr && m_ndf >0 && m_ntracks >=2 ;
495
496
497 ATH_MSG_DEBUG( " xAOD::Vertex : " << ( myxAODVertex != nullptr ? 1 : 0 )
498 << ", #dof = " << m_ndf << ", #tracks (weight>0.01) = " << m_ntracks );
499
500
501// below is while-okay
502 if (!goodVertex)
503 {
504
505 ATH_MSG_DEBUG( " Going to new iteration with: " << seedTracks.size() << " seed tracks after BAD VERTEX. " );
506
507
508 if (myxAODVertex)
509 {
510// DUMMY !! Only for validation study
511 removeCompatibleTracks( myxAODVertex, perigeesToFit, seedTracks);
512#ifdef MONITORTUNES
513 myxAODVertex->setVertexType( xAOD::VxType::KinkVtx );
514 theVertexContainer->push_back( myxAODVertex );
515 mDecor_mass( *myxAODVertex ) = m_v0mass ;
516 mDecor_energy( *myxAODVertex ) = m_v0ee ;
517 mDecor_NumTrk( *myxAODVertex ) = m_ntracks ;
518 mDecor_HitsFilter( *myxAODVertex ) = m_hif ;
519 mDecor_direction( *myxAODVertex ) = m_dir ;
520 mDecor_intrk( *myxAODVertex ) = numberOfTracks ;
521#else
522 delete myxAODVertex;
523 myxAODVertex=0;
524#endif
525 }
526 else
527 {
528 removeAllFrom( perigeesToFit, seedTracks );
529#ifdef MONITORTUNES
530 xAOD::Vertex * seededxAODVertex = new xAOD::Vertex;
531 theVertexContainer->push_back( seededxAODVertex ); // have to add vertex to container here first so it can use its aux store
532
533 seededxAODVertex->setPosition( seedVertex );
534 Amg::MatrixX dummyCovariance(3,3);
535 dummyCovariance.setIdentity();
536 seededxAODVertex->setCovariancePosition( dummyCovariance );
537 seededxAODVertex->vxTrackAtVertex() = std::vector<Trk::VxTrackAtVertex>();
538 seededxAODVertex->setFitQuality( 99.9, m_ndf ) ;
539 seededxAODVertex->setVertexType( xAOD::VxType::NoVtx );
540
541 mDecor_mass( *seededxAODVertex ) = m_v0mass ;
542 mDecor_energy( *seededxAODVertex ) = m_v0ee ;
543 mDecor_NumTrk( *seededxAODVertex ) = m_ntracks ;
544 mDecor_HitsFilter( *seededxAODVertex ) = m_hif ;
545 mDecor_direction( *seededxAODVertex ) = m_dir ;
546 mDecor_intrk( *seededxAODVertex ) = numberOfTracks ;
547#endif
548 }
549
550 continue ; // try next seed
551
552 }
553
554// here is a boundary XXXXXXXXXXXXXX, at least a "good" vertex is found
555 // Now the goodVertex could be not so good ...
556 mDecor_nrobbed( *myxAODVertex ) = 0 ;
558 {
559
560
561 ATH_MSG_VERBOSE( " N tracks used for fit before reallocating: " << perigeesToFit.size() );
562
563 //now you HAVE a good vertex
564 //but you want to add the tracks which you missed...
565
566 int numberOfAddedTracks=0;
567 const AmgSymMatrix(3) covariance = (&(*myxAODVertex))->covariancePosition() ;
568 const Amg::Vector3D position = (&(*myxAODVertex))->position() ;
569
570 //iterate on remaining vertices and cross-check if tracks need to be attached
571 //to new vertex
572 xAOD::VertexContainer::iterator vxBegin=theVertexContainer->begin();
573 xAOD::VertexContainer::iterator vxEnd=theVertexContainer->end();
574
575 for (xAOD::VertexContainer::iterator vxIter=vxBegin;vxIter!=vxEnd;++vxIter)
576 {
577 // A vertex should not rob tracks from itself
578#ifdef MONITORTUNES
579 if ( (*vxIter)->vertexType() == xAOD::VxType::NoVtx ) continue ;
580#endif
581 std::vector<Trk::VxTrackAtVertex>* myVxTracksAtVtx=(&(*vxIter)->vxTrackAtVertex());
582
583 if ( ! myVxTracksAtVtx ) continue;
584
585 int nrobbed = 0 ;
586
587 const AmgSymMatrix(3) oldcovariance = (*vxIter)->covariancePosition() ;
588 const Amg::Vector3D oldposition = (*vxIter)->position() ;
589
590 //now iterate on tracks at vertex
591 std::vector<Trk::VxTrackAtVertex>::const_iterator tracksBegin=myVxTracksAtVtx->begin();
592 std::vector<Trk::VxTrackAtVertex>::const_iterator tracksEnd=myVxTracksAtVtx->end();
593
594
595 ATH_MSG_VERBOSE( " Iterating over new vertex to look for tracks to reallocate... " );
596
597
598 for (std::vector<Trk::VxTrackAtVertex>::const_iterator tracksIter=tracksBegin;
599 tracksIter!=tracksEnd;)
600 {
601
602 //only try with tracks which are not too tightly assigned to another vertex
603 if ((*tracksIter).weight() > m_minWghtAtVtx )
604 {
605 ++tracksIter;
606 continue;
607 }
608
609 const Trk::TrackParameters* trackPerigee=(*tracksIter).initialPerigee();
610
611 if (trackPerigee==nullptr)
612 {
613 ATH_MSG_ERROR( " Cast to perigee gives 0 pointer, cannot continue " );
614 return invalidResponse;
615;
616 }
617
618 double chi2_newvtx=compatibility(*trackPerigee, covariance, position );
619 double chi2_oldvtx=compatibility(*trackPerigee, oldcovariance, oldposition );
620
621 double minDist = VrtVrtDist( myxAODVertex, *vxIter ) ;
622
623
624 ATH_MSG_VERBOSE( "Compatibility to old vtx is : " << chi2_oldvtx <<
625 " to new vtx is: " << chi2_newvtx );
626
627
628 bool remove=false;
629
630 if ( chi2_newvtx < chi2_oldvtx
632 )
633 {
634
635
636 ATH_MSG_DEBUG( " Found track of old vertex (chi2= " << chi2_oldvtx <<
637 ") more compatible to new one (chi2= " << chi2_newvtx << ")" );
638
639
640 perigeesToFit.push_back(trackPerigee);
641 //but you need to re-add it to the seedTracks too...
642
643 bool isFound=false;
644
645 std::vector<Trk::ITrackLink*>::iterator origBegin=origTracks.begin();
646 std::vector<Trk::ITrackLink*>::iterator origEnd=origTracks.end();
647
648 for (std::vector<Trk::ITrackLink*>::iterator origIter=origBegin;
649 origIter!=origEnd;++origIter)
650 {
651 if ( (*origIter)->parameters()==trackPerigee )
652 {
653
654 ATH_MSG_VERBOSE( " found the old perigee to be re-added to seedTracks in order to be deleted again!" );
655
656 isFound=true;
657 seedTracks.push_back(*origIter);
658 break;
659 }
660 }
661
662 if (!isFound)
663 {
664 ATH_MSG_WARNING( " Cannot find old perigee to re-add back to seed tracks... " );
665 }
666
667 numberOfAddedTracks+=1;
668
669 remove=true;
670 }
671
672 if (remove)
673 {
674 //now you have to delete the track from the old vertex...
675 //easy???
676 nrobbed ++ ;
677 tracksIter=myVxTracksAtVtx->erase(tracksIter);
678 tracksBegin=myVxTracksAtVtx->begin();
679 tracksEnd=myVxTracksAtVtx->end();
680 }
681 else
682 {
683 ++tracksIter;
684 }
685 }//end of iterating on tracks at previous vertices
686
687 if ( nrobbed > 0 ) mDecor_nrobbed( *(*vxIter) ) ++;
688
689 }//end of iterating on already found vertices in event
690
691
692 ATH_MSG_VERBOSE( " N tracks used for fit after reallocating: " << perigeesToFit.size() );
693
694
695 //now you have to delete the previous xAOD::Vertex, do a new fit, i
696 // then check if you still have a good vertex
697
698 if ( numberOfAddedTracks > 0 )
699 {
700
701 ATH_MSG_DEBUG( " refit with additional " << numberOfAddedTracks
702 << " from other vertices " );
703
704 Amg::Vector3D fitposition = myxAODVertex->position() ;
705 delete myxAODVertex;
706
707 myxAODVertex=nullptr;
708
709 if ( perigeesToFit.size()>1)
710 {
711 myxAODVertex=m_iVertexFitter->fit( perigeesToFit, fitposition ) ;
712 mDecor_nrobbed( *myxAODVertex ) = 0 ; // robbing but also refit
713 }
714
715 countTracksAndNdf( myxAODVertex, m_ndf, m_ntracks);
716
717 goodVertex = myxAODVertex != nullptr && m_ndf >0 && m_ntracks >=2 ;
718
719
720 ATH_MSG_DEBUG( " Refitted xAODVertex is pointer: " << myxAODVertex <<
721 " #dof = " << m_ndf << " #tracks (with weight>0.01) " << m_ntracks );
722
723
724 if ( ! goodVertex )
725 {
726// it was ever good vertex at least, it become nothing after eating new a track ...
727 if (myxAODVertex)
728 {
729 removeCompatibleTracks( myxAODVertex, perigeesToFit, seedTracks);
730#ifdef MONITORTUNES
731 theVertexContainer->push_back(myxAODVertex);
732 myxAODVertex->setVertexType( xAOD::VxType::KinkVtx );
733
734 mDecor_mass( *myxAODVertex ) = m_v0mass ;
735 mDecor_energy( *myxAODVertex ) = m_v0ee ;
736 mDecor_HitsFilter( *myxAODVertex ) = m_hif ;
737 mDecor_direction( *myxAODVertex ) = m_dir ;
738 mDecor_NumTrk( *myxAODVertex ) = m_ntracks ;
739 mDecor_intrk( *myxAODVertex ) = numberOfTracks ;
740#else
741 delete myxAODVertex;
742 myxAODVertex=0;
743#endif
744 }
745 else
746 {
747 removeAllFrom(perigeesToFit,seedTracks);
748
749 ATH_MSG_DEBUG( " Adding tracks resulted in an invalid vertex. Should be rare... " );
750 ATH_MSG_DEBUG( " Going to new iteration with " << seedTracks.size()
751 << " seed tracks after BAD VERTEX. " );
752
753#ifdef MONITORTUNES
754
755 xAOD::Vertex * seededxAODVertex = new xAOD::Vertex;
756 theVertexContainer->push_back( seededxAODVertex );
757 seededxAODVertex->setPosition( position );
758 seededxAODVertex->setCovariancePosition( covariance );
759 seededxAODVertex->vxTrackAtVertex() = std::vector<Trk::VxTrackAtVertex>();
760
761 seededxAODVertex->setVertexType( xAOD::VxType::NoVtx );
762 seededxAODVertex->setFitQuality( 99.9, m_ndf ) ;
763 mDecor_NumTrk( *seededxAODVertex ) = m_ntracks ;
764 mDecor_mass( *seededxAODVertex ) = m_v0mass ;
765 mDecor_energy( *seededxAODVertex ) = m_v0ee ;
766 mDecor_HitsFilter( *seededxAODVertex ) = m_hif ;
767 mDecor_direction( *seededxAODVertex ) = m_dir ;
768 mDecor_intrk( *seededxAODVertex ) = numberOfTracks ;
769#endif
770 }
771 continue ;
772
773 }
774 }// end if tracks were added...
775
776 }//end reassign tracks from previous vertices and refitting if needed
777
778
779// here is a boundary YYYYYYYYYYYYYYYYY, a "good" vertex is anyway found
780 m_v0mass = -99.9 ;
781
782 myxAODVertex->setVertexType( xAOD::VxType::SecVtx );
783
784 countTracksAndNdf( myxAODVertex, m_ndf, m_ntracks );
785 bool isv0 = ( m_ntracks < 2 ? false
786 : V0kine( getVertexMomenta( myxAODVertex ) ,
787 (&(*myxAODVertex))->position(), m_v0mass, m_dir ) ) ;
788
789#ifdef MONITORTUNES
790 mDecor_HitsFilter( *myxAODVertex ) = m_hif ;
791 mDecor_direction( *myxAODVertex ) = m_dir ;
792 mDecor_NumTrk( *myxAODVertex ) = m_ntracks ;
793#else
794 if ( ( m_filterLevel >= 3 && isv0 )
795 || ( m_filterLevel >= 2 && m_dir < m_privtxRef )
796 )
797 {
798 removeCompatibleTracks( myxAODVertex, perigeesToFit, seedTracks);
799 delete myxAODVertex;
800 myxAODVertex=0;
801 continue ; // fakes, just skip and have its tracks removed from seeds
802 }
803#endif
804
805// now we arrived at a V0 tagging :
806 mDecor_mass( *myxAODVertex ) = m_v0mass ;
807 mDecor_energy( *myxAODVertex ) = m_v0ee ;
808 mDecor_intrk( *myxAODVertex ) = numberOfTracks ;
809 if ( isv0 ) myxAODVertex->setVertexType(xAOD::VxType::V0Vtx);
810
811// Now check the hits inner and outer the seed position
812 m_hif = 0. ;
813#ifdef MONITORTUNES
814 if ( m_filterLevel > 4 && myxAODVertex->vertexType() != xAOD::VxType::V0Vtx )
815#else
816 if ( m_filterLevel > 4 && m_ntracks == 2 && myxAODVertex->vertexType() != xAOD::VxType::V0Vtx )
817#endif
818 {
819
820 xAOD::Vertex * groomed = new xAOD::Vertex() ;
821 * groomed = * myxAODVertex ;
822 m_hif = removeTracksInBadSeed( groomed, perigeesToFit ) ;
823
824 ATH_MSG_DEBUG( " radiiPattern Filter : " << m_hif <<" "<< m_ntracks );
825
826#ifndef MONITORTUNES
827 if ( m_hif > m_CutHitsFilter )
828 {
829 removeCompatibleTracks( groomed, perigeesToFit, seedTracks);
830 delete groomed ;
831 groomed = 0 ;
832 delete myxAODVertex ;
833 myxAODVertex=0;
834 continue ;
835 }
836#endif
837
838 if ( m_hif > 0 )
839 {
840 bool goodgroom = false ;
841 int ngroom = 0 ;
842 if ( perigeesToFit.size() >= 2 )
843 {
844 groomed = nullptr;
845 groomed=m_iVertexFitter->fit( perigeesToFit, myxAODVertex->position() ) ;
846
847 countTracksAndNdf( groomed, m_ndf, ngroom );
848
849 goodgroom = ( groomed != nullptr && m_ndf > 0 && ngroom >= 2 ) ;
850 }
851
852 if ( perigeesToFit.size() < 2 || ! goodgroom )
853 {
854#ifdef MONITORTUNES
855 theVertexContainer->push_back(myxAODVertex);
856
857 if ( perigeesToFit.size() < 2 || ngroom < 1 )
858 myxAODVertex->setVertexType( xAOD::VxType::NoVtx );
859 else
860 myxAODVertex->setVertexType( xAOD::VxType::KinkVtx );
861
862 mDecor_HitsFilter( *myxAODVertex ) = m_hif ;
863 mDecor_NumTrk( *myxAODVertex ) = m_ntracks ;
864#endif
865 if ( groomed )
866 {
867// removeCompatibleTracks( groomed, perigeesToFit, seedTracks);
868 delete groomed ;
869 groomed=nullptr;
870 }
871
872 removeCompatibleTracks( myxAODVertex, perigeesToFit, seedTracks);
873 continue ; // a vertex went bad and skipped
874 }
875
876// found a good vertex after groom and refit, re-check V0 with reduced tracks
877 * myxAODVertex = * groomed ;
878// myxAODVertex->setTrackParticleLinks( groomed->trackParticleLinks() ) ;
879
880 delete groomed ;
881 groomed=nullptr;
882
883 ATH_MSG_DEBUG( " new vertex after grooming with reminded tracks : " << ngroom
884 << " with Chi2/dof " << myxAODVertex->chiSquared()/myxAODVertex->numberDoF() );
885
886 isv0 = false ;
887 if ( ngroom == 2 )
888 isv0 = V0kine( getVertexMomenta( myxAODVertex ) , (&(*myxAODVertex))->position(), m_v0mass, m_dir ) ;
889
890#ifndef MONITORTUNES
891 if ( isv0 || m_dir < m_privtxRef )
892 {
893 removeCompatibleTracks( myxAODVertex, perigeesToFit, seedTracks);
894 delete myxAODVertex;
895 myxAODVertex=0;
896 continue ; // fakes, just skip and have its tracks removed from seeds
897 }
898#endif
899 if ( isv0 )
900 myxAODVertex->setVertexType( xAOD::VxType::V0Vtx );
901 else
902 myxAODVertex->setVertexType( xAOD::VxType::SecVtx );
903
904 mDecor_mass( *myxAODVertex ) = m_v0mass ;
905 mDecor_energy( *myxAODVertex ) = m_v0ee ;
906 mDecor_intrk( *myxAODVertex ) = numberOfTracks ;
907#ifdef MONITORTUNES
908 mDecor_HitsFilter( *myxAODVertex ) = m_hif ;
909 mDecor_direction( *myxAODVertex ) = m_dir ;
910 mDecor_NumTrk( *myxAODVertex ) = ngroom ;
911#endif
912
913 ATH_MSG_DEBUG( " radiiPattern filter worked with VxType : " << myxAODVertex->vertexType() );
914
915 } else
916 {
917 ATH_MSG_DEBUG( " radiiPattern filter hasn't affect in tracks : " << m_ntracks
918 << " with VxType : " << myxAODVertex->vertexType() );
919 delete groomed ;
920 groomed = nullptr ;
921#ifdef MONITORTUNES
922 mDecor_HitsFilter( *myxAODVertex ) = m_hif ;
923#endif
924
925 }
926
927 }
928
929 removeCompatibleTracks( myxAODVertex, perigeesToFit, seedTracks);
930 // SeedTracks are updated so that the remained could be seeded/fit in next iteration
931
932
933 ATH_MSG_VERBOSE( "Number of seeds after removal of outliers: " << seedTracks.size() );
934
935
936 if ( ! m_createSplitVertices )
937 {
938 theVertexContainer->push_back(myxAODVertex);
939 ATH_MSG_DEBUG( " A GOOD vertex is stored in its container with tracks : " << m_ntracks );
940 }
941 // there is no necessary to split for 2ndVtx... lets forget it.
942
943 } while ( seedTracks.size() > 1 && m_iterations < m_maxVertices ) ;
944
945 ATH_MSG_DEBUG( " Iterations done " ) ;
946
947 //unfortunately you have still a problem with the track to links!!!
948
949 xAOD::VertexContainer::iterator vxBegin=theVertexContainer->begin();
950 xAOD::VertexContainer::iterator vxEnd=theVertexContainer->end();
951
952 //prepare iterators for tracks only necessary for seeding
953 std::vector<Trk::ITrackLink*>::iterator origtrkbegin=origTracks.begin();
954 std::vector<Trk::ITrackLink*>::iterator origtrkend=origTracks.end();
955
956 // refit a vertex if there were tracks eventualy ever removed/robbed from it
957// for (xAOD::VertexContainer::const_iterator vxIter=vxBegin;vxIter!=vxEnd;++vxIter)
958 for (xAOD::VertexContainer::iterator vxIter=vxBegin;vxIter!=vxEnd;)
959 {
960#ifdef MONITORTUNES
961 if ( (*vxIter)->vertexType() == xAOD::VxType::NoVtx || (*vxIter)->vertexType() == xAOD::VxType::KinkVtx )
962 {
963 ++vxIter ;
964 continue ;
965 }
966#endif
967 std::vector<Trk::VxTrackAtVertex>* myVxTracksAtVtx = &((*vxIter)->vxTrackAtVertex());
968 if ( !myVxTracksAtVtx )
969 {
970 ++vxIter ;
971 continue ;
972 }
973
974 int nrobbed = mDecor_nrobbed (**vxIter);
975 if ( nrobbed < 1 )
976 {
977 ++vxIter ;
978 continue ;
979 }
980
981 std::vector<Trk::VxTrackAtVertex>::const_iterator tracksBegin=myVxTracksAtVtx->begin();
982 std::vector<Trk::VxTrackAtVertex>::const_iterator tracksEnd=myVxTracksAtVtx->end();
983
984 std::vector<const Trk::TrackParameters*> perigeesToFit ;
985 std::vector<Trk::ITrackLink*> nullseedTracks ;
986 int ntrk = 0 ;
987 for (std::vector<Trk::VxTrackAtVertex>::const_iterator tracksIter=tracksBegin;
988 tracksIter!=tracksEnd;++tracksIter)
989 {
990 for (std::vector<Trk::ITrackLink*>::const_iterator origtrkiter=origtrkbegin;
991 origtrkiter!=origtrkend;++origtrkiter)
992 {
993 if ((*origtrkiter)->parameters()==(*tracksIter).initialPerigee())
994 {
995 if ( (*tracksIter).weight() > m_minWghtAtVtx ) ntrk ++ ;
996 nullseedTracks.push_back( *origtrkiter ) ;
997 perigeesToFit.push_back( (*origtrkiter)->parameters() ) ;
998 }
999 }
1000 }
1001// redo the fit
1002 xAOD::Vertex * QxAODVertex = new xAOD::Vertex() ;
1003 QxAODVertex = m_iVertexFitter->fit( perigeesToFit, (*vxIter)->position() );
1004 if ( QxAODVertex == nullptr )
1005 {
1006 ++vxIter ;
1007 continue ;
1008 }
1009
1010 removeCompatibleTracks( QxAODVertex, perigeesToFit, nullseedTracks);
1011
1012 float chi2dof1 = (*vxIter)->chiSquared()/(*vxIter)->numberDoF() ;
1013 float chi2dof2 = QxAODVertex->chiSquared()/QxAODVertex->numberDoF() ;
1014
1015 float oldhf = 0. ;
1016 if ( mDecor_HitsFilter.isAvailable (**vxIter) )
1017 oldhf = mDecor_HitsFilter (**vxIter);
1018 if ( chi2dof2 >= chi2dof1 )
1019 {
1020 ++vxIter ;
1021 continue ;
1022 }
1023
1024 int nit = mDecor_intrk (**vxIter);
1025 *(*vxIter) = *QxAODVertex ;
1026
1027 bool isv0 = V0kine( getVertexMomenta( QxAODVertex ), (&(*QxAODVertex))->position(), m_v0mass, m_dir ) ;
1028 mDecor_mass( *(*vxIter) ) = m_v0mass ;
1029 mDecor_energy( *(*vxIter) ) = m_v0ee ;
1030 mDecor_intrk( *(*vxIter) ) = nit ;
1031
1032 if ( isv0 && ntrk == 2 ) (*vxIter)->setVertexType(xAOD::VxType::V0Vtx);
1033
1034#ifdef MONITORTUNES
1035 mDecor_HitsFilter( *(*vxIter) ) = oldhf ;
1036 mDecor_direction( *(*vxIter) ) = m_dir ;
1037 mDecor_NumTrk( *(*vxIter) ) = ntrk ;
1038#else
1039 if ( ( m_filterLevel >= 3 && m_ntracks == 2 && isv0 )
1040 || ( m_filterLevel > 2 && m_dir < m_privtxRef )
1041 )
1042 {
1043 vxIter = theVertexContainer->erase( vxIter ) ;
1044 vxBegin=theVertexContainer->begin();
1045 vxEnd=theVertexContainer->end();
1046 }
1047#endif
1048
1049 ++ vxIter ;
1050
1051
1052 delete QxAODVertex ;
1053 }
1054
1055 // update vertices, while fix the TrackParticles linked to the vertex
1056 int nv = 0 ;
1057
1058 for (xAOD::VertexContainer::iterator vxIter=vxBegin;vxIter!=vxEnd;++vxIter, nv++)
1059 {
1060
1061 ATH_MSG_DEBUG( " filled " << nv << " 'th vertex : x= " << (*vxIter)->position().x() <<" , y= "
1062 << (*vxIter)->position().y() <<" , z= " << (*vxIter)->position().z()
1063 << " vxType = " << (*vxIter)->vertexType() ) ;
1064
1065 std::vector<float> trkWght ;
1066
1067 std::vector<Trk::VxTrackAtVertex>* myVxTracksAtVtx = &((*vxIter)->vxTrackAtVertex());
1068
1069 if ( ! myVxTracksAtVtx )
1070 {
1071#ifdef MONITORTUNES
1072 mDecor_trkDOE( *(*vxIter) ) = trkWght ;
1073 mDecor_trkWght( *(*vxIter) ) = trkWght ;
1074 mDecor_sumPt2( *(*vxIter) ) = -99.9 ;
1075 mDecor_nrobbed( *(*vxIter) ) = 0 ;
1076#endif
1077 continue ;
1078 }
1079
1080 float pt2 = 0.0 ;
1081 std::vector<float> xdoe ;
1082
1083 std::vector<Trk::VxTrackAtVertex>::iterator tracksBegin=myVxTracksAtVtx->begin();
1084 std::vector<Trk::VxTrackAtVertex>::iterator tracksEnd=myVxTracksAtVtx->end();
1085
1086 for (std::vector<Trk::VxTrackAtVertex>::iterator tracksIter=tracksBegin;
1087 tracksIter!=tracksEnd;++tracksIter)
1088 {
1089 bool found = false ;
1090
1091 trkWght.push_back( (*tracksIter).weight() ) ;
1092
1093 //now look for corresponding ITrackLink
1094 for (std::vector<Trk::ITrackLink*>::iterator origtrkiter=origtrkbegin;
1095 origtrkiter!=origtrkend;++origtrkiter)
1096 {
1097
1098 if ((*origtrkiter)->parameters()==(*tracksIter).initialPerigee())
1099 {
1100
1101 found=true;
1102
1103 const Trk::TrackParameters * svperigee = (*tracksIter).perigeeAtVertex() ;
1104
1105 if ( svperigee )
1106 {
1107 ATH_MSG_VERBOSE( " svperigee retrieved " );
1108
1109 float pt = svperigee->parameters()[Trk::theta]/std::abs( svperigee->parameters()[Trk::qOverP] ) ;
1110 ATH_MSG_VERBOSE( " track perigee parameters retrieved " );
1111 pt2 += pt*pt*0.000001 ;
1112
1113#ifdef MONITORTUNES
1114 double distance=0.;
1115 try
1116 {
1117 std::unique_ptr<Trk::PlaneSurface> mySurface=
1118 m_ImpactPoint3dEstimator->Estimate3dIP( (*tracksIter).initialPerigee(), &((*vxIter)->position()), distance );
1119 }
1120 catch (error::ImpactPoint3dEstimatorProblem err)
1121 {
1122 ATH_MSG_WARNING( " ImpactPoint3dEstimator failed " << err.p );
1123 }
1124
1125 double error= 1. ;
1126
1127 if( svperigee->covariance() )
1128 error = std::sqrt( (*svperigee->covariance())(Trk::d0,Trk::d0)
1129 + (*svperigee->covariance())(Trk::z0,Trk::z0)
1130 );
1131 xdoe.push_back( std::abs( distance/error ) ) ;
1132#endif
1133 }
1134
1135 //assigning the input track to the fitted vertex through VxTrackAtVertices vector
1136 (*tracksIter).setOrigTrack (*origtrkiter );
1137
1138 // See if the trklink is to an xAOD::TrackParticle
1139 Trk::LinkToXAODTrackParticle* linkToXAODTP = dynamic_cast<Trk::LinkToXAODTrackParticle*>( *origtrkiter );
1140
1141 // If track is an xAOD::TrackParticle, set directly in xAOD::Vertex
1142 if (linkToXAODTP)
1143 {
1144 (*vxIter)->addTrackAtVertex(*linkToXAODTP, (*tracksIter).weight());
1145
1146 } //TODO: else write in a warning? (if tracks were Trk::Tracks or Trk::TrackParticleBase)
1147
1148 origTracks.erase(origtrkiter);
1149 origtrkbegin=origTracks.begin();
1150 origtrkend=origTracks.end();
1151 break;
1152 } // end of matched
1153 } // end of matching loop over original trakcs
1154
1155 if (!found)
1156 {
1157 ATH_MSG_ERROR( " Cannot find vector element to fix links (step 4)! " );
1158 }
1159
1160 }//end iterate on tracks at vtx
1161
1162// update the decorations
1163
1164 const xAOD::Vertex * vtx = (*vxIter);
1165
1166#ifdef MONITORTUNES
1167 mDecor_trkWght( *vtx ) = trkWght ;
1168 mDecor_trkDOE( *vtx ) = xdoe ;
1169#endif
1170
1171 mDecor_sumPt2( *vtx ) = pt2 ;
1172
1173 }//end iterate on vertices
1174
1175 ATH_MSG_DEBUG(" #Vtx "<< theVertexContainer->size() <<" with track-vertex association fixed " );
1176
1177
1178 for (unsigned int i = 0 ; i < theVertexContainer->size() ; i++)
1179 {
1180
1181 ATH_MSG_VERBOSE( " Vtx: " << i <<
1182 " x= " << (*theVertexContainer)[i]->position().x() <<
1183 " y= " << (*theVertexContainer)[i]->position().y() <<
1184 " z= " << (*theVertexContainer)[i]->position().z() <<
1185 " ntracks= " << (*theVertexContainer)[i]->vxTrackAtVertex().size() <<
1186 " chi2= " << (*theVertexContainer)[i]->chiSquared() <<
1187 " #dof = " << (*theVertexContainer)[i]->numberDoF() );
1188 }
1189
1190 for (std::vector<Trk::ITrackLink*>::iterator origtrkiter=origtrkbegin;
1191 origtrkiter!=origtrkend;++origtrkiter)
1192 {
1193 if ((*origtrkiter)!=0) {
1194 delete *origtrkiter;
1195 *origtrkiter=0;
1196 }
1197 }
1198
1199#ifdef MONITORTUNES
1200 m_OTree->Fill() ;
1201#endif
1202
1203 return std::make_pair(theVertexContainer, theVertexAuxContainer);
1204}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
Eigen::Matrix< double, 3, 1 > Vector3D
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
#define y
#define x
#define z
value_type push_back(value_type pElem)
Add an element to the end of the collection.
DataModel_detail::iterator< DataVector > iterator
Definition DataVector.h:842
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
iterator erase(iterator position)
Remove element at a given position.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
size_type size() const noexcept
Returns the number of elements in the collection.
void countTracksAndNdf(xAOD::Vertex *myxAODVertex, float &ndf, int &ntracks) const
double compatibility(const Trk::TrackParameters &measPerigee, const AmgSymMatrix(3) &covariancePosition, const Amg::Vector3D &position) const
ToolHandle< Trk::IImpactPoint3dEstimator > m_ImpactPoint3dEstimator
ToolHandle< Trk::AdaptiveVertexFitter > m_iVertexFitter
bool V0kine(const std::vector< Amg::Vector3D > &, const Amg::Vector3D &position, float &, float &) const
ToolHandle< Trk::IVertexSeedFinder > m_SeedFinder
Gaudi::Property< int > m_splitVerticesTrkInvFraction
Integer: 1./fraction of tracks to be assigned to the tag split vertex.
void removeCompatibleTracks(xAOD::Vertex *myxAODVertex, std::vector< const Trk::TrackParameters * > &perigeesToFit, std::vector< Trk::ITrackLink * > &seedTracks) const
const std::vector< Amg::Vector3D > getVertexMomenta(xAOD::Vertex *myxAODVertex) const
static double VrtVrtDist(xAOD::Vertex *v1, xAOD::Vertex *v2)
void removeAllFrom(std::vector< const Trk::TrackParameters * > &perigeesToFit, std::vector< Trk::ITrackLink * > &seedTracks) const
float removeTracksInBadSeed(xAOD::Vertex *myxAODVertex, std::vector< const Trk::TrackParameters * > &) const
SG::Decorator< T, ALLOC > Decorator
Definition AuxElement.h:575
void setCovariancePosition(const AmgSymMatrix(3)&covariancePosition)
Sets the vertex covariance matrix.
void setVertexType(VxType::VertexType vType)
Set the type of the vertex.
VxType::VertexType vertexType() const
The type of the vertex.
float chiSquared() const
Returns the of the vertex fit as float.
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.
float distance(const Amg::Vector3D &p1, const Amg::Vector3D &p2)
calculates the distance between two point in 3D space
float chiSquared(const U &p)
@ theta
Definition ParamDefs.h:66
@ qOverP
perigee
Definition ParamDefs.h:67
@ d0
Definition ParamDefs.h:63
@ z0
Definition ParamDefs.h:64
ParametersBase< TrackParametersDim, Charged > TrackParameters
@ 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".

◆ findVertex() [3/5]

std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > InDet::InDetIterativeSecVtxFinderTool::findVertex ( const TrackCollection * trackTES)
override

Definition at line 1577 of file InDetIterativeSecVtxFinderTool.cxx.

1578{
1579
1580 ATH_MSG_DEBUG( " Number of input tracks before track selection: " << trackTES->size() );
1581
1582 std::vector<Trk::ITrackLink*> selectedTracks;
1583
1584 using TrackDataVecIter = DataVector<Trk::Track>::const_iterator;
1585
1586 bool selectionPassed;
1587 for (TrackDataVecIter itr = (*trackTES).begin(); itr != (*trackTES).end(); ++itr) {
1588 Trk::Vertex null(Amg::Vector3D(0,0,0));
1589 selectionPassed=static_cast<bool>(m_trkFilter->accept(**itr,&null));
1590 if ( selectionPassed ) selectionPassed=static_cast<bool>(m_SVtrkFilter->accept(**itr,&null));
1591 if (selectionPassed)
1592 {
1593 ElementLink<TrackCollection> link;
1594 link.setElement(*itr);
1595 Trk::LinkToTrack * linkTT = new Trk::LinkToTrack(link);
1596 linkTT->setStorableObject(*trackTES);
1597 selectedTracks.push_back(linkTT);
1598 }
1599 }
1600
1601 ATH_MSG_DEBUG( "Of " << trackTES->size() << " tracks "
1602 << selectedTracks.size() << " survived the preselection." );
1603
1604
1605 std::pair<xAOD::VertexContainer*, xAOD::VertexAuxContainer*> returnContainers
1606 = findVertex( selectedTracks );
1607
1608 return returnContainers;
1609
1610}
DataModel_detail::const_iterator< DataVector > const_iterator
Standard const_iterator.
Definition DataVector.h:838

◆ findVertex() [4/5]

std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > InDet::InDetIterativeSecVtxFinderTool::findVertex ( const Trk::TrackParticleBaseCollection * trackTES)
override

Definition at line 1612 of file InDetIterativeSecVtxFinderTool.cxx.

1612 {
1613
1614 ATH_MSG_DEBUG( " Number of input tracks before track selection: " << trackTES->size() );
1615
1616 std::vector<Trk::ITrackLink*> selectedTracks;
1617
1618
1619 using TrackParticleDataVecIter = DataVector<Trk::TrackParticleBase>::const_iterator;
1620
1621 bool selectionPassed;
1622 for (TrackParticleDataVecIter itr = (*trackTES).begin(); itr != (*trackTES).end(); ++itr) {
1623
1624 Trk::Vertex null(Amg::Vector3D(0,0,0));
1625 selectionPassed=static_cast<bool>(m_trkFilter->accept( *((*itr)->originalTrack()), &null));
1626 if ( selectionPassed ) selectionPassed =static_cast<bool>( m_SVtrkFilter->accept( *((*itr)->originalTrack()), &null));
1627
1628 if (selectionPassed)
1629 {
1630 ElementLink<Trk::TrackParticleBaseCollection> link;
1631 link.setElement(*itr);
1632 Trk::LinkToTrackParticleBase * linkTT = new Trk::LinkToTrackParticleBase(link);
1633 linkTT->setStorableObject(*trackTES);
1634 selectedTracks.push_back(linkTT);
1635 }
1636 }
1637
1638 ATH_MSG_DEBUG( "Of " << trackTES->size() << " tracks "
1639 << selectedTracks.size() << " survived the preselection." );
1640
1641 std::pair<xAOD::VertexContainer*, xAOD::VertexAuxContainer*> returnContainers
1642 = findVertex( selectedTracks );
1643
1644 return returnContainers;
1645
1646}

◆ findVertex() [5/5]

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

Finding method.

Has as input a track collection and as output a VxContainer.

Definition at line 76 of file InDetIterativeSecVtxFinderTool.cxx.

77{
78 ATH_MSG_DEBUG(" Number of input tracks before track selection: " << trackParticles->size());
79
80 m_evtNum ++ ;
81
82 std::vector<Trk::ITrackLink*> selectedTracks;
83
84
85 using TrackParticleDataVecIter = DataVector<xAOD::TrackParticle>::const_iterator;
86
87 bool selectionPassed;
88 m_trkdefiPars.clear() ;
89 xAOD::Vertex null;
90 null.makePrivateStore();
91 null.setPosition(Amg::Vector3D(0,0,0));
92 AmgSymMatrix(3) vertexError;
93 vertexError.setZero();
94 null.setCovariancePosition(vertexError);
95 for (TrackParticleDataVecIter itr = trackParticles->begin(); itr != trackParticles->end(); ++itr) {
96
97 selectionPassed=static_cast<bool>(m_trkFilter->accept(**itr,&null));
98 if ( selectionPassed ) selectionPassed =static_cast<bool>(m_SVtrkFilter->accept(**itr,&null));
99 if (selectionPassed)
100 {
101 Amg::VectorX par = (*itr)->definingParameters();
102 par[0] = 1.0*(*itr)->hitPattern() ;
103 m_trkdefiPars.push_back( par ) ;
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
117 std::pair<xAOD::VertexContainer*, xAOD::VertexAuxContainer*> returnContainers
118 =findVertex( selectedTracks );
119
120 return returnContainers;
121
122}

◆ getModes1d()

int InDet::InDetIterativeSecVtxFinderTool::getModes1d ( std::vector< int > * ,
std::vector< int > * ,
std::vector< int > *  ) const

◆ getVertexMomenta()

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

Definition at line 1879 of file InDetIterativeSecVtxFinderTool.cxx.

1880{
1881 std::vector< Amg::Vector3D > TrkAtVtxMomenta ;
1882
1883 std::vector<Trk::VxTrackAtVertex>* tracksAtVertex= &(myxAODVertex->vxTrackAtVertex());
1884
1885 ATH_MSG_DEBUG( " getVertexMomenta ... #Tracks associated at vertex : " << tracksAtVertex->size() );
1886
1887 std::vector<Trk::VxTrackAtVertex>::const_iterator tracksAtVertexBegin=tracksAtVertex->begin();
1888 std::vector<Trk::VxTrackAtVertex>::const_iterator tracksAtVertexEnd=tracksAtVertex->end();
1889
1890 for (std::vector<Trk::VxTrackAtVertex>::const_iterator tracksAtVertexIter=tracksAtVertexBegin;
1891 tracksAtVertexIter!=tracksAtVertexEnd;++tracksAtVertexIter)
1892 {
1893 if ((*tracksAtVertexIter).weight() <= m_minWghtAtVtx ) continue ;
1894 {
1895 const Trk::TrackParameters* sv_perigee = (*tracksAtVertexIter).perigeeAtVertex() ;
1896
1897 double qp = 1/std::abs( sv_perigee->parameters()[Trk::qOverP] ) ;
1898 double theta = sv_perigee->parameters()[Trk::theta];
1899 double phi = sv_perigee->parameters()[Trk::phi];
1900
1901 m_TrkAtVtxMomenta.emplace_back( qp*sin(theta)*cos(phi), qp*sin(theta)*sin(phi), qp*cos(theta) ) ;
1902
1903 }
1904
1905 }
1906
1907 return m_TrkAtVtxMomenta ;
1908}
Scalar phi() const
phi method
Scalar theta() const
theta method
@ phi
Definition ParamDefs.h:75

◆ initialize()

StatusCode InDet::InDetIterativeSecVtxFinderTool::initialize ( )
override

Definition at line 1508 of file InDetIterativeSecVtxFinderTool.cxx.

1509{
1510 StatusCode sc;
1511
1512 m_evtNum = 0 ;
1513 ATH_CHECK(m_iVertexFitter.retrieve());
1514 ATH_CHECK(m_SeedFinder.retrieve());
1517
1518 ATH_CHECK(m_trkFilter.retrieve());
1519 ATH_CHECK(m_SVtrkFilter.retrieve());
1520
1521 // since some parameters special to an inherited class this method
1522 // will be overloaded by the inherited class
1524
1525#ifdef MONITORTUNES
1526 SmartIF<ITHistSvc> hist_root{Gaudi::svcLocator()->service("THistSvc")};
1527 ATH_CHECK(hist_root.isValid());
1528
1529 m_leastmodes = new std::vector<int>() ;
1530 m_sdFsmwX = new std::vector< std::vector < float > >() ;
1531 m_sdFsmwY = new std::vector< std::vector < float > >() ;
1532 m_sdFsmwZ = new std::vector< std::vector < float > >() ;
1533 m_sdcrsWght = new std::vector< std::vector < float > >() ;
1534
1535 m_nperiseed = new std::vector<int>() ;
1536 m_seedX = new std::vector < float >() ;
1537 m_seedY = new std::vector < float >() ;
1538 m_seedZ = new std::vector < float >() ;
1539 m_seedXYdist = new std::vector< float >() ;
1540 m_seedZdist = new std::vector< float >() ;
1541 m_seedac = new std::vector < int >() ;
1542
1543 m_OTree = new TTree( "IncSecVtxUnder", "TTree of underlying/upstream techinfo for InclusiveSecVtx" ) ;
1544
1545 // the number of event processed ( calling of this tool )
1546 m_OTree->Branch("EvtNum", &m_evtNum, "EvtNum/I" ) ;
1547
1548 m_OTree->Branch("IteraNum", &m_iterations, "IteraNum/I" ) ;
1549
1550 // the minumal number of FsmwModes ( candidates ) for each seed
1551 // for each searching one has one or more modes, data type of vector<vecvtot<float>> :
1552 m_OTree->Branch("ModeNum", &m_leastmodes ) ;
1553 m_OTree->Branch("sdFsmwX", &m_sdFsmwX ) ;
1554 m_OTree->Branch("sdFsmwY", &m_sdFsmwY ) ;
1555 m_OTree->Branch("sdFsmwZ", &m_sdFsmwZ ) ;
1556 m_OTree->Branch("crossingWght", &m_sdcrsWght ) ;
1557 // seeds as candidate ( number of trials ) of SecVtx
1558 m_OTree->Branch("nperigeeseed", &m_nperiseed ) ;
1559 m_OTree->Branch("seedX", &m_seedX ) ;
1560 m_OTree->Branch("seedY", &m_seedY ) ;
1561 m_OTree->Branch("seedZ", &m_seedZ ) ;
1562 m_OTree->Branch("seedXYdist", &m_seedXYdist ) ;
1563 m_OTree->Branch("seedZdist", &m_seedZdist ) ;
1564 m_OTree->Branch("seedaccepted", &m_seedac ) ;
1565
1566 // seeds filter
1567
1568 // seeds compatibility ( after fit )
1569
1570 ATH_CHECK( hist_root->regTree("/AANT/SecVtxIncunder", m_OTree ) ) ;
1571#endif
1572
1573 ATH_MSG_DEBUG( "Initialization successful" );
1574 return StatusCode::SUCCESS;
1575}
#define ATH_CHECK
Evaluate an expression and check for errors.
static Double_t sc
::StatusCode StatusCode
StatusCode definition for legacy code.

◆ passHitsFilter()

bool InDet::InDetIterativeSecVtxFinderTool::passHitsFilter ( const Trk::TrackParameters * perigee,
float vtxR,
float absvz ) const
private

Borrowed from Reconstruction/VKalVrt/VrtSecInclusive

Definition at line 1299 of file InDetIterativeSecVtxFinderTool.cxx.

1301{
1302 bool pass = true ;
1303
1304 // PrimaryVtx information is used, who was assumpted be always be there
1305 ATH_MSG_DEBUG( " seed input : " << rad <<" "<< absz );
1306
1308
1309 //
1310 // rough guesses for active layers:
1311 // BeamPipe: 23.5-24.3
1312 // IBL: 31.0-38.4
1313 // Pix0 (BLayer): 47.7-54.4, Pix1: 85.5-92.2, Pix2: 119.3-126.1
1314 // Sct0: 290-315, Sct1: 360-390, Sct2: 430-460, Sct3:500-530
1315 //
1316
1317 enum vertexArea {
1318 insideBeamPipe,
1319
1320 insidePixelBarrel0,
1321 aroundPixelBarrel0,
1322
1323 outsidePixelBarrel0_and_insidePixelBarrel1,
1324 aroundPixelBarrel1,
1325
1326 outsidePixelBarrel1_and_insidePixelBarrel2,
1327 aroundPixelBarrel2,
1328
1329 outsidePixelBarrel2_and_insidePixelBarrel3,
1330 aroundPixelBarrel3,
1331
1332 outsidePixelBarrel3_and_insideSctBarrel0,
1333 aroundSctBarrel0,
1334
1335 outsideSctBarrel0_and_insideSctBarrel1,
1336 aroundSctBarrel1,
1337
1338 insideSilicon
1339 };
1340
1341 int vertex_pattern = 0;
1342 if( rad < 23.50 ) {
1343 vertex_pattern = insideBeamPipe;
1344 } else if( rad < 31.0 && absz < 331.5 ) {
1345 vertex_pattern = insidePixelBarrel0;
1346 } else if( rad < 38.4 && absz < 331.5 ) {
1347 vertex_pattern = aroundPixelBarrel0;
1348 } else if( rad < 47.7 && absz < 400.5 ) {
1349 vertex_pattern = outsidePixelBarrel0_and_insidePixelBarrel1;
1350 } else if( rad < 54.4 && absz < 400.5 ) {
1351 vertex_pattern = aroundPixelBarrel1;
1352 } else if( rad < 85.5 && absz < 400.5 ) {
1353 vertex_pattern = outsidePixelBarrel1_and_insidePixelBarrel2;
1354 } else if( rad < 92.2 && absz < 400.5 ) {
1355 vertex_pattern = aroundPixelBarrel2;
1356 } else if( rad < 119.3 && absz < 400.5 ) {
1357 vertex_pattern = outsidePixelBarrel2_and_insidePixelBarrel3;
1358 } else if( rad < 126.1 && absz < 400.5 ) {
1359 vertex_pattern = aroundPixelBarrel3;
1360 } else if( rad < 290 && absz < 749.0 ) {
1361 vertex_pattern = outsidePixelBarrel3_and_insideSctBarrel0;
1362 } else if( rad < 315 && absz < 749.0 ) {
1363 vertex_pattern = aroundSctBarrel0;
1364 } else if( rad < 360 && absz < 749.0 ) {
1365 vertex_pattern = outsideSctBarrel0_and_insideSctBarrel1;
1366 } else if( rad < 390 && absz < 749.0 ) {
1367 vertex_pattern = aroundSctBarrel1;
1368 } else {
1369 vertex_pattern = insideSilicon ;
1370 }
1371
1372
1373
1374 for ( unsigned int p = 0 ; p < m_trkdefiPars.size() ; p++ )
1375 {
1376
1377 Amg::VectorX tpperigee = m_trkdefiPars.at(p) ;
1378
1379 if ( tpperigee[4] != perigee->parameters()[Trk::qOverP]
1380 || tpperigee[3] != perigee->parameters()[Trk::theta]
1381 || tpperigee[2] != perigee->parameters()[Trk::phi]
1382 || tpperigee[1] != perigee->parameters()[Trk::z0]
1383
1384 ) continue ;
1385// for each perigee only single element in trkdefiPars is expected
1386
1387 uint32_t HitPattern= tpperigee[0] ;
1388 ATH_MSG_VERBOSE( " Track HitPattern : " << HitPattern );
1389
1390
1391
1392// shift 1 left with bits of Trk::InnerDet_LayerX_DiscY to define a mask
1393// int mask = 1<<Trk::InnerDet_LayerX_DiscY ;
1394// check whether this bit is set in the pattern
1395// bool set = HitPattern&mask ;
1396
1397 if( vertex_pattern == insideBeamPipe ) {
1398
1399 if( ! (HitPattern & (1<<Trk::pixelBarrel0)) ) return false;
1400
1401 } else if( vertex_pattern == insidePixelBarrel0 ) {
1402
1403 if( ! (HitPattern & (1<<Trk::pixelBarrel0)) ) return false;
1404 }
1405
1406 else if( vertex_pattern == aroundPixelBarrel0 ) {
1407 // require nothing for PixelBarrel0
1408 if( ! (HitPattern & (1<<Trk::pixelBarrel1)) ) return false;
1409 }
1410
1411 else if( vertex_pattern == outsidePixelBarrel0_and_insidePixelBarrel1 ) {
1412
1413 if( (HitPattern & (1<<Trk::pixelBarrel0)) ) return false;
1414 if( ! (HitPattern & (1<<Trk::pixelBarrel1)) ) return false;
1415 }
1416
1417 else if( vertex_pattern == aroundPixelBarrel1 ) {
1418
1419 if( (HitPattern & (1<<Trk::pixelBarrel0)) ) return false;
1420 // require nothing for PixelBarrel1
1421 if( ! (HitPattern & (1<<Trk::pixelBarrel2)) ) return false;
1422 }
1423
1424 else if( vertex_pattern == outsidePixelBarrel1_and_insidePixelBarrel2 ) {
1425
1426 if( (HitPattern & (1<<Trk::pixelBarrel0)) ) return false;
1427 if( (HitPattern & (1<<Trk::pixelBarrel1)) ) return false;
1428 if( ! (HitPattern & (1<<Trk::pixelBarrel2)) ) return false;
1429 }
1430
1431 else if( vertex_pattern == aroundPixelBarrel2 ) {
1432
1433 if( (HitPattern & (1<<Trk::pixelBarrel0)) ) return false;
1434 if( (HitPattern & (1<<Trk::pixelBarrel1)) ) return false;
1435 // require nothing for PixelBarrel2
1436 if( ! (HitPattern & (1<<Trk::pixelBarrel3)) ) return false;
1437 }
1438
1439 else if( vertex_pattern == outsidePixelBarrel2_and_insidePixelBarrel3 ) {
1440
1441 if( (HitPattern & (1<<Trk::pixelBarrel0)) ) return false;
1442 if( (HitPattern & (1<<Trk::pixelBarrel1)) ) return false;
1443 if( (HitPattern & (1<<Trk::pixelBarrel2)) ) return false;
1444 if( ! (HitPattern & (1<<Trk::pixelBarrel3)) ) return false;
1445 }
1446
1447 else if( vertex_pattern == aroundPixelBarrel3 ) {
1448
1449 if( (HitPattern & (1<<Trk::pixelBarrel0)) ) return false;
1450 if( (HitPattern & (1<<Trk::pixelBarrel1)) ) return false;
1451 if( (HitPattern & (1<<Trk::pixelBarrel2)) ) return false;
1452 // require nothing for PixelBarrel3
1453 if( ! (HitPattern & (1<<Trk::sctBarrel0)) ) return false;
1454 }
1455
1456 else if( vertex_pattern == outsidePixelBarrel3_and_insideSctBarrel0 ) {
1457
1458 if( (HitPattern & (1<<Trk::pixelBarrel0)) ) return false;
1459 if( (HitPattern & (1<<Trk::pixelBarrel1)) ) return false;
1460 if( (HitPattern & (1<<Trk::pixelBarrel2)) ) return false;
1461 if( (HitPattern & (1<<Trk::pixelBarrel3)) ) return false;
1462 if( ! (HitPattern & (1<<Trk::sctBarrel0)) ) return false;
1463 }
1464
1465 else if( vertex_pattern == aroundSctBarrel0 ) {
1466
1467 if( (HitPattern & (1<<Trk::pixelBarrel0)) ) return false;
1468 if( (HitPattern & (1<<Trk::pixelBarrel1)) ) return false;
1469 if( (HitPattern & (1<<Trk::pixelBarrel2)) ) return false;
1470 if( (HitPattern & (1<<Trk::pixelBarrel3)) ) return false;
1471 // require nothing for SctBarrel0
1472 if( ! (HitPattern & (1<<Trk::sctBarrel1)) ) return false;
1473 }
1474
1475 else if( vertex_pattern == outsideSctBarrel0_and_insideSctBarrel1 ) {
1476
1477 if( (HitPattern & (1<<Trk::pixelBarrel0)) ) return false;
1478 if( (HitPattern & (1<<Trk::pixelBarrel1)) ) return false;
1479 if( (HitPattern & (1<<Trk::pixelBarrel2)) ) return false;
1480 if( (HitPattern & (1<<Trk::pixelBarrel3)) ) return false;
1481 if( (HitPattern & (1<<Trk::sctBarrel0)) ) return false;
1482 if( ! (HitPattern & (1<<Trk::sctBarrel1)) ) return false;
1483 }
1484
1485 else if( vertex_pattern == aroundSctBarrel1 ) {
1486 if( (HitPattern & (1<<Trk::pixelBarrel0)) ) return false;
1487 if( (HitPattern & (1<<Trk::pixelBarrel1)) ) return false;
1488 if( (HitPattern & (1<<Trk::pixelBarrel2)) ) return false;
1489 if( (HitPattern & (1<<Trk::pixelBarrel3)) ) return false;
1490 if( (HitPattern & (1<<Trk::sctBarrel0)) ) return false;
1491 // require nothing for SctBarrel1
1492 if( ! (HitPattern & (1<<Trk::sctBarrel2)) ) return false;
1493 } else
1494 pass = true ;
1495 }
1496
1497 ATH_MSG_DEBUG( " passHitsFilter : " << pass );
1498
1499 return pass ;
1500}
@ pixelBarrel0
there are three or four pixel barrel layers (R1/R2)
setEventNumber uint32_t

◆ printParameterSettings()

void InDet::InDetIterativeSecVtxFinderTool::printParameterSettings ( )
privatevirtual

Definition at line 1685 of file InDetIterativeSecVtxFinderTool.cxx.

1686{
1687 ATH_MSG_DEBUG( "InDetIterativeSecVtxFinderTool initialize(): Parametersettings " );
1688 ATH_MSG_DEBUG( "VertexFitter " << m_iVertexFitter );
1689}

◆ removeAllFrom()

void InDet::InDetIterativeSecVtxFinderTool::removeAllFrom ( std::vector< const Trk::TrackParameters * > & perigeesToFit,
std::vector< Trk::ITrackLink * > & seedTracks ) const
private

Definition at line 1723 of file InDetIterativeSecVtxFinderTool.cxx.

1725{
1726 //remove all perigeesToFit and go on...
1727
1728 std::vector<Trk::ITrackLink*>::iterator seedBegin=seedTracks.begin();
1729 std::vector<Trk::ITrackLink*>::iterator seedEnd=seedTracks.end();
1730
1731 std::vector<const Trk::TrackParameters*>::iterator perigeesToFitBegin=perigeesToFit.begin();
1732 std::vector<const Trk::TrackParameters*>::iterator perigeesToFitEnd=perigeesToFit.end();
1733
1734 for (std::vector<const Trk::TrackParameters*>::iterator perigeesToFitIter=perigeesToFitBegin;
1735 perigeesToFitIter!=perigeesToFitEnd;++perigeesToFitIter)
1736 {
1737
1738 ATH_MSG_VERBOSE( " Iterating on new track in original perigeesToFit list of BAD VERTEX..." );
1739
1740
1741 bool found=false;
1742
1743 for (std::vector<Trk::ITrackLink*>::iterator seedIter=seedTracks.begin();
1744 seedIter!=seedEnd;++seedIter)
1745 {
1746 if ((*seedIter)->parameters()==*perigeesToFitIter)
1747 {
1748
1749 ATH_MSG_VERBOSE( " found and deleted from seeds!" );
1750
1751 found=true;
1752 seedTracks.erase(seedIter);
1753 seedBegin=seedTracks.begin();
1754 seedEnd=seedTracks.end();
1755 break;
1756 }
1757 }
1758
1759 if (!found)
1760 {
1761 ATH_MSG_ERROR( " Cannot find vector element to delete when removing BAD vertex! " );
1762 }
1763
1764 }//for cycle
1765
1766}

◆ removeCompatibleTracks()

void InDet::InDetIterativeSecVtxFinderTool::removeCompatibleTracks ( xAOD::Vertex * myxAODVertex,
std::vector< const Trk::TrackParameters * > & perigeesToFit,
std::vector< Trk::ITrackLink * > & seedTracks ) const
private

Definition at line 1910 of file InDetIterativeSecVtxFinderTool.cxx.

1913{
1914 //now you have your new vertex with its new tracks
1915 //now you have to get the compatibility also of all tracks which DIDN'T BELONG to the vertex!
1916 std::vector<Trk::VxTrackAtVertex>* tracksAtVertex= &(myxAODVertex->vxTrackAtVertex());
1917
1918 ATH_MSG_VERBOSE( " removeCompatibleTracks ... #Tracks associated at vertex : " << tracksAtVertex->size() );
1919
1920 std::vector<Trk::VxTrackAtVertex>::const_iterator tracksAtVertexBegin=tracksAtVertex->begin();
1921 std::vector<Trk::VxTrackAtVertex>::const_iterator tracksAtVertexEnd=tracksAtVertex->end();
1922
1923 std::vector<Trk::ITrackLink*>::iterator seedBegin=seedTracks.begin();
1924 std::vector<Trk::ITrackLink*>::iterator seedEnd=seedTracks.end();
1925
1926 std::vector<const Trk::TrackParameters*>::iterator perigeesToFitBegin=perigeesToFit.begin();
1927 std::vector<const Trk::TrackParameters*>::iterator perigeesToFitEnd=perigeesToFit.end();
1928
1929 const AmgSymMatrix(3) covariance = myxAODVertex->covariancePosition() ;
1930 const Amg::Vector3D position = myxAODVertex->position() ;
1931
1932 for (std::vector<Trk::VxTrackAtVertex>::const_iterator tracksAtVertexIter=tracksAtVertexBegin;
1933 tracksAtVertexIter!=tracksAtVertexEnd;++tracksAtVertexIter)
1934 {
1935
1936
1937 ATH_MSG_VERBOSE( " new track..." );
1938
1939
1940 bool found=false;
1941
1942 for (std::vector<Trk::ITrackLink*>::iterator seedIter=seedBegin;
1943 seedIter!=seedEnd;++seedIter)
1944 {
1945 if ((*seedIter)->parameters()==(*tracksAtVertexIter).initialPerigee() )
1946 {
1947 found=true;
1948 if ((*tracksAtVertexIter).weight()> m_minWghtAtVtx )
1949 {
1950
1951 ATH_MSG_VERBOSE( " found and deleted from seeds!" );
1952
1953 seedTracks.erase(seedIter);
1954 seedBegin=seedTracks.begin();
1955 seedEnd=seedTracks.end();
1956 }
1957 break;
1958 }
1959 }
1960
1961 if (!found)
1962 {
1963 ATH_MSG_ERROR( " Cannot find vector element to delete (step 1)! " );
1964 }
1965
1966 found=false;
1967 for (std::vector<const Trk::TrackParameters*>::iterator perigeesToFitIter=perigeesToFitBegin;
1968 perigeesToFitIter!=perigeesToFitEnd;++perigeesToFitIter)
1969 {
1970 if (*perigeesToFitIter==(*tracksAtVertexIter).initialPerigee())
1971 {
1972 found=true;
1973
1974 if ((*tracksAtVertexIter).weight()> m_minWghtAtVtx )
1975 {
1976
1977 ATH_MSG_VERBOSE( " found and deleted from perigeesToFit!" );
1978
1979 perigeesToFit.erase(perigeesToFitIter);
1980 perigeesToFitBegin=perigeesToFit.begin();
1981 perigeesToFitEnd=perigeesToFit.end();
1982
1983 }
1984 break;
1985 }
1986 }
1987
1988#ifndef MONITORTUNES
1989 if (!found)
1990 {
1991 ATH_MSG_WARNING( " Cannot find vector element to delete (step 2)! " );
1992 }
1993#endif
1994
1995 }//finishing iterating on tracks at vertex
1996
1997
1998 ATH_MSG_DEBUG( " Outliers still to be considered: " << perigeesToFit.size() );
1999
2000 ATH_MSG_VERBOSE( "Number of seedtracks after removal of inliers: " << seedTracks.size() );
2001
2002
2003
2004 std::vector<Trk::VxTrackAtVertex>* myVxTracksAtVertex= &(myxAODVertex->vxTrackAtVertex());
2005
2006 std::vector<Trk::VxTrackAtVertex>::iterator tracksBegin=myVxTracksAtVertex->begin();
2007 std::vector<Trk::VxTrackAtVertex>::iterator tracksEnd=myVxTracksAtVertex->end();
2008
2009 for (std::vector<const Trk::TrackParameters*>::iterator perigeesToFitIter=perigeesToFitBegin;
2010 perigeesToFitIter!=perigeesToFitEnd;++perigeesToFitIter)
2011 {
2012
2013 bool found=false;
2014
2015 //compute the chi2 with respect to the last fitted vertex!
2016 //(easy since track was NOT used in the last vertex fit...)
2017
2018 const Trk::TrackParameters* myPerigee=(*perigeesToFitIter);
2019
2020 if (myPerigee==nullptr)
2021 {
2022 ATH_MSG_ERROR( " Cast to perigee gives 0 pointer " );
2023 return;
2024
2025 }
2026
2027 double chi2=compatibility(*myPerigee, covariance, position );
2028
2029
2030 ATH_MSG_VERBOSE( "Compatibility is : " << chi2 );
2031
2032
2033 //check if still sufficiently compatible to previous vertex
2034 //(CAN BE VERY LOOSE TO BE CONSERVATIVE AGAINST FAR OUTLIERS)
2036 {
2037
2038 ATH_MSG_VERBOSE( " Found track with compatibility: " << chi2 <<
2039 " to be removed from the seeds... " );
2040
2041
2042// These seed tracks did NOT participate the fitting of vertex, but they seem compatible...
2043 for (std::vector<Trk::ITrackLink*>::iterator seedIter=seedTracks.begin();
2044 seedIter!=seedEnd;++seedIter)
2045 {
2046 if ((*seedIter)->parameters()==*perigeesToFitIter)
2047 {
2048
2049 ATH_MSG_VERBOSE( " found and deleted from seeds!" );
2050
2051 found=true;
2052 seedTracks.erase(seedIter);
2053 seedBegin=seedTracks.begin();
2054 seedEnd=seedTracks.end();
2055 break;
2056 }
2057 }
2058
2059 if (!found)
2060 {
2061 ATH_MSG_ERROR( " Cannot find vector element to delete (step 3)! " );
2062 }
2063 }
2064 else
2065 {
2066 //look if the track is attached to the vertex. If this is the case you should
2067 //delete the track from the vertex!
2068
2069
2070 ATH_MSG_VERBOSE( " Found track with compatibility: " << chi2 <<
2071 " to be further considered and thus to be removed from previous vertex if it was there... " );
2072
2073
2074 bool found=false;
2075
2076// some tracks seem INcompatible/INconsistet even though they succeed the fit
2077 for (std::vector<Trk::VxTrackAtVertex>::iterator tracksIter=tracksBegin;
2078 tracksIter!=tracksEnd;++tracksIter)
2079 {
2080 if ((*tracksIter).initialPerigee()==*perigeesToFitIter)
2081 {
2082
2083
2084 ATH_MSG_VERBOSE( " OK, removing track with compatibility:" << (*tracksIter).trackQuality().chiSquared() <<
2085 " or vtx compatibility" << (*tracksIter).vtxCompatibility() << " which was found attached to the vertex... " );
2086
2087 // this delete is no longer needed because objects in myVxTracksAtVertex are no longer pointers - memory deletion of this VxTrackAtVertex
2088 // was already taken care of inside fitter the moment the VxTrackAtVertex was added to the vector stored in xAOD::Vertex
2089 // -David S.
2090 //
2091 //delete *tracksIter; // delete has to happen BEFORE the erase (because the iterator will point to the next object in the vector AFTER the erase!)
2092 myVxTracksAtVertex->erase(tracksIter);
2093 tracksBegin=myVxTracksAtVertex->begin();
2094 tracksEnd=myVxTracksAtVertex->end();
2095 found=true;
2096 break;
2097 }
2098
2099 }
2100
2101 if (!found)
2102 {
2103
2104 ATH_MSG_VERBOSE( "Track not found: probably it was already not attached to the vertex" );
2105
2106 }
2107 }
2108 }//iterate on all perigeesToFit
2109
2110 ATH_MSG_VERBOSE( " #CompatibleTracks associated at vertex : " << myVxTracksAtVertex->size() );
2111
2112}
double chi2(TH1 *h0, TH1 *h1)

◆ removeTracksInBadSeed()

float InDet::InDetIterativeSecVtxFinderTool::removeTracksInBadSeed ( xAOD::Vertex * myxAODVertex,
std::vector< const Trk::TrackParameters * > & perigeesToFit ) const
private

Definition at line 1793 of file InDetIterativeSecVtxFinderTool.cxx.

1795{
1796 int removed = 0 ;
1797
1798 std::vector<Trk::VxTrackAtVertex>* tracksAtVertex= &(myxAODVertex->vxTrackAtVertex());
1799
1800 int tot = tracksAtVertex->size() ;
1801 ATH_MSG_DEBUG( " removeTracksInBadSeed ... #Tracks associated at vertex : " << tot );
1802
1803 if ( tot == 0 ) return 0. ;
1804
1805 std::vector<Trk::VxTrackAtVertex>::const_iterator tracksBegin=tracksAtVertex->begin();
1806 std::vector<Trk::VxTrackAtVertex>::const_iterator tracksEnd=tracksAtVertex->end();
1807
1808 std::vector<const Trk::TrackParameters*>::iterator perigeesToFitBegin=perigeesToFit.begin();
1809 std::vector<const Trk::TrackParameters*>::iterator perigeesToFitEnd=perigeesToFit.end();
1810
1811 const Amg::Vector3D position = myxAODVertex->position() ;
1812
1813 tot = 0 ;
1814 float pt_tot = 0. , pt_hf = 0. ;
1815
1816 std::vector<const Trk::TrackParameters*> perigees_deleted ;
1817 for (std::vector<Trk::VxTrackAtVertex>::const_iterator tracksAtVertexIter=tracksBegin;
1818 tracksAtVertexIter!=tracksEnd;++tracksAtVertexIter)
1819 {
1820 if ( (*tracksAtVertexIter).weight() <= m_minWghtAtVtx ) continue ;
1821
1822 const Trk::TrackParameters* measPerigee = (*tracksAtVertexIter).initialPerigee() ;
1823
1824 if ( measPerigee == nullptr ) continue ;
1825
1826 tot ++ ;
1827 bool hf = passHitsFilter( measPerigee, position.perp(), std::abs( position.z() ) ) ;
1828
1829 float pt = 0. ;
1830 const Trk::TrackParameters * svperigee = (*tracksAtVertexIter).perigeeAtVertex() ;
1831 if ( svperigee )
1832 pt = svperigee->parameters()[Trk::theta]/std::abs( svperigee->parameters()[Trk::qOverP] ) ;
1833
1834 pt_tot += pt ;
1835
1836 if ( hf ) continue ;
1837
1838 ATH_MSG_DEBUG( " found and will delete from perigeesToFit !" );
1839
1840 auto found=std::find(perigeesToFitBegin,perigeesToFitEnd,measPerigee);
1841 if(found!=perigeesToFitEnd){
1842 perigees_deleted.push_back(*found);
1843 perigeesToFit.erase(found);
1844 perigeesToFitBegin=perigeesToFit.begin();
1845 perigeesToFitEnd=perigeesToFit.end();
1846 pt_hf += pt ;
1847 removed++;
1848 }
1849 }
1850
1851 if ( tot == 0 || pt_tot == 0. ) return 0. ;
1852
1853//#ifndef MONITORTUNES
1854 for (std::vector<const Trk::TrackParameters*>::iterator perigeesToFitIter=perigees_deleted.begin();
1855 perigeesToFitIter!=perigees_deleted.end();++perigeesToFitIter)
1856 {
1857 for (std::vector<Trk::VxTrackAtVertex>::const_iterator tracksAtVertexIter=tracksBegin;
1858 tracksAtVertexIter!=tracksEnd;++tracksAtVertexIter)
1859 {
1860 const Trk::TrackParameters* measPerigee = (*tracksAtVertexIter).initialPerigee() ;
1861
1862 if (*perigeesToFitIter==measPerigee)
1863 {
1864 tracksAtVertex->erase( tracksAtVertexIter );
1865 tracksBegin=tracksAtVertex->begin();
1866 tracksEnd=tracksAtVertex->end();
1867 ATH_MSG_DEBUG( " deleted tracksAtVertex !" );
1868 break ;
1869 }
1870 }
1871 }
1872//#endif
1873
1874 ATH_MSG_DEBUG( " Bad tracks removed from seed : " << removed );
1875 return 0.4*removed/( 1.0*tot ) + 0.6*pt_hf/pt_tot ;
1876
1877}
bool passHitsFilter(const Trk::TrackParameters *, float vtxR, float absvz) const

◆ setPriVtxPosition()

void InDet::InDetIterativeSecVtxFinderTool::setPriVtxPosition ( double vx,
double vy,
double vz )
override

Definition at line 71 of file InDetIterativeSecVtxFinderTool.cxx.

72{
73 m_privtx = Amg::Vector3D( vx, vy, vz ) ;
74}

◆ SGError()

void InDet::InDetIterativeSecVtxFinderTool::SGError ( const std::string & errService)
private

Definition at line 1691 of file InDetIterativeSecVtxFinderTool.cxx.

1692{
1693 ATH_MSG_FATAL( errService << " not found. Exiting !" );
1694}
#define ATH_MSG_FATAL(x)

◆ V0kine()

bool InDet::InDetIterativeSecVtxFinderTool::V0kine ( const std::vector< Amg::Vector3D > & momenta,
const Amg::Vector3D & position,
float & mass,
float & modir ) const
private

Definition at line 1206 of file InDetIterativeSecVtxFinderTool.cxx.

1208{
1209 mass = -99.9 ;
1210 modir = -99999.9 ;
1211
1212 int ntrk = momenta.size() ;
1213
1214 if ( ntrk < 2 )
1215 {
1216 ATH_MSG_DEBUG( " ntrk < 2 , Meaningless to test mass " );
1217 return false ;
1218 }
1219
1220 std::vector<double> Pv (ntrk);
1221 double vx = 0., vy = 0., vz = 0., eK0 = 0. ;
1223
1224 for ( int t = 0 ; t < ntrk ; t ++ )
1225 {
1226 Amg::Vector3D trk = momenta[t] ;
1227
1228 vz += trk.z() ;
1229 vx += trk.x() ;
1230 vy += trk.y() ;
1231 Pv[t] = trk.x()*trk.x() + trk.y()*trk.y() + trk.z()*trk.z() ;
1232 eK0 += std::sqrt( Pv[t] + pi2 ) ;
1233 }
1234
1235 double mnt2 = vx*vx + vy*vy + vz*vz ;
1236 mass = eK0*eK0 - mnt2 ;
1237 mass = 0.001*( mass >= 0 ? std::sqrt( mass ) : std::sqrt( -mass ) ) ;
1238
1239 Amg::Vector3D vdif = posi - m_privtx ;
1240 Amg::Vector3D vmoment = Amg::Vector3D( vx, vy, vz ) ;
1241
1242 modir = vmoment.dot( vdif )/std::sqrt( mnt2 ) ;
1243
1244 // borrowed from InnerDetector/InDetRecAlgs/InDetV0Finder/InDetV0FinderTool
1245 double a0z = ( vdif + vmoment*vmoment.dot( vdif )/( mnt2 + 0.00001 ) ).z() ;
1246 double Rxy = vdif.perp() ;
1247
1248 ATH_MSG_DEBUG( " V0kine : a0z = " << a0z << " Rxy = " << Rxy <<" direction "<< modir );
1249
1250 if ( ntrk != 2 )
1251 {
1252 ATH_MSG_VERBOSE( " ntrk != 2 , Meaningless to test V0 " );
1253 return false ;
1254 }
1255
1256 if ( a0z > 15. || Rxy > 500. ) { return false ; }
1257
1258 // 1 eV^(-1) of time = hbar / eV = 6.582173*10^(-16) second, for energy-time in natural unit
1259// double planck = 6.582173 ;
1260
1262 double mGam = eGam*eGam - mnt2 ;
1263
1265 double eLam = Pv[0] > Pv[1] ? std::sqrt( Pv[0] + prtn2 ) + std::sqrt( Pv[1] + pi2 ) :
1266 sqrt( Pv[0] + pi2 ) + std::sqrt( Pv[1] + prtn2 ) ;
1267 double mLam = eLam*eLam - mnt2 ;
1268
1269 ATH_MSG_DEBUG( " V0 masses : " << mass
1270 <<" "<< ( mGam >= 0 ? std::sqrt( mGam ) : std::sqrt( -mGam ) )
1271 <<" "<< ( mLam >= 0 ? std::sqrt( mLam ) : std::sqrt( -mLam ) ) );
1272
1273 if ( ( std::abs( mass - ParticleConstants::KZeroMassInMeV ) < 100. ) // K short
1274 || ( mGam > 0 && sqrt( mGam ) < 40. ) // gamma conversion ;
1275 || ( mLam > 0 && std::abs( sqrt( mLam ) - ParticleConstants::lambdaMassInMeV ) < 200. ) // Lambda
1276 ) return true ;
1277
1278 return false ;
1279}
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)

◆ VrtVrtDist()

double InDet::InDetIterativeSecVtxFinderTool::VrtVrtDist ( xAOD::Vertex * v1,
xAOD::Vertex * v2 )
staticprivate

Definition at line 1283 of file InDetIterativeSecVtxFinderTool.cxx.

1284{
1285 double dist = 50 ;
1286
1287 Amg::Vector3D vdif = v1->position() - v2->position() ;
1288 AmgSymMatrix(3) vErrs = v1->covariancePosition() + v2->covariancePosition() ;
1289
1290 vErrs = vErrs.inverse().eval();
1291
1292 dist = vdif.dot( vErrs * vdif ) ;
1293
1294 return dist ;
1295
1296}

Member Data Documentation

◆ m_createSplitVertices

Gaudi::Property<bool> InDet::InDetIterativeSecVtxFinderTool::m_createSplitVertices {this, "createSplitVertices", false, ""}
private

Definition at line 165 of file InDetIterativeSecVtxFinderTool.h.

165{this, "createSplitVertices", false, ""};

◆ m_CutHitsFilter

Gaudi::Property<double> InDet::InDetIterativeSecVtxFinderTool::m_CutHitsFilter {this, "TrackInnerOuterFraction",0.95, ""}
private

Definition at line 162 of file InDetIterativeSecVtxFinderTool.h.

162{this, "TrackInnerOuterFraction",0.95, ""};

◆ m_dir

float InDet::InDetIterativeSecVtxFinderTool::m_dir {}
private

Definition at line 196 of file InDetIterativeSecVtxFinderTool.h.

196{}, m_v0ee{}, m_dir{}, m_ndf{}, m_hif{} ;

◆ m_doMaxTracksCut

Gaudi::Property<bool> InDet::InDetIterativeSecVtxFinderTool::m_doMaxTracksCut {this, "doMaxTracksCut", false, ""}
private

Definition at line 170 of file InDetIterativeSecVtxFinderTool.h.

170{this, "doMaxTracksCut", false, ""};

◆ m_evtNum

long int InDet::InDetIterativeSecVtxFinderTool::m_evtNum {}
private

Definition at line 179 of file InDetIterativeSecVtxFinderTool.h.

179{} ;

◆ m_filterLevel

Gaudi::Property<int> InDet::InDetIterativeSecVtxFinderTool::m_filterLevel {this,"VertexFilterLevel",0,""}
private

Definition at line 157 of file InDetIterativeSecVtxFinderTool.h.

157{this,"VertexFilterLevel",0,""} ;

◆ m_hif

float InDet::InDetIterativeSecVtxFinderTool::m_hif {}
private

Definition at line 196 of file InDetIterativeSecVtxFinderTool.h.

196{}, m_v0ee{}, m_dir{}, m_ndf{}, m_hif{} ;

◆ m_ImpactPoint3dEstimator

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

Definition at line 153 of file InDetIterativeSecVtxFinderTool.h.

153{this, "ImpactPoint3dEstimator", "Trk::ImpactPoint3dEstimator", "impact point estimator"};

◆ m_iterations

int InDet::InDetIterativeSecVtxFinderTool::m_iterations {}
private

Definition at line 180 of file InDetIterativeSecVtxFinderTool.h.

180{} ;

◆ m_iVertexFitter

ToolHandle< Trk::AdaptiveVertexFitter > InDet::InDetIterativeSecVtxFinderTool::m_iVertexFitter {this, "VertexFitterTool", "Trk::AdaptiveVertexFitter","Vertex Fitter"}
private

Definition at line 149 of file InDetIterativeSecVtxFinderTool.h.

149{this, "VertexFitterTool", "Trk::AdaptiveVertexFitter","Vertex Fitter"};

◆ m_leastmodes

std::vector<int>* InDet::InDetIterativeSecVtxFinderTool::m_leastmodes {}
private

Definition at line 182 of file InDetIterativeSecVtxFinderTool.h.

182{} ;

◆ m_LinearizedTrackFactory

ToolHandle< Trk::IVertexLinearizedTrackFactory > InDet::InDetIterativeSecVtxFinderTool::m_LinearizedTrackFactory {this, "LinearizedTrackFactory", "Trk::FullLinearizedTrackFactory", "linearized track factory"}
private

Definition at line 154 of file InDetIterativeSecVtxFinderTool.h.

154{this, "LinearizedTrackFactory", "Trk::FullLinearizedTrackFactory", "linearized track factory"};

◆ m_maximumChi2cutForSeeding

Gaudi::Property<double> InDet::InDetIterativeSecVtxFinderTool::m_maximumChi2cutForSeeding {this, "maxCompatibilityCutSeeding",18., ""}
private

Definition at line 159 of file InDetIterativeSecVtxFinderTool.h.

159{this, "maxCompatibilityCutSeeding",18., ""};

◆ m_maxTracks

Gaudi::Property<unsigned int> InDet::InDetIterativeSecVtxFinderTool::m_maxTracks {this, "MaxTracks", 5000,""}
private

Definition at line 171 of file InDetIterativeSecVtxFinderTool.h.

171{this, "MaxTracks", 5000,""};

◆ m_maxVertices

Gaudi::Property<double> InDet::InDetIterativeSecVtxFinderTool::m_maxVertices {this, "maxVertices",20, ""}
private

Definition at line 161 of file InDetIterativeSecVtxFinderTool.h.

161{this, "maxVertices",20, ""};

◆ m_minVtxDist

Gaudi::Property<float> InDet::InDetIterativeSecVtxFinderTool::m_minVtxDist {this, "SeedsMinimumDistance",0.1,""}
private

Definition at line 164 of file InDetIterativeSecVtxFinderTool.h.

164{this, "SeedsMinimumDistance",0.1,""};

◆ m_minWghtAtVtx

Gaudi::Property<double> InDet::InDetIterativeSecVtxFinderTool::m_minWghtAtVtx {this, "minTrackWeightAtVtx",0.02, ""}
private

Definition at line 160 of file InDetIterativeSecVtxFinderTool.h.

160{this, "minTrackWeightAtVtx",0.02, ""};

◆ m_ndf

float InDet::InDetIterativeSecVtxFinderTool::m_ndf {}
private

Definition at line 196 of file InDetIterativeSecVtxFinderTool.h.

196{}, m_v0ee{}, m_dir{}, m_ndf{}, m_hif{} ;

◆ m_nperiseed

std::vector<int>* InDet::InDetIterativeSecVtxFinderTool::m_nperiseed {}
private

Definition at line 188 of file InDetIterativeSecVtxFinderTool.h.

188{} ;

◆ m_ntracks

int InDet::InDetIterativeSecVtxFinderTool::m_ntracks {}
mutableprivate

Definition at line 197 of file InDetIterativeSecVtxFinderTool.h.

197{} ;

◆ m_OTree

TTree* InDet::InDetIterativeSecVtxFinderTool::m_OTree {}
private

Definition at line 177 of file InDetIterativeSecVtxFinderTool.h.

177{} ;

◆ m_privtx

Amg::Vector3D InDet::InDetIterativeSecVtxFinderTool::m_privtx
private

Definition at line 128 of file InDetIterativeSecVtxFinderTool.h.

◆ m_privtxRef

Gaudi::Property<float> InDet::InDetIterativeSecVtxFinderTool::m_privtxRef {this, "MomentumProjectionOnDirection",-999.9,""}
private

Definition at line 163 of file InDetIterativeSecVtxFinderTool.h.

163{this, "MomentumProjectionOnDirection",-999.9,""};

◆ m_reassignTracksAfterFirstFit

Gaudi::Property<bool> InDet::InDetIterativeSecVtxFinderTool::m_reassignTracksAfterFirstFit {this, "reassignTracksAfterFirstFit", true, ""}
private

Definition at line 168 of file InDetIterativeSecVtxFinderTool.h.

168{this, "reassignTracksAfterFirstFit", true, ""};

◆ m_sdcrsWght

std::vector< std::vector < float > >* InDet::InDetIterativeSecVtxFinderTool::m_sdcrsWght {}
private

Definition at line 186 of file InDetIterativeSecVtxFinderTool.h.

186{};

◆ m_sdFsmwX

std::vector< std::vector < float > >* InDet::InDetIterativeSecVtxFinderTool::m_sdFsmwX {}
private

Definition at line 183 of file InDetIterativeSecVtxFinderTool.h.

183{};

◆ m_sdFsmwY

std::vector< std::vector < float > >* InDet::InDetIterativeSecVtxFinderTool::m_sdFsmwY {}
private

Definition at line 184 of file InDetIterativeSecVtxFinderTool.h.

184{};

◆ m_sdFsmwZ

std::vector< std::vector < float > >* InDet::InDetIterativeSecVtxFinderTool::m_sdFsmwZ {}
private

Definition at line 185 of file InDetIterativeSecVtxFinderTool.h.

185{};

◆ m_seedac

std::vector< int >* InDet::InDetIterativeSecVtxFinderTool::m_seedac {}
private

Definition at line 194 of file InDetIterativeSecVtxFinderTool.h.

194{} ;

◆ m_SeedFinder

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

Definition at line 152 of file InDetIterativeSecVtxFinderTool.h.

152{this, "SeedFinder", "Trk::IndexedCrossDistancesSeedFinder", "seed finder"};

◆ m_seedX

std::vector< float >* InDet::InDetIterativeSecVtxFinderTool::m_seedX {}
private

Definition at line 189 of file InDetIterativeSecVtxFinderTool.h.

189{} ;

◆ m_seedXYdist

std::vector< float >* InDet::InDetIterativeSecVtxFinderTool::m_seedXYdist {}
private

Definition at line 192 of file InDetIterativeSecVtxFinderTool.h.

192{} ;

◆ m_seedY

std::vector< float >* InDet::InDetIterativeSecVtxFinderTool::m_seedY {}
private

Definition at line 190 of file InDetIterativeSecVtxFinderTool.h.

190{} ;

◆ m_seedZ

std::vector< float >* InDet::InDetIterativeSecVtxFinderTool::m_seedZ {}
private

Definition at line 191 of file InDetIterativeSecVtxFinderTool.h.

191{} ;

◆ m_seedZdist

std::vector< float >* InDet::InDetIterativeSecVtxFinderTool::m_seedZdist {}
private

Definition at line 193 of file InDetIterativeSecVtxFinderTool.h.

193{} ;

◆ m_significanceCutSeeding

Gaudi::Property<double> InDet::InDetIterativeSecVtxFinderTool::m_significanceCutSeeding {this, "significanceCutSeeding", 9., ""}
private

Definition at line 158 of file InDetIterativeSecVtxFinderTool.h.

158{this, "significanceCutSeeding", 9., ""};

◆ m_splitVerticesTrkInvFraction

Gaudi::Property<int> InDet::InDetIterativeSecVtxFinderTool::m_splitVerticesTrkInvFraction {this, "splitVerticesTrkInvFraction",2,""}
private

Integer: 1./fraction of tracks to be assigned to the tag split vertex.

Definition at line 166 of file InDetIterativeSecVtxFinderTool.h.

166{this, "splitVerticesTrkInvFraction",2,""};

◆ m_SVtrkFilter

ToolHandle< InDet::IInDetTrackSelectionTool > InDet::InDetIterativeSecVtxFinderTool::m_SVtrkFilter {this, "SecVtxTrackSelector", "InDet::InDetSecVtxTrackSelection", "SV track selector"}
private

Definition at line 151 of file InDetIterativeSecVtxFinderTool.h.

151{this, "SecVtxTrackSelector", "InDet::InDetSecVtxTrackSelection", "SV track selector"};

◆ m_TrkAtVtxMomenta

std::vector< Amg::Vector3D > InDet::InDetIterativeSecVtxFinderTool::m_TrkAtVtxMomenta
mutableprivate

Definition at line 129 of file InDetIterativeSecVtxFinderTool.h.

◆ m_trkdefiPars

std::vector< Amg::VectorX > InDet::InDetIterativeSecVtxFinderTool::m_trkdefiPars
mutableprivate

Definition at line 198 of file InDetIterativeSecVtxFinderTool.h.

◆ m_trkFilter

ToolHandle< InDet::IInDetTrackSelectionTool > InDet::InDetIterativeSecVtxFinderTool::m_trkFilter {this, "BaseTrackSelector", "InDet::InDetTrackSelection", "base track selector"}
private

Definition at line 150 of file InDetIterativeSecVtxFinderTool.h.

150{this, "BaseTrackSelector", "InDet::InDetTrackSelection", "base track selector"};

◆ m_v0ee

float InDet::InDetIterativeSecVtxFinderTool::m_v0ee {}
private

Definition at line 196 of file InDetIterativeSecVtxFinderTool.h.

196{}, m_v0ee{}, m_dir{}, m_ndf{}, m_hif{} ;

◆ m_v0mass

float InDet::InDetIterativeSecVtxFinderTool::m_v0mass {}
mutableprivate

Definition at line 196 of file InDetIterativeSecVtxFinderTool.h.

196{}, m_v0ee{}, m_dir{}, m_ndf{}, m_hif{} ;

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