ATLAS Offline Software
Loading...
Searching...
No Matches
InDetSecVtxTruthMatchTool.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef InDetSecVtxTruthMatchTool_h
6#define InDetSecVtxTruthMatchTool_h
7
8// Framework include(s):
9#include "AsgTools/AsgTool.h"
11
12// EDM include(s):
20
21
22// standard includes
23#include<vector>
24#include <bitset>
25
27
28 //Namespace for useful analysis things on the truth matching decorations applied to the VertexContainer
29
30 //How the matching info is stored; link to the truth vertex, a float with its contribution to the relative track weight, and a float with its contribution to the track sumpt2 of the truth vertex
31 typedef std::tuple<ElementLink<xAOD::TruthVertexContainer>, float, float> VertexTruthMatchInfo;
32
33 //type codes for vertex matching on all vertices
35 Matched=0, // > threshold (default 50%) from one truth interaction
36 Merged, // not matched
37 Split, // highest weight truth interaction contributes to >1 vtx (vtx with highest fraction of sumpT2 remains matched/merged)
38 Fake, // highest contribution is fake (if pile-up MC info not present those tracks end up as "fakes")
40 };
41
42 //type codes for optional SM origin matching
43 //used to determine the content of "other" type vertices above for in-depth studies
44 //should only be used where extended truth info is available -- e.g when LLP1 has been run with saveFullTruth=True
45 //NOTE: types are NOT exclusive -- vertex can be both from B and D decay, signal and hadronic interation, etc
65
66 //type codes for truth vertices
67 //NOTE: types are subsets of subsequent types
69 Reconstructable=0, // fiducial cuts, >= 2 charged daughters
70 Accepted, // >= 2 reco tracks
71 Seeded, // tracks pass cuts
72 Reconstructed, // matched to reco vtx
73 ReconstructedSplit // matched to multiple vtx
74 };
75
76 inline bool isMatched(int matchInfo) {
77 if (matchInfo & (0x1 << Matched)) return true;
78 return false;
79 }
80 inline bool isMerged(int matchInfo) {
81 if (matchInfo & (0x1 << Merged)) return true;
82 return false;
83 }
84 inline bool isSplit(int matchInfo) {
85 if (matchInfo & (0x1 << Split)) return true;
86 return false;
87 }
88 inline bool isFake(int matchInfo) {
89 if (matchInfo & (0x1 << Fake)) return true;
90 return false;
91 }
92 inline bool isOther(int matchInfo) {
93 if (matchInfo & (0x1 << Other)) return true;
94 return false;
95 }
96
97 inline bool isReconstructable(int matchInfo) {
98 if (matchInfo & (0x1 << Reconstructable)) return true;
99 return false;
100 }
101 inline bool isAccepted(int matchInfo) {
102 if (matchInfo & (0x1 << Accepted)) return true;
103 return false;
104 }
105 inline bool isSeeded(int matchInfo) {
106 if (matchInfo & (0x1 << Seeded)) return true;
107 return false;
108 }
109 inline bool isReconstructed(int matchInfo) {
110 if (matchInfo & (0x1 << Reconstructed)) return true;
111 return false;
112 }
113 inline bool isReconstructedSplit(int matchInfo) {
114 if (matchInfo & (0x1 << ReconstructedSplit)) return true;
115 return false;
116 }
117
118 inline bool isOriginType(int matchInfo, VertexMatchOriginType type) {
119 if (type < 0 || type > Signal) { // Signal is the last valid value
120 return false;
121 }
122 return matchInfo & (0x1 << type);
123 }
124}
125
131 public asg::AsgTool {
132
134
135 public:
136
137 InDetSecVtxTruthMatchTool( const std::string & name );
138
139 virtual StatusCode initialize() override final;
140
141 // take collection of vertices, match them, and decorate with matching info
142 virtual StatusCode matchVertices( std::vector<const xAOD::Vertex*> recoVerticesToMatch, std::vector<const xAOD::TruthVertex*> truthVerticesToMatch, const xAOD::TrackParticleContainer* trackParticles ) override;
143
144 private:
145
146 Gaudi::Property<float> m_trkMatchProb{this, "trackMatchProb", 0.5, "Required MC match probability to consider track a good match" };
147 Gaudi::Property<float> m_vxMatchWeight{this, "vertexMatchWeight", 0.5, "Relative weight threshold to consider vertex matched"};
148 Gaudi::Property<float> m_trkPtCut{this, "trackPtCut", 1000., "pt cut to apply on tracks"};
149 Gaudi::Property<std::string> m_selectedTrackFlag{this, "selectedTrackFlag", "is_selected", "Aux decoration on tracks for seeding efficiencies"};
150 Gaudi::Property<bool> m_doMuSA{this, "doMuSA", false, "Combination flag for special MuSA logic" };
151 Gaudi::Property<bool> m_doSMOrigin{this, "doSMOrigin", false, "Enable decoration of SM origin types"};
152
153 ToolHandle<InDet::IInDetTrackTruthOriginTool> m_trackTruthOriginTool{this, "TrackTruthOriginTool", "InDet::InDetTrackTruthOriginTool/TrackTruthOriginTool"};
154 Gaudi::Property<std::string> m_muonContainerName{this, "MuonContainer", "StdWithLRTMuons", "Primary muon container name used in MuSA mode"};
155 Gaudi::Property<std::string> m_muonFallbackContainerName{this, "FallbackMuonContainer", "Muons", "Fallback muon container if primary is unavailable"};
156
157 //private methods to check if particles are good to use
158 //returns barcode of LLP production truth vertex
159 int checkProduction( const xAOD::TruthParticle& truthPart, std::vector<const xAOD::TruthVertex*> truthVerticesToMatch ) const;
161 std::vector<const xAOD::TruthParticle*>& set, int counter) const;
162 std::vector<int> checkParticle( const xAOD::TruthParticle& part, const xAOD::TrackParticleContainer* tkCont, const xAOD::MuonContainer* muonCont ) const;
163 bool isFrom(const xAOD::TruthParticle& truth, int flav) const;
164 int checkSMProduction( const xAOD::TruthParticle& truthPart) const;
166};
167
168#endif
#define ASG_TOOL_CLASS(CLASSNAME, INT1)
xAOD::MuonContainer * muonContainer
Class for vertex truth matching.
ToolHandle< InDet::IInDetTrackTruthOriginTool > m_trackTruthOriginTool
Gaudi::Property< std::string > m_selectedTrackFlag
void countReconstructibleDescendentParticles(const xAOD::TruthVertex &signalTruthVertex, std::vector< const xAOD::TruthParticle * > &set, int counter) const
virtual StatusCode initialize() override final
Dummy implementation of the initialisation function.
virtual StatusCode matchVertices(std::vector< const xAOD::Vertex * > recoVerticesToMatch, std::vector< const xAOD::TruthVertex * > truthVerticesToMatch, const xAOD::TrackParticleContainer *trackParticles) override
Gaudi::Property< bool > m_doSMOrigin
InDetSecVtxTruthMatchTool(const std::string &name)
Gaudi::Property< std::string > m_muonFallbackContainerName
int checkProduction(const xAOD::TruthParticle &truthPart, std::vector< const xAOD::TruthVertex * > truthVerticesToMatch) const
int checkSMProduction(const xAOD::TruthParticle &truthPart) const
const xAOD::Muon * findStandAloneMuon(const xAOD::TrackParticle &mstp, const xAOD::MuonContainer *muonContainer) const
std::vector< int > checkParticle(const xAOD::TruthParticle &part, const xAOD::TrackParticleContainer *tkCont, const xAOD::MuonContainer *muonCont) const
Gaudi::Property< std::string > m_muonContainerName
Gaudi::Property< float > m_trkPtCut
Gaudi::Property< float > m_vxMatchWeight
Gaudi::Property< float > m_trkMatchProb
bool isFrom(const xAOD::TruthParticle &truth, int flav) const
Base class for the dual-use tool implementation classes.
Definition AsgTool.h:47
STL class.
std::tuple< ElementLink< xAOD::TruthVertexContainer >, float, float > VertexTruthMatchInfo
bool isOriginType(int matchInfo, VertexMatchOriginType type)
TruthVertex_v1 TruthVertex
Typedef to implementation.
Definition TruthVertex.h:15
TrackParticle_v1 TrackParticle
Reference the current persistent version:
TruthParticle_v1 TruthParticle
Typedef to implementation.
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".
Muon_v1 Muon
Reference the current persistent version:
MuonContainer_v1 MuonContainer
Definition of the current "Muon container version".