ATLAS Offline Software
xAODtoHepMCTool.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 #ifndef XAOD_STANDALONE
8 #endif
10 
12  : asg::AsgTool(name)
13 {
14 }
15 
17 {
18  ATH_MSG_INFO("Initializing xAODtoHepMCTool " << name() << "...");
19  ATH_MSG_INFO("SignalOnly = " << m_signalOnly);
20  return StatusCode::SUCCESS;
21 }
22 
24 {
25  ATH_MSG_INFO("==============================================================");
26  ATH_MSG_INFO("========== xAOD -> HepMC Tool :: Run Summary ==========");
27  ATH_MSG_INFO("==============================================================");
28  if (m_noProdVtx)
29  {
30  ATH_MSG_INFO("Number of events that have a missing production vertex = " << m_noProdVtx);
31  }
32  else
33  {
34  ATH_MSG_INFO("No missing production vertices");
35  }
36  if (m_badBeams)
37  {
38  ATH_MSG_INFO("Number events with improperly defined beamparticles = " << m_badBeams);
39  ATH_MSG_INFO("Inconsistencies with the beam particles storage in the event record!");
40  ATH_MSG_INFO("... check xAOD::TruthEvent record ...");
41  }
42  else
43  {
44  ATH_MSG_INFO("No events with undefined beams.");
45  }
46  ATH_MSG_INFO("==============================================================");
47  ATH_MSG_INFO("=================== End Run Summary ===================");
48  ATH_MSG_INFO("==============================================================");
49  return StatusCode::SUCCESS;
50 }
51 
52 std::vector<HepMC::GenEvent> xAODtoHepMCTool ::getHepMCEvents(const xAOD::TruthEventContainer *xTruthEventContainer, const xAOD::EventInfo *eventInfo) const
53 {
54  ++m_evtCount;
55  bool doPrint = m_evtCount < m_maxCount;
56  // CREATE HEPMC EVENT COLLECTION
57  std::vector<HepMC::GenEvent> mcEventCollection;
58  // Loop over xAOD truth container
59  // Signal event is always first, followed by pileup events
60  ATH_MSG_DEBUG("Start event loop");
61  for (const xAOD::TruthEvent *xAODEvent : *xTruthEventContainer)
62  {
63  if (doPrint)
64  ATH_MSG_DEBUG("XXX Printing xAOD Event");
65  if (doPrint)
66  printxAODEvent(xAODEvent, eventInfo);
67  // Create GenEvent for each xAOD truth event
68  ATH_MSG_DEBUG("Create new GenEvent");
69  HepMC::GenEvent&& hepmcEvent = createHepMCEvent(xAODEvent, eventInfo);
70  #ifdef HEPMC3
71  std::shared_ptr<HepMC3::GenRunInfo> runinfo = std::make_shared<HepMC3::GenRunInfo>(*(hepmcEvent.run_info().get()));
72  #endif
73  // Insert into McEventCollection
74  mcEventCollection.push_back(std::move(hepmcEvent));
75  #ifdef HEPMC3
76  mcEventCollection[mcEventCollection.size()-1].set_run_info(runinfo);
77  #endif
78  if (doPrint)
79  ATH_MSG_DEBUG("XXX Printing HepMC Event");
80  if (doPrint)
81  HepMC::Print::line(std::cout, hepmcEvent);
82  // Quit if signal only
83  if (m_signalOnly)
84  break;
85  }
86  return mcEventCollection;
87 }
88 
89 HepMC::GenEvent xAODtoHepMCTool::createHepMCEvent(const xAOD::TruthEvent *xEvt, const xAOD::EventInfo *eventInfo) const
90 {
91 
93  HepMC::GenEvent genEvt;
94  #ifdef HEPMC3
95  genEvt.set_units(HepMC3::Units::MEV, HepMC3::Units::MM);
96  #endif
97 
98  long long int evtNum = eventInfo->eventNumber();
99  genEvt.set_event_number(evtNum);
100  ATH_MSG_DEBUG("Start createHepMCEvent for event " << evtNum);
101 
102  //Weights
103  #ifndef XAOD_STANDALONE
104  const std::vector<float> weights = xEvt->weights();
105  std::map<std::string, int> weightNameMap;
106  if (AthAnalysisHelper::retrieveMetadata("/Generation/Parameters","HepMCWeightNames", weightNameMap).isFailure()) {
107  ATH_MSG_DEBUG("Couldn't find meta-data for weight names.");
108  }
109  #ifdef HEPMC3
110  std::shared_ptr<HepMC3::GenRunInfo> runinfo = std::make_shared<HepMC3::GenRunInfo>();
111  genEvt.set_run_info(runinfo);
112  std::vector<std::string> wnames;
113  wnames.reserve(weights.size());
114  for (int idx = 0; idx < int(weights.size()); ++idx) {
115  for (const auto& it : weightNameMap) {
116  if (it.second == idx) {
117  wnames.push_back(it.first);
118  break;
119  }
120  }
121  }
122  genEvt.run_info()->set_weight_names(wnames);
123  for ( std::vector<float>::const_iterator wgt = weights.begin(); wgt != weights.end(); ++wgt ) {
124  genEvt.weights().push_back(*wgt);
125  }
126  #else
127  if (weightNameMap.size()) {
128  HepMC::WeightContainer& wc = genEvt.weights();
129  wc.clear();
130  for (int idx = 0; idx < int(weights.size()); ++idx) {
131  for (const auto& it : weightNameMap) {
132  if (it.second == idx) {
133  wc[ it.first ] = weights[idx];
134  break;
135  }
136  }
137  }
138  }
139  else {
140  for ( std::vector<float>::const_iterator wgt = weights.begin(); wgt != weights.end(); ++wgt ) {
141  genEvt.weights().push_back(*wgt);
142  }
143  }
144  #endif
145  #endif
146 
147  // PARTICLES AND VERTICES
148  // Map of existing vertices - needed for the tree linking
149  std::map<const xAOD::TruthVertex *, HepMC::GenVertexPtr> vertexMap;
150 
151  // Loop over all of the particles in the event, call particle builder
152  for (auto tlink : xEvt->truthParticleLinks()) {
153  if (!tlink.isValid()) { continue; }
154  const xAOD::TruthParticle *xPart = *tlink;
155 
156  // sanity check
157  if (xPart == nullptr) {
158  ATH_MSG_WARNING("xAOD TruthParticle is equal to NULL. This should not happen!");
159  continue;
160  }
161 
162  if (!xPart->hasProdVtx() && !xPart->hasDecayVtx()) {
163  ATH_MSG_WARNING("xAOD particle with no vertices, uniqueID = " << HepMC::uniqueID(xPart));
164  continue;
165  }
166 
167  // skip particles which are Geant4 secondaries
168  if (HepMC::is_simulation_particle(xPart)) { continue; }
169 
170  // Create GenParticle
171  // presumably the GenEvent takes ownership of this, but creating a unique_ptr here as that will only happen if there's an associated vertex
172 #ifdef HEPMC3
173  auto hepmcParticle = createHepMCParticle(xPart);
174 #else
175  std::unique_ptr<HepMC::GenParticle> hepmcParticle(createHepMCParticle(xPart));
176 #endif
177  int bcpart = HepMC::barcode(xPart);
178 
179  // Get the production and decay vertices
180  if (xPart->hasProdVtx()) {
181  const xAOD::TruthVertex *xAODProdVtx = xPart->prodVtx();
182  // skip production vertices which are Geant4 secondaries
183  if (HepMC::is_simulation_vertex(xAODProdVtx)) { continue; }
184  bool prodVtxSeenBefore(false); // is this new?
185  auto hepmcProdVtx = vertexHelper(xAODProdVtx, vertexMap, prodVtxSeenBefore);
186  // Set the decay/production links
187 #ifdef HEPMC3
188  hepmcProdVtx->add_particle_out(hepmcParticle);
189 #else
190  hepmcProdVtx->add_particle_out(hepmcParticle.get());
191 #endif
192  // Insert into Event
193  if (!prodVtxSeenBefore) {
194  genEvt.add_vertex(hepmcProdVtx);
195  if (!HepMC::suggest_barcode(hepmcProdVtx, HepMC::barcode(xAODProdVtx))) {
196  ATH_MSG_WARNING("suggest_barcode failed for vertex " << HepMC::barcode(xAODProdVtx));
197  ++m_badSuggest;
198  }
199  }
200  if (!HepMC::suggest_barcode(hepmcParticle, bcpart)) {
201  ATH_MSG_VERBOSE("suggest_barcode failed for particle " << bcpart);
202  ++m_badSuggest;
203  }
204  bcpart = 0;
205  }
206  else {
207  ATH_MSG_VERBOSE("No production vertex found for particle " << HepMC::uniqueID(xPart));
208  }
209 
210  if (xPart->hasDecayVtx()) {
211  const xAOD::TruthVertex *xAODDecayVtx = xPart->decayVtx();
212  // skip decay vertices which are Geant4 secondaries
213  if (HepMC::is_simulation_vertex(xAODDecayVtx)) {
215 #ifndef HEPMC3
216  if (xPart->hasProdVtx()) { (void)hepmcParticle.release(); }
217 #endif
218  continue;
219  }
220  bool decayVtxSeenBefore(false); // is this new?
221  auto hepmcDecayVtx = vertexHelper(xAODDecayVtx, vertexMap, decayVtxSeenBefore);
222  // Set the decay/production links
223 #ifdef HEPMC3
224  hepmcDecayVtx->add_particle_in(hepmcParticle);
225 #else
226  hepmcDecayVtx->add_particle_in(hepmcParticle.get());
227 #endif
228  // Insert into Event
229  if (!decayVtxSeenBefore) {
230  genEvt.add_vertex(hepmcDecayVtx);
231  if (!HepMC::suggest_barcode(hepmcDecayVtx, HepMC::barcode(xAODDecayVtx))) {
232  ATH_MSG_WARNING("suggest_barcode failed for vertex "
233  << HepMC::barcode(xAODDecayVtx));
234  ++m_badSuggest;
235  }
236  }
237  if (bcpart != 0) {
238  if (!HepMC::suggest_barcode(hepmcParticle, bcpart)) {
239  ATH_MSG_DEBUG("suggest_barcode failed for particle " << bcpart);
240  ++m_badSuggest;
241  }
242  bcpart = 0;
243  }
244  }
245 #ifndef HEPMC3
246  (void)hepmcParticle.release();
247 #endif
248 
249  } // end of particle loop
250 
251  return genEvt;
252 }
253 
254 // Helper to check whether a vertex exists or not using a map;
255 // calls createHepMCVertex if not
257  std::map<const xAOD::TruthVertex *, HepMC::GenVertexPtr> &vertexMap,
258  bool &seenBefore) const
259 {
260 
261  HepMC::GenVertexPtr hepmcVertex;
263  vMapItr = vertexMap.find(xaodVertex);
264  // Vertex seen before?
265  if (vMapItr != vertexMap.end()) {
266  // YES: use the HepMC::Vertex already in the map
267  hepmcVertex = (*vMapItr).second;
268  seenBefore = true;
269  }
270  else {
271  // NO: create a new HepMC::Vertex
272  vertexMap[xaodVertex] = createHepMCVertex(xaodVertex);
273  hepmcVertex = vertexMap[xaodVertex];
274  seenBefore = false;
275  }
276  return hepmcVertex;
277 }
278 
279 // Create the HepMC GenParticle
281 {
282  ATH_MSG_VERBOSE("Creating GenParticle for uniqueID " << HepMC::uniqueID(particle));
283  const HepMC::FourVector fourVec(m_momFac * particle->px(), m_momFac * particle->py(), m_momFac * particle->pz(), m_momFac * particle->e());
284  auto hepmcParticle = HepMC::newGenParticlePtr(fourVec, particle->pdgId(), particle->status());
285  hepmcParticle->set_generated_mass(m_momFac * particle->m());
286  return hepmcParticle;
287 }
288 
289 // Create the HepMC GenVertex
291 {
292  ATH_MSG_VERBOSE("Creating GenVertex for uniqueID " << HepMC::uniqueID(vertex));
293  HepMC::FourVector prod_pos(m_lenFac * vertex->x(), m_lenFac * vertex->y(), m_lenFac * vertex->z(), m_lenFac * vertex->t());
294  auto genVertex = HepMC::newGenVertexPtr(prod_pos);
295  return genVertex;
296 }
297 
298 // Print xAODTruth Event. The printout is particle oriented, unlike the HepMC particle/vertex printout. Regenerated Geant (i.e. ones surving an interaction) and pileup particles (e.g. suppressed pileup barcodes), with barcode > 1000000 are omitted. Thus we output only the original particles created in Geant4.
299 
301 {
302 
303  std::vector<int> uidPars;
304  std::vector<int> uidKids;
305 
306  long long int evtNum = eventInfo->eventNumber();
307 
308  std::cout << "======================================================================================" << std::endl;
309  std::cout << "xAODTruth Event " << evtNum << std::endl;
310  std::cout << " UniqueID PDG Id Status px(GeV) py(GeV) pz(GeV) E(GeV) Parent: Decay" << std::endl;
311  std::cout << " -----------------------------------------------------------------------------------" << std::endl;
312 
313  int nPart = event->nTruthParticles();
314  for (int i = 0; i < nPart; ++i)
315  {
316  const xAOD::TruthParticle *part = event->truthParticle(i);
317  if (part == nullptr) continue;
318  if (HepMC::generations(part) > 0) continue;
319  int id = part->pdgId();
320  if (id != 25) continue;
321  int stat = part->status();
322  float px = part->px() / 1000.;
323  float py = part->py() / 1000.;
324  float pz = part->pz() / 1000.;
325  float e = part->e() / 1000.;
326  uidPars.clear();
327  uidKids.clear();
328 
329  if (part->hasProdVtx())
330  {
331  const xAOD::TruthVertex *pvtx = part->prodVtx();
332  if (pvtx)
333  uidPars.push_back(HepMC::uniqueID(pvtx));
334  }
335 
336  if (part->hasDecayVtx())
337  {
338  const xAOD::TruthVertex *dvtx = part->decayVtx();
339  if (dvtx)
340  uidKids.push_back(HepMC::uniqueID(dvtx));
341  }
342 
343  std::cout << std::setw(10) << HepMC::uniqueID(part) << std::setw(12) << id
344  << std::setw(8) << stat
345  << std::setprecision(2) << std::fixed
346  << std::setw(10) << px << std::setw(10) << py
347  << std::setw(10) << pz << std::setw(10) << e << " ";
348  std::cout << "P: ";
349  for (unsigned int k = 0; k < uidPars.size(); ++k)
350  {
351  std::cout << uidPars[k] << " ";
352  }
353  std::cout << " D: ";
354  for (unsigned int k = 0; k < uidKids.size(); ++k)
355  {
356  std::cout << uidKids[k] << " ";
357  }
358  std::cout << std::endl;
359  }
360  std::cout << "======================================================================================" << std::endl;
361 }
HepMC::GenVertexPtr
HepMC::GenVertex * GenVertexPtr
Definition: GenVertex.h:59
LArG4FSStartPointFilter.part
part
Definition: LArG4FSStartPointFilter.py:21
xAOD::iterator
JetConstituentVector::iterator iterator
Definition: JetConstituentVector.cxx:68
AllowedVariables::e
e
Definition: AsgElectronSelectorTool.cxx:37
HepMC::suggest_barcode
bool suggest_barcode(T &p, int i)
Definition: GenEvent.h:670
Trk::ParticleSwitcher::particle
constexpr ParticleHypothesis particle[PARTICLEHYPOTHESES]
the array of masses
Definition: ParticleHypothesis.h:76
xAODtoHepMCTool::m_signalOnly
Gaudi::Property< bool > m_signalOnly
Definition: xAODtoHepMCTool.h:49
xAODtoHepMCTool::xAODtoHepMCTool
xAODtoHepMCTool(const std::string &name)
Definition: xAODtoHepMCTool.cxx:11
test_pyathena.px
px
Definition: test_pyathena.py:18
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
xAOD::EventInfo_v1::eventNumber
uint64_t eventNumber() const
The current event's event number.
CaloCellPos2Ntuple.int
int
Definition: CaloCellPos2Ntuple.py:24
HepMC::GenParticlePtr
GenParticle * GenParticlePtr
Definition: GenParticle.h:37
skel.it
it
Definition: skel.GENtoEVGEN.py:396
asg
Definition: DataHandleTestTool.h:28
MM
@ MM
Definition: RegSelEnums.h:38
ATH_MSG_VERBOSE
#define ATH_MSG_VERBOSE(x)
Definition: AthMsgStreamMacros.h:28
HepMC::Print::line
void line(std::ostream &os, const GenEvent &e)
Definition: GenEvent.h:676
DeMoUpdate.runinfo
def runinfo
Definition: DeMoUpdate.py:542
xAODtoHepMCTool::m_badBeams
int m_badBeams
Definition: xAODtoHepMCTool.h:60
xAODtoHepMCTool::printxAODEvent
void printxAODEvent(const xAOD::TruthEvent *event, const xAOD::EventInfo *eventInfo) const
Definition: xAODtoHepMCTool.cxx:300
xAODtoHepMCTool::m_evtCount
std::atomic< int > m_evtCount
Counters.
Definition: xAODtoHepMCTool.h:57
xAODtoHepMCTool::m_noProdVtx
int m_noProdVtx
Definition: xAODtoHepMCTool.h:59
xAODtoHepMCTool.h
xAODtoHepMCTool::m_maxCount
Gaudi::Property< int > m_maxCount
Definition: xAODtoHepMCTool.h:52
xAODtoHepMCTool::createHepMCEvent
HepMC::GenEvent createHepMCEvent(const xAOD::TruthEvent *xEvt, const xAOD::EventInfo *eventInfo) const
Definition: xAODtoHepMCTool.cxx:89
xAOD::TruthParticle_v1::hasDecayVtx
bool hasDecayVtx() const
Check for a decay vertex on this particle.
xAODtoHepMCTool::finalize
StatusCode finalize() override
Definition: xAODtoHepMCTool.cxx:23
event
POOL::TEvent event(POOL::TEvent::kClassAccess)
HepMC::is_simulation_particle
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...
Definition: MagicNumbers.h:355
HepMC::newGenVertexPtr
GenVertexPtr newGenVertexPtr(const HepMC::FourVector &pos=HepMC::FourVector(0.0, 0.0, 0.0, 0.0), const int i=0)
Definition: GenVertex.h:64
lumiFormat.i
int i
Definition: lumiFormat.py:85
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
HepMC::barcode
int barcode(const T *p)
Definition: Barcode.h:16
xAOD::TruthParticle_v1
Class describing a truth particle in the MC record.
Definition: TruthParticle_v1.h:37
xAODtoHepMCTool::m_lenFac
Gaudi::Property< float > m_lenFac
Definition: xAODtoHepMCTool.h:46
xAOD::TruthParticle_v1::hasProdVtx
bool hasProdVtx() const
Check for a production vertex on this particle.
Definition: TruthParticle_v1.cxx:74
HepMC::uniqueID
int uniqueID(const T &p)
Definition: MagicNumbers.h:116
Amg::pz
@ pz
Definition: GeoPrimitives.h:40
xAOD::TruthEvent_v1
Class describing a signal truth event in the MC record.
Definition: TruthEvent_v1.h:35
HepMC::is_simulation_vertex
bool is_simulation_vertex(const T &v)
Method to establish if the vertex was created during simulation (TODO migrate to be based on status).
Definition: MagicNumbers.h:361
DataVector
Derived DataVector<T>.
Definition: DataVector.h:794
beamspotman.stat
stat
Definition: beamspotman.py:266
xAOD::TruthParticle_v1::decayVtx
const TruthVertex_v1 * decayVtx() const
The decay vertex of this particle.
xAOD::TruthEvent_v1::weights
const std::vector< float > & weights() const
Const access to the weights vector.
xAOD::TruthParticle_v1::prodVtx
const TruthVertex_v1 * prodVtx() const
The production vertex of this particle.
Definition: TruthParticle_v1.cxx:80
xAOD::TruthVertex_v1
Class describing a truth vertex in the MC record.
Definition: TruthVertex_v1.h:37
xAODtoHepMCTool::createHepMCParticle
HepMC::GenParticlePtr createHepMCParticle(const xAOD::TruthParticle *) const
Definition: xAODtoHepMCTool.cxx:280
Amg::py
@ py
Definition: GeoPrimitives.h:39
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:228
MagicNumbers.h
weights
Definition: herwig7_interface.h:38
xAOD::EventInfo_v1
Class describing the basic event information.
Definition: EventInfo_v1.h:43
AthAnalysisHelper.h
Trk::vertex
@ vertex
Definition: MeasurementType.h:21
HepMC::newGenParticlePtr
GenParticlePtr newGenParticlePtr(const HepMC::FourVector &mom=HepMC::FourVector(0.0, 0.0, 0.0, 0.0), int pid=0, int status=0)
Definition: GenParticle.h:39
xAODtoHepMCTool::initialize
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
Definition: xAODtoHepMCTool.cxx:16
xAODtoHepMCTool::vertexHelper
HepMC::GenVertexPtr vertexHelper(const xAOD::TruthVertex *, std::map< const xAOD::TruthVertex *, HepMC::GenVertexPtr > &, bool &) const
Definition: xAODtoHepMCTool.cxx:256
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
LArNewCalib_DelayDump_OFC_Cali.idx
idx
Definition: LArNewCalib_DelayDump_OFC_Cali.py:69
xAODtoHepMCTool::createHepMCVertex
HepMC::GenVertexPtr createHepMCVertex(const xAOD::TruthVertex *) const
Definition: xAODtoHepMCTool.cxx:290
xAODtoHepMCTool::getHepMCEvents
std::vector< HepMC::GenEvent > getHepMCEvents(const xAOD::TruthEventContainer *xTruthEventContainer, const xAOD::EventInfo *eventInfo) const override
Definition: xAODtoHepMCTool.cxx:52
xAODtoHepMCTool::m_momFac
Gaudi::Property< float > m_momFac
Input container key (job property)
Definition: xAODtoHepMCTool.h:43
xAOD::TruthEventBase_v1::truthParticleLinks
const TruthParticleLinks_t & truthParticleLinks() const
Get all the truth particles.
HepMC::generations
int generations(const T &p)
Method to return how many interactions a particle has undergone during simulation (TODO migrate to be...
Definition: MagicNumbers.h:358
xAODtoHepMCTool::m_badSuggest
std::atomic< int > m_badSuggest
Definition: xAODtoHepMCTool.h:58
fitman.k
k
Definition: fitman.py:528
AthAnalysisHelper::retrieveMetadata
static std::string retrieveMetadata(const std::string &folder, const std::string &key, const ServiceHandle< StoreGateSvc > &inputMetaStore)
method that always returns as a string you can use from, e.g, pyROOT with evt = ROOT....
Definition: AthAnalysisHelper.h:254