ATLAS Offline Software
Loading...
Searching...
No Matches
VrtSecInclusive.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// VKalVrt.h
6//
7#ifndef VRTSECINCLUSIVE_VRTSECINCLUSIVE_H
8#define VRTSECINCLUSIVE_VRTSECINCLUSIVE_H
9
10
11#include "VrtSecInclusive/Constants.h"
12
15
16// Gaudi includes
17#include "GaudiKernel/ToolHandle.h"
18#include "GaudiKernel/ITHistSvc.h"
19//
21
22// for truth
24
32#include "GaudiKernel/ServiceHandle.h"
36
37// xAOD Classes
48#include "xAODMuon/Muon.h"
50#include "xAODEgamma/Electron.h"
51
52// Normal STL and physical vectors
53#include <vector>
54#include <deque>
55#include <functional>
56#include <optional>
57
58
60
61class TH1;
62
63namespace Trk {
64 class ITruthToTrack;
65 class ITrackToVertex;
67 class IExtrapolator;
68 class IVertexMapper;
69 struct MappedVertex;
70}
71
72namespace VKalVrtAthena {
73
74 namespace GeoModel { enum GeoModel { Run1=1, Run2=2 }; }
75
76 class IntersectionPos;
79
80 class NtupleVars;
81}
82
83
84namespace VKalVrtAthena {
85
87 public:
89 VrtSecInclusive(const std::string& name, ISvcLocator* pSvcLocator);
90
93
94 virtual StatusCode initialize();
95 virtual StatusCode finalize();
96 virtual StatusCode execute();
97 virtual StatusCode initEvent();
98
99 private:
100
102 //
103 // Member Variables
104 //
105
107 // JO: GeoModel
109
110 std::string TrackLocation;
111 std::string MuonLocation;
112 std::string ElectronLocation;
113 std::string PrimVrtLocation;
116 std::string augVerString;
118
121
122 // common feature flags
129 bool extrapPV; //extrapolate reco and prim truth trks to PV (for testing only)
130
132
133 bool doFastMode; // flag for running in rapid finder mode instead of using graph
134
135 // track selection conditions
136 unsigned int SelTrkMaxCutoff;
138
139 /* impact parameters */
147
158
159 /* pT anc chi2 */
161 double TrkPtCut;
162
163 /* hit requirements */
164 bool doTRTPixCut; // Kazuki
170 int CutTRTHits; // Kazuki
173
174 /* track extrpolator; 1==VKalGetImpact, 2==m_trackToVertexTool*/
176
177 // Vertex reconstruction
192 double VertexMergeFinalDistCut; // Kazuki
196
200
201 // When truncateWrkVertices set to true, maxWrkVertices is the maximum
202 // number of potential vertices to process, used to truncate the list
203 // in rare caes where several thousands are found, to avoid algorithm
204 // timeout issues
207
213
218
220
221 // vertexing using muons (test implementation)
227
228 // vertexing using disapperaing track
232
233 // When doSelectTracksWithLRTCuts is set to true, the addtional track cuts
234 // be applied to the selected tracks to reduce the number of fake tracks in
235 // the selected track collected. These cuts are inspired by the improvments that
236 // were implmented for LRT Run 3.
238
239 // Additional dressing option
240 bool doAugmentDVimpactParametersToMuons; // potentially useful for DV + muon search
241 bool doAugmentDVimpactParametersToElectrons; // potentially useful for analyses involving electrons
242
243 // MC truth
246
247 };
248
250
251 // Indicates give-up modes during vertexing
252 // 0 if no errors occured
253 // 1 if too few selected tracks
254 // 2 if too many selected tracks
255 // 3 if wrkVertices container is truncated
256 // 4 if any uncaught exception is raised at the top level of VSI execute()
258
259 // xAOD Accessors
262 std::vector<const xAOD::TrackParticle*> m_selectedTracks;
263 std::vector<const xAOD::TrackParticle*> m_associatedTracks;
264 std::vector<const xAOD::TrackParticle*> m_leptonicTracks;
265 std::vector<double> m_BeamPosition;
266
268 //
269 // Athena JobOption Properties
270 //
271
272 ToolHandle <Trk::ITrkVKalVrtFitter> m_fitSvc; // VKalVrtFitter tool
273 ToolHandle <Trk::ITruthToTrack> m_truthToTrack; // tool to create trkParam from genPart
274
276 ToolHandle< Reco::ITrackToVertex > m_trackToVertexTool;
277 ToolHandle<Trk::ITrackToVertexIPEstimator> m_trackToVertexIPEstimatorTool;
278 ToolHandle<Trk::IExtrapolator> m_extrapolator;
279 ToolHandle<Trk::IVertexMapper> m_vertexMapper;
280
282 ToolHandle<IInDetConditionsTool> m_pixelCondSummaryTool{this, "PixelConditionsSummaryTool", "PixelConditionsSummaryTool", "Tool to retrieve Pixel Conditions summary"};
283 ToolHandle<IInDetConditionsTool> m_sctCondSummaryTool{this, "InDetSCT_ConditionsSummaryTool", "SCT_ConditionsSummaryTool/InDetSCT_ConditionsSummaryTool", "Tool to retrieve SCT conditions summary"};
284
285 const AtlasDetectorID* m_atlasId = nullptr;
286 const PixelID* m_pixelId = nullptr;
287 const SCT_ID* m_sctId = nullptr;
288
290 using PatternStrategyFunc = bool (VrtSecInclusive::*) ( const xAOD::TrackParticle *trk, const Amg::Vector3D& vertex );
291 std::map<std::string, PatternStrategyFunc> m_patternStrategyFuncs;
292
293 // AuxElement decorators
294 std::optional< SG::Decorator< char > > m_decor_isSelected;
295 std::optional< SG::Decorator< char > > m_decor_isAssociated;
296 std::optional< SG::Decorator< char > > m_decor_is_svtrk_final;
297 std::map< unsigned, SG::Decorator<float> > m_trkDecors;
298
300 SG::ReadHandleKey<xAOD::EventInfo> m_eventInfoKey{this,"EventInfoKey", "EventInfo", "EventInfo name"};
302
304 std::vector< IPDecoratorType > m_ipDecors;
305
307 std::optional< VertexELType > m_decor_svLink;
308
310 //
311 // define ntuple variables here
312 //
313
314 // The standard AANT, CollectionTree, is bare bones
316 std::unique_ptr<NtupleVars> m_ntupleVars;
317
318 // Histograms for stats
319 std::map<std::string, TH1*> m_hists;
320
321
323 //
324 // Private member functions
325 //
326
327 // for event info to new ntuple (used to go by default in CollectionTree)
328 void declareProperties();
329
330 StatusCode addEventInfo();
331 StatusCode setupNtupleVariables();
332 StatusCode setupNtuple();
333 StatusCode clearNtupleVariables();
334 StatusCode deleteNtupleVariables();
335 StatusCode processPrimaryVertices();
336 StatusCode fillAANT_SelectedBaseTracks();
338
339 //
340 struct WrkVrt {
341 bool isGood = false;
342 std::deque<long int> selectedTrackIndices;
343 std::deque<long int> associatedTrackIndices;
345 TLorentzVector vertexMom;
346 std::vector<double> vertexCov;
347 double Chi2 = 0;
348 double Chi2_core = 0;
349 std::vector<double> Chi2PerTrk;
350 long int Charge = 0;
351 std::vector< std::vector<double> > TrkAtVrt;
352 unsigned long closestWrkVrtIndex = 0;
354
355 inline double ndof() const { return 2.0*( selectedTrackIndices.size() + associatedTrackIndices.size() ) - 3.0; }
356 inline double ndof_core() const { return 2.0*( selectedTrackIndices.size() ) - 3.0; }
357 inline unsigned nTracksTotal() const { return selectedTrackIndices.size() + associatedTrackIndices.size(); }
358 inline double fitQuality() const { return Chi2 / ndof(); }
359 };
360
361
362 using Detector = int;
363 using Bec = int;
364 using Layer = int;
365 using Flag = int;
366 using ExtrapolatedPoint = std::tuple<const TVector3, Detector, Bec, Layer, Flag>;
367 using ExtrapolatedPattern = std::vector< ExtrapolatedPoint >;
368 using PatternBank = std::map<const xAOD::TrackParticle*, std::pair< std::unique_ptr<ExtrapolatedPattern>, std::unique_ptr<ExtrapolatedPattern> > >;
369
371
372 std::vector< std::pair<int, int> > m_incomp;
373
374 // the map used by printWrkSet
375 std::map<const xAOD::TruthVertex*, bool> m_matchMap;
376
381
383 void selectTrack( const xAOD::TrackParticle* );
384 StatusCode selectTracksInDet();
385 StatusCode selectTracksFromMuons();
386 StatusCode selectTracksFromElectrons();
387 StatusCode selectInDetAndGSFTracks();
388
389 using TrackSelectionAlg = StatusCode (VrtSecInclusive::*)();
390 std::vector<TrackSelectionAlg> m_trackSelectionAlgs;
391
393 using CutFunc = bool (VrtSecInclusive::*) ( const xAOD::TrackParticle* ) const;
394 std::vector<CutFunc> m_trackSelectionFuncs;
395
398 bool selectTrack_pTCut ( const xAOD::TrackParticle* ) const;
399 bool selectTrack_chi2Cut ( const xAOD::TrackParticle* ) const;
400 bool selectTrack_hitPattern ( const xAOD::TrackParticle* ) const;
402 bool selectTrack_d0Cut ( const xAOD::TrackParticle* ) const;
403 bool selectTrack_z0Cut ( const xAOD::TrackParticle* ) const;
404 bool selectTrack_d0errCut ( const xAOD::TrackParticle* ) const;
405 bool selectTrack_z0errCut ( const xAOD::TrackParticle* ) const;
406 static bool selectTrack_d0signifCut ( const xAOD::TrackParticle* ) ;
407 static bool selectTrack_z0signifCut ( const xAOD::TrackParticle* ) ;
408 bool selectTrack_LRTR3Cut ( const xAOD::TrackParticle* ) const;
409
411 StatusCode extractIncompatibleTrackPairs( std::vector<WrkVrt>* );
412 StatusCode findNtrackVertices(std::vector<WrkVrt>* );
413 StatusCode rearrangeTracks( std::vector<WrkVrt>* );
414
416 StatusCode reassembleVertices( std::vector<WrkVrt>* );
417
420 StatusCode mergeByShuffling( std::vector<WrkVrt>* );
421
423 StatusCode mergeFinalVertices( std::vector<WrkVrt>* ); // Kazuki
424
426 StatusCode associateNonSelectedTracks( std::vector<WrkVrt>* );
427
429 StatusCode refitAndSelectGoodQualityVertices( std::vector<WrkVrt>* );
430
432 bool getSVImpactParameters(const xAOD::TrackParticle* trk, const Amg::Vector3D& vertex, std::vector<double>& impactParameters, std::vector<double>& impactParErrors);
433
434 enum TrkParameter { k_d0=0, k_z0=1, k_theta=2, k_phi=3, k_qOverP=4 ,k_nTP=5 };
436
437 using vertexingAlg = StatusCode (VrtSecInclusive::*)( std::vector<WrkVrt>* );
438 std::vector< std::pair<std::string, vertexingAlg> > m_vertexingAlgorithms;
440
441
443 //
444 // Supporting utility functions
445
447 void printWrkSet(const std::vector<WrkVrt> *WrkVrtSet, const std::string& name);
448
450 StatusCode refitVertex( WrkVrt& );
451 StatusCode refitVertex( WrkVrt&, Trk::IVKalState& istate );
452
454 StatusCode refitVertexWithSuggestion( WrkVrt&, const Amg::Vector3D& );
455 StatusCode refitVertexWithSuggestion( WrkVrt&, const Amg::Vector3D&, Trk::IVKalState& istate );
456
459 double improveVertexChi2( WrkVrt& );
460
461 static void removeTrackFromVertex(std::vector<WrkVrt>*,
462 std::vector< std::deque<long int> > *,
463 const long int & ,const long int & );
464
465 StatusCode disassembleVertex(std::vector<WrkVrt> *, const unsigned& vertexIndex );
466
467 void trackClassification(std::vector< WrkVrt >* , std::map< long int, std::vector<long int> >& );
468
469 double findWorstChi2ofMaximallySharedTrack(std::vector<WrkVrt>*, std::map< long int, std::vector<long int> >&, long int & ,long int & );
470
472 static size_t nTrkCommon( std::vector<WrkVrt> *WrkVrtSet, const std::pair<unsigned, unsigned>& pairIndex ) ;
473
475 double significanceBetweenVertices( const WrkVrt&, const WrkVrt& ) const;
476
478 double distanceBetweenVertices( const WrkVrt&, const WrkVrt& ) const;
479
480 using AlgForVerticesPair = double (VrtSecInclusive::*)( const WrkVrt&, const WrkVrt& ) const;
481
483 double findMinVerticesPair( std::vector<WrkVrt>*, std::pair<unsigned, unsigned>&, const AlgForVerticesPair& );
484
486 static double findMinVerticesNextPair( std::vector<WrkVrt>*, std::pair<unsigned, unsigned>& );
487
489 StatusCode mergeVertices( WrkVrt& destination, WrkVrt& source );
490
492
514
516 static void fillTrackSummary( track_summary& summary, const xAOD::TrackParticle *trk );
517
519
520 bool patternCheck ( const uint32_t& pattern, const Amg::Vector3D& vertex );
521 static bool patternCheckRun1( const uint32_t& pattern, const Amg::Vector3D& vertex );
522 static bool patternCheckRun2( const uint32_t& pattern, const Amg::Vector3D& vertex );
523
524 bool patternCheckOuterOnly ( const uint32_t& pattern, const Amg::Vector3D& vertex );
525 static bool patternCheckRun1OuterOnly( const uint32_t& pattern, const Amg::Vector3D& vertex );
526 static bool patternCheckRun2OuterOnly( const uint32_t& pattern, const Amg::Vector3D& vertex );
527
529 bool checkTrackHitPatternToVertex( const xAOD::TrackParticle *trk, const Amg::Vector3D& vertex );
530
533
536
539
541 bool passedFakeReject( const Amg::Vector3D& FitVertex, const xAOD::TrackParticle *itrk, const xAOD::TrackParticle *jtrk );
542
544 void removeInconsistentTracks( WrkVrt& );
545
546 template<class Track> void getIntersection(Track *trk, std::vector<IntersectionPos*>& layers, const Trk::Perigee* per);
547 template<class Track> void setIntersection(Track *trk, IntersectionPos *bec, const Trk::Perigee* per);
548
550 StatusCode monitorVertexingAlgorithmStep( std::vector<WrkVrt>*, const std::string& name, bool final = false );
551
553 //
554 // Truth Information Algorithms Member Functions
555 //
556 //
557
559
560 StatusCode categorizeVertexTruthTopology( xAOD::Vertex *vertex );
561
563
564 std::vector<const xAOD::TruthVertex*> m_tracingTruthVertices;
565
567 //
568 // Additional augmentation
569 //
570 //
571
572 template<class LeptonFlavor>
573 StatusCode augmentDVimpactParametersToLeptons( const std::string& containerName );
574
576 void lockTrackDecorations( const xAOD::TrackParticle* trk, bool onlySelection ) const;
577 void lockLeptonDecorations( const SG::AuxVectorData* cont ) const;
578 StatusCode lockTrackDecorations( bool onlySelection ) const;
579
580 };
581
582} // end of namespace bracket
583
584
585// This header file contains the definition of member templates
586#include "details/Utilities.h"
587
588
589#endif /* VRTSECINCLUSIVE_VRTSECINCLUSIVE_H */
Helper class to provide type-safe access to aux data.
Extrapolation for HepMC particles.
This is an Identifier helper class for the Pixel subdetector.
This is an Identifier helper class for the SCT subdetector.
AthAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters:
This class provides an interface to generate or decode an identifier for the upper levels of the dete...
This is an Identifier helper class for the Pixel subdetector.
Definition PixelID.h:67
This is an Identifier helper class for the SCT subdetector.
Definition SCT_ID.h:68
SG::Decorator< T, ALLOC > Decorator
Definition AuxElement.h:575
Manage lookup of vectors of auxiliary data.
Property holding a SG store/key/clid from which a ReadHandle is made.
Property holding a SG store/key/clid/attr name from which a WriteDecorHandle is made.
Interface class for the extrapolation AlgTool, it inherits from IAlgTool Detailed information about p...
ITruthToTrack is an interface to create Trk::TrackParameters from a HepMC::GenParticle.
Interface class IVertexMapper.
bool patternCheck(const uint32_t &pattern, const Amg::Vector3D &vertex)
SG::AuxElement::Decorator< std::vector< std::vector< float > > > IPDecoratorType
StatusCode associateNonSelectedTracks(std::vector< WrkVrt > *)
in addition to selected tracks, associate as much tracks as possible
VrtSecInclusive(const std::string &name, ISvcLocator *pSvcLocator)
Standard Athena-Algorithm Constructor.
bool selectTrack_pTCut(const xAOD::TrackParticle *) const
std::map< std::string, PatternStrategyFunc > m_patternStrategyFuncs
static const xAOD::TruthParticle * getTrkGenParticle(const xAOD::TrackParticle *)
Definition TruthAlgs.cxx:26
StatusCode refitVertexWithSuggestion(WrkVrt &, const Amg::Vector3D &)
refit the vertex with suggestion
bool passedFakeReject(const Amg::Vector3D &FitVertex, const xAOD::TrackParticle *itrk, const xAOD::TrackParticle *jtrk)
Flag false if the consistituent tracks are not consistent with the vertex position.
bool checkTrackHitPatternToVertex(const xAOD::TrackParticle *trk, const Amg::Vector3D &vertex)
A classical method with hard-coded geometry.
static bool patternCheckRun1(const uint32_t &pattern, const Amg::Vector3D &vertex)
bool checkTrackHitPatternToVertexByExtrapolationAssist(const xAOD::TrackParticle *trk, const Amg::Vector3D &vertex)
New method with track extrapolation.
void printWrkSet(const std::vector< WrkVrt > *WrkVrtSet, const std::string &name)
print the contents of reconstructed vertices
static bool patternCheckRun2(const uint32_t &pattern, const Amg::Vector3D &vertex)
bool checkTrackHitPatternToVertexByExtrapolation(const xAOD::TrackParticle *trk, const Amg::Vector3D &vertex)
New method with track extrapolation.
ToolHandle< IInDetConditionsTool > m_pixelCondSummaryTool
Condition service.
StatusCode mergeFinalVertices(std::vector< WrkVrt > *)
attempt to merge vertices by lookng at the distance between two vertices
StatusCode refitAndSelectGoodQualityVertices(std::vector< WrkVrt > *)
finalization of the vertex and store to xAOD::VertexContainer
bool selectTrack_chi2Cut(const xAOD::TrackParticle *) const
bool selectTrack_notPVassociated(const xAOD::TrackParticle *) const
track-by-track selection strategies
void lockTrackDecorations(const xAOD::TrackParticle *trk, bool onlySelection) const
lock decorations at the end of the algorithm
void trackClassification(std::vector< WrkVrt > *, std::map< long int, std::vector< long int > > &)
std::map< unsigned, SG::Decorator< float > > m_trkDecors
ExtrapolatedPattern * extrapolatedPattern(const xAOD::TrackParticle *, enum Trk::PropDirection)
void selectTrack(const xAOD::TrackParticle *)
Vertexing Algorithm Member Functions.
StatusCode(VrtSecInclusive::*)(std::vector< WrkVrt > *) vertexingAlg
void lockLeptonDecorations(const SG::AuxVectorData *cont) const
bool selectTrack_hitPatternTight(const xAOD::TrackParticle *) const
double findWorstChi2ofMaximallySharedTrack(std::vector< WrkVrt > *, std::map< long int, std::vector< long int > > &, long int &, long int &)
StatusCode reassembleVertices(std::vector< WrkVrt > *)
attempt to merge vertices when all tracks of a vertex A is close to vertex B in terms of impact param...
static bool selectTrack_d0signifCut(const xAOD::TrackParticle *)
bool selectTrack_d0errCut(const xAOD::TrackParticle *) const
std::unique_ptr< NtupleVars > m_ntupleVars
static bool patternCheckRun1OuterOnly(const uint32_t &pattern, const Amg::Vector3D &vertex)
bool getSVImpactParameters(const xAOD::TrackParticle *trk, const Amg::Vector3D &vertex, std::vector< double > &impactParameters, std::vector< double > &impactParErrors)
get secondary vertex impact parameters
StatusCode categorizeVertexTruthTopology(xAOD::Vertex *vertex)
Definition TruthAlgs.cxx:44
void removeInconsistentTracks(WrkVrt &)
Remove inconsistent tracks from vertices.
StatusCode rearrangeTracks(std::vector< WrkVrt > *)
ToolHandle< IInDetConditionsTool > m_sctCondSummaryTool
bool selectTrack_hitPattern(const xAOD::TrackParticle *) const
ToolHandle< Trk::IVertexMapper > m_vertexMapper
StatusCode extractIncompatibleTrackPairs(std::vector< WrkVrt > *)
related to the graph method and verte finding
bool checkTrackHitPatternToVertexOuterOnly(const xAOD::TrackParticle *trk, const Amg::Vector3D &vertex)
A classical method with hard-coded geometry.
std::vector< CutFunc > m_trackSelectionFuncs
std::optional< SG::Decorator< char > > m_decor_isSelected
void getIntersection(Track *trk, std::vector< IntersectionPos * > &layers, const Trk::Perigee *per)
ToolHandle< Trk::IExtrapolator > m_extrapolator
ToolHandle< Trk::ITrackToVertexIPEstimator > m_trackToVertexIPEstimatorTool
std::map< std::string, TH1 * > m_hists
static void fillTrackSummary(track_summary &summary, const xAOD::TrackParticle *trk)
retrieve the track hit information
double improveVertexChi2(WrkVrt &)
attempt to improve the vertex chi2 by removing the most-outlier track one by one until the vertex chi...
SG::ReadHandleKey< xAOD::EventInfo > m_eventInfoKey
Read/Write Handle Keys.
double distanceBetweenVertices(const WrkVrt &, const WrkVrt &) const
calculate the physical distance
StatusCode(VrtSecInclusive::*)() TrackSelectionAlg
StatusCode findNtrackVertices(std::vector< WrkVrt > *)
StatusCode monitorVertexingAlgorithmStep(std::vector< WrkVrt > *, const std::string &name, bool final=false)
monitor the intermediate status of vertexing
static size_t nTrkCommon(std::vector< WrkVrt > *WrkVrtSet, const std::pair< unsigned, unsigned > &pairIndex)
returns the number of tracks commonly present in both vertices
static double findMinVerticesNextPair(std::vector< WrkVrt > *, std::pair< unsigned, unsigned > &)
returns the next pair of vertices that give next-to-minimum distance significance
std::tuple< const TVector3, Detector, Bec, Layer, Flag > ExtrapolatedPoint
bool(VrtSecInclusive::*)(const xAOD::TrackParticle *) const CutFunc
track selection
bool selectTrack_LRTR3Cut(const xAOD::TrackParticle *) const
std::optional< SG::Decorator< char > > m_decor_is_svtrk_final
double findMinVerticesPair(std::vector< WrkVrt > *, std::pair< unsigned, unsigned > &, const AlgForVerticesPair &)
returns the pair of vertices that give minimum in terms of some observable (e.g.
std::optional< VertexELType > m_decor_svLink
std::map< const xAOD::TrackParticle *, std::pair< std::unique_ptr< ExtrapolatedPattern >, std::unique_ptr< ExtrapolatedPattern > > > PatternBank
const AtlasDetectorID * m_atlasId
static void removeTrackFromVertex(std::vector< WrkVrt > *, std::vector< std::deque< long int > > *, const long int &, const long int &)
ToolHandle< Reco::ITrackToVertex > m_trackToVertexTool
get a handle on the Track to Vertex tool
bool selectTrack_d0Cut(const xAOD::TrackParticle *) const
bool selectTrack_z0Cut(const xAOD::TrackParticle *) const
std::vector< double > m_BeamPosition
SG::AuxElement::Decorator< std::vector< ElementLink< xAOD::VertexContainer > > > VertexELType
std::vector< IPDecoratorType > m_ipDecors
std::vector< ExtrapolatedPoint > ExtrapolatedPattern
const xAOD::VertexContainer * m_primaryVertices
StatusCode mergeVertices(WrkVrt &destination, WrkVrt &source)
the 2nd vertex is merged into the 1st vertex.
std::vector< const xAOD::TrackParticle * > m_associatedTracks
bool selectTrack_z0errCut(const xAOD::TrackParticle *) const
ToolHandle< Trk::ITrkVKalVrtFitter > m_fitSvc
std::vector< const xAOD::TrackParticle * > m_selectedTracks
StatusCode disassembleVertex(std::vector< WrkVrt > *, const unsigned &vertexIndex)
StatusCode fillAANT_SecondaryVertices(xAOD::VertexContainer *)
double(VrtSecInclusive::*)(const WrkVrt &, const WrkVrt &) const AlgForVerticesPair
StatusCode mergeByShuffling(std::vector< WrkVrt > *)
attempt to merge splitted vertices when they are significantly distant due to the long-tail behavior ...
bool patternCheckOuterOnly(const uint32_t &pattern, const Amg::Vector3D &vertex)
static bool selectTrack_z0signifCut(const xAOD::TrackParticle *)
std::vector< std::pair< std::string, vertexingAlg > > m_vertexingAlgorithms
~VrtSecInclusive()
Default Destructor.
struct VKalVrtAthena::VrtSecInclusive::track_summary_properties track_summary
SG::WriteDecorHandleKey< xAOD::EventInfo > m_vertexingStatusKey
void setIntersection(Track *trk, IntersectionPos *bec, const Trk::Perigee *per)
std::vector< const xAOD::TrackParticle * > m_leptonicTracks
std::vector< std::pair< int, int > > m_incomp
ToolHandle< Trk::ITruthToTrack > m_truthToTrack
static bool patternCheckRun2OuterOnly(const uint32_t &pattern, const Amg::Vector3D &vertex)
std::optional< SG::Decorator< char > > m_decor_isAssociated
std::vector< const xAOD::TruthVertex * > m_tracingTruthVertices
bool(VrtSecInclusive::*)(const xAOD::TrackParticle *trk, const Amg::Vector3D &vertex) PatternStrategyFunc
double significanceBetweenVertices(const WrkVrt &, const WrkVrt &) const
calculate the significance (Mahalanobis distance) between two reconstructed vertices
std::vector< TrackSelectionAlg > m_trackSelectionAlgs
std::map< const xAOD::TruthVertex *, bool > m_matchMap
Eigen::Matrix< double, 3, 1 > Vector3D
USAGE: openCoraCool.exe "COOLONL_SCT/COMP200".
Ensure that the ATLAS eigen extensions are properly loaded.
PropDirection
PropDirection, enum for direction of the propagation.
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
TrackParticle_v1 TrackParticle
Reference the current persistent version:
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.
TruthParticle_v1 TruthParticle
Typedef to implementation.
double Chi2
VKalVrt fit covariance.
std::vector< double > vertexCov
VKalVrt fit vertex 4-momentum.
std::vector< double > Chi2PerTrk
VKalVrt fit chi2 result.
std::deque< long int > selectedTrackIndices
flagged true for good vertex candidates
std::deque< long int > associatedTrackIndices
list if indices in TrackParticleContainer for selectedBaseTracks
long int Charge
list of VKalVrt fit chi2 for each track
double Chi2_core
VKalVrt fit chi2 result.
Amg::Vector3D vertex
list if indices in TrackParticleContainer for associatedTracks
TLorentzVector vertexMom
VKalVrt fit vertex position.
double closestWrkVrtValue
stores the index of the closest WrkVrt in std::vector<WrkVrt>
double ndof() const
stores the value of some observable to the closest WrkVrt ( observable = e.g. significance )
std::vector< std::vector< double > > TrkAtVrt
total charge of the vertex
unsigned long closestWrkVrtIndex
list of track parameters wrt the reconstructed vertex