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

Algorithm to skim the xAOD truth particle container for xAOD MET filter. More...

#include <xAODTruthParticleSlimmerMET.h>

Inheritance diagram for xAODTruthParticleSlimmerMET:
Collaboration diagram for xAODTruthParticleSlimmerMET:

Public Member Functions

 xAODTruthParticleSlimmerMET (const std::string &name, ISvcLocator *svcLoc)
 Regular algorithm constructor.
virtual StatusCode initialize ()
 Function initialising the algorithm.
virtual StatusCode execute ()
 Function executing the algorithm.
virtual StatusCode sysInitialize () override
 Override sysInitialize.
virtual const DataObjIDColl & extraOutputDeps () const override
 Return the list of extra output dependencies.
ServiceHandle< StoreGateSvc > & evtStore ()
 The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
const ServiceHandle< StoreGateSvc > & detStore () const
 The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
virtual StatusCode sysStart () override
 Handle START transition.
virtual std::vector< Gaudi::DataHandle * > inputHandles () const override
 Return this algorithm's input handles.
virtual std::vector< Gaudi::DataHandle * > outputHandles () const override
 Return this algorithm's output handles.
Gaudi::Details::PropertyBase & declareProperty (Gaudi::Property< T, V, H > &t)
void updateVHKA (Gaudi::Details::PropertyBase &)
MsgStream & msg () const
bool msgLvl (const MSG::Level lvl) const

Protected Member Functions

void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce (T &h)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

SG::ReadHandleKey< xAOD::TruthEventContainerm_xaodTruthEventContainerName {this, "xAODTruthEventContainerName", "TruthEvents"}
SG::WriteHandleKey< xAOD::TruthParticleContainerm_xaodTruthParticleContainerNameMET {this, "xAODTruthParticleContainerNameMET","TruthMET"}
 The key for the output xAOD truth containers.
PublicToolHandle< IMCTruthClassifierm_classif {this, "MCTruthClassifier", "MCTruthClassifier/DFCommonTruthClassifier"}
DataObjIDColl m_extendedExtraObjects
StoreGateSvc_t m_evtStore
 Pointer to StoreGate (event store by default)
StoreGateSvc_t m_detStore
 Pointer to StoreGate (detector store by default)
std::vector< SG::VarHandleKeyArray * > m_vhka
bool m_varHandleArraysDeclared

Detailed Description

Algorithm to skim the xAOD truth particle container for xAOD MET filter.

This algorithm is used to copy and skim the particles from the xAOD TruthParticles container, keeping just relevant MET particles from the event. The design of this class heavily mirrors the DerivationFramework::TruthCollectionMaker.

Author
Jeff Dandoy Jeff..nosp@m.Dand.nosp@m.oy@ce.nosp@m.rn.c.nosp@m.h

Definition at line 22 of file xAODTruthParticleSlimmerMET.h.

Member Typedef Documentation

◆ StoreGateSvc_t

typedef ServiceHandle<StoreGateSvc> AthCommonDataStore< AthCommonMsg< Algorithm > >::StoreGateSvc_t
privateinherited

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ xAODTruthParticleSlimmerMET()

xAODTruthParticleSlimmerMET::xAODTruthParticleSlimmerMET ( const std::string & name,
ISvcLocator * svcLoc )

Regular algorithm constructor.

Definition at line 25 of file xAODTruthParticleSlimmerMET.cxx.

26 : AthAlgorithm(name, svcLoc)
27{
28}
AthAlgorithm()
Default constructor:

Member Function Documentation

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Algorithm > >::declareGaudiProperty ( Gaudi::Property< T, V, H > & hndl,
const SG::VarHandleKeyType &  )
inlineprivateinherited

specialization for handling Gaudi::Property<SG::VarHandleKey>

Definition at line 156 of file AthCommonDataStore.h.

158 {
160 hndl.value(),
161 hndl.documentation());
162
163 }
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)

◆ declareProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Algorithm > >::declareProperty ( Gaudi::Property< T, V, H > & t)
inlineinherited

Definition at line 145 of file AthCommonDataStore.h.

145 {
146 typedef typename SG::HandleClassifier<T>::type htype;
148 }
Gaudi::Details::PropertyBase & declareGaudiProperty(Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
specialization for handling Gaudi::Property<SG::VarHandleKey>

◆ detStore()

const ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< Algorithm > >::detStore ( ) const
inlineinherited

The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 95 of file AthCommonDataStore.h.

◆ evtStore()

ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< Algorithm > >::evtStore ( )
inlineinherited

The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 85 of file AthCommonDataStore.h.

◆ execute()

StatusCode xAODTruthParticleSlimmerMET::execute ( )
virtual

Function executing the algorithm.

Definition at line 40 of file xAODTruthParticleSlimmerMET.cxx.

41{
42 // If the containers already exists then assume that nothing needs to be done
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
50 SG::WriteHandle<xAOD::TruthParticleContainer> xTruthParticleContainerMET(m_xaodTruthParticleContainerNameMET);
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
55 SG::ReadHandle<xAOD::TruthEventContainer> xTruthEventContainer{m_xaodTruthEventContainerName};
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}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
virtual bool isValid() override final
Can the handle be successfully dereferenced?
SG::WriteHandleKey< xAOD::TruthParticleContainer > m_xaodTruthParticleContainerNameMET
The key for the output xAOD truth containers.
PublicToolHandle< IMCTruthClassifier > m_classif
SG::ReadHandleKey< xAOD::TruthEventContainer > m_xaodTruthEventContainerName
bool isGenStable() const
Check if this is generator stable particle.
bool contains(const std::string &s, const std::string &regx)
does a string contain the substring
Definition hcg.cxx:114
bool prompt(const xAOD::TruthParticle *part, ToolHandle< IMCTruthClassifier > &m_classif)
Definition Common.cxx:7
int uniqueID(const T &p)
constexpr int UNDEFINED_ID
int isPrompt(const unsigned int classify, bool allow_prompt_tau_decays=true)
bool isInteracting(const T &p)
Identify if the particle with given PDG ID would not interact with the detector, i....
TruthParticle_v1 TruthParticle
Typedef to implementation.

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::extraDeps_update_handler ( Gaudi::Details::PropertyBase & ExtraDeps)
protectedinherited

Add StoreName to extra input/output deps as needed.

use the logic of the VarHandleKey to parse the DataObjID keys supplied via the ExtraInputs and ExtraOuputs Properties to add the StoreName if it's not explicitly given

◆ extraOutputDeps()

const DataObjIDColl & AthAlgorithm::extraOutputDeps ( ) const
overridevirtualinherited

Return the list of extra output dependencies.

This list is extended to include symlinks implied by inheritance relations.

Definition at line 50 of file AthAlgorithm.cxx.

51{
52 // If we didn't find any symlinks to add, just return the collection
53 // from the base class. Otherwise, return the extended collection.
54 if (!m_extendedExtraObjects.empty()) {
56 }
57 return Algorithm::extraOutputDeps();
58}
DataObjIDColl m_extendedExtraObjects

◆ initialize()

StatusCode xAODTruthParticleSlimmerMET::initialize ( )
virtual

Function initialising the algorithm.

Definition at line 30 of file xAODTruthParticleSlimmerMET.cxx.

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}

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Algorithm > >::inputHandles ( ) const
overridevirtualinherited

Return this algorithm's input handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ msg()

MsgStream & AthCommonMsg< Algorithm >::msg ( ) const
inlineinherited

Definition at line 24 of file AthCommonMsg.h.

24 {
25 return this->msgStream();
26 }

◆ msgLvl()

bool AthCommonMsg< Algorithm >::msgLvl ( const MSG::Level lvl) const
inlineinherited

Definition at line 30 of file AthCommonMsg.h.

30 {
31 return this->msgLevel(lvl);
32 }

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Algorithm > >::outputHandles ( ) const
overridevirtualinherited

Return this algorithm's output handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ renounce()

std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > AthCommonDataStore< AthCommonMsg< Algorithm > >::renounce ( T & h)
inlineprotectedinherited

Definition at line 380 of file AthCommonDataStore.h.

381 {
382 h.renounce();
384 }
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)

◆ renounceArray()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ sysInitialize()

StatusCode AthAlgorithm::sysInitialize ( )
overridevirtualinherited

Override sysInitialize.

Override sysInitialize from the base class.

Loop through all output handles, and if they're WriteCondHandles, automatically register them and this Algorithm with the CondSvc

Scan through all outputHandles, and if they're WriteCondHandles, register them with the CondSvc

Reimplemented from AthCommonDataStore< AthCommonMsg< Algorithm > >.

Reimplemented in AthAnalysisAlgorithm, AthFilterAlgorithm, AthHistogramAlgorithm, and PyAthena::Alg.

Definition at line 66 of file AthAlgorithm.cxx.

66 {
68
69 if (sc.isFailure()) {
70 return sc;
71 }
72 ServiceHandle<ICondSvc> cs("CondSvc",name());
73 for (auto h : outputHandles()) {
74 if (h->isCondition() && h->mode() == Gaudi::DataHandle::Writer) {
75 // do this inside the loop so we don't create the CondSvc until needed
76 if ( cs.retrieve().isFailure() ) {
77 ATH_MSG_WARNING("no CondSvc found: won't autoreg WriteCondHandles");
78 return StatusCode::SUCCESS;
79 }
80 if (cs->regHandle(this,*h).isFailure()) {
81 sc = StatusCode::FAILURE;
82 ATH_MSG_ERROR("unable to register WriteCondHandle " << h->fullKey()
83 << " with CondSvc");
84 }
85 }
86 }
87 return sc;
88}
static Double_t sc
virtual StatusCode sysInitialize() override
Override sysInitialize.
AthCommonDataStore(const std::string &name, T... args)
virtual std::vector< Gaudi::DataHandle * > outputHandles() const override
::StatusCode StatusCode
StatusCode definition for legacy code.

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< Algorithm > >::sysStart ( )
overridevirtualinherited

Handle START transition.

We override this in order to make sure that conditions handle keys can cache a pointer to the conditions container.

◆ updateVHKA()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::updateVHKA ( Gaudi::Details::PropertyBase & )
inlineinherited

Definition at line 308 of file AthCommonDataStore.h.

308 {
309 // debug() << "updateVHKA for property " << p.name() << " " << p.toString()
310 // << " size: " << m_vhka.size() << endmsg;
311 for (auto &a : m_vhka) {
313 for (auto k : keys) {
314 k->setOwner(this);
315 }
316 }
317 }
std::vector< SG::VarHandleKeyArray * > m_vhka

Member Data Documentation

◆ m_classif

PublicToolHandle<IMCTruthClassifier> xAODTruthParticleSlimmerMET::m_classif {this, "MCTruthClassifier", "MCTruthClassifier/DFCommonTruthClassifier"}
private

Definition at line 39 of file xAODTruthParticleSlimmerMET.h.

39{this, "MCTruthClassifier", "MCTruthClassifier/DFCommonTruthClassifier"};

◆ m_detStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< Algorithm > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_evtStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< Algorithm > >::m_evtStore
privateinherited

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_extendedExtraObjects

DataObjIDColl AthAlgorithm::m_extendedExtraObjects
privateinherited

Definition at line 79 of file AthAlgorithm.h.

◆ m_varHandleArraysDeclared

bool AthCommonDataStore< AthCommonMsg< Algorithm > >::m_varHandleArraysDeclared
privateinherited

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

std::vector<SG::VarHandleKeyArray*> AthCommonDataStore< AthCommonMsg< Algorithm > >::m_vhka
privateinherited

Definition at line 398 of file AthCommonDataStore.h.

◆ m_xaodTruthEventContainerName

SG::ReadHandleKey<xAOD::TruthEventContainer> xAODTruthParticleSlimmerMET::m_xaodTruthEventContainerName {this, "xAODTruthEventContainerName", "TruthEvents"}
private

Definition at line 33 of file xAODTruthParticleSlimmerMET.h.

34{this, "xAODTruthEventContainerName", "TruthEvents"};

◆ m_xaodTruthParticleContainerNameMET

SG::WriteHandleKey<xAOD::TruthParticleContainer> xAODTruthParticleSlimmerMET::m_xaodTruthParticleContainerNameMET {this, "xAODTruthParticleContainerNameMET","TruthMET"}
private

The key for the output xAOD truth containers.

Definition at line 36 of file xAODTruthParticleSlimmerMET.h.

37{this, "xAODTruthParticleContainerNameMET","TruthMET"};

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