Loading [MathJax]/extensions/tex2jax.js
ATLAS Offline Software
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
xAODTruthParticleSlimmerMET.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 #include "AthLinks/ElementLink.h"
7 
9 
10 #include "GaudiKernel/MsgStream.h"
11 #include "GaudiKernel/DataSvc.h"
12 #include "GaudiKernel/PhysicalConstants.h"
13 
17 
19 
22 
24 
25 xAODTruthParticleSlimmerMET::xAODTruthParticleSlimmerMET(const std::string &name, ISvcLocator *svcLoc)
26  : AthAlgorithm(name, svcLoc)
27 {
28 }
29 
31 {
33  ATH_MSG_INFO("xAOD output TruthParticleContainerMET name = " << m_xaodTruthParticleContainerNameMET.key());
35  ATH_MSG_INFO("xAOD input xAODTruthEventContainerName name = " << m_xaodTruthEventContainerName.key());
36  ATH_CHECK(m_classif.retrieve());
37  return StatusCode::SUCCESS;
38 }
39 
41 {
42  // If the containers already exists then assume that nothing needs to be done
43  if (evtStore()->contains<xAOD::TruthParticleContainer>(m_xaodTruthParticleContainerNameMET.key()))
44  {
45  ATH_MSG_WARNING("xAOD MET Truth Particles are already available in the event");
46  return StatusCode::SUCCESS;
47  }
48 
49  // Create new output container
51  ATH_CHECK(xTruthParticleContainerMET.record(std::make_unique<xAOD::TruthParticleContainer>(), std::make_unique<xAOD::TruthParticleAuxContainer>()));
52  ATH_MSG_INFO("Recorded TruthParticleContainerMET with key: " << m_xaodTruthParticleContainerNameMET.key());
53 
54  // Retrieve full TruthEventContainer container
56  if ( !xTruthEventContainer.isValid() )
57  {
58  ATH_MSG_ERROR("No TruthEvent collection with name " << m_xaodTruthEventContainerName.key() << " found in StoreGate!");
59  return StatusCode::FAILURE;
60  }
61 
62  // Set up decorators if needed
63  const static SG::AuxElement::Decorator<bool> isPrompt("isPrompt");
64 
65  // Loop over full TruthParticle container
67  for (itr = xTruthEventContainer->begin(); itr!=xTruthEventContainer->end(); ++itr) {
68 
69  unsigned int nPart = (*itr)->nTruthParticles();
70  std::vector<int> uniqueID_list;
71  int zero_uniqueID=0;
72  int dup_uniqueID=0;
73 
74  for (unsigned int iPart = 0; iPart < nPart; ++iPart) {
75  const xAOD::TruthParticle* theParticle = (*itr)->truthParticle(iPart);
76 
77  int my_uniqueID = HepMC::uniqueID(theParticle);
78  if ( my_uniqueID == HepMC::UNDEFINED_ID ) {
79  zero_uniqueID++;
80  continue;
81  }
82  bool found = false;
83  if (uniqueID_list.size() > 0){
84  found = (std::find(uniqueID_list.begin(), uniqueID_list.end(), my_uniqueID) != uniqueID_list.end());
85  if(found) {
86  dup_uniqueID++;
87  continue;}
88  }
89  uniqueID_list.push_back(my_uniqueID);
90 
91 
92 
93  // stable and non-interacting, implemented from DerivationFramework
94  //https://gitlab.cern.ch/atlas/athena/-/blob/master/PhysicsAnalysis/DerivationFramework/DerivationFrameworkMCTruth/python/MCTruthCommon.py#L183
95  // which in turn use the implementation from Reconstruction
96  //https://gitlab.cern.ch/atlas/athena/blob/21.0/Reconstruction/MET/METReconstruction/Root/METTruthTool.cxx#L143
97  if (!theParticle->isGenStable()) continue;
98  if (MC::isInteracting(theParticle)) continue;
99 
100 
101  xAOD::TruthParticle *xTruthParticle = new xAOD::TruthParticle();
102  xTruthParticleContainerMET->push_back( xTruthParticle );
103 
104  // Fill with numerical content
105  *xTruthParticle=*theParticle;
106 
107  //Decorate
108  isPrompt(*xTruthParticle) = Common::prompt(theParticle,m_classif);
109  }
110  if(zero_uniqueID != 0 || dup_uniqueID !=0 ) ATH_MSG_INFO("Found " << zero_uniqueID << " uniqueID=0 particles and " <<dup_uniqueID <<"duplicated");
111  }
112 
113  return StatusCode::SUCCESS;
114 }
115 
116 
DataModel_detail::const_iterator
Const iterator class for DataVector/DataList.
Definition: DVLIterator.h:82
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
find
std::string find(const std::string &s)
return a remapped string
Definition: hcg.cxx:135
SG::ReadHandle
Definition: StoreGate/StoreGate/ReadHandle.h:67
TruthParticleContainer.h
xAOD::TruthParticle_v1::isGenStable
bool isGenStable() const
Check if this is generator stable particle.
Definition: TruthParticle_v1.cxx:316
Common::prompt
bool prompt(const xAOD::TruthParticle *part, ToolHandle< IMCTruthClassifier > &m_classif)
Definition: Common.cxx:7
SG::VarHandleKey::key
const std::string & key() const
Return the StoreGate ID for the referenced object.
Definition: AthToolSupport/AsgDataHandles/Root/VarHandleKey.cxx:141
TruthParticleAuxContainer.h
AthCommonDataStore< AthCommonMsg< Algorithm > >::evtStore
ServiceHandle< StoreGateSvc > & evtStore()
The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
Definition: AthCommonDataStore.h:85
xAODTruthParticleSlimmerMET.h
IMCTruthClassifier.h
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
SG::Decorator
Helper class to provide type-safe access to aux data.
Definition: Decorator.h:59
xAODTruthParticleSlimmerMET::initialize
virtual StatusCode initialize()
Function initialising the algorithm.
Definition: xAODTruthParticleSlimmerMET.cxx:30
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
xAOD::TruthParticle_v1
Class describing a truth particle in the MC record.
Definition: TruthParticle_v1.h:37
xAOD::TruthParticle
TruthParticle_v1 TruthParticle
Typedef to implementation.
Definition: Event/xAOD/xAODTruth/xAODTruth/TruthParticle.h:15
HepMC::uniqueID
int uniqueID(const T &p)
Definition: MagicNumbers.h:116
ATH_CHECK
#define ATH_CHECK
Definition: AthCheckMacros.h:40
xAODTruthParticleSlimmerMET::execute
virtual StatusCode execute()
Function executing the algorithm.
Definition: xAODTruthParticleSlimmerMET.cxx:40
SG::VarHandleKey::initialize
StatusCode initialize(bool used=true)
If this object is used as a property, then this should be called during the initialize phase.
Definition: AthToolSupport/AsgDataHandles/Root/VarHandleKey.cxx:103
AthAlgorithm
Definition: AthAlgorithm.h:47
HepMC::UNDEFINED_ID
constexpr int UNDEFINED_ID
Definition: MagicNumbers.h:56
Common.h
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:240
errorcheck.h
Helpers for checking error return status codes and reporting errors.
DataVector::push_back
value_type push_back(value_type pElem)
Add an element to the end of the collection.
checkTriggerxAOD.found
found
Definition: checkTriggerxAOD.py:328
xAODTruthParticleSlimmerMET::m_xaodTruthParticleContainerNameMET
SG::WriteHandleKey< xAOD::TruthParticleContainer > m_xaodTruthParticleContainerNameMET
The key for the output xAOD truth containers.
Definition: xAODTruthParticleSlimmerMET.h:37
SG::WriteHandle
Definition: StoreGate/StoreGate/WriteHandle.h:73
xAODTruthParticleSlimmerMET::m_xaodTruthEventContainerName
SG::ReadHandleKey< xAOD::TruthEventContainer > m_xaodTruthEventContainerName
Definition: xAODTruthParticleSlimmerMET.h:34
SG::WriteHandle::record
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
MCTruthPartClassifier::isPrompt
int isPrompt(const unsigned int classify, bool allow_prompt_tau_decays=true)
Definition: TruthClassifiers.h:180
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
xAODTruthParticleSlimmerMET::xAODTruthParticleSlimmerMET
xAODTruthParticleSlimmerMET(const std::string &name, ISvcLocator *svcLoc)
Regular algorithm constructor.
Definition: xAODTruthParticleSlimmerMET.cxx:25
xAODTruthParticleSlimmerMET::m_classif
PublicToolHandle< IMCTruthClassifier > m_classif
Definition: xAODTruthParticleSlimmerMET.h:39
MC::isInteracting
bool isInteracting(const T &p)
Identify if the particle with given PDG ID would not interact with the detector, i....
Definition: HepMCHelpers.h:33
TruthParticle.h
HepMCHelpers.h