ATLAS Offline Software
Loading...
Searching...
No Matches
TrigTauPrecisionDiKaonHypoTool Class Reference

Precision step hypothesis tool for applying meson kinematic cuts (meson chains) More...

#include <TrigTauPrecisionDiKaonHypoTool.h>

Inheritance diagram for TrigTauPrecisionDiKaonHypoTool:
Collaboration diagram for TrigTauPrecisionDiKaonHypoTool:

Public Member Functions

 TrigTauPrecisionDiKaonHypoTool (const std::string &type, const std::string &name, const IInterface *parent)
virtual StatusCode initialize () override
virtual StatusCode decide (std::vector< ITrigTauJetHypoTool::ToolInfo > &input) const override
virtual bool decide (const ITrigTauJetHypoTool::ToolInfo &i) const override

Private Attributes

HLT::Identifier m_decisionId
Gaudi::Property< float > m_ptMin {this, "PtMin", 0, "Minimum tau pT value"}
Gaudi::Property< int > m_numTrackMax {this, "NTracksMax", 2, "Maximum number of Tracks"}
Gaudi::Property< int > m_numTrackMin {this, "NTracksMin", 1, "Minimum number of Tracks"}
Gaudi::Property< int > m_numIsoTrackMax {this, "NIsoTracksMax", 1, "Maximum number of wide Tracks"}
Gaudi::Property< float > m_massTrkSysMin {this, "massTrkSysMin", 0, "Minimum DiTrack mass value"}
Gaudi::Property< float > m_massTrkSysMax {this, "massTrkSysMax", 1000000000, "Maximum DiTrack mass value"}
Gaudi::Property< float > m_massTrkSysKaonMin {this, "massTrkSysKaonMin", 0, "Minimum DiKaon mass value"}
Gaudi::Property< float > m_massTrkSysKaonMax {this, "massTrkSysKaonMax", 1000000000, "Maximum DiKaon mass value"}
Gaudi::Property< float > m_massTrkSysKaonPiMin {this, "massTrkSysKaonPiMin", 0, "Minimum KaonPi mass value"}
Gaudi::Property< float > m_massTrkSysKaonPiMax {this, "massTrkSysKaonPiMax", 1000000000, "Maximum KaonPi mass value"}
Gaudi::Property< float > m_targetMassTrkSysKaonPi {this, "targetMassTrkSysKaonPi", 0, "Target KaonPi mass value (parameter)"}
Gaudi::Property< float > m_leadTrkPtMin {this, "leadTrkPtMin", 0, "Minimum Pt of Lead Track"}
Gaudi::Property< float > m_dRmaxMax {this, "dRmaxMax", 10, "Maximum dRmax value"}
Gaudi::Property< float > m_etOverPtLeadTrkMax {this, "etOverPtLeadTrkMax", 10, "Maximum et/pt(lead track)"}
Gaudi::Property< float > m_etOverPtLeadTrkMin {this, "etOverPtLeadTrkMin", 0, "Minimum et/pt(lead track)"}
Gaudi::Property< float > m_EMPOverTrkSysPMax {this, "EMPOverTrkSysPMax", 5, "Maximum Cluster pt over ditrack pt"}
Gaudi::Property< bool > m_acceptAll {this, "AcceptAll", false, "Ignore selection"}
ToolHandle< GenericMonitoringToolm_monTool {this, "MonTool", "", "Monitoring tool"}

Detailed Description

Precision step hypothesis tool for applying meson kinematic cuts (meson chains)

Definition at line 29 of file TrigTauPrecisionDiKaonHypoTool.h.

Constructor & Destructor Documentation

◆ TrigTauPrecisionDiKaonHypoTool()

TrigTauPrecisionDiKaonHypoTool::TrigTauPrecisionDiKaonHypoTool ( const std::string & type,
const std::string & name,
const IInterface * parent )

Definition at line 27 of file TrigTauPrecisionDiKaonHypoTool.cxx.

30 : base_class(type, name, parent), m_decisionId(HLT::Identifier::fromToolName(name))
31{
32
33}
static HLT::Identifier fromToolName(const std::string &tname)

Member Function Documentation

◆ decide() [1/2]

bool TrigTauPrecisionDiKaonHypoTool::decide ( const ITrigTauJetHypoTool::ToolInfo & i) const
overridevirtual

Definition at line 70 of file TrigTauPrecisionDiKaonHypoTool.cxx.

71{
72 ATH_MSG_DEBUG(name() << ": in execute()");
73
74 auto NInputTaus = Monitored::Scalar<int>("NInputTaus", -1);
75 auto passedCuts = Monitored::Scalar<int>("CutCounter", 0);
76 auto PtAccepted = Monitored::Scalar<float>("PtAccepted", -1);
77 auto NTracksAccepted = Monitored::Scalar<int>("NTracksAccepted", -1);
78 auto NIsoTracksAccepted = Monitored::Scalar<int>("NIsoTracksAccepted", -1);
79 auto massTrkSysAccepted = Monitored::Scalar<float>("massTrkSysAccepted", -10);
80 auto massTrkSysKaonAccepted = Monitored::Scalar<float>("massTrkSysKaonAccepted", -10);
81 auto massTrkSysKaonPiAccepted = Monitored::Scalar<float>("massTrkSysKaonPiAccepted", -10);
82 auto leadTrkPtAccepted = Monitored::Scalar<float>("leadTrkPtAccepted", -10);
83 auto dRAccepted = Monitored::Scalar<float>("dRAccepted", -1);
84 auto etOverPtLeadTrkAccepted = Monitored::Scalar<float>("etOverPtLeadTrkAccepted", -10);
85 auto EMOverTrkSysPAccepted = Monitored::Scalar<float>("EMOverTrkSysPAccepted", -10);
86
87 auto monitorIt = Monitored::Group(m_monTool,
88 NInputTaus, passedCuts,
89 PtAccepted, NTracksAccepted, NIsoTracksAccepted,
90 massTrkSysAccepted, massTrkSysKaonAccepted, massTrkSysKaonPiAccepted,
91 leadTrkPtAccepted, dRAccepted, etOverPtLeadTrkAccepted, EMOverTrkSysPAccepted);
92
93
94 // Tau pass flag
95 bool pass = false;
96
97 if(m_acceptAll) {
98 pass = true;
99 ATH_MSG_DEBUG("AcceptAll property is set: taking all events");
100 }
101
102 // Debugging location of the TauJet RoI
103 ATH_MSG_DEBUG("Input RoI eta: " << input.roi->eta() << ", phi: " << input.roi->phi() << ", z: " << input.roi->zed());
104
105 const xAOD::TauJetContainer* TauContainer = input.tauContainer;
106 NInputTaus = TauContainer->size();
107 // There should only be a single TauJet in the TauJetContainer; just in case we still run the loop
108 for(const xAOD::TauJet* Tau : *TauContainer) {
109 ATH_MSG_DEBUG(" New HLT TauJet candidate:");
110 passedCuts++;
111
112
113 //---------------------------------------------------------------
114 // Calibrated tau pT cut ('idperf' step)
115 //---------------------------------------------------------------
116 float pT = Tau->pt();
117 ATH_MSG_DEBUG(" pT: " << pT / Gaudi::Units::GeV);
118
119 if(!(pT > m_ptMin)) continue;
120 passedCuts++;
121 PtAccepted = pT / Gaudi::Units::GeV;
122
123
124 //---------------------------------------------------------------
125 // Track counting ('perf' step)
126 //---------------------------------------------------------------
127 int numTrack = Tau->nTracks();
128 int numIsoTrack = Tau->nTracksIsolation();
129
130 ATH_MSG_DEBUG(" N Tracks: " << numTrack);
131 ATH_MSG_DEBUG(" N Iso Tracks: " << numIsoTrack);
132
133 // Apply NTrackMin/Max and NIsoTrackMax cuts:
134 if(!(numTrack >= m_numTrackMin && numTrack <= m_numTrackMax)) continue;
135 if(!(numIsoTrack <= m_numIsoTrackMax)) continue;
136 passedCuts++;
137 NTracksAccepted = numTrack;
138 NIsoTracksAccepted = numIsoTrack;
139
140
141 //---------------------------------------------------------------
142 // Leading track pT cut
143 //---------------------------------------------------------------
144 float leadTrkPt = -1;
145 Tau->detail(xAOD::TauJetParameters::leadTrkPt, leadTrkPt);
146 ATH_MSG_DEBUG(" leadTrkPt: " << leadTrkPt / Gaudi::Units::GeV);
147 if(!(leadTrkPt > m_leadTrkPtMin)) continue;
148 passedCuts++;
149 leadTrkPtAccepted = leadTrkPt / Gaudi::Units::GeV;
150
151
152 //---------------------------------------------------------------
153 // Cuts on the track-system mass, for different meson hypothesis
154 //---------------------------------------------------------------
155 float massTrkSys = -1.;
156 Tau->detail(xAOD::TauJetParameters::massTrkSys, massTrkSys);
157 ATH_MSG_DEBUG(" massTrkSys: " << massTrkSys / Gaudi::Units::GeV);
158
159 // For the dikaon mass hypothesis, compute the invariant mass with the kaon mass
160 // Also cache the p4 of tracks, to use for other mass hypothesis
161 TLorentzVector my_kaons;
162 std::vector<TLorentzVector> my_trks;
163 for(unsigned int i = 0; i < Tau->nTracks(); ++i) {
164 const xAOD::TrackParticle* trk = nullptr;
165 TLorentzVector tmpKaon;;
166
167 try {
168 trk = Tau->track(i)->track();
169 } catch(const std::exception& e) {
170 ATH_MSG_WARNING(" Failed to get tau track link!");
171 }
172
173 if(trk) {
174 tmpKaon.SetPtEtaPhiM(trk->pt(), trk->eta(), trk->phi(), ParticleConstants::chargedKaonMassInMeV);
175 my_trks.push_back(trk->p4());
176 }
177
178 my_kaons = my_kaons + tmpKaon;
179 }
180
181 float massTrkSysKaon = my_kaons.M();
182 ATH_MSG_DEBUG(" massTrkSys with kaon mass hypo: " << massTrkSysKaon / Gaudi::Units::GeV);
183
184 // kaon+pi mass hypo
185 float finalKPiMass = 0;
186 if(my_trks.size() == 2) {
187 TLorentzVector tmpKaon;
188
189 tmpKaon.SetPtEtaPhiM(my_trks.at(0).Pt(), my_trks.at(0).Eta(), my_trks.at(0).Phi(), ParticleConstants::chargedKaonMassInMeV);
190 TLorentzVector tmpPion = my_trks.at(1);
191 float kPiMass1 = (tmpKaon+tmpPion).M();
192
193 tmpKaon.SetPtEtaPhiM(my_trks.at(1).Pt(), my_trks.at(1).Eta(), my_trks.at(1).Phi(), ParticleConstants::chargedKaonMassInMeV);
194 tmpPion = my_trks.at(0);
195 float kPiMass2 = (tmpKaon+tmpPion).M();
196
197 if(std::abs(kPiMass1 - m_targetMassTrkSysKaonPi) < std::abs(kPiMass2 - m_targetMassTrkSysKaonPi)) {
198 finalKPiMass = kPiMass1;
199 } else {
200 finalKPiMass = kPiMass2;
201 }
202 }
203 float massTrkSysKaonPi = finalKPiMass;
204 ATH_MSG_DEBUG(" massTrkSys with kaon+pi mass hypo: " << massTrkSysKaonPi / Gaudi::Units::GeV);
205
206 if(!(massTrkSys > m_massTrkSysMin && massTrkSys < m_massTrkSysMax)) continue;
207 passedCuts++;
208 massTrkSysAccepted = massTrkSys / Gaudi::Units::GeV;
209
210 if(!(massTrkSysKaon > m_massTrkSysKaonMin && massTrkSysKaon < m_massTrkSysKaonMax)) continue;
211 passedCuts++;
212 massTrkSysKaonAccepted = massTrkSysKaon / Gaudi::Units::GeV;
213
214 // Use '>=' here, otherwise the singlepion chain would fail!
215 if(!(massTrkSysKaonPi >= m_massTrkSysKaonPiMin && massTrkSysKaonPi < m_massTrkSysKaonPiMax)) continue;
216 passedCuts++;
217 massTrkSysKaonPiAccepted = massTrkSysKaonPi / Gaudi::Units::GeV;
218
219
220 //---------------------------------------------------------------
221 // Cut on E_{T}^{EM} / p_T^{trk sys}
222 //---------------------------------------------------------------
223 float EMPOverTrkSysP = -1;
224 Tau->detail(xAOD::TauJetParameters::EMPOverTrkSysP, EMPOverTrkSysP);
225 ATH_MSG_DEBUG(" EMPOverTrkSysP: " << EMPOverTrkSysP);
226 if(!(EMPOverTrkSysP < m_EMPOverTrkSysPMax)) continue;
227 passedCuts++;
228 EMOverTrkSysPAccepted = EMPOverTrkSysP;
229
230
231 //---------------------------------------------------------------
232 // Cut on E_{T} / p_T^{lead trk}
233 //---------------------------------------------------------------
234 float etOverPtLeadTrk = -1;
235 Tau->detail(xAOD::TauJetParameters::etOverPtLeadTrk, etOverPtLeadTrk);
236 ATH_MSG_DEBUG(" etOverPtLeadTrk: " << etOverPtLeadTrk);
237 if(!(etOverPtLeadTrk > m_etOverPtLeadTrkMin && etOverPtLeadTrk < m_etOverPtLeadTrkMax)) continue;
238 passedCuts++;
239 etOverPtLeadTrkAccepted = etOverPtLeadTrk;
240
241
242 //---------------------------------------------------------------
243 // Cut on Max dR(track, tau)
244 //---------------------------------------------------------------
245 float dRmax = -1;
246 Tau->detail(xAOD::TauJetParameters::dRmax, dRmax);
247 ATH_MSG_DEBUG(" dRmax: " << dRmax);
248 if(!(dRmax < m_dRmaxMax)) continue;
249 passedCuts++;
250 dRAccepted = dRmax;
251
252
253 //---------------------------------------------------------
254 // At least one Tau passed all the cuts. Accept the event!
255 //---------------------------------------------------------
256 pass = true;
257
258 ATH_MSG_DEBUG(" Pass hypo tool: " << pass);
259 }
260
261 return pass;
262}
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
size_type size() const noexcept
Returns the number of elements in the collection.
ToolHandle< GenericMonitoringTool > m_monTool
virtual FourMom_t p4() const override final
The full 4-momentum of the particle.
virtual double phi() const override final
The azimuthal angle ( ) of the particle (has range to .)
virtual double pt() const override final
The transverse momentum ( ) of the particle.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
constexpr double chargedKaonMassInMeV
the mass of the charged kaon (in MeV)
@ Tau
The object is a tau (jet)
Definition ObjectType.h:49
@ dRmax
Get maximal dR of tracks associated to calo-seeded tau.
Definition TauDefs.h:226
TrackParticle_v1 TrackParticle
Reference the current persistent version:
TauJet_v3 TauJet
Definition of the current "tau version".
TauJetContainer_v3 TauJetContainer
Definition of the current "taujet container version".

◆ decide() [2/2]

StatusCode TrigTauPrecisionDiKaonHypoTool::decide ( std::vector< ITrigTauJetHypoTool::ToolInfo > & input) const
overridevirtual

Definition at line 264 of file TrigTauPrecisionDiKaonHypoTool.cxx.

265{
266 for(ITrigTauJetHypoTool::ToolInfo& i : input) {
267 if(passed(m_decisionId.numeric(), i.previousDecisionIDs)) {
268 if(decide(i)) {
269 addDecisionID(m_decisionId, i.decision);
270 }
271 }
272 }
273
274 return StatusCode::SUCCESS;
275}
virtual StatusCode decide(std::vector< ITrigTauJetHypoTool::ToolInfo > &input) const override
bool passed(DecisionID id, const DecisionIDContainer &idSet)
checks if required decision ID is in the set of IDs in the container
void addDecisionID(DecisionID id, Decision *d)
Appends the decision (given as ID) to the decision object.

◆ initialize()

StatusCode TrigTauPrecisionDiKaonHypoTool::initialize ( )
overridevirtual

Definition at line 36 of file TrigTauPrecisionDiKaonHypoTool.cxx.

37{
38 ATH_MSG_DEBUG(name() << ": in initialize()");
39
40 ATH_MSG_DEBUG("TrigTauPrecisionDiKaonHypoTool will cut on:");
41 ATH_MSG_DEBUG(" - PtMin: " << m_ptMin.value());
42 ATH_MSG_DEBUG(" - NTracksMin: " << m_numTrackMin.value());
43 ATH_MSG_DEBUG(" - NTracksMax: " << m_numTrackMax.value());
44 ATH_MSG_DEBUG(" - NIsoTracksMax: " << m_numIsoTrackMax.value());
45
46 ATH_MSG_DEBUG(" - massTrkSysMin: " << m_massTrkSysMin.value());
47 ATH_MSG_DEBUG(" - massTrkSysMax: " << m_massTrkSysMax.value());
48 ATH_MSG_DEBUG(" - massTrkSysKaonMin: " << m_massTrkSysKaonMin.value());
49 ATH_MSG_DEBUG(" - massTrkSysKaonMax: " << m_massTrkSysKaonMax.value());
50 ATH_MSG_DEBUG(" - massTrkSysKaonPiMin: " << m_massTrkSysKaonPiMin.value());
51 ATH_MSG_DEBUG(" - massTrkSysKaonPiMax: " << m_massTrkSysKaonPiMax.value());
52 ATH_MSG_DEBUG(" - targetMassTrkSysKaonPi: " << m_targetMassTrkSysKaonPi.value());
53 ATH_MSG_DEBUG(" - leadTrkPtMin: " << m_leadTrkPtMin.value());
54 ATH_MSG_DEBUG(" - dRmaxMax: " << m_dRmaxMax.value());
55 ATH_MSG_DEBUG(" - etOverPtLeadTrkMin: " << m_etOverPtLeadTrkMin.value());
56 ATH_MSG_DEBUG(" - etOverPtLeadTrkMax: " << m_etOverPtLeadTrkMax.value());
57 ATH_MSG_DEBUG(" - EMPOverTrkSysPMax: " << m_EMPOverTrkSysPMax.value());
58
62 ATH_MSG_ERROR("Invalid tool configuration!");
63 return StatusCode::FAILURE;
64 }
65
66 return StatusCode::SUCCESS;
67}
#define ATH_MSG_ERROR(x)

Member Data Documentation

◆ m_acceptAll

Gaudi::Property<bool> TrigTauPrecisionDiKaonHypoTool::m_acceptAll {this, "AcceptAll", false, "Ignore selection"}
private

Definition at line 64 of file TrigTauPrecisionDiKaonHypoTool.h.

64{this, "AcceptAll", false, "Ignore selection"};

◆ m_decisionId

HLT::Identifier TrigTauPrecisionDiKaonHypoTool::m_decisionId
private

Definition at line 40 of file TrigTauPrecisionDiKaonHypoTool.h.

◆ m_dRmaxMax

Gaudi::Property<float> TrigTauPrecisionDiKaonHypoTool::m_dRmaxMax {this, "dRmaxMax", 10, "Maximum dRmax value"}
private

Definition at line 59 of file TrigTauPrecisionDiKaonHypoTool.h.

59{this, "dRmaxMax", 10, "Maximum dRmax value"};

◆ m_EMPOverTrkSysPMax

Gaudi::Property<float> TrigTauPrecisionDiKaonHypoTool::m_EMPOverTrkSysPMax {this, "EMPOverTrkSysPMax", 5, "Maximum Cluster pt over ditrack pt"}
private

Definition at line 62 of file TrigTauPrecisionDiKaonHypoTool.h.

62{this, "EMPOverTrkSysPMax", 5, "Maximum Cluster pt over ditrack pt"};

◆ m_etOverPtLeadTrkMax

Gaudi::Property<float> TrigTauPrecisionDiKaonHypoTool::m_etOverPtLeadTrkMax {this, "etOverPtLeadTrkMax", 10, "Maximum et/pt(lead track)"}
private

Definition at line 60 of file TrigTauPrecisionDiKaonHypoTool.h.

60{this, "etOverPtLeadTrkMax", 10, "Maximum et/pt(lead track)"};

◆ m_etOverPtLeadTrkMin

Gaudi::Property<float> TrigTauPrecisionDiKaonHypoTool::m_etOverPtLeadTrkMin {this, "etOverPtLeadTrkMin", 0, "Minimum et/pt(lead track)"}
private

Definition at line 61 of file TrigTauPrecisionDiKaonHypoTool.h.

61{this, "etOverPtLeadTrkMin", 0, "Minimum et/pt(lead track)"};

◆ m_leadTrkPtMin

Gaudi::Property<float> TrigTauPrecisionDiKaonHypoTool::m_leadTrkPtMin {this, "leadTrkPtMin", 0, "Minimum Pt of Lead Track"}
private

Definition at line 58 of file TrigTauPrecisionDiKaonHypoTool.h.

58{this, "leadTrkPtMin", 0, "Minimum Pt of Lead Track"};

◆ m_massTrkSysKaonMax

Gaudi::Property<float> TrigTauPrecisionDiKaonHypoTool::m_massTrkSysKaonMax {this, "massTrkSysKaonMax", 1000000000, "Maximum DiKaon mass value"}
private

Definition at line 52 of file TrigTauPrecisionDiKaonHypoTool.h.

52{this, "massTrkSysKaonMax", 1000000000, "Maximum DiKaon mass value"};

◆ m_massTrkSysKaonMin

Gaudi::Property<float> TrigTauPrecisionDiKaonHypoTool::m_massTrkSysKaonMin {this, "massTrkSysKaonMin", 0, "Minimum DiKaon mass value"}
private

Definition at line 51 of file TrigTauPrecisionDiKaonHypoTool.h.

51{this, "massTrkSysKaonMin", 0, "Minimum DiKaon mass value"};

◆ m_massTrkSysKaonPiMax

Gaudi::Property<float> TrigTauPrecisionDiKaonHypoTool::m_massTrkSysKaonPiMax {this, "massTrkSysKaonPiMax", 1000000000, "Maximum KaonPi mass value"}
private

Definition at line 55 of file TrigTauPrecisionDiKaonHypoTool.h.

55{this, "massTrkSysKaonPiMax", 1000000000, "Maximum KaonPi mass value"};

◆ m_massTrkSysKaonPiMin

Gaudi::Property<float> TrigTauPrecisionDiKaonHypoTool::m_massTrkSysKaonPiMin {this, "massTrkSysKaonPiMin", 0, "Minimum KaonPi mass value"}
private

Definition at line 54 of file TrigTauPrecisionDiKaonHypoTool.h.

54{this, "massTrkSysKaonPiMin", 0, "Minimum KaonPi mass value"};

◆ m_massTrkSysMax

Gaudi::Property<float> TrigTauPrecisionDiKaonHypoTool::m_massTrkSysMax {this, "massTrkSysMax", 1000000000, "Maximum DiTrack mass value"}
private

Definition at line 49 of file TrigTauPrecisionDiKaonHypoTool.h.

49{this, "massTrkSysMax", 1000000000, "Maximum DiTrack mass value"};

◆ m_massTrkSysMin

Gaudi::Property<float> TrigTauPrecisionDiKaonHypoTool::m_massTrkSysMin {this, "massTrkSysMin", 0, "Minimum DiTrack mass value"}
private

Definition at line 48 of file TrigTauPrecisionDiKaonHypoTool.h.

48{this, "massTrkSysMin", 0, "Minimum DiTrack mass value"};

◆ m_monTool

ToolHandle<GenericMonitoringTool> TrigTauPrecisionDiKaonHypoTool::m_monTool {this, "MonTool", "", "Monitoring tool"}
private

Definition at line 66 of file TrigTauPrecisionDiKaonHypoTool.h.

66{this, "MonTool", "", "Monitoring tool"};

◆ m_numIsoTrackMax

Gaudi::Property<int> TrigTauPrecisionDiKaonHypoTool::m_numIsoTrackMax {this, "NIsoTracksMax", 1, "Maximum number of wide Tracks"}
private

Definition at line 45 of file TrigTauPrecisionDiKaonHypoTool.h.

45{this, "NIsoTracksMax", 1, "Maximum number of wide Tracks"};

◆ m_numTrackMax

Gaudi::Property<int> TrigTauPrecisionDiKaonHypoTool::m_numTrackMax {this, "NTracksMax", 2, "Maximum number of Tracks"}
private

Definition at line 43 of file TrigTauPrecisionDiKaonHypoTool.h.

43{this, "NTracksMax", 2, "Maximum number of Tracks"};

◆ m_numTrackMin

Gaudi::Property<int> TrigTauPrecisionDiKaonHypoTool::m_numTrackMin {this, "NTracksMin", 1, "Minimum number of Tracks"}
private

Definition at line 44 of file TrigTauPrecisionDiKaonHypoTool.h.

44{this, "NTracksMin", 1, "Minimum number of Tracks"};

◆ m_ptMin

Gaudi::Property<float> TrigTauPrecisionDiKaonHypoTool::m_ptMin {this, "PtMin", 0, "Minimum tau pT value"}
private

Definition at line 42 of file TrigTauPrecisionDiKaonHypoTool.h.

42{this, "PtMin", 0, "Minimum tau pT value"};

◆ m_targetMassTrkSysKaonPi

Gaudi::Property<float> TrigTauPrecisionDiKaonHypoTool::m_targetMassTrkSysKaonPi {this, "targetMassTrkSysKaonPi", 0, "Target KaonPi mass value (parameter)"}
private

Definition at line 56 of file TrigTauPrecisionDiKaonHypoTool.h.

56{this, "targetMassTrkSysKaonPi", 0, "Target KaonPi mass value (parameter)"};

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