ATLAS Offline Software
Loading...
Searching...
No Matches
BuildTruthTaus.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// Local include(s)
11
12// Core include(s):
13#include "AthLinks/ElementLink.h"
15
16// EDM include(s):
19
21
22using namespace TauAnalysisTools;
23
24//=================================PUBLIC-PART==================================
25//______________________________________________________________________________
26BuildTruthTaus::BuildTruthTaus( const std::string& name )
27 : AsgMetadataTool(name)
28{
29}
30
31//______________________________________________________________________________
33{
35 ATH_MSG_INFO( "Initializing BuildTruthTaus, will generate " << m_truthTauOutputContainer.key() << " from " << m_truthParticleContainer.key() << " container" );
36 }
37 else {
38 ATH_MSG_INFO( "Initializing BuildTruthTaus in truth matching mode, using input container " << m_truthTauInputContainer.key() );
39 }
40
41 // input containers
47
48 // in truth matching mode, truth electron/muon/jet containers are optional but expected
50 if (m_truthElectronContainer.empty()) {
51 ATH_MSG_WARNING("Truth electron container is not available, won't perform matching to truth electrons");
52 }
53 if (m_truthMuonContainer.empty()) {
54 ATH_MSG_WARNING("Truth muon container is not available, won't perform matching to truth muons");
55 }
56 if (m_truthJetContainer.empty()) {
57 ATH_MSG_WARNING("Truth jet container is not available, won't perform matching to truth jets");
58 }
59 }
60
61 // output container
63
64 // input decorations
69 // output decorations
75
76 // drop at earliest occasion
78
79 return StatusCode::SUCCESS;
80}
81
82
84{
85 const EventContext& ctx = Gaudi::Hive::currentContext();
87}
88
89
90StatusCode BuildTruthTaus::retrieveTruthTaus(ITruthTausEvent& truthTausEvent, const EventContext& ctx) const
91{
92 return retrieveTruthTaus (dynamic_cast<TruthTausEvent&> (truthTausEvent), ctx);
93}
94
95
96StatusCode BuildTruthTaus::retrieveTruthTaus(TruthTausEvent& truthTausEvent, const EventContext& ctx) const
97{
98 // truth matching mode
100 if (!m_truthElectronContainer.empty()) {
102 if (!truthElectronsHandle.isValid()) {
103 ATH_MSG_ERROR ("Could not retrieve " << truthElectronsHandle.key());
104 return StatusCode::FAILURE;
105 }
106 truthTausEvent.m_xTruthElectronContainerConst = truthElectronsHandle.cptr();
107 }
108
109 if (!m_truthMuonContainer.empty()) {
111 if (!truthMuonsHandle.isValid()) {
112 ATH_MSG_ERROR ("Could not retrieve " << truthMuonsHandle.key());
113 return StatusCode::FAILURE;
114 }
115 truthTausEvent.m_xTruthMuonContainerConst = truthMuonsHandle.cptr();
116 }
117
118 if (!m_truthJetContainer.empty()) {
120 if (!truthJetsHandle.isValid()) {
121 ATH_MSG_ERROR ("Could not retrieve " << truthJetsHandle.key());
122 return StatusCode::FAILURE;
123 }
124 truthTausEvent.m_xTruthJetContainerConst = truthJetsHandle.cptr();
125 }
126
128 if (!truthTausHandle.isValid()) {
129 ATH_MSG_ERROR ("Could not retrieve " << truthTausHandle.key());
130 return StatusCode::FAILURE;
131 }
132 truthTausEvent.m_xTruthTauContainerConst = truthTausHandle.cptr();
133 }
134 // truth tau building mode
135 else {
137 if (!truthParticlesHandle.isValid()) {
138 ATH_MSG_ERROR ("Could not retrieve " << truthParticlesHandle.key());
139 return StatusCode::FAILURE;
140 }
141 truthTausEvent.m_xTruthParticleContainer = truthParticlesHandle.cptr();
142
143 auto truthTausOutput = std::make_unique<xAOD::TruthParticleContainer>();
144 auto truthTausOutputAux = std::make_unique<xAOD::TruthParticleAuxContainer>();
145 truthTausOutput->setStore(truthTausOutputAux.get());
146 truthTausEvent.m_xTruthTauContainer = truthTausOutput.get();
147
148 auto writeHandle = SG::makeHandle(m_truthTauOutputContainer, ctx);
149 ATH_CHECK(writeHandle.record(std::move(truthTausOutput), std::move(truthTausOutputAux)));
150
151 ATH_CHECK( buildTruthTausFromTruthParticles(truthTausEvent, ctx) );
152 }
153
154 return StatusCode::SUCCESS;
155}
156
157//=================================PRIVATE-PART=================================
158//______________________________________________________________________________
159//______________________________________________________________________________
160StatusCode
161BuildTruthTaus::buildTruthTausFromTruthParticles(TruthTausEvent& truthTausEvent, const EventContext& ctx) const
162{
167
173
174 for (auto xTruthParticle : *truthTausEvent.m_xTruthParticleContainer)
175 {
176 if ( xTruthParticle->isTau() )
177 {
178 auto xTruthTau = std::make_unique<xAOD::TruthParticle>();
179 xTruthTau->makePrivateStore( *xTruthParticle );
180
181 if ( examineTruthTau(*xTruthTau).isFailure() )
182 {
183 continue;
184 }
185
186 // adding decorations onto the truth tau requires that it is in the final Aux store
187 truthTausEvent.m_xTruthTauContainer->push_back(std::move(xTruthTau));
188 xAOD::TruthParticle* truthTau = truthTausEvent.m_xTruthTauContainer->back();
189
190 // propagate MCTruthClassifier decorations
191 typeDecorator(*truthTau) = typeReadDecor.isAvailable() ? typeReadDecor(*xTruthParticle) : -1234;
192 originDecorator(*truthTau) = originReadDecor.isAvailable() ? originReadDecor(*xTruthParticle) : -1234;
193 outcomeDecorator(*truthTau) = outcomeReadDecor.isAvailable() ? outcomeReadDecor(*xTruthParticle) : -1234;
194 classificationDecorator(*truthTau) = classificationReadDecor.isAvailable() ? classificationReadDecor(*xTruthParticle) : -1234;
195
196 // create link to the original TruthParticle
197 ElementLink < xAOD::TruthParticleContainer > lTruthParticleLink(xTruthParticle, *truthTausEvent.m_xTruthParticleContainer);
198 linkDecorator(*truthTau) = lTruthParticleLink;
199 }
200 }
201 return StatusCode::SUCCESS;
202}
203
204//______________________________________________________________________________
205StatusCode BuildTruthTaus::examineTruthTau(const xAOD::TruthParticle& xTruthParticle) const
206{
207 // skip this tau if it has no decay vertex, should not happen
208 if ( !xTruthParticle.hasDecayVtx() )
209 return StatusCode::FAILURE;
210
211 ATH_MSG_VERBOSE("looking for charged daughters of a truth tau");
212
213 TauTruthInfo truthInfo;
214
215 const xAOD::TruthVertex* xDecayVertex = xTruthParticle.decayVtx();
216 if (xDecayVertex == nullptr)
217 return StatusCode::FAILURE;
218 for ( size_t iOutgoingParticle = 0; iOutgoingParticle < xDecayVertex->nOutgoingParticles(); ++iOutgoingParticle )
219 {
220 const xAOD::TruthParticle* xTruthDaughter = xDecayVertex->outgoingParticle(iOutgoingParticle);
221 if (xTruthDaughter == nullptr)
222 {
223 ATH_MSG_ERROR("Truth daughter of tau decay was not found in "<< m_truthParticleContainer.key() <<" container. Please ensure that this container has the full tau decay information or produce the TruthTaus container in AtlasDerivation.");
224 return StatusCode::FAILURE;
225 }
226
227 // if tau decays into tau this is not a proper tau decay
228 if ( xTruthDaughter->isTau() )
229 {
230 ATH_MSG_VERBOSE("Tau decays into a tau itself -> skip this decay");
231 return StatusCode::FAILURE;
232 }
233 }
234
235 examineTruthTauDecay(xTruthParticle, truthInfo).ignore();
236
237 if (truthInfo.m_bIsHadronicTau)
238 ATH_MSG_VERBOSE(truthInfo.m_iNChargedDaughters << " prong hadronic truth tau was found with uniqueID "<<HepMC::uniqueID(xTruthParticle));
239 else
240 ATH_MSG_VERBOSE(truthInfo.m_iNChargedDaughters << " prong leptonic truth tau was found with uniqueID "<<HepMC::uniqueID(xTruthParticle));
241 if ( truthInfo.m_iNChargedDaughters%2 == 0 )
242 {
243 ATH_MSG_WARNING("found tau with even multiplicity: " << truthInfo.m_iNChargedDaughters);
244 printDecay(xTruthParticle);
245 }
246
247 static const SG::Decorator<double> decPtVis("pt_vis");
248 static const SG::Decorator<double> decEtaVis("eta_vis");
249 static const SG::Decorator<double> decPhiVis("phi_vis");
250 static const SG::Decorator<double> decMVis("m_vis");
251
252 static const SG::Decorator<size_t> decNumCharged("numCharged");
253 static const SG::Decorator<size_t> decNumChargedPion("numChargedPion");
254 static const SG::Decorator<size_t> decNumNeutral("numNeutral");
255 static const SG::Decorator<size_t> decNumNeutralPion("numNeutralPion");
256
257 decPtVis(xTruthParticle) = truthInfo.m_vTruthVisTLV.Pt();
258 decEtaVis(xTruthParticle) = truthInfo.m_vTruthVisTLV.Eta();
259 decPhiVis(xTruthParticle) = truthInfo.m_vTruthVisTLV.Phi();
260 decMVis(xTruthParticle) = truthInfo.m_vTruthVisTLV.M();
261
262 decNumCharged(xTruthParticle) = truthInfo.m_iNChargedDaughters;
263 decNumChargedPion(xTruthParticle) = truthInfo.m_iNChargedPions;
264 decNumNeutral(xTruthParticle) = truthInfo.m_iNNeutralPions+truthInfo.m_iNNeutralOthers;
265 decNumNeutralPion(xTruthParticle) = truthInfo.m_iNNeutralPions;
266
267 static const SG::Decorator<char> decIsHadronicTau("IsHadronicTau");
268 decIsHadronicTau(xTruthParticle) = (char)truthInfo.m_bIsHadronicTau;
269
271 {
272 TLorentzVector vTruthInvisTLV = xTruthParticle.p4() - truthInfo.m_vTruthVisTLV;
273 static const SG::Decorator<double> decPtInvis("pt_invis");
274 static const SG::Decorator<double> decEtaInvis("eta_invis");
275 static const SG::Decorator<double> decPhiInvis("phi_invis");
276 static const SG::Decorator<double> decMInvis("m_invis");
277 decPtInvis(xTruthParticle) = vTruthInvisTLV.Pt();
278 decEtaInvis(xTruthParticle) = vTruthInvisTLV.Eta();
279 decPhiInvis(xTruthParticle) = vTruthInvisTLV.Phi();
280 decMInvis(xTruthParticle) = vTruthInvisTLV.M();
281 }
282
284 {
285 static const SG::Decorator<double> decPtVisCharged("pt_vis_charged");
286 static const SG::Decorator<double> decEtaVisCharged("eta_vis_charged");
287 static const SG::Decorator<double> decPhiVisCharged("phi_vis_charged");
288 static const SG::Decorator<double> decMVisCharged("m_vis_charged");
289 decPtVisCharged(xTruthParticle) = truthInfo.m_vTruthVisTLVCharged.Pt();
290 decEtaVisCharged(xTruthParticle) = truthInfo.m_vTruthVisTLVCharged.Eta();
291 decPhiVisCharged(xTruthParticle) = truthInfo.m_vTruthVisTLVCharged.Phi();
292 decMVisCharged(xTruthParticle) = truthInfo.m_vTruthVisTLVCharged.M();
293 }
294
296 {
297 static const SG::Decorator<double> decPtVisNeutral("pt_vis_neutral");
298 static const SG::Decorator<double> decEtaVisNeutral("eta_vis_neutral");
299 static const SG::Decorator<double> decPhiVisNeutral("phi_vis_neutral");
300 static const SG::Decorator<double> decMVisNeutral("m_vis_neutral");
301 decPtVisNeutral(xTruthParticle) = truthInfo.m_vTruthVisTLVNeutral.Pt();
302 decEtaVisNeutral(xTruthParticle) = truthInfo.m_vTruthVisTLVNeutral.Eta();
303 decPhiVisNeutral(xTruthParticle) = truthInfo.m_vTruthVisTLVNeutral.Phi();
304 decMVisNeutral(xTruthParticle) = truthInfo.m_vTruthVisTLVNeutral.M();
305 }
306
308 {
309 static const SG::Decorator<std::vector<int> > decDecayModeVector("DecayModeVector");
310 decDecayModeVector(xTruthParticle) = truthInfo.m_vDecayMode;
311 }
312
313 if ( m_bWriteVertices )
314 {
315 // tau decay vertex
316 static const SG::Decorator<float> decDecayVertexX("decay_vertex_x");
317 static const SG::Decorator<float> decDecayVertexY("decay_vertex_y");
318 static const SG::Decorator<float> decDecayVertexZ("decay_vertex_z");
319
320 decDecayVertexX(xTruthParticle) = truthInfo.m_vDecayVertex.X();
321 decDecayVertexY(xTruthParticle) = truthInfo.m_vDecayVertex.Y();
322 decDecayVertexZ(xTruthParticle) = truthInfo.m_vDecayVertex.Z();
323
324 // tau production vertex
325 static const SG::Decorator<float> decProdVertexX("prod_vertex_x");
326 static const SG::Decorator<float> decProdVertexY("prod_vertex_y");
327 static const SG::Decorator<float> decProdVertexZ("prod_vertex_z");
328
329 decProdVertexX(xTruthParticle) = truthInfo.m_vProdVertex.X();
330 decProdVertexY(xTruthParticle) = truthInfo.m_vProdVertex.Y();
331 decProdVertexZ(xTruthParticle) = truthInfo.m_vProdVertex.Z();
332 }
333
334 return StatusCode::SUCCESS;
335}
336
337//______________________________________________________________________________
339 TauTruthInfo& truthInfo) const
340{
341 // get vertex and check if it exists
342 const xAOD::TruthVertex* xDecayVertex = xTruthParticle.decayVtx();
343 if (!xDecayVertex)
344 return StatusCode::SUCCESS;
345
346 truthInfo.m_vDecayVertex.SetXYZ(xDecayVertex->x(),xDecayVertex->y(),xDecayVertex->z());
347
348 if (xTruthParticle.hasProdVtx() ) {
349 const xAOD::TruthVertex* xProdVertex = xTruthParticle.prodVtx();
350 truthInfo.m_vProdVertex.SetXYZ(xProdVertex->x(),xProdVertex->y(),xProdVertex->z());
351 } else {
352 truthInfo.m_vProdVertex.SetXYZ(-1234,-1234,-1234);
353 }
354
355 for ( size_t iOutgoingParticle = 0; iOutgoingParticle < xDecayVertex->nOutgoingParticles(); ++iOutgoingParticle )
356 {
357 const xAOD::TruthParticle* xTruthDaughter = xDecayVertex->outgoingParticle(iOutgoingParticle);
358 if (xTruthDaughter == nullptr)
359 {
360 ATH_MSG_ERROR("Truth daughter of tau decay was not found in "<< m_truthParticleContainer.key() <<" container. Please ensure that this container has the full tau decay information or produce the TruthTaus container in AtlasDerivation.");
361 return StatusCode::FAILURE;
362 }
363
364 int iAbsPdgId = xTruthDaughter->absPdgId();
365 int iPdgId = xTruthDaughter->pdgId();
366
367 // look at decay of unstable particles
368 if (MC::isDecayed(xTruthDaughter) || !MC::isPhysical(xTruthDaughter))
369 {
370 if ( iAbsPdgId != 111 && iAbsPdgId != 311 && iAbsPdgId != 310 && iAbsPdgId != 130 )
371 {
372 examineTruthTauDecay(*xTruthDaughter, truthInfo).ignore();
373 continue;
374 }
375 }
376
377 // only process stable particles
378 if (!MC::isStable(xTruthDaughter) && !MC::isDecayed(xTruthDaughter))
379 continue;
380
381 // add pdgID to vector for decay mode classification
382 truthInfo.m_vDecayMode.push_back(iPdgId);
383
384 // if tau decays leptonically, indicated by an electron/muon neutrino then
385 // it is not a hadronic decay
386 if ( xTruthDaughter->isHadron() )
387 truthInfo.m_bIsHadronicTau = true;
388
389 // ignore neutrinos for further progress
390 if ( xTruthDaughter->isNeutrino() )
391 {
392 ATH_MSG_VERBOSE("found neutrino decay particle with PdgId "<<iPdgId);
393 continue;
394 }
395
396 // add momentum of non-neutrino particle to visible momentum
397 truthInfo.m_vTruthVisTLV += xTruthDaughter->p4();
399 if ( xTruthDaughter->isCharged() )
400 truthInfo.m_vTruthVisTLVCharged += xTruthDaughter->p4();
402 if ( xTruthDaughter->isNeutral() )
403 truthInfo.m_vTruthVisTLVNeutral += xTruthDaughter->p4();
404
405 // only count charged decay particles
406 if ( xTruthDaughter->isCharged() )
407 {
408 ATH_MSG_VERBOSE("found charged decay particle with PdgId "<<iPdgId);
409 truthInfo.m_iNChargedDaughters++;
410 // count charged pions
411 if (iAbsPdgId==211) truthInfo.m_iNChargedPions++;
412 else truthInfo.m_iNChargedOthers++;
413 }
414 else
415 {
416 ATH_MSG_VERBOSE("found neutral decay particle with PdgId "<<iPdgId);
417 // count neutral pions
418 if (iAbsPdgId==111) truthInfo.m_iNNeutralPions++;
419 else truthInfo.m_iNNeutralOthers++;
420 }
421 }
422 return StatusCode::SUCCESS;
423}
424
425void BuildTruthTaus::printDecay(const xAOD::TruthParticle& xTruthParticle, int depth) const
426{
427 // loop over all decay particles, print their kinematic and other properties
428
429 const xAOD::TruthVertex* xDecayVertex = xTruthParticle.decayVtx();
430 if (xDecayVertex == nullptr)
431 return;
432
433 for ( size_t iOutgoingParticle = 0; iOutgoingParticle < xDecayVertex->nOutgoingParticles(); ++iOutgoingParticle )
434 {
435 const xAOD::TruthParticle* xTruthDaughter = xDecayVertex->outgoingParticle(iOutgoingParticle);
436 if (xTruthDaughter == nullptr)
437 {
438 ATH_MSG_WARNING("Truth daughter of tau decay was not found in "<< m_truthParticleContainer.key() <<" container. Please ensure that this container has the full tau decay information or produce the TruthTaus container in AtlasDerivation.");
439 return;
440 }
441 ATH_MSG_WARNING("depth "<<depth
442 <<" e "<<xTruthDaughter->e()
443 <<" eta "<<xTruthDaughter->p4().Eta()
444 <<" phi "<<xTruthDaughter->p4().Phi()
445 <<" pdgid "<<xTruthDaughter->pdgId()
446 <<" status "<<xTruthDaughter->status()
447 <<" uniqueID "<<HepMC::uniqueID(xTruthDaughter));
448 printDecay(*xTruthDaughter, depth+1);
449 }
450}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
Handle class for reading a decoration on an object.
Handle class for reading from StoreGate.
Handle class for adding a decoration to an object.
Handle class for recording to StoreGate.
Helper class to provide type-safe access to aux data.
Helper class to provide type-safe access to aux data.
Definition Decorator.h:59
Handle class for reading a decoration on an object.
bool isAvailable()
Test to see if this variable exists in the store, for the referenced object.
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type cptr()
Dereference the pointer.
virtual const std::string & key() const override final
Return the StoreGate ID for the referenced object.
Handle class for adding a decoration to an object.
Gaudi::Property< bool > m_bWriteDecayModeVector
SG::ReadHandleKey< xAOD::TruthParticleContainer > m_truthElectronContainer
SG::ReadHandleKey< xAOD::TruthParticleContainer > m_truthMuonContainer
void printDecay(const xAOD::TruthParticle &xTruthParticle, int depth=0) const
SG::ReadDecorHandleKey< xAOD::TruthParticleContainer > m_outcomeReadDecorKey
SG::ReadHandleKey< xAOD::TruthParticleContainer > m_truthParticleContainer
SG::ReadDecorHandleKey< xAOD::TruthParticleContainer > m_originReadDecorKey
SG::ReadHandleKey< xAOD::JetContainer > m_truthJetContainer
Gaudi::Property< bool > m_bWriteVisibleChargedFourMomentum
ASG_TOOL_CLASS(BuildTruthTaus, TauAnalysisTools::IBuildTruthTaus) public BuildTruthTaus(const std::string &name)
Create a proper constructor for Athena.
Gaudi::Property< bool > m_bWriteVertices
SG::WriteDecorHandleKey< xAOD::TruthParticleContainer > m_originDecoratorKey
SG::ReadDecorHandleKey< xAOD::TruthParticleContainer > m_classificationReadDecorKey
SG::ReadDecorHandleKey< xAOD::TruthParticleContainer > m_typeReadDecorKey
virtual StatusCode retrieveTruthTaus() override
SG::WriteDecorHandleKey< xAOD::TruthParticleContainer > m_linkDecoratorKey
SG::WriteDecorHandleKey< xAOD::TruthParticleContainer > m_typeDecoratorKey
Gaudi::Property< bool > m_bWriteInvisibleFourMomentum
Gaudi::Property< bool > m_bWriteVisibleNeutralFourMomentum
StatusCode examineTruthTauDecay(const xAOD::TruthParticle &xTruthParticle, TauTruthInfo &truthInfo) const
SG::WriteHandleKey< xAOD::TruthParticleContainer > m_truthTauOutputContainer
virtual StatusCode initialize() override
StatusCode examineTruthTau(const xAOD::TruthParticle &xTruthParticle) const
StatusCode buildTruthTausFromTruthParticles(TruthTausEvent &truthTausEvent, const EventContext &ctx) const
SG::ReadHandleKey< xAOD::TruthParticleContainer > m_truthTauInputContainer
SG::WriteDecorHandleKey< xAOD::TruthParticleContainer > m_classificationDecoratorKey
SG::WriteDecorHandleKey< xAOD::TruthParticleContainer > m_outcomeDecoratorKey
Declare the interface that the class provides.
AsgMetadataTool(const std::string &name)
Normal ASG tool constructor with a name.
bool isNeutrino() const
Whether the particle is a neutrino (or antineutrino)
int status() const
Status code.
bool isHadron() const
Whether the particle is a hadron.
int pdgId() const
PDG ID code.
const TruthVertex_v1 * decayVtx() const
The decay vertex of this particle.
bool isCharged() const
Whether the particle is electrically charged.
int absPdgId() const
Absolute PDG ID code (often useful)
bool hasProdVtx() const
Check for a production vertex on this particle.
bool hasDecayVtx() const
Check for a decay vertex on this particle.
virtual double e() const override final
The total energy of the particle.
bool isNeutral() const
Whether the particle is electrically neutral.
const TruthVertex_v1 * prodVtx() const
The production vertex of this particle.
virtual FourMom_t p4() const override final
The full 4-momentum of the particle.
bool isTau() const
Whether the particle is a tau (or antitau)
const TruthParticle_v1 * outgoingParticle(size_t index) const
Get one of the outgoing particles.
float z() const
Vertex longitudinal distance along the beam line form the origin.
float y() const
Vertex y displacement.
size_t nOutgoingParticles() const
Get the number of outgoing particles.
float x() const
Vertex x displacement.
std::string depth
tag string for intendation
Definition fastadd.cxx:46
int uniqueID(const T &p)
bool isStable(const T &p)
Identify if the particle is stable, i.e. has not decayed.
bool isDecayed(const T &p)
Identify if the particle decayed.
bool isPhysical(const T &p)
Identify if the particle is physical, i.e. is stable or decayed.
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
TruthVertex_v1 TruthVertex
Typedef to implementation.
Definition TruthVertex.h:15
TruthParticle_v1 TruthParticle
Typedef to implementation.