ATLAS Offline Software
Loading...
Searching...
No Matches
ParticleSortingTool.cxx
Go to the documentation of this file.
1
2
3/*
4 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
5*/
6
7// ParticleSortingTool.cxx
8// Implementation file for class ParticleSortingTool
9// Author: Karsten Koeneke <karsten.koeneke@cern.ch>
11
12// EventUtils includes
13#include "ParticleSortingTool.h"
14
15// EDM includes
16#include "xAODBase/IParticle.h"
31
32// Constructors
35 const std::string& name,
36 const IInterface* parent ) :
37 ::AthAlgTool ( type, name, parent )
38{
39 declareInterface< DerivationFramework::IAugmentationTool >(this);
40}
41
42
43// Destructor
47
48
49
50// Athena algtool's Hooks
53{
54 ATH_MSG_DEBUG ("Initializing " << name() << "...");
55
56 // Print out the used configuration
57 ATH_MSG_DEBUG ( " using = " << m_inCollKey );
58 ATH_MSG_DEBUG ( " using = " << m_outCollKey );
59
60 // initialize the counters
61 m_sortID = 0;
63
64 // Figure out how to sort
65 if ( m_sortVar.value() == "pt" ) { m_sortID = 1; }
66 else if ( m_sortVar.value() == "eta" ) { m_sortID = 2; }
67 else if ( m_sortVar.value() == "phi" ) { m_sortID = 3; }
68 else if ( m_sortVar.value() == "m" ) { m_sortID = 4; }
69 else if ( m_sortVar.value() == "e" ) { m_sortID = 5; }
70 else if ( m_sortVar.value() == "rapidity" ) { m_sortID = 6; }
71 else {
72 ATH_MSG_INFO("Didn't find a valid value for SortVariable=" << m_sortVar.value() << "."
73 << " Assuming it's an auxdata member");
74 m_sortID = 7;
75 }
76 if ( m_sortDescending.value() ) { m_sortID *= -1; }
77
78 return StatusCode::SUCCESS;
79}
80
81
82
83
85{
86 ATH_MSG_DEBUG ("Finalizing " << name() << "...");
87
88 return StatusCode::SUCCESS;
89}
90
91
92
93// Declare a short pre-processor macro to deal with the different container types
94#define COPY_AND_SORT_CONTAINER( CONTAINERTYPE ) \
95else if ( evtStore()->contains<CONTAINERTYPE>( m_inCollKey.value() ) ) { \
96 ATH_MSG_DEBUG("Trying to copy, sort, and record container of type "#CONTAINERTYPE ); \
97 const CONTAINERTYPE* inCont; \
98 ATH_CHECK( evtStore()->retrieve( inCont, m_inCollKey.value() ) ); \
99 CONTAINERTYPE* outCont = new CONTAINERTYPE( SG::VIEW_ELEMENTS ); \
100 *outCont = *inCont; \
101 ATH_CHECK( evtStore()->record ( outCont, m_outCollKey.value() ) ); \
102 ATH_CHECK( this->doSort(outCont) ); \
103}
104
105
106// Declare a short pre-processor macro to deal with the different container types
107#define OVERWRITE_AND_SORT_CONTAINER( CONTAINERTYPE ) \
108else if ( evtStore()->contains<CONTAINERTYPE>( m_inCollKey.value() ) ) { \
109 ATH_MSG_DEBUG("Trying to copy, sort, and overwrite container of type "#CONTAINERTYPE ); \
110 const CONTAINERTYPE* inCont; \
111 ATH_CHECK( evtStore()->retrieve( inCont, m_inCollKey.value() ) ); \
112 ConstDataVector<CONTAINERTYPE>* outCont = new ConstDataVector<CONTAINERTYPE>( SG::VIEW_ELEMENTS ); \
113 for ( const CONTAINERTYPE::base_value_type* inPart : *inCont ){ \
114 outCont->push_back(inPart); \
115 } \
116 ATH_CHECK( evtStore()->overwrite( outCont, m_inCollKey.value() ) ); \
117 ATH_CHECK( this->doSortConst<CONTAINERTYPE>(outCont) ); \
118}
119
120
121
122StatusCode ParticleSortingTool::addBranches(const EventContext&) const
123{
124 // Increase the event counter
126
127 // Simple status message at the beginning of each event execute,
128 ATH_MSG_DEBUG ( "==> addBranches " << name() << " on " << m_nEventsProcessed << ". event..." );
129
130 if ( m_outCollKey.value().empty() ) {
131 // Try to get the input container as non-const
132 ATH_MSG_DEBUG("Got an empty 'OutputCollection' property. "
133 << "Trying to retrieve a non-const version of the 'InputContainer'...");
134 xAOD::IParticleContainer* inCont = evtStore()->tryRetrieve<xAOD::IParticleContainer>( m_inCollKey.value() ); // FIXME Use Handles
135 if (inCont){ ATH_CHECK( this->doSort(inCont) ); }
136 else {
137 ATH_MSG_DEBUG("We couldn't retrieve a non-const version of the input container... try const.");
138 const xAOD::IParticleContainer* inCont2 = nullptr;
139 ATH_CHECK( evtStore()->retrieve( inCont2, m_inCollKey.value()) ); // FIXME Use Handles
140 // Now, do the copy and sorting and overwriting of all known container types
141 if (false) {
142 }
155 else {
156 ATH_MSG_ERROR("Couln't find the provided intput container in store gate for later overwriting");
157 return StatusCode::FAILURE;
158 }
159 }
160 }
161 else {
162 ATH_MSG_DEBUG("Got a non-empty 'OutputCollection' property. "
163 << "Trying to retrieve a const version of the 'InputContainer'...");
164
165 // Now, do the copy and sorting of all known container types
166 if (false) {
167 }
180 else {
181 ATH_MSG_ERROR("Couln't find the provided intput container in store gate");
182 return StatusCode::FAILURE;
183 }
184
185 }
186
187 return StatusCode::SUCCESS;
188}
189
190
191
193{
194 if ( !cont ) {
195 ATH_MSG_ERROR("No container to be sorted");
196 return StatusCode::FAILURE;
197 }
198 // Actually do the sorting, using a C++11 lambda function construct
199 // to be able to use the member function here
200 if ( std::abs(m_sortID) == 1 ) {
201 cont->sort( [this](const xAOD::IParticle* a, const xAOD::IParticle* b) {
202 return this->comparePt(a,b);
203 } );
204 }
205 else if ( std::abs(m_sortID) == 2 ) {
206 cont->sort( [this](const xAOD::IParticle* a, const xAOD::IParticle* b) {
207 return this->compareEta(a,b);
208 } );
209 }
210 else if ( std::abs(m_sortID) == 3 ) {
211 cont->sort( [this](const xAOD::IParticle* a, const xAOD::IParticle* b) {
212 return this->comparePhi(a,b);
213 } );
214 }
215 else if ( std::abs(m_sortID) == 4 ) {
216 cont->sort( [this](const xAOD::IParticle* a, const xAOD::IParticle* b) {
217 return this->compareMass(a,b);
218 } );
219 }
220 else if ( std::abs(m_sortID) == 5 ) {
221 cont->sort( [this](const xAOD::IParticle* a, const xAOD::IParticle* b) {
222 return this->compareEnergy(a,b);
223 } );
224 }
225 else if ( std::abs(m_sortID) == 6 ) {
226 cont->sort( [this](const xAOD::IParticle* a, const xAOD::IParticle* b) {
227 return this->compareRapidity(a,b);
228 } );
229 }
230 else if ( std::abs(m_sortID) == 7 ) {
231 cont->sort( [this](const xAOD::IParticle* a, const xAOD::IParticle* b) {
232 return this->compareAuxData(a,b);
233 } );
234 }
235
236 return StatusCode::SUCCESS;
237}
238
239
241 const xAOD::IParticle* partB ) const
242{
243 const double a = partA->pt();
244 const double b = partB->pt();
245 return this->compareDouble(a,b);
246}
247
248
250 const xAOD::IParticle* partB ) const
251{
252 const double a = partA->eta();
253 const double b = partB->eta();
254 return this->compareDouble(a,b);
255}
256
257
259 const xAOD::IParticle* partB ) const
260{
261 const double a = partA->phi();
262 const double b = partB->phi();
263 return this->compareDouble(a,b);
264}
265
266
268 const xAOD::IParticle* partB ) const
269{
270 const double a = partA->m();
271 const double b = partB->m();
272 return this->compareDouble(a,b);
273}
274
275
277 const xAOD::IParticle* partB ) const
278{
279 const double a = partA->e();
280 const double b = partB->e();
281 return this->compareDouble(a,b);
282}
283
284
286 const xAOD::IParticle* partB ) const
287{
288 const double a = partA->rapidity();
289 const double b = partB->rapidity();
290 return this->compareDouble(a,b);
291}
292
294 const xAOD::IParticle* partB ) const
295{
296 SG::ConstAccessor<float> acc( this->m_sortVar.value() );
297 return this->compareDouble(acc(*partA),acc(*partB));
298}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
DataVector adapter that acts like it holds const pointers.
static Double_t a
#define OVERWRITE_AND_SORT_CONTAINER(CONTAINERTYPE)
#define COPY_AND_SORT_CONTAINER(CONTAINERTYPE)
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
ServiceHandle< StoreGateSvc > & evtStore()
void sort()
Sort the container.
StatusCode doSort(xAOD::IParticleContainer *cont) const
Helper method that implements the call to the right sort function.
StringProperty m_inCollKey
Input container name.
bool compareEnergy(const xAOD::IParticle *partA, const xAOD::IParticle *partB) const
The method to compare the particle's energy.
virtual ~ParticleSortingTool()
Destructor:
bool compareAuxData(const xAOD::IParticle *partA, const xAOD::IParticle *partB) const
The method to compare an auxdata member of the particle.
BooleanProperty m_sortDescending
Define if the container should be sorted in a descending order (default=true)
int m_sortID
Internal identifier for the type of sorting.
std::atomic< unsigned long > m_nEventsProcessed
Internal event counter.
StringProperty m_sortVar
Define by what parameter to sort (default: 'pt')
bool comparePt(const xAOD::IParticle *partA, const xAOD::IParticle *partB) const
The method to compare the particle's pt.
ParticleSortingTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
virtual StatusCode addBranches(const EventContext &ctx) const final override
Implement the method from the ISkimmingTool interface.
bool compareRapidity(const xAOD::IParticle *partA, const xAOD::IParticle *partB) const
The method to compare the particle's rapidity.
bool compareEta(const xAOD::IParticle *partA, const xAOD::IParticle *partB) const
The method to compare the particle's eta.
StringProperty m_outCollKey
The name of the output container (with SG::VIEW_ELEMENTS) with the sorted copy of input objects.
virtual StatusCode initialize() override
Athena algtool's initialize.
bool compareDouble(double a, double b) const
Method to compare two doubles.
bool compareMass(const xAOD::IParticle *partA, const xAOD::IParticle *partB) const
The method to compare the particle's mass.
virtual StatusCode finalize() override
Athena algtool's finalize.
bool comparePhi(const xAOD::IParticle *partA, const xAOD::IParticle *partB) const
The method to compare the particle's phi.
Helper class to provide constant type-safe access to aux data.
Class providing the definition of the 4-vector interface.
virtual double eta() const =0
The pseudorapidity ( ) of the particle.
virtual double pt() const =0
The transverse momentum ( ) of the particle.
virtual double m() const =0
The invariant mass of the particle.
virtual double rapidity() const =0
The true rapidity (y) of the particle.
virtual double phi() const =0
The azimuthal angle ( ) of the particle.
virtual double e() const =0
The total energy of the particle.
PhotonContainer_v1 PhotonContainer
Definition of the current "photon container version".
ElectronContainer_v1 ElectronContainer
Definition of the current "electron container version".
PFOContainer_v1 PFOContainer
Definition of the current "pfo container version".
NeutralParticleContainer_v1 NeutralParticleContainer
Definition of the current "NeutralParticle container version".
ParticleContainer_v1 ParticleContainer
Define the latest version of the particle class.
CompositeParticleContainer_v1 CompositeParticleContainer
Define the latest version of the CompositeParticle class.
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".
CaloClusterContainer_v1 CaloClusterContainer
Define the latest version of the calorimeter cluster container.
JetContainer_v1 JetContainer
Definition of the current "jet container version".
TauJetContainer_v3 TauJetContainer
Definition of the current "taujet container version".
MuonContainer_v1 MuonContainer
Definition of the current "Muon container version".
TruthParticleContainer_v1 TruthParticleContainer
Declare the latest version of the truth particle container.
DataVector< IParticle > IParticleContainer
Simple convenience declaration of IParticleContainer.