ATLAS Offline Software
Loading...
Searching...
No Matches
HardTruthThinning.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// HardTruthThinning.cxx, (c) ATLAS Detector software
7// Author: Frank Paige
8// Based on MenuTruthThinning and should perhaps be merged with it.
9// Intended for use with CompactHardTruth.
10// Preserves graph only for decays of selected particles.
11//
12// No treatment of Geant particles (yet).
13//
15
18#include "xAODBase/IParticle.h"
20#include "GaudiKernel/SystemOfUnits.h"
22#include "GaudiKernel/ThreadLocalContext.h"
23
25#include <vector>
26#include <string>
27#include <format>
28
29using Gaudi::Units::GeV;
30
32// Constructor
34
35// Required parameters:
36// EventInfo, TruthParticles, TruthVertices, HardParticles
37// Parameters that turn on features:
38// JetName If not empty, save constituents with cuts
39// KeepIds If pdgId list not empty, save particles and decay chains
40// IsolRadius If positive, save stable particles in isolation cones
41
43 const std::string& t, const std::string& n, const IInterface* p ) :
44 base_class(t,n,p),
45 m_jetPtCut(0.),
46 m_jetEtaCut(5.0),
49 m_isolR(-1),
50 m_isolPtCut(0),
51 m_maxCount(0)
52{
53 declareProperty("JetPtCut", m_jetPtCut,
54 "truth jet minumum pt");
55
56 declareProperty("JetEtaCut", m_jetEtaCut,
57 "truth jet maximum abs(eta)");
58
59 declareProperty("JetConstPtCut", m_jetConstPtCut,
60 "jet constituent minimum pt");
61
62 declareProperty("JetPhotonPtCut", m_jetPhotonPtCut,
63 "jet constituent photon minimum pt");
64
65 declareProperty("KeepIds", m_keepIds,
66 "list of abs(pdgID) to keep");
67
68 declareProperty("IsolRadius", m_isolR,
69 "isolation radius for leptons and photons");
70
71 declareProperty("IsolPtCut", m_isolPtCut,
72 "isolation particle minimum pt");
73
74 declareProperty("MaxCount",
76 "maximum number of events to print");
77}
78
79
81// Destructor
83
86
87
89// Athena initialize and finalize
91
93{
94 ATH_CHECK( m_eventInfoKey.initialize() );
95 ATH_CHECK( m_hardParticleKey.initialize() );
99
100 m_evtCount = -1;
101 m_errCount = 0;
102
103 // Check for unset photon cut
105
106 return StatusCode::SUCCESS;
107}
108
109
111{
112 ATH_MSG_INFO("finalize() ...");
113 if( m_errCount > 0 ){
114 ATH_MSG_WARNING("TruthHard/TruthEvent pdgId mismatches " <<m_errCount);
115 } else {
116 ATH_MSG_INFO("No TruthHard/TruthEvent pdgId mismatches");
117 }
118
119 return StatusCode::SUCCESS;
120}
121
122
124// doThinning
126
128{
129 const EventContext& ctx = Gaudi::Hive::currentContext();
130
131 ++m_evtCount;
132 bool doPrint = m_evtCount < m_maxCount;
133
134 // Retrieve truth particles and vertices
136 (m_truthParticleName, ctx);
138 (m_truthVertexName, ctx);
139
141 if (!inHardParts.isValid()) {
142 ATH_MSG_ERROR("Cannot retrieve TruthParticleContainer " << m_hardParticleKey);
143 return StatusCode::FAILURE;
144 }
145
146 // Hard particles are different objects but have matching unique IDs.
147 // Find hard particles with status==1.
148 // Save 4vectors of leptons and photons for matching if requested.
149 // Do we need a photon pt cut for soft photons from quarks?
150
151 std::vector<const xAOD::TruthParticle*> hardPart;
152 std::vector<TLorentzVector> pLepGam;
153
154 for (const auto* pItr : *inHardParts) {
155 if( MC::isStable(pItr) ){
156 hardPart.push_back(pItr);
157 if( m_isolR > 0 ){
158 int ida = pItr->absPdgId();
159 if( MC::isElectron(ida) || MC::isMuon(ida) || MC::isTau(ida) || MC::isPhoton(ida) ){
160 pLepGam.push_back( pItr->p4() );
161 }
162 }
163 }
164 }
165
166 // Print full input event
167 long long int evtNum{};
168 if( doPrint ){
170 if(!evt.isValid()) {
171 ATH_MSG_ERROR("Failed to retrieve EventInfo");
172 return StatusCode::FAILURE;
173 }
174 evtNum = evt->eventNumber();
175
176 printxAODTruth(evtNum, inTruthParts.cptr());
177 }
178
179
180 // Set up masks for particles/vertices
181 // Default is false for all
182
183 std::vector<bool> partMask;
184 std::vector<bool> vertMask;
185 unsigned int inPartNum = inTruthParts->size();
186 unsigned int inVertNum = inTruthVerts->size();
187 partMask.assign(inPartNum, false);
188 vertMask.assign(inVertNum, false);
189
190 ATH_MSG_DEBUG("Initial particles/vertices = " <<inPartNum <<" " <<inVertNum);
191
192
193 // Particles from hard event are stable => keep no parent/child vertices.
194 // No descendants (except Geant).
195 // Particles from keep list need descendants and their parent vertices.
196 // There could be overlap, e.g. for taus.
197 // Expect no pdgId mismatches.
198
199 std::vector<const xAOD::TruthParticle*> kids;
200
201 for (const auto* pItr : *inTruthParts) {
202 // Stable hard particle matches
203 int uid = HepMC::uniqueID(pItr);
204 bool isHard = false;
205 for(unsigned int i=0; i<hardPart.size(); ++i){
206 if( uid == HepMC::uniqueID(hardPart[i]) ){
207 isHard = true;
208 if( pItr->pdgId() != (hardPart[i])->pdgId() ){
209 ATH_MSG_WARNING("pdgID mismatch, TruthParticle uid/pdgid "
210 <<HepMC::uniqueID(pItr) <<" " <<pItr->pdgId()
211 <<" Hard uid/pdgid " <<HepMC::uniqueID(hardPart[i])
212 <<" " <<(hardPart[i])->pdgId());
213 ++m_errCount;
214 break;
215 }
216 }
217 }
218 if( isHard ){
219 if( doPrint ) ATH_MSG_DEBUG("ParticleMask isHard " <<uid);
220 partMask[ pItr->index() ] = true;
221 }
222
223 // Keep particles
224 int ida = pItr->absPdgId();
225 bool isKeep = false;
226 for(unsigned int i=0; i<m_keepIds.size(); ++i){
227 if( ida == m_keepIds[i] ){
228 isKeep = true;
229 break;
230 }
231 }
232 if( isKeep ){
233 if( doPrint ) ATH_MSG_DEBUG("ParticleMask isKeep " <<uid);
234 partMask[pItr->index()] = true;
235 int nkids = getDescendants( pItr, kids );
236 for(int i=0; i<nkids; ++i){
237 if( doPrint ) ATH_MSG_DEBUG("ParticleMask isKeep kid "
238 <<uid <<" " <<HepMC::uniqueID(kids[i]));
239 partMask[ (kids[i])->index() ] = true;
240 const xAOD::TruthVertex* v = (kids[i])->prodVtx();
241 if( v ) vertMask[ v->index() ] = true;
242 }
243 }
244
245 // Delta(R) matches to hard truth leptons/photons
246 if( !pLepGam.empty() && pItr->pt() > m_isolPtCut ){
247 TLorentzVector pp4 = pItr->p4();
248 for(unsigned int lep=0; lep<pLepGam.size(); ++lep){
249 double r = pp4.DeltaR( pLepGam[lep] );
250 if( r < m_isolR ){
251 if( doPrint ) ATH_MSG_DEBUG("ParticleMask isol " <<uid);
252 partMask[ pItr->index() ] = true;
253 }
254 }
255 }
256 }
257
258
259 // Retrieve optional jets
260 // Add particles that are constituents of selected jets using unique IDs.
261 // Is index() for JetConstituentVector or TruthParticleContainer or??
262
263 if( !m_truthJetsKey.empty() ){
264
266 if (!inJets.isValid()) {
267 ATH_MSG_ERROR("Cannot retrieve JetContainer " << m_truthJetsKey);
268 return StatusCode::FAILURE;
269 }
270
271 std::vector<int> uidJetConst;
272
273 for(const auto* ajet : *inJets){
274 if( ajet->pt() < m_jetPtCut ) continue;
275 if( std::abs(ajet->eta()) > m_jetEtaCut ) continue;
276
277 xAOD::JetConstituentVector aconst = ajet->getConstituents();
280 for( ; aItr != aItrE; ++aItr){
281 const xAOD::JetConstituent* aip = (*aItr);
282 const xAOD::IParticle* aipraw = aip->rawConstituent();
283 const xAOD::TruthParticle* pp =
284 dynamic_cast<const xAOD::TruthParticle*>(aipraw);
285 if( pp ) {
286 if( pp->pt()>m_jetConstPtCut && !MC::isPhoton(pp) ){
287 uidJetConst.push_back( HepMC::uniqueID(pp) );
288 }
289 if( pp->pt()>m_jetPhotonPtCut && !MC::isPhoton(pp) ){
290 uidJetConst.push_back( HepMC::uniqueID(pp) );
291 }
292 } else {
293 ATH_MSG_WARNING("Bad cast for particle in jet " <<ajet->index());
294 }
295 }
296 }
297
298 for (const auto* pItr : *inTruthParts) {
299 int uid = HepMC::uniqueID(pItr);
300 bool isJet = false;
301 for(unsigned int i=0; i<uidJetConst.size(); ++i){
302 if( uid == uidJetConst[i] ){
303 isJet = true;
304 break;
305 }
306 }
307 if( isJet ){
308 if( doPrint ) ATH_MSG_DEBUG("ParticleMask isJet " <<uid);
309 partMask[ pItr->index() ] = true;
310 }
311 }
312
313 } //end optional jets
314
315
316 // Execute the thinning service based on the mask. Finish.
317 inTruthParts.keep (partMask);
318 inTruthVerts.keep (vertMask);
319
320 // Final statistics
321 int outPartNum = 0;
322 for(unsigned int i=0; i<partMask.size(); ++i){
323 if( partMask[i] ) ++outPartNum;
324 }
325 int outVertNum = 0;
326 for(unsigned int i=0; i<vertMask.size(); ++i){
327 if( vertMask[i] ) ++outVertNum;
328 }
329
330 ATH_MSG_DEBUG("Final particles/vertices = " <<outPartNum <<" " <<outVertNum);
331
332 if( doPrint ){
333 std::cout <<"======================================================================================" <<std::endl;
334 std::cout <<"HardTruthThinning complete for event " <<evtNum <<std::endl;
335 std::cout <<"Saved " <<outPartNum <<" particles" <<std::endl;
336 std::cout <<"Particle unique IDs = ";
337 for(unsigned int i=0; i<partMask.size(); ++i){
338 if( partMask[i] ) std::cout << HepMC::uniqueID((*inTruthParts)[i]) <<" ";
339 }
340 std::cout <<std::endl;
341
342 std::cout <<"Saved " <<outVertNum <<" vertices" <<std::endl;
343 std::cout <<"Vertex unique IDs = ";
344 for(unsigned int i=0; i<vertMask.size(); ++i){
345 if( vertMask[i] ) std::cout << HepMC::uniqueID((*inTruthVerts)[i]) <<" ";
346 }
347 std::cout <<std::endl;
348 std::cout <<"======================================================================================" <<std::endl;
349 }
350
351 return StatusCode::SUCCESS;
352
353}
354
355
357// Utility functions
359
360// Emulate HepMC descendant iterator
361// Multiple particles can give same descendant (string/cluster)
362// Remove Geant descendants
363// MUST check ElementLink validity with isValid() for thinned samples
364
366 const xAOD::TruthParticle* p,
367 std::vector<const xAOD::TruthParticle*>& descendants ) {
368 descendants.clear();
369 if( ! (p->hasDecayVtx()) ) return 0;
370 const xAOD::TruthVertex* dvtx = p->decayVtx();
371 if( !dvtx ) return 0;
372 if( dvtx->nOutgoingParticles() == 0 ) return 0;
373 if(HepMC::is_simulation_vertex(dvtx)) return 0;
374 const std::vector< ElementLink< xAOD::TruthParticleContainer > >& outPart =
375 dvtx->outgoingParticleLinks();
376 for(unsigned int k=0; k<outPart.size(); ++k){
377 if( ! (outPart[k]).isValid() ) continue;
378 const xAOD::TruthParticle* kid = *(outPart[k]);
379 descendants.push_back(kid);
380 }
381
382 int nstart = 0;
383 int nstop = descendants.size();
384
385 while( nstop > nstart ){
386 for(int i=nstart; i<nstop; ++i){
387 const xAOD::TruthParticle* pp = descendants[i];
388 if( ! (pp->hasDecayVtx()) ) continue;
389 const xAOD::TruthVertex* vpp = pp->decayVtx();
390 if( !vpp ) continue;
391 if( vpp->nOutgoingParticles() == 0 ) continue;
392 if( HepMC::is_simulation_vertex(vpp)) continue;
393 const std::vector< ElementLink< xAOD::TruthParticleContainer > >&
394 outPart2 = vpp->outgoingParticleLinks();
395 for(unsigned int k=0; k<outPart2.size(); ++k){
396 if( ! (outPart2[k]).isValid() ) continue;
397 const xAOD::TruthParticle* kpp = *(outPart2[k]);
398 if( HepMC::is_simulation_particle(kpp)) continue;
399 bool isIn = false;
400 for(unsigned int kk=0; kk<descendants.size(); ++kk){
401 if(kpp==descendants[kk]) isIn = true;
402 }
403 if( !isIn ) descendants.push_back(kpp);
404 }
405 }
406 nstart = nstop;
407 nstop = descendants.size();
408 }
409
410 return nstop;
411}
412
413// Print xAODTruth Event. The printout is particle oriented, unlike
414// the HepMC particle/vertex printout. Geant and pileup particles are
415// omitted.
416
418 const xAOD::TruthParticleContainer* truths) {
419
420 std::vector<int> uidPars;
421 std::vector<int> uidKids;
422
423 std::cout <<"======================================================================================\n" ;
424 std::cout <<"xAODTruth Event " <<evnum <<"\n";
425 std::cout <<" Unique ID PDG Id Status px(GeV) py(GeV) pz(GeV) E(GeV) Parent: Decay\n" ;
426 std::cout <<" -----------------------------------------------------------------------------------\n" ;
427
428 for (const auto* tpItr : *truths) {
429 if (HepMC::is_simulation_particle(tpItr)) continue;
430 int uid = HepMC::uniqueID(tpItr);
431 int id = tpItr->pdgId();
432 int stat = tpItr->status();
433 float px = tpItr->px()/GeV;
434 float py = tpItr->py()/GeV;
435 float pz = tpItr->pz()/GeV;
436 float e = tpItr->e()/GeV;
437 uidPars.clear();
438 uidKids.clear();
439
440 if( tpItr->hasProdVtx() ){
441 const xAOD::TruthVertex* pvtx = tpItr->prodVtx();
442 if( pvtx ){
443 const std::vector< ElementLink< xAOD::TruthParticleContainer > >& pars =
444 pvtx->incomingParticleLinks();
445 for(unsigned int k=0; k<pars.size(); ++k){
446 if( ! (pars[k]).isValid() ) continue;
447 const xAOD::TruthParticle* par = *(pars[k]);
448 uidPars.push_back(HepMC::uniqueID(par));
449 }
450 }
451 }
452 if( tpItr->hasDecayVtx() ){
453 const xAOD::TruthVertex* dvtx = tpItr->decayVtx();
454 if( dvtx ){
455 const std::vector< ElementLink< xAOD::TruthParticleContainer > >& kids =
456 dvtx->outgoingParticleLinks();
457 for(unsigned int k=0; k<kids.size(); ++k){
458 if( ! (kids[k]).isValid() ) continue;
459 const xAOD::TruthParticle* kid = *(kids[k]);
460 uidKids.push_back(HepMC::uniqueID(kid));
461 }
462 }
463 }
464
465 std::cout << std::format("{:>10}{:>12}{:>8}{:>10.2f}{:>10.2f}{:>10.2f}{:>10.2f} P: ",
466 uid, id, stat, px, py, pz, e);
467 for(unsigned int k=0; k<uidPars.size(); ++k){
468 std::cout <<uidPars[k] <<" ";
469 }
470 std::cout <<" D: ";
471 for(unsigned int k=0; k<uidKids.size(); ++k){
472 std::cout <<uidKids[k] <<" ";
473 }
474 std::cout <<"\n";
475 }
476 std::cout <<"======================================================================================" <<std::endl;
477}
#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)
#define ATH_MSG_DEBUG(x)
bool isValid(const T &p)
Av: we implement here an ATLAS-sepcific convention: all particles which are 99xxxxx are fine.
Definition AtlasPID.h:878
Helpers for checking error return status codes and reporting errors.
ATLAS-specific HepMC functions.
Handle for requesting thinning for a data object.
virtual StatusCode doThinning() const override
SG::ThinningHandleKey< xAOD::TruthVertexContainer > m_truthVertexName
virtual StatusCode initialize() override
static void printxAODTruth(long long evnum, const xAOD::TruthParticleContainer *truths)
HardTruthThinning(const std::string &t, const std::string &n, const IInterface *p)
SG::ThinningHandleKey< xAOD::TruthParticleContainer > m_truthParticleName
virtual StatusCode finalize() override
SG::ReadHandleKey< xAOD::EventInfo > m_eventInfoKey
SG::ReadHandleKey< xAOD::TruthParticleContainer > m_hardParticleKey
SG::ReadHandleKey< xAOD::JetContainer > m_truthJetsKey
static int getDescendants(const xAOD::TruthParticle *p, std::vector< const xAOD::TruthParticle * > &d)
virtual bool isValid() override final
Can the handle be successfully dereferenced?
void keep(size_t ndx)
Mark that index ndx in the container should be kept (not thinned away).
Handle for requesting thinning for a data object.
Class providing the definition of the 4-vector interface.
A vector of jet constituents at the scale used during jet finding.
iterator begin() const
iterator on the first constituent
iterator end() const
iterator after the last constituent
4-vector of jet constituent at the scale used during jet finding.
const IParticle * rawConstituent() const
Access the real underlying IParticle.
const TruthVertex_v1 * decayVtx() const
The decay vertex of this particle.
bool hasDecayVtx() const
Check for a decay vertex on this particle.
virtual double pt() const override final
The transverse momentum ( ) of the particle.
const TPLinks_t & outgoingParticleLinks() const
Get all the outgoing particles.
size_t nOutgoingParticles() const
Get the number of outgoing particles.
const TPLinks_t & incomingParticleLinks() const
Get all the incoming particles.
int r
Definition globals.cxx:22
bool is_simulation_vertex(const T &v)
Method to establish if the vertex was created during simulation (TODO migrate to be based on status).
int uniqueID(const T &p)
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 isPhoton(const T &p)
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)
bool isTau(const T &p)
Definition index.py:1
TruthVertex_v1 TruthVertex
Typedef to implementation.
Definition TruthVertex.h:15
TruthParticle_v1 TruthParticle
Typedef to implementation.
TruthParticleContainer_v1 TruthParticleContainer
Declare the latest version of the truth particle container.