ATLAS Offline Software
Loading...
Searching...
No Matches
DerivationFramework::GenFilterTool Class Reference

#include <GenFilterTool.h>

Inheritance diagram for DerivationFramework::GenFilterTool:
Collaboration diagram for DerivationFramework::GenFilterTool:

Public Member Functions

virtual StatusCode initialize () override final
virtual StatusCode addBranches (const EventContext &ctx) const override final

Private Member Functions

StatusCode getGenFiltVars (const EventContext &ctx, float &genFiltHT, float &genFiltHTinclNu, float &genFiltMET, float &genFiltPTZ, float &genFiltFatJ) const
bool isPrompt (const xAOD::TruthParticle *tp) const

Private Attributes

SG::ReadHandleKey< xAOD::EventInfom_eventInfoKey {this,"EventInfoName" , "EventInfo"}
SG::ReadHandleKey< xAOD::TruthParticleContainerm_mcKey {this, "MCCollectionName", "TruthParticles"}
SG::ReadHandleKey< xAOD::JetContainerm_truthJetsKey {this, "TruthJetCollectionName", "AntiKt4TruthWZJets"}
SG::ReadHandleKey< xAOD::JetContainerm_truthFatJetsKey {this, "TruthFatJetCollectionName", "AntiKt10TruthJets"}
SG::ReadDecorHandleKey< xAOD::TruthParticleContainerm_mcReadDecor {this, "TruthClassKey",m_mcKey,"classifierParticleOrigin"}
SG::WriteDecorHandleKey< xAOD::EventInfom_dec_genFiltHTKey { this, "GenFiltHTKey", m_eventInfoKey, "GenFiltHT", "GenFiltHT EventInfo decoration name" }
SG::WriteDecorHandleKey< xAOD::EventInfom_dec_genFiltHTinclNuKey { this, "GenFiltHTinclNuKey", m_eventInfoKey, "GenFiltHTinclNu", "GenFiltHTinclNu EventInfo decoration name" }
SG::WriteDecorHandleKey< xAOD::EventInfom_dec_genFiltMETKey { this, "GenFiltMET", m_eventInfoKey, "GenFiltMET", "GenFiltMET EventInfo decoration name" }
SG::WriteDecorHandleKey< xAOD::EventInfom_dec_genFiltPTZKey { this, "GenFiltPTZ", m_eventInfoKey, "GenFiltPTZ", "GenFiltPTZ EventInfo decoration name" }
SG::WriteDecorHandleKey< xAOD::EventInfom_dec_genFiltFatJKey { this, "GenFiltFatJ", m_eventInfoKey, "GenFiltFatJ", "GenFiltFatJ EventInfo decoration name" }
Gaudi::Property< float > m_MinJetPt {this, "MinJetPt", 35.* Gaudi::Units::GeV}
 Min pT for the truth jets.
Gaudi::Property< float > m_MaxJetEta {this, "MaxJetEta", 2.5}
 Max eta for the truth jets.
Gaudi::Property< float > m_MinLepPt {this,"MinLeptonPt", 25.*Gaudi::Units::GeV}
 Min pT for the truth leptons.
Gaudi::Property< float > m_MaxLepEta {this, "MaxLeptonEta", 2.5}
 Max eta for the truth leptons.
PublicToolHandle< IMCTruthClassifierm_classif {this, "TruthClassifier", ""}

Detailed Description

Definition at line 34 of file GenFilterTool.h.

Member Function Documentation

◆ addBranches()

StatusCode DerivationFramework::GenFilterTool::addBranches ( const EventContext & ctx) const
finaloverridevirtual

Definition at line 76 of file GenFilterTool.cxx.

76 {
77 ATH_MSG_VERBOSE("GenFilterTool::addBranches(const EventContext& ctx)");
78 SG::ReadHandle<xAOD::EventInfo> eventInfo{m_eventInfoKey, ctx};
79 if (!eventInfo.isValid()) {
80 ATH_MSG_ERROR("could not retrieve event info " <<m_eventInfoKey.fullKey());
81 return StatusCode::FAILURE;
82 }
83
84
85 float genFiltHT{0.f}, genFiltHTinclNu{0.f}, genFiltMET{0.f}, genFiltPTZ{0.f}, genFiltFatJ{0.f};
86 ATH_CHECK( getGenFiltVars(ctx, genFiltHT, genFiltHTinclNu, genFiltMET, genFiltPTZ, genFiltFatJ) );
87
88 ATH_MSG_DEBUG("Computed generator filter quantities: HT " << genFiltHT/1e3 << ", HTinclNu " << genFiltHTinclNu/1e3 << ", MET " << genFiltMET/1e3 << ", PTZ " << genFiltPTZ/1e3 << ", FatJ " << genFiltFatJ/1e3 );
89
90 SG::makeHandle<float>(m_dec_genFiltHTKey, ctx)(*eventInfo) = genFiltHT;
91 SG::makeHandle<float>(m_dec_genFiltHTinclNuKey, ctx)(*eventInfo) = genFiltHTinclNu;
92 SG::makeHandle<float>(m_dec_genFiltMETKey, ctx)(*eventInfo) = genFiltMET;
93 SG::makeHandle<float>(m_dec_genFiltPTZKey, ctx)(*eventInfo) = genFiltPTZ;
94 SG::makeHandle<float>(m_dec_genFiltFatJKey, ctx)(*eventInfo) = genFiltFatJ;
95
96 return StatusCode::SUCCESS;
97 }
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_DEBUG(x)
SG::WriteDecorHandleKey< xAOD::EventInfo > m_dec_genFiltHTinclNuKey
StatusCode getGenFiltVars(const EventContext &ctx, float &genFiltHT, float &genFiltHTinclNu, float &genFiltMET, float &genFiltPTZ, float &genFiltFatJ) const
SG::WriteDecorHandleKey< xAOD::EventInfo > m_dec_genFiltPTZKey
SG::WriteDecorHandleKey< xAOD::EventInfo > m_dec_genFiltMETKey
SG::WriteDecorHandleKey< xAOD::EventInfo > m_dec_genFiltHTKey
SG::WriteDecorHandleKey< xAOD::EventInfo > m_dec_genFiltFatJKey
SG::ReadHandleKey< xAOD::EventInfo > m_eventInfoKey
virtual bool isValid() override final
Can the handle be successfully dereferenced?
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())

◆ getGenFiltVars()

StatusCode DerivationFramework::GenFilterTool::getGenFiltVars ( const EventContext & ctx,
float & genFiltHT,
float & genFiltHTinclNu,
float & genFiltMET,
float & genFiltPTZ,
float & genFiltFatJ ) const
private

Definition at line 99 of file GenFilterTool.cxx.

99 {
100 // Get jet container out
101
102 SG::ReadDecorHandle<xAOD::TruthParticleContainer, unsigned int> mcParticleOrigin{m_mcReadDecor, ctx} ;
103 if (!mcParticleOrigin.isValid()) {
104 ATH_MSG_ERROR("WARNING could not retrieve TruthParticleContainer " <<m_mcKey.fullKey());
105 return StatusCode::FAILURE;
106 }
107
108 SG::ReadHandle<xAOD::JetContainer> truthjets{m_truthJetsKey, ctx};
109 if (!truthjets.isValid()){
110 ATH_MSG_ERROR( "No xAOD::JetContainer found in StoreGate with key " << m_truthJetsKey );
111 return StatusCode::FAILURE;
112 }
113
114 // Get HT
115 genFiltHT = 0.;
116 genFiltHTinclNu = 0.; // this HT definition includes neutrinos from W/Z/tau
117 for (const auto *const tj : *truthjets) {
118 if ( tj->pt()>m_MinJetPt && std::abs(tj->eta())<m_MaxJetEta ) {
119 ATH_MSG_VERBOSE("Adding truth jet with pt " << tj->pt()
120 << ", eta " << tj->eta()
121 << ", phi " << tj->phi()
122 << ", nconst = " << tj->numConstituents());
123 genFiltHT += tj->pt();
124 genFiltHTinclNu += tj->pt();
125 }
126 }
127
128 // Get MET and add leptons to HT
129 float MEx(0.), MEy(0.);
130 for (const auto *const tp : *mcParticleOrigin){
131 int pdgid = tp->pdgId();
132 if (HepMC::is_simulation_particle(tp)) continue; // Particle is from G4
133 if (MC::isZeroEnergyPhoton(tp)) continue; // Work around for an old generator bug
134 if ( !MC::isStable(tp)) continue; // Stable!
135
136 if ((MC::isElectron(pdgid) || MC::isMuon(pdgid)) && tp->pt()>m_MinLepPt && std::abs(tp->eta())<m_MaxLepEta) {
137 if( isPrompt(tp) ) {
138 ATH_MSG_VERBOSE("Adding prompt lepton " << tp);
139 genFiltHT += tp->pt();
140 genFiltHTinclNu += tp->pt();
141 }
142 }
143
144 // include neutrinos from W/Z/Tau in one of the HT definitions
145 // this corresponds to a subset of HT-sliced samples where the HT filter
146 // was configured to include these particles
147 unsigned int orig = mcParticleOrigin (*tp);
148 if (tp->isNeutrino() && isFromWZTau(orig)) {
149 ATH_MSG_VERBOSE("Adding neutrino from W/Z/Tau with pt " << tp->pt()
150 << ", eta " << tp->eta()
151 << ", phi " << tp->phi()
152 << ", status " << tp->status()
153 << ", pdgId " << pdgid);
154 genFiltHTinclNu += tp->pt();
155 }
156
157 if (MC::isSpecialNonInteracting(tp) && isPrompt(tp) ) {
158 ATH_MSG_VERBOSE("Found prompt nonInteracting particle with pt " << tp->pt()
159 << ", eta " << tp->eta()
160 << ", phi " << tp->phi()
161 << ", status " << tp->status()
162 << ", pdgId " << pdgid);
163 MEx += tp->px();
164 MEy += tp->py();
165 }
166 }
167 genFiltMET = std::sqrt(MEx*MEx+MEy*MEy);
168
169 // Get PTZ
170 float PtZ(.0);
171 float MinPt_PTZ(5000.), MaxEta_PTZ(5.0), MinMass_PTZ(20000.), MaxMass_PTZ(14000000.);
172 bool AllowElecMu_PTZ = false;
173 bool AllowSameCharge_PTZ = false;
174 for (const xAOD::TruthParticle* pitr1 : *mcParticleOrigin){
175 const int pdgId1 = pitr1->pdgId();
176 if (HepMC::is_simulation_particle(pitr1)) continue;
177 if (!MC::isStable(pitr1)) continue;
178 // Pick electrons or muons with Pt > MinPt_PTZ and |eta| < m_maxEta
179 if (MC::isElectron(pdgId1) || MC::isMuon(pdgId1)) {
180 if (pitr1->pt() >= MinPt_PTZ && std::abs(pitr1->eta()) <= MaxEta_PTZ){
181 for (const xAOD::TruthParticle* pitr2 : *mcParticleOrigin){
182 if (pitr2==pitr1) continue;
183 if (HepMC::is_simulation_particle(pitr2)) continue;
184 if (!MC::isStable(pitr2)) continue;
185 const int pdgId2 = pitr2->pdgId();
186 // Pick electrons or muons with Pt > MinPt_PTZ and |eta| < MaxEta_PTZ
187 // If AllowSameCharge_PTZ is not true only pick those with opposite charge to the first particle
188 // If AllowElecMu_PTZ is true allow also Z -> emu compinations (with charge requirements as above)
189 if (!(MC::isElectron(pdgId2) || MC::isMuon(pdgId2))) continue;
190 if ( ( AllowSameCharge_PTZ && ( AllowElecMu_PTZ || std::abs(pdgId2) == std::abs(pdgId1) ) ) ||
191 ( !AllowSameCharge_PTZ && ( pdgId2 == -1*pdgId1 || (AllowElecMu_PTZ && MC::charge3(pdgId1) == MC::charge3(pdgId2) ) ) )
192 ) {
193 if (pitr2->pt() >= MinPt_PTZ && std::abs(pitr2->eta()) <= MaxEta_PTZ){
194 double invMass = (pitr1->p4()+pitr2->p4()).M();
195 double dilepPt = (pitr1->p4()+pitr2->p4()).Pt();
196 // Only consider pair that fall in the mass window
197 if (MinMass_PTZ < invMass && invMass < MaxMass_PTZ) {
198 if (dilepPt > PtZ) PtZ = dilepPt;
199 }
200 }
201 }
202 }
203 }
204 }
205 }
206 genFiltPTZ = PtZ;
207
208 //Get FatJ
209 // Get correct jet container
210 SG::ReadHandle<xAOD::JetContainer> truthjets10{m_truthFatJetsKey, ctx} ;
211 if ( !truthjets10.isValid()){
212 ATH_MSG_ERROR( "No xAOD::JetContainer found in StoreGate with key "<<m_truthFatJetsKey.fullKey() );
213 return StatusCode::FAILURE;
214 }
215 genFiltFatJ=0.;
216 for (const auto *const j : *truthjets10) {
217 if (j->pt()>genFiltFatJ) genFiltFatJ=j->pt();
218 }
219
220
221 return StatusCode::SUCCESS;
222 }
Gaudi::Property< float > m_MinLepPt
Min pT for the truth leptons.
SG::ReadHandleKey< xAOD::JetContainer > m_truthFatJetsKey
SG::ReadDecorHandleKey< xAOD::TruthParticleContainer > m_mcReadDecor
SG::ReadHandleKey< xAOD::TruthParticleContainer > m_mcKey
Gaudi::Property< float > m_MaxJetEta
Max eta for the truth jets.
SG::ReadHandleKey< xAOD::JetContainer > m_truthJetsKey
Gaudi::Property< float > m_MinJetPt
Min pT for the truth jets.
Gaudi::Property< float > m_MaxLepEta
Max eta for the truth leptons.
bool isPrompt(const xAOD::TruthParticle *tp) const
static bool isFromWZTau(unsigned int orig)
bool is_simulation_particle(const T &p)
Method to establish if a particle (or barcode) was created during the simulation (TODO update to be s...
bool isZeroEnergyPhoton(const T &p)
Identify a photon with zero energy. Probably a workaround for a generator bug.
int charge3(const T &p)
bool isSpecialNonInteracting(const T &p)
Identify a special non-interacting particles.
bool isElectron(const T &p)
bool isStable(const T &p)
Identify if the particle is stable, i.e. has not decayed.
bool isMuon(const T &p)
double invMass(const I4Momentum &pA, const I4Momentum &pB)
invariant mass from two I4momentum references
Definition P4Helpers.h:252
TruthParticle_v1 TruthParticle
Typedef to implementation.

◆ initialize()

StatusCode DerivationFramework::GenFilterTool::initialize ( )
finaloverridevirtual

Definition at line 63 of file GenFilterTool.cxx.

63 {
64 ATH_CHECK(m_eventInfoKey.initialize());
65 ATH_CHECK(m_mcKey.initialize());
66 ATH_CHECK(m_truthJetsKey.initialize());
67 ATH_CHECK(m_truthFatJetsKey.initialize());
68 ATH_CHECK(m_dec_genFiltHTKey.initialize());
70 ATH_CHECK(m_dec_genFiltMETKey.initialize());
71 ATH_CHECK(m_dec_genFiltPTZKey.initialize());
72 ATH_CHECK(m_dec_genFiltFatJKey.initialize());
73 ATH_CHECK(m_mcReadDecor.initialize());
74 return StatusCode::SUCCESS;
75 }

◆ isPrompt()

bool DerivationFramework::GenFilterTool::isPrompt ( const xAOD::TruthParticle * tp) const
private

Definition at line 33 of file GenFilterTool.cxx.

34 {
35 ParticleOrigin orig = m_classif->particleTruthClassifier( tp ).second;
36 ATH_MSG_VERBOSE("Particle has origin " << orig);
37
38 switch(orig) {
39 case NonDefined:
40 case PhotonConv:
41 case DalitzDec:
42 case ElMagProc:
43 case Mu:
44 case LightMeson:
45 case StrangeMeson:
46 case CharmedMeson:
47 case BottomMeson:
48 case CCbarMeson:
49 case JPsi:
50 case BBbarMeson:
51 case LightBaryon:
52 case StrangeBaryon:
53 case CharmedBaryon:
54 case BottomBaryon:
55 case PionDecay:
56 case KaonDecay:
57 return false;
58 default:
59 break;
60 }
61 return true;
62 }
PublicToolHandle< IMCTruthClassifier > m_classif

Member Data Documentation

◆ m_classif

PublicToolHandle<IMCTruthClassifier> DerivationFramework::GenFilterTool::m_classif {this, "TruthClassifier", ""}
private

Definition at line 72 of file GenFilterTool.h.

72{this, "TruthClassifier", ""};

◆ m_dec_genFiltFatJKey

SG::WriteDecorHandleKey<xAOD::EventInfo> DerivationFramework::GenFilterTool::m_dec_genFiltFatJKey { this, "GenFiltFatJ", m_eventInfoKey, "GenFiltFatJ", "GenFiltFatJ EventInfo decoration name" }
private

Definition at line 64 of file GenFilterTool.h.

65{ this, "GenFiltFatJ", m_eventInfoKey, "GenFiltFatJ", "GenFiltFatJ EventInfo decoration name" };

◆ m_dec_genFiltHTinclNuKey

SG::WriteDecorHandleKey<xAOD::EventInfo> DerivationFramework::GenFilterTool::m_dec_genFiltHTinclNuKey { this, "GenFiltHTinclNuKey", m_eventInfoKey, "GenFiltHTinclNu", "GenFiltHTinclNu EventInfo decoration name" }
private

Definition at line 58 of file GenFilterTool.h.

59{ this, "GenFiltHTinclNuKey", m_eventInfoKey, "GenFiltHTinclNu", "GenFiltHTinclNu EventInfo decoration name" };

◆ m_dec_genFiltHTKey

SG::WriteDecorHandleKey<xAOD::EventInfo> DerivationFramework::GenFilterTool::m_dec_genFiltHTKey { this, "GenFiltHTKey", m_eventInfoKey, "GenFiltHT", "GenFiltHT EventInfo decoration name" }
private

Definition at line 56 of file GenFilterTool.h.

57{ this, "GenFiltHTKey", m_eventInfoKey, "GenFiltHT", "GenFiltHT EventInfo decoration name" };

◆ m_dec_genFiltMETKey

SG::WriteDecorHandleKey<xAOD::EventInfo> DerivationFramework::GenFilterTool::m_dec_genFiltMETKey { this, "GenFiltMET", m_eventInfoKey, "GenFiltMET", "GenFiltMET EventInfo decoration name" }
private

Definition at line 60 of file GenFilterTool.h.

61{ this, "GenFiltMET", m_eventInfoKey, "GenFiltMET", "GenFiltMET EventInfo decoration name" };

◆ m_dec_genFiltPTZKey

SG::WriteDecorHandleKey<xAOD::EventInfo> DerivationFramework::GenFilterTool::m_dec_genFiltPTZKey { this, "GenFiltPTZ", m_eventInfoKey, "GenFiltPTZ", "GenFiltPTZ EventInfo decoration name" }
private

Definition at line 62 of file GenFilterTool.h.

63{ this, "GenFiltPTZ", m_eventInfoKey, "GenFiltPTZ", "GenFiltPTZ EventInfo decoration name" };

◆ m_eventInfoKey

SG::ReadHandleKey<xAOD::EventInfo> DerivationFramework::GenFilterTool::m_eventInfoKey {this,"EventInfoName" , "EventInfo"}
private

Definition at line 49 of file GenFilterTool.h.

49{this,"EventInfoName" , "EventInfo"};

◆ m_MaxJetEta

Gaudi::Property<float> DerivationFramework::GenFilterTool::m_MaxJetEta {this, "MaxJetEta", 2.5}
private

Max eta for the truth jets.

Definition at line 68 of file GenFilterTool.h.

68{this, "MaxJetEta", 2.5};

◆ m_MaxLepEta

Gaudi::Property<float> DerivationFramework::GenFilterTool::m_MaxLepEta {this, "MaxLeptonEta", 2.5}
private

Max eta for the truth leptons.

Definition at line 70 of file GenFilterTool.h.

70{this, "MaxLeptonEta", 2.5};

◆ m_mcKey

SG::ReadHandleKey<xAOD::TruthParticleContainer> DerivationFramework::GenFilterTool::m_mcKey {this, "MCCollectionName", "TruthParticles"}
private

Definition at line 50 of file GenFilterTool.h.

50{this, "MCCollectionName", "TruthParticles"};

◆ m_mcReadDecor

SG::ReadDecorHandleKey<xAOD::TruthParticleContainer> DerivationFramework::GenFilterTool::m_mcReadDecor {this, "TruthClassKey",m_mcKey,"classifierParticleOrigin"}
private

Definition at line 54 of file GenFilterTool.h.

54{this, "TruthClassKey",m_mcKey,"classifierParticleOrigin"};

◆ m_MinJetPt

Gaudi::Property<float> DerivationFramework::GenFilterTool::m_MinJetPt {this, "MinJetPt", 35.* Gaudi::Units::GeV}
private

Min pT for the truth jets.

Definition at line 67 of file GenFilterTool.h.

67{this, "MinJetPt", 35.* Gaudi::Units::GeV};

◆ m_MinLepPt

Gaudi::Property<float> DerivationFramework::GenFilterTool::m_MinLepPt {this,"MinLeptonPt", 25.*Gaudi::Units::GeV}
private

Min pT for the truth leptons.

Definition at line 69 of file GenFilterTool.h.

69{this,"MinLeptonPt", 25.*Gaudi::Units::GeV};

◆ m_truthFatJetsKey

SG::ReadHandleKey<xAOD::JetContainer> DerivationFramework::GenFilterTool::m_truthFatJetsKey {this, "TruthFatJetCollectionName", "AntiKt10TruthJets"}
private

Definition at line 52 of file GenFilterTool.h.

52{this, "TruthFatJetCollectionName", "AntiKt10TruthJets"};

◆ m_truthJetsKey

SG::ReadHandleKey<xAOD::JetContainer> DerivationFramework::GenFilterTool::m_truthJetsKey {this, "TruthJetCollectionName", "AntiKt4TruthWZJets"}
private

Definition at line 51 of file GenFilterTool.h.

51{this, "TruthJetCollectionName", "AntiKt4TruthWZJets"};

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