ATLAS Offline Software
Generators
GeneratorFilters
src
xAODChargedTracksFilter.cxx
Go to the documentation of this file.
1
/*
2
Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3
*/
4
5
#include "
GeneratorFilters/xAODChargedTracksFilter.h
"
6
#include "
TruthUtils/HepMCHelpers.h
"
7
8
9
xAODChargedTracksFilter::xAODChargedTracksFilter
(
const
std::string&
name
, ISvcLocator* pSvcLocator)
10
:
GenFilter
(
name
, pSvcLocator)
11
{
13
declareProperty
(
"Ptcut"
,
m_Ptmin
= 50.0);
14
declareProperty
(
"Etacut"
,
m_EtaRange
= 2.5);
15
declareProperty
(
"NTracks"
,
m_NTracks
= 40);
16
}
17
18
19
StatusCode
xAODChargedTracksFilter::filterEvent
() {
20
21
// Retrieve TruthGen container from xAOD Gen slimmer, contains all particles witout barcode_zero and
22
// duplicated barcode ones
23
const
xAOD::TruthParticleContainer
* xTruthParticleContainer;
24
if
(
evtStore
()->
retrieve
(xTruthParticleContainer,
"TruthGen"
).isFailure()) {
25
ATH_MSG_ERROR
(
"No TruthParticle collection with name "
<<
"TruthGen"
<<
" found in StoreGate!"
);
26
return
StatusCode::FAILURE;
27
}
28
29
30
int
nChargedTracks
= 0;
31
unsigned
int
nPart = xTruthParticleContainer->
size
();
32
for
(
unsigned
int
iPart = 0; iPart < nPart; ++iPart) {
33
const
xAOD::TruthParticle
*
part
= (*xTruthParticleContainer)[iPart];
34
// We only care about stable particles
35
if
(!
part
->isGenStable())
continue
;
36
37
// Particle's charge
38
int
pID =
part
->pdgId();
39
double
pCharge =
MC::charge
(pID);
40
41
// Count tracks in specified acceptance
42
const
double
pT
=
part
->pt();
43
const
double
eta
=
part
->eta();
44
if
(
pT
>=
m_Ptmin
&& std::abs(
eta
) <=
m_EtaRange
&& pCharge != 0) {
45
ATH_MSG_DEBUG
(
"Found particle, "
<<
46
" pT = "
<<
pT
<<
47
" eta = "
<<
eta
<<
48
" pdg id = "
<< pID <<
49
" charge = "
<< pCharge <<
" in acceptance"
);
50
nChargedTracks
+= 1;
51
}
52
53
}
//end loop on particles
54
55
// Summarise event
56
ATH_MSG_DEBUG
(
"# of tracks "
<<
nChargedTracks
<<
57
" with pT >= "
<<
m_Ptmin
<<
58
" |eta| < "
<<
m_EtaRange
<<
59
" minNTracks = "
<<
m_NTracks
);
60
61
// Record passed status
62
setFilterPassed(
nChargedTracks
>
m_NTracks
);
63
return
StatusCode::SUCCESS;
64
}
LArG4FSStartPointFilter.part
part
Definition:
LArG4FSStartPointFilter.py:21
python.PyKernel.retrieve
def retrieve(aClass, aKey=None)
Definition:
PyKernel.py:110
CalculateHighPtTerm.pT
pT
Definition:
ICHEP2016/CalculateHighPtTerm.py:57
eta
Scalar eta() const
pseudorapidity method
Definition:
AmgMatrixBasePlugin.h:79
AthCommonDataStore< AthCommonMsg< Algorithm > >::declareProperty
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T > &t)
Definition:
AthCommonDataStore.h:145
xAODChargedTracksFilter::m_NTracks
double m_NTracks
Definition:
xAODChargedTracksFilter.h:32
xAODChargedTracksFilter::xAODChargedTracksFilter
xAODChargedTracksFilter(const std::string &name, ISvcLocator *pSvcLocator)
Definition:
xAODChargedTracksFilter.cxx:9
AthCommonDataStore< AthCommonMsg< Algorithm > >::evtStore
ServiceHandle< StoreGateSvc > & evtStore()
The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
Definition:
AthCommonDataStore.h:85
GenFilter
Base class for event generator filtering modules.
Definition:
GenFilter.h:30
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition:
AthMsgStreamMacros.h:33
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition:
PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition:
AthMsgStreamMacros.h:29
xAOD::TruthParticle_v1
Class describing a truth particle in the MC record.
Definition:
TruthParticle_v1.h:41
xAOD::TauJetParameters::nChargedTracks
@ nChargedTracks
Definition:
TauDefs.h:322
DataVector
Derived DataVector<T>.
Definition:
DataVector.h:581
xAODChargedTracksFilter::m_Ptmin
double m_Ptmin
Definition:
xAODChargedTracksFilter.h:26
name
std::string name
Definition:
Control/AthContainers/Root/debug.cxx:195
charge
double charge(const T &p)
Definition:
AtlasPID.h:494
xAODChargedTracksFilter.h
xAODChargedTracksFilter::m_EtaRange
double m_EtaRange
Definition:
xAODChargedTracksFilter.h:29
DataVector::size
size_type size() const noexcept
Returns the number of elements in the collection.
HepMCHelpers.h
xAODChargedTracksFilter::filterEvent
StatusCode filterEvent()
Definition:
xAODChargedTracksFilter.cxx:19
Generated on Fri Jul 5 2024 21:33:00 for ATLAS Offline Software by
1.8.18