ATLAS Offline Software
TruthDressingTool.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 // TruthDressingTool.cxx
7 // Author: Kevin Finelli (kevin.finelli@cern.ch)
8 // Create dressed (i.e. including FSR photons) 4-vectors of truth objects
9 
17 
19 
20 #include "fastjet/PseudoJet.hh"
21 #include "fastjet/ClusterSequence.hh"
22 #include <vector>
23 #include <string>
24 #include <algorithm>
25 #include <memory>
26 namespace {
27  static const SG::ConstAccessor<unsigned int> acc_origin("Classification");
28 }
29 // Constructor
31  const std::string& n,
32  const IInterface* p )
33  : AthAlgTool(t,n,p)
34 {
35  declareInterface<DerivationFramework::IAugmentationTool>(this);
36 }
37 
38 // Destructor
40 }
41 
42 // Athena initialize and finalize
44 {
45  // Initialise handle keys
46  ATH_CHECK(m_particlesKey.initialize());
47  ATH_CHECK(m_dressParticlesKey.initialize());
48  m_decorator_eKey = m_dressParticlesKey.key() + ".e_dressed";
49  ATH_CHECK(m_decorator_eKey.initialize());
50  m_decorator_ptKey = m_dressParticlesKey.key() + ".pt_dressed";
51  ATH_CHECK(m_decorator_ptKey.initialize());
52  m_decorator_etaKey= m_dressParticlesKey.key() + ".eta_dressed";
53  ATH_CHECK(m_decorator_etaKey.initialize());
54  m_decorator_phiKey= m_dressParticlesKey.key() + ".phi_dressed";
55  ATH_CHECK(m_decorator_phiKey.initialize());
56  m_decorator_pt_visKey = m_dressParticlesKey.key() + ".pt_vis_dressed";
57  ATH_CHECK(m_decorator_pt_visKey.initialize());
58  m_decorator_eta_visKey= m_dressParticlesKey.key() + ".eta_vis_dressed";
59  ATH_CHECK(m_decorator_eta_visKey.initialize());
60  m_decorator_phi_visKey= m_dressParticlesKey.key() + ".phi_vis_dressed";
61  ATH_CHECK(m_decorator_phi_visKey.initialize());
62  m_decorator_m_visKey= m_dressParticlesKey.key() + ".m_vis_dressed";
63  ATH_CHECK(m_decorator_m_visKey.initialize());
64  m_decorator_nphotonKey = m_dressParticlesKey.key() + ".nPhotons_dressed";
65  ATH_CHECK(m_decorator_nphotonKey.initialize());
66  if (!m_decorationName.empty()) {m_decorationKey = m_particlesKey.key()+"."+m_decorationName;}
67  else {m_decorationKey = m_particlesKey.key()+".unusedPhotonDecoration";}
68  ATH_CHECK(m_decorationKey.initialize());
69  m_truthClassKey = m_dressParticlesKey.key() + "." + SG::AuxTypeRegistry::instance().getName(acc_origin.auxid());
70  ATH_CHECK(m_truthClassKey.initialize());
71  return StatusCode::SUCCESS;
72 }
73 
74 // Function to do dressing, implements interface in IAugmentationTool
76 {
77  // Get the event context
78  const EventContext& ctx = Gaudi::Hive::currentContext();
79 
80  // Retrieve the truth collections
81  SG::ReadHandle<xAOD::TruthParticleContainer> truthParticles(m_particlesKey,ctx);
82  if (!truthParticles.isValid()) {
83  ATH_MSG_ERROR("Couldn't retrieve TruthParticle collection with name " << m_particlesKey);
84  return StatusCode::FAILURE;
85  }
86 
87  SG::ReadHandle<xAOD::TruthParticleContainer> dressTruthParticles(m_dressParticlesKey,ctx);
88  if (!dressTruthParticles.isValid()) {
89  ATH_MSG_ERROR("Couldn't retrieve TruthParticle collection with name " << m_dressParticlesKey);
90  return StatusCode::FAILURE;
91  }
92 
93  // Decorators
94  SG::WriteDecorHandle< xAOD::TruthParticleContainer,float > decorator_e(m_decorator_eKey, ctx);
95  SG::WriteDecorHandle< xAOD::TruthParticleContainer,float > decorator_pt(m_decorator_ptKey, ctx);
96  SG::WriteDecorHandle< xAOD::TruthParticleContainer,float > decorator_eta(m_decorator_etaKey, ctx);
97  SG::WriteDecorHandle< xAOD::TruthParticleContainer,float > decorator_phi(m_decorator_phiKey, ctx);
98  // for truth taus, use 'vis' in the decoration name to avoid prompt/visible tau momentum ambiguity
99  // use (pt,eta,phi,m) for taus, for consistency with other TauAnalysisTools decorations
100  SG::WriteDecorHandle< xAOD::TruthParticleContainer,float > decorator_pt_vis(m_decorator_pt_visKey, ctx);
101  SG::WriteDecorHandle< xAOD::TruthParticleContainer,float > decorator_eta_vis(m_decorator_eta_visKey, ctx);
102  SG::WriteDecorHandle< xAOD::TruthParticleContainer,float > decorator_phi_vis(m_decorator_phi_visKey, ctx);
103  SG::WriteDecorHandle< xAOD::TruthParticleContainer,float > decorator_m_vis(m_decorator_m_visKey, ctx);
104  SG::WriteDecorHandle< xAOD::TruthParticleContainer,int > decorator_nphoton(m_decorator_nphotonKey, ctx);
105  // One for the photons as well
106  // Can't use a handle here, as this decoration gets touched by
107  // multiple algorithms. Need to explicitly schedule a LockDecorations
108  // algorithm to lock it after all modifications.
109  // FIXME: This is not MT-safe.
110  SG::Decorator< char > dressDec (SG::decorKeyFromKey (m_decorationKey.key()));
111  // If we want to decorate, then we need to decorate everything with false to begin with
112  if (!m_decorationName.empty()){
113  if (!dressDec.isAvailable(*truthParticles)) {
114  for (const auto * particle : *truthParticles){
115  dressDec(*particle);
116  }
117  } // Loop over particles
118  } // We are using the decoration
119 
120  //get struct of helper functions
122 
123  std::vector<const xAOD::TruthParticle*> listOfParticlesToDress;
124  std::vector<xAOD::TruthParticle::FourMom_t> listOfDressedParticles;
125  std::vector<int> dressedParticlesNPhot;
126 
127  if(m_listOfPIDs.size()==1 && abs(m_listOfPIDs[0])==15) {
128  // when dressing only truth taus, it is assumed that the truth tau container has
129  // been built beforehand and is used as input
130  for (auto *pItr : *dressTruthParticles) {
131  listOfParticlesToDress.push_back(pItr);
132  }
133  }
134  else {
135  // non-prompt particles are still included here to ensure all particles
136  // will get the decoration; however further down only the prompt particles
137  // are actually dressed depending on the value of m_useLeptonsFromHadrons
138  decayHelper.constructListOfFinalParticles(dressTruthParticles.ptr(), listOfParticlesToDress, m_listOfPIDs, true);
139  }
140 
141  //initialize list of dressed particles
142  for (const auto& part : listOfParticlesToDress) {
143  listOfDressedParticles.push_back(part->p4());
144  dressedParticlesNPhot.push_back(0);
145  }
146 
147  //fill the photon list
148  std::vector<const xAOD::TruthParticle*> photonsFSRList;
149  std::vector<int> photonPID{22};
150  const bool pass = decayHelper.constructListOfFinalParticles(truthParticles.ptr(), photonsFSRList,
151  photonPID, m_usePhotonsFromHadrons);
152  if (!pass) {
153  ATH_MSG_WARNING("Cannot construct the list of final state particles "<<m_truthClassKey.fullKey());
154  }
155 
156  static const SG::ConstAccessor<double> pt_visAcc("pt_vis");
157  static const SG::ConstAccessor<double> eta_visAcc("eta_vis");
158  static const SG::ConstAccessor<double> phi_visAcc("phi_vis");
159  static const SG::ConstAccessor<double> mvisAcc("m_vis");
160 
161  // Do dR-based photon dressing (default)
162  if (!m_useAntiKt){
163  //loop over photons, uniquely associate each to nearest bare particle
164  for (const auto& phot : photonsFSRList ) {
165  double dRmin = m_coneSize;
166  int idx = -1;
167 
168  for (size_t i = 0; i < listOfParticlesToDress.size(); ++i) {
169  if (!m_useLeptonsFromHadrons) {
170  if (!acc_origin.isAvailable(*listOfParticlesToDress[i])) {
171  ATH_MSG_WARNING("MCTruthClassifier "<<m_truthClassKey.fullKey() <<" not available, cannot apply notFromHadron veto!");
172  }
173  unsigned int result = acc_origin(*listOfParticlesToDress[i]);
174  const bool isPrompt = MCTruthPartClassifier::isPrompt(result, true);
175  if (!isPrompt) continue;
176  }
178  if(listOfParticlesToDress[i]->isTau()) {
179 
180  if( !pt_visAcc.isAvailable(*listOfParticlesToDress[i]) ||
181  !eta_visAcc.isAvailable(*listOfParticlesToDress[i]) ||
182  !phi_visAcc.isAvailable(*listOfParticlesToDress[i]) ||
183  !mvisAcc.isAvailable(*listOfParticlesToDress[i])) {
184  ATH_MSG_ERROR("Visible momentum not available for truth taus, cannot perform dressing!");
185  return StatusCode::FAILURE;
186  }
187 
188  bare_part.SetPtEtaPhiM(pt_visAcc(*listOfParticlesToDress[i]),
189  eta_visAcc(*listOfParticlesToDress[i]),
190  phi_visAcc(*listOfParticlesToDress[i]),
191  mvisAcc(*listOfParticlesToDress[i]));
192  }
193  else {
194  bare_part = listOfParticlesToDress[i]->p4();
195  }
196 
197  double dR = bare_part.DeltaR(phot->p4());
198  if (dR < dRmin) {
199  dRmin = dR;
200  idx = i;
201  }
202  }
203 
204  if(idx > -1) {
205  listOfDressedParticles[idx] += phot->p4();
206  dressedParticlesNPhot[idx]++;
207  if (!m_decorationName.empty()){
208  dressDec(*phot) = 1;
209  }
210  }
211  }
212 
213  //loop over particles and add decorators
214  //for (const auto& part : listOfDressedParticles) {
215  for (size_t i = 0; i < listOfParticlesToDress.size(); ++i) {
216  const xAOD::TruthParticle* part = listOfParticlesToDress[i];
217  xAOD::TruthParticle::FourMom_t& dressedVec = listOfDressedParticles[i];
218 
219  if(part->isTau()) {
220  decorator_pt_vis(*part) = dressedVec.Pt();
221  decorator_eta_vis(*part) = dressedVec.Eta();
222  decorator_phi_vis(*part) = dressedVec.Phi();
223  decorator_m_vis(*part) = dressedVec.M();
224  }
225  else {
226  decorator_e(*part) = dressedVec.E();
227  decorator_pt(*part) = dressedVec.Pt();
228  decorator_eta(*part) = dressedVec.Eta();
229  decorator_phi(*part) = dressedVec.Phi();
230  }
231  decorator_nphoton(*part) = dressedParticlesNPhot[i];
232  }
233  } // end of the dR matching part
234 
235  //build the anti-kt jet list
236  if (m_useAntiKt) {
237  std::vector<fastjet::PseudoJet> sorted_jets;
238  std::vector<fastjet::PseudoJet> fj_particles;
239  for (const auto& part : listOfParticlesToDress) {
240 
241  if(part->isTau()) {
242  if(!pt_visAcc.isAvailable(*part) || !eta_visAcc.isAvailable(*part)
243  || !phi_visAcc.isAvailable(*part) || !mvisAcc.isAvailable(*part)) {
244  ATH_MSG_ERROR("Visible momentum not available for truth taus, cannot perform dressing!");
245  return StatusCode::FAILURE;
246  }
247 
248  TLorentzVector tauvis;
249  tauvis.SetPtEtaPhiM(pt_visAcc(*part),
250  eta_visAcc(*part),
251  phi_visAcc(*part),
252  mvisAcc(*part));
253  fj_particles.emplace_back(tauvis.Px(), tauvis.Py(), tauvis.Pz(), tauvis.E());
254  }
255  else {
256  fj_particles.emplace_back(part->px(), part->py(), part->pz(), part->e());
257  }
258 
259  fj_particles.back().set_user_index(HepMC::uniqueID(part));
260  }
261  for (const auto& part : photonsFSRList) {
262  fj_particles.emplace_back(part->px(), part->py(), part->pz(), part->e());
263  fj_particles.back().set_user_index(HepMC::uniqueID(part));
264  }
265 
266  //run the clustering
267  fastjet::JetAlgorithm alg=fastjet::antikt_algorithm;
268  const fastjet::JetDefinition jet_def(alg, m_coneSize);
269  fastjet::ClusterSequence cseq(fj_particles, jet_def);
270  sorted_jets = sorted_by_pt(cseq.inclusive_jets(0));
271  //associate clustered jets back to bare particles
272  std::vector<int> photon_uniqueIDs(50);
273  photon_uniqueIDs.clear();
274  for (const auto& part : listOfParticlesToDress) {
275  //loop over fastjet pseudojets and associate one with this particle
276  bool found=false;
277  auto pjItr=sorted_jets.begin();
278  while(!found && pjItr!=sorted_jets.end()) {
279  std::vector<fastjet::PseudoJet> constituents = pjItr->constituents();
280  for(const auto& constit : constituents) {
281  if (HepMC::uniqueID(part)==constit.user_index()) {
282 
283  // shall we count the number of photons among the pseudojet constituents
284  // to decorate leptons with the number of dressing photons found by the anti-kt algorithm?
285 
286  if(part->isTau()) {
287  decorator_pt_vis(*part) = pjItr->pt();
288  decorator_eta_vis(*part) = pjItr->pseudorapidity();
289  decorator_phi_vis(*part) = pjItr->phi_std(); //returns phi in [-pi,pi]
290  decorator_m_vis(*part) = pjItr->m();
291  }
292  else {
293  decorator_e(*part) = pjItr->e();
294  decorator_pt(*part) = pjItr->pt();
295  decorator_eta(*part) = pjItr->pseudorapidity();
296  decorator_phi(*part) = pjItr->phi_std(); //returns phi in [-pi,pi]
297  }
298  found=true;
299  break;
300  } // Found the matching unique ID
301  } // Loop over the jet constituents
302  if (found){
303  for(const auto& constit : constituents) {
304  photon_uniqueIDs.push_back(constit.user_index());
305  } // Loop over the constituents
306  } // Found one of the key leptons in this jet
307  ++pjItr;
308  }
309  if (!found) {
310  if(part->isTau()) {
311  decorator_pt_vis(*part) = 0.;
312  decorator_eta_vis(*part) = 0.;
313  decorator_phi_vis(*part) = 0.;
314  decorator_m_vis(*part) = 0.;
315  }
316  else {
317  decorator_e(*part) = 0;
318  decorator_pt(*part) = 0;
319  decorator_eta(*part) = 0;
320  decorator_phi(*part) = 0;
321  }
322  ATH_MSG_WARNING("Bare particle not found in constituents ");
323  }
324  }
325  // Check if we wanted to decorate photons used for dressing
326  if (!m_decorationName.empty()){
327  //loop over photons, uniquely associate each to nearest bare particle
328  for (const auto& phot : photonsFSRList ) {
329  bool found = std::find(photon_uniqueIDs.begin(), photon_uniqueIDs.end(), HepMC::uniqueID(phot)) != photon_uniqueIDs.end();
330  if (found) {
331  dressDec(*phot) = 1;
332  }
333  } // End of loop over photons
334  } // End of decoration of photons used in dressing
335  } // End of anti-kT dressing
336 
337  return StatusCode::SUCCESS;
338 }
LArG4FSStartPointFilter.part
part
Definition: LArG4FSStartPointFilter.py:21
DerivationFramework::TruthDressingTool::initialize
StatusCode initialize()
Definition: TruthDressingTool.cxx:43
Trk::ParticleSwitcher::particle
constexpr ParticleHypothesis particle[PARTICLEHYPOTHESES]
the array of masses
Definition: ParticleHypothesis.h:76
SGout2dot.alg
alg
Definition: SGout2dot.py:243
get_generator_info.result
result
Definition: get_generator_info.py:21
find
std::string find(const std::string &s)
return a remapped string
Definition: hcg.cxx:135
DerivationFramework::TruthDressingTool::addBranches
virtual StatusCode addBranches() const
Pass the thinning service
Definition: TruthDressingTool.cxx:75
SG::AuxTypeRegistry::instance
static AuxTypeRegistry & instance()
Return the singleton registry instance.
Definition: AuxTypeRegistry.cxx:639
SG::ReadHandle
Definition: StoreGate/StoreGate/ReadHandle.h:70
SG::decorKeyFromKey
std::string decorKeyFromKey(const std::string &key)
Extract the decoration part of key.
Definition: DecorKeyHelpers.cxx:41
SG::AuxTypeRegistry::getName
std::string getName(SG::auxid_t auxid) const
Return the name of an aux data item.
Definition: AuxTypeRegistry.cxx:881
SG::ConstAccessor< unsigned int >
MCTruthClassifier.h
DerivationFramework::DecayGraphHelper::constructListOfFinalParticles
bool constructListOfFinalParticles(const xAOD::TruthParticleContainer *allParticles, std::vector< const xAOD::TruthParticle * > &selectedlist, const std::vector< int > &pdgId, bool allowFromHadron=false, bool chargedOnly=false) const
Definition: DecayGraphHelper.h:218
read_hist_ntuple.t
t
Definition: read_hist_ntuple.py:5
DecorKeyHelpers.h
Some common helper functions used by decoration handles.
python.utils.AtlRunQueryDQUtils.p
p
Definition: AtlRunQueryDQUtils.py:210
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
SG::Decorator< char >
jet::ClusterSequence
fastjet::ClusterSequence ClusterSequence
Definition: ClusterSequence.h:21
lumiFormat.i
int i
Definition: lumiFormat.py:85
beamspotman.n
n
Definition: beamspotman.py:731
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
xAOD::TruthParticle_v1
Class describing a truth particle in the MC record.
Definition: TruthParticle_v1.h:37
SG::WriteDecorHandle
Handle class for adding a decoration to an object.
Definition: StoreGate/StoreGate/WriteDecorHandle.h:100
WriteDecorHandle.h
Handle class for adding a decoration to an object.
HepMC::uniqueID
int uniqueID(const T &p)
Definition: MagicNumbers.h:116
ATH_CHECK
#define ATH_CHECK
Definition: AthCheckMacros.h:40
TruthDressingTool.h
DerivationFramework::TruthDressingTool::TruthDressingTool
TruthDressingTool(const std::string &t, const std::string &n, const IInterface *p)
Definition: TruthDressingTool.cxx:30
SG::ReadHandle::isValid
virtual bool isValid() override final
Can the handle be successfully dereferenced?
isTau
bool isTau(const T &p)
Definition: AtlasPID.h:173
MagicNumbers.h
DerivationFramework::DecayGraphHelper
Definition: DecayGraphHelper.h:26
xAOD::JetAlgorithmType::antikt_algorithm
@ antikt_algorithm
Definition: JetContainerInfo.h:33
SG::ReadHandle::ptr
const_pointer_type ptr()
Dereference the pointer.
DerivationFramework::TruthDressingTool::~TruthDressingTool
~TruthDressingTool()
Definition: TruthDressingTool.cxx:39
CondAlgsOpts.found
int found
Definition: CondAlgsOpts.py:101
MCTruthPartClassifier::isPrompt
int isPrompt(const unsigned int classify, bool allow_prompt_tau_decays=true)
Definition: TruthClassifiers.h:180
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
SG::Decorator::isAvailable
bool isAvailable(const ELT &e) const
Test to see if this variable exists in the store.
LArNewCalib_DelayDump_OFC_Cali.idx
idx
Definition: LArNewCalib_DelayDump_OFC_Cali.py:69
ReadDecorHandle.h
Handle class for reading a decoration on an object.
SG::ConstAccessor::isAvailable
bool isAvailable(const ELT &e) const
Test to see if this variable exists in the store.
ConstAccessor.h
Helper class to provide constant type-safe access to aux data.
AthAlgTool
Definition: AthAlgTool.h:26
TruthEventContainer.h
xAOD::TruthParticle_v1::FourMom_t
IParticle::FourMom_t FourMom_t
Definition of the 4-momentum type.
Definition: TruthParticle_v1.h:147