ATLAS Offline Software
Loading...
Searching...
No Matches
xAODTruthCnvAlg.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4#include "xAODTruthCnvAlg.h"
5
7#include "AthLinks/ElementLink.h"
8
9#include "GaudiKernel/MsgStream.h"
10#include "GaudiKernel/DataSvc.h"
11#include "GaudiKernel/PhysicalConstants.h"
14
18
22
26
30
33
36
40
41bool isSeparatorGenEvent(const HepMC::GenEvent *genEvt) {
42 // Separator defined by pid==0 and eventNumber==-1 as per
43 // https://twiki.cern.ch/twiki/bin/viewauth/AtlasComputing/PileupDigitization#Arrangement_of_Truth_Information
44 const int pid = HepMC::signal_process_id(genEvt);
45 const int eventNumber = genEvt->event_number();
46 return (pid==0 && eventNumber==-1);
47}
48
49using namespace std;
50
51namespace xAODMaker {
52
53
54 xAODTruthCnvAlg::xAODTruthCnvAlg( const string& name, ISvcLocator* svcLoc )
55 : AthReentrantAlgorithm( name, svcLoc )
56 , m_metaStore( "StoreGateSvc/MetaDataStore", name )
57 , m_firstBeginRun(true)
58 {
59 // leaving metadata alone for now--to be updated
60 declareProperty( "MetaObjectName", m_metaName = "TruthMetaData" );
61 declareProperty( "MetaDataStore", m_metaStore );
62 }
63
64
67 ATH_MSG_FATAL( "Contradictory xAOD truth pile-up setting: all pile-up AND in-time alone requested simultaneously. Check settings." );
68 return StatusCode::FAILURE;
69 }
70
71 if (m_writeMetaData) {
72 ATH_CHECK( m_meta.initialize (m_metaStore, m_metaName) );
73 }
74
75 // initialize handles
77
78 // only if doing full truth
79 ATH_CHECK(m_aodContainerKey.initialize());
86 }
87
88 ATH_CHECK(m_evtInfo.initialize());
89
90 ServiceHandle<IIncidentSvc> incSvc("IncidentSvc", name());
91 ATH_CHECK(incSvc.retrieve());
92 incSvc->addListener( this, "BeginRun", 10);
93
94 ATH_MSG_DEBUG("AODContainerName = " << m_aodContainerKey.key() );
95 ATH_MSG_DEBUG("xAOD TruthEventContainer name = " << m_xaodTruthEventContainerKey.key() );
96 ATH_MSG_DEBUG("xAOD TruthPileupEventContainer name = " << m_xaodTruthPUEventContainerKey.key() );
97 ATH_MSG_DEBUG("xAOD TruthParticleContainer name = " << m_xaodTruthParticleContainerKey.key() );
98 ATH_MSG_DEBUG("xAOD TruthVertexContainer name = " << m_xaodTruthVertexContainerKey.key() );
100 ATH_MSG_DEBUG("xAOD TruthLHEParticleContainer name = " << m_lheTruthParticleContainerKey.key() );
101 }
102
103 if (m_doAllPileUp) ATH_MSG_INFO( "All pile-up truth (including out-of-time) will be written" );
104 if (m_doInTimePileUp) ATH_MSG_INFO( "In-time pile-up truth (but not out-of-time) will be written" );
105 if (!m_doAllPileUp && !m_doInTimePileUp) ATH_MSG_INFO( "No pile-up truth will be written" );
106
107 return StatusCode::SUCCESS;
108 }
109
110
111 StatusCode xAODTruthCnvAlg::execute (const EventContext& ctx) const {
112
114 ATH_CHECK(truthLinkVec.record(std::make_unique<xAODTruthParticleLinkVector>()));
115
116 // Retrieve the HepMC truth:
118 // validity check is only really needed for serial running. Remove when MT is only way.
119 if (!mcColl.isValid()) {
120 ATH_MSG_ERROR("Could not retrieve HepMC with key:" << m_aodContainerKey.key());
121 return StatusCode::FAILURE;
122 } else {
123 ATH_MSG_DEBUG( "Retrieved HepMC with key: " << m_aodContainerKey.key() );
124 }
125
126 // **************************************************************
127 // Create the xAOD containers and their auxiliary stores:
128 // **************************************************************
129 // Signal event
131 ATH_CHECK(xTruthEventContainer.record(std::make_unique<xAOD::TruthEventContainer>(),
132 std::make_unique<xAOD::TruthEventAuxContainer>()));
133 ATH_MSG_DEBUG( "Recorded TruthEventContainer with key: " << m_xaodTruthEventContainerKey.key() );
134
135 // Pile-up events
136 SG::WriteHandle<xAOD::TruthPileupEventContainer> xTruthPileupEventContainer;
139 ATH_CHECK(xTruthPileupEventContainer.record(std::make_unique<xAOD::TruthPileupEventContainer>(),
140 std::make_unique<xAOD::TruthPileupEventAuxContainer>()));
141 ATH_MSG_DEBUG( "Recorded TruthPileupEventContainer with key: " << m_xaodTruthPUEventContainerKey.key() );
142 }
143
144 // Particles
146 ATH_CHECK(xTruthParticleContainer.record(std::make_unique<xAOD::TruthParticleContainer>(),
147 std::make_unique<xAOD::TruthParticleAuxContainer>()));
148 ATH_MSG_DEBUG( "Recorded TruthParticleContainer with key: " << m_xaodTruthParticleContainerKey.key() );
149 // Vertices
151 ATH_CHECK(xTruthVertexContainer.record(std::make_unique<xAOD::TruthVertexContainer>(),
152 std::make_unique<xAOD::TruthVertexAuxContainer>()));
153 ATH_MSG_DEBUG( "Recorded TruthVertexContainer with key: " << m_xaodTruthVertexContainerKey.key() );
154
155 // To keep track of whether we wrote an LHE event already or not
156#ifdef HEPMC3
157 bool hadLHERecord = false;
158#endif
159
160 // ***********************************************************************************
161 // Create the xAOD objects
162 // This consists of three parts:
163 // (1) For each Athena event, loop over the GenEvents and build TruthEvent collections
164 // In principle there can be more than one GenEvent per event
165 // (2) For each GenEvent, loop over the GenParticles. For each GenParticle:
166 // (a) Create a TruthParticle.
167 // (b) Call fillParticle
168 // (c) Add the TruthParticle to the TruthParticle container, and add
169 // an EL to this TruthParticle to the truthParticles in TruthEvent
170 // (call this EL eltp)
171 // (d) For the GenVertex * that's this particle's production vertex,
172 // (i) see if it is in tempMap. If not, add it.
173 // (ii) add a copy of eltp (the EL to this truth particle) to map[gv].second
174 // (e) For the GenVertex * that's this particle's decay vertex,
175 // (i) see if it is in tempMap. If not, add it.
176 // (ii) add a copy of eltp (the EL to this truth particle) to map[gv].first
177 // (3) Iterate over tempMap. For each GenVertex *:
178 // (a) Create a TruthVertex
179 // (b) Call fillVertex
180 // (c) Add the TruthVertex to the TruthTruth container, and add an EL to this TruthVertex
181 // to the truthVertices in TruthEvent. (call this EL eltv)
182 // (d) call tv->setIncomingParticles(mapiter.second.first) <- I think mapiter.second.first is the first of the pair
183 // (e) call tv->setOutgoingParticles(mapiter.second.second)
184 // (f) Iterate through the incomingParticles, and set the decay vertex EL as eltv.
185 // (g) Iterate through the outgoingParticles, and set the incoming vertex EL as eltv.
186 //
187 // Comment lines below follow this recipe
188 // ************************************************************************************
189
190 // (1) Build TruthEvents
191 ATH_MSG_DEBUG("Number of GenEvents in this Athena event = " << mcColl->size());
192#ifdef HEPMC3
193 bool newAttributesPresent(false);
194#endif
195 for (unsigned int cntr = 0; cntr < mcColl->size(); ++cntr) {
196 const HepMC::GenEvent* genEvt = (*mcColl)[cntr];
197 bool isSignalProcess(false);
198 if (cntr==0) {
199 isSignalProcess=true;
200#ifdef HEPMC3
201 auto bunchCrossingTime = genEvt->attribute<HepMC3::IntAttribute>("BunchCrossingTime");
202 if (bunchCrossingTime) {
203 newAttributesPresent = true;
204 ATH_MSG_VERBOSE("New attributes present.");
205 }
206 else {
207 ATH_MSG_VERBOSE("New attributes missing.");
208 }
209#else
210 ATH_MSG_VERBOSE("New attributes missing.");
211#endif
212 }
213 if (cntr>0) {
214 // Handle pile-up events
215 if (!m_doInTimePileUp && !m_doAllPileUp) break;
216 isSignalProcess=false;
217#ifdef HEPMC3
218 auto bunchCrossingTime = genEvt->attribute<HepMC3::IntAttribute>("BunchCrossingTime");
219 if (bunchCrossingTime) {
220 // New approach based on checking the bunch crossing
221 // time directly.
222 if (m_doInTimePileUp && std::abs(bunchCrossingTime->value()) > 0) {
223 // Skip out-of-time pile-up events
224 continue;
225 }
226 }
227 else {
228 // Old approach based on McEventCollection structure. If
229 // in-time pileup only is requested, loop stops when the
230 // separator GenEvent between out-of-time and in-time is
231 // reached
232 if (m_doInTimePileUp && isSeparatorGenEvent(genEvt)) {
233 if (newAttributesPresent) {
234 // Old structure with new Attributes?! Check all
235 // GenEvents just in case.
236 ATH_MSG_VERBOSE("New-style, but seeing separator GenEvents");
237 continue;
238 }
239 // Old structure - stop at the first separator
240 // GenEvent.
241 break;
242 }
243 }
244#else
245 // Old approach based on McEventCollection structure. If
246 // in-time pileup only is requested, loop stops when the
247 // separator GenEvent between out-of-time and in-time is
248 // reached
249 if (m_doInTimePileUp && isSeparatorGenEvent(genEvt)) {
250 // Old structure - stop at the first separator
251 // GenEvent.
252 break;
253 }
254#endif
255 }
256
257 xAOD::TruthEvent* xTruthEvent = nullptr;
258 xAOD::TruthPileupEvent* xTruthPileupEvent = nullptr;
259
260
261 if (isSignalProcess) {
262 xTruthEvent = xTruthEventContainer->push_back( std::make_unique<xAOD::TruthEvent>() );
263 // Cross-section
264 auto crossSection = genEvt->cross_section();
265#ifdef HEPMC3
266 xTruthEvent->setCrossSection(crossSection ? (float)crossSection->xsec() : -1);
267 xTruthEvent->setCrossSectionError(crossSection ? (float)crossSection->xsec_err() : -1);
268#else
269 xTruthEvent->setCrossSection(crossSection ? (float)crossSection->cross_section() : -1);
270 xTruthEvent->setCrossSectionError(crossSection ? (float)crossSection->cross_section_error() : -1);
271#endif
272
273 if (m_writeMetaData) {
274 //The mcChannelNumber is used as a unique identifier for which truth meta data belongs to
275 uint32_t mcChannelNumber = 0;
277 if (evtInfo.isValid()) {
278 mcChannelNumber = evtInfo->mcChannelNumber();
279 if (mcChannelNumber==0) mcChannelNumber = evtInfo->runNumber();
280 }
281 else {
282 ATH_MSG_FATAL("Faied to retrieve EventInfo");
283 return StatusCode::FAILURE;
284 }
285
286 ATH_CHECK( m_meta.maybeWrite (mcChannelNumber, *genEvt, m_metaFields) );
287 }
288 // Event weights
289 vector<float> weights;
290 for (const double& w : genEvt->weights()) weights.push_back((float)(w));
291 //AV This to be decided. It is always a good idea to have a default weight 1.0.
292 //if (weights.empty()) weights.push_back(1.0);
293 xTruthEvent->setWeights(weights);
294
295 // Heavy ion info
296 auto const hiInfo = genEvt->heavy_ion();
297 if (hiInfo) {
298#ifdef HEPMC3
299 /* Please note HepMC3 as well as more recent HePMC2 versions have more Hi parameters */
300 xTruthEvent->setHeavyIonParameter(hiInfo->Ncoll_hard, xAOD::TruthEvent::NCOLLHARD);
301 xTruthEvent->setHeavyIonParameter(hiInfo->Npart_proj, xAOD::TruthEvent::NPARTPROJ);
302 xTruthEvent->setHeavyIonParameter(hiInfo->Npart_targ, xAOD::TruthEvent::NPARTTARG);
303 xTruthEvent->setHeavyIonParameter(hiInfo->Ncoll, xAOD::TruthEvent::NCOLL);
304 xTruthEvent->setHeavyIonParameter(hiInfo->spectator_neutrons, xAOD::TruthEvent::SPECTATORNEUTRONS);
305 xTruthEvent->setHeavyIonParameter(hiInfo->spectator_protons, xAOD::TruthEvent::SPECTATORPROTONS);
306 xTruthEvent->setHeavyIonParameter(hiInfo->N_Nwounded_collisions, xAOD::TruthEvent::NNWOUNDEDCOLLISIONS);
307 xTruthEvent->setHeavyIonParameter(hiInfo->Nwounded_N_collisions, xAOD::TruthEvent::NWOUNDEDNCOLLISIONS);
308 xTruthEvent->setHeavyIonParameter(hiInfo->Nwounded_Nwounded_collisions, xAOD::TruthEvent::NWOUNDEDNWOUNDEDCOLLISIONS);
309 xTruthEvent->setHeavyIonParameter((float)hiInfo->impact_parameter, xAOD::TruthEvent::IMPACTPARAMETER);
310 xTruthEvent->setHeavyIonParameter((float)hiInfo->event_plane_angle, xAOD::TruthEvent::EVENTPLANEANGLE);
311 xTruthEvent->setHeavyIonParameter((float)hiInfo->eccentricity, xAOD::TruthEvent::ECCENTRICITY);
312 xTruthEvent->setHeavyIonParameter((float)hiInfo->sigma_inel_NN, xAOD::TruthEvent::SIGMAINELNN);
313#else
314 xTruthEvent->setHeavyIonParameter(hiInfo->Ncoll_hard(), xAOD::TruthEvent::NCOLLHARD);
315 xTruthEvent->setHeavyIonParameter(hiInfo->Npart_proj(), xAOD::TruthEvent::NPARTPROJ);
316 xTruthEvent->setHeavyIonParameter(hiInfo->Npart_targ(), xAOD::TruthEvent::NPARTTARG);
317 xTruthEvent->setHeavyIonParameter(hiInfo->Ncoll(), xAOD::TruthEvent::NCOLL);
318 xTruthEvent->setHeavyIonParameter(hiInfo->spectator_neutrons(), xAOD::TruthEvent::SPECTATORNEUTRONS);
319 xTruthEvent->setHeavyIonParameter(hiInfo->spectator_protons(), xAOD::TruthEvent::SPECTATORPROTONS);
320 xTruthEvent->setHeavyIonParameter(hiInfo->N_Nwounded_collisions(), xAOD::TruthEvent::NNWOUNDEDCOLLISIONS);
321 xTruthEvent->setHeavyIonParameter(hiInfo->Nwounded_N_collisions(), xAOD::TruthEvent::NWOUNDEDNCOLLISIONS);
322 xTruthEvent->setHeavyIonParameter(hiInfo->Nwounded_Nwounded_collisions(), xAOD::TruthEvent::NWOUNDEDNWOUNDEDCOLLISIONS);
323 xTruthEvent->setHeavyIonParameter(hiInfo->impact_parameter(), xAOD::TruthEvent::IMPACTPARAMETER);
324 xTruthEvent->setHeavyIonParameter(hiInfo->event_plane_angle(), xAOD::TruthEvent::EVENTPLANEANGLE);
325 xTruthEvent->setHeavyIonParameter(hiInfo->eccentricity(), xAOD::TruthEvent::ECCENTRICITY);
326 xTruthEvent->setHeavyIonParameter(hiInfo->sigma_inel_NN(), xAOD::TruthEvent::SIGMAINELNN);
327#endif
328 // This doesn't yet exist in our version of HepMC
329 // xTruthEvent->setHeavyIonParameter(hiInfo->centrality(),xAOD::TruthEvent::CENTRALITY);
330 }
331
332 // Parton density info
333 // This will exist 99% of the time, except for e.g. cosmic or particle gun simulation
334 auto const pdfInfo = genEvt->pdf_info();
335 if (pdfInfo) {
336#ifdef HEPMC3
337 xTruthEvent->setPdfInfoParameter(pdfInfo->parton_id[0], xAOD::TruthEvent::PDGID1);
338 xTruthEvent->setPdfInfoParameter(pdfInfo->parton_id[1], xAOD::TruthEvent::PDGID2);
339 xTruthEvent->setPdfInfoParameter(pdfInfo->pdf_id[1], xAOD::TruthEvent::PDFID1);
340 xTruthEvent->setPdfInfoParameter(pdfInfo->pdf_id[1], xAOD::TruthEvent::PDFID2);
341
342 xTruthEvent->setPdfInfoParameter((float)pdfInfo->x[0], xAOD::TruthEvent::X1);
343 xTruthEvent->setPdfInfoParameter((float)pdfInfo->x[1], xAOD::TruthEvent::X2);
344 xTruthEvent->setPdfInfoParameter((float)pdfInfo->scale, xAOD::TruthEvent::Q);
345 xTruthEvent->setPdfInfoParameter((float)pdfInfo->xf[0], xAOD::TruthEvent::XF1);
346 xTruthEvent->setPdfInfoParameter((float)pdfInfo->xf[1], xAOD::TruthEvent::XF2);
347#else
348 xTruthEvent->setPdfInfoParameter(pdfInfo->id1(), xAOD::TruthEvent::PDGID1);
349 xTruthEvent->setPdfInfoParameter(pdfInfo->id2(), xAOD::TruthEvent::PDGID2);
350 xTruthEvent->setPdfInfoParameter(pdfInfo->pdf_id1(), xAOD::TruthEvent::PDFID1);
351 xTruthEvent->setPdfInfoParameter(pdfInfo->pdf_id2(), xAOD::TruthEvent::PDFID2);
352
353 xTruthEvent->setPdfInfoParameter((float)pdfInfo->x1(), xAOD::TruthEvent::X1);
354 xTruthEvent->setPdfInfoParameter((float)pdfInfo->x2(), xAOD::TruthEvent::X2);
355 xTruthEvent->setPdfInfoParameter((float)pdfInfo->scalePDF(), xAOD::TruthEvent::Q);
356 xTruthEvent->setPdfInfoParameter((float)pdfInfo->pdf1(), xAOD::TruthEvent::XF1);
357 xTruthEvent->setPdfInfoParameter((float)pdfInfo->pdf2(), xAOD::TruthEvent::XF2);
358#endif
359 }
360
361 // Handle LHE particles, only supported for HEPMC3
362#ifdef HEPMC3
363 auto lhe_record_attribute = genEvt->attribute<HepMC::ShortEventAttribute>("LHERecord");
364
365 if (lhe_record_attribute && !hadLHERecord && !m_lheTruthParticleContainerKey.empty()){
366 hadLHERecord=true;
367 // The event had an LHE record, so let's record it. This will only happen once per event.
369 ATH_CHECK(xTruthLHEParticleContainer.record(std::make_unique<xAOD::TruthParticleContainer>(),
370 std::make_unique<xAOD::TruthParticleAuxContainer>()));
371 ATH_MSG_DEBUG( "Recorded TruthLHEParticleContainer with key: " << m_lheTruthParticleContainerKey.key() );
372 // The LHE record is stored in a struct with old-style LHE format, so we have to re-encode it
373 for (int nPart=0;nPart<lhe_record_attribute->NUP;++nPart) {
374 // Create TruthParticle
375 xAOD::TruthParticle* xTruthParticle = new xAOD::TruthParticle();
376 // Put particle into container;
377 xTruthLHEParticleContainer->push_back( xTruthParticle );
378 // Copy LHE info into the new particle; good description is in https://arxiv.org/abs/hep-ph/0609017
379 xTruthParticle->setPdgId( lhe_record_attribute->IDUP[nPart] );
380 xTruthParticle->setUid( nPart+1 );
381 xTruthParticle->setStatus( lhe_record_attribute->ISTUP[nPart] );
382 xTruthParticle->setPx( lhe_record_attribute->PUP[nPart][0] );
383 xTruthParticle->setPy( lhe_record_attribute->PUP[nPart][1] );
384 xTruthParticle->setPz( lhe_record_attribute->PUP[nPart][2] );
385 xTruthParticle->setE( lhe_record_attribute->PUP[nPart][3] );
386 xTruthParticle->setM( lhe_record_attribute->PUP[nPart][4] );
387 } // End of loop over particles
388 } // End of if we found the LHE record attribute
389 else if (hadLHERecord){
390 ATH_MSG_WARNING("Truth record appeared to have two LHE records; this should not be possible");
391 }
392#else
393 if (!m_lheTruthParticleContainerKey.empty()){
394 ATH_MSG_WARNING("HEPMC2 does not support LHE truth record storage. Skipping.");
395 }
396#endif
397 }else{//not isSignalProcess
398 xTruthPileupEvent = xTruthPileupEventContainer->push_back( std::make_unique<xAOD::TruthPileupEvent>() );
399 }
400
401 // (2) Build particles and vertices
402 // Map for building associations between particles and vertices
403 // The pair in the map is the (incomingParticles . outgoingParticles) of the given vertex
404 // If signal process vertex is a disconnected vertex (no incoming/outgoing particles), add it manually
405 VertexMap vertexMap;
406 VertexMap::iterator mapItr;
407 vector<HepMC::ConstGenVertexPtr> vertices;
408
409 // Check signal process vertex
410 // If this is a disconnected vertex, add it manually or won't be added from the loop over particles below.
411 auto disconnectedSignalProcessVtx = HepMC::signal_process_vertex(genEvt); // Get the signal process vertex
412 if (disconnectedSignalProcessVtx) {
413 if (disconnectedSignalProcessVtx->particles_in_size() == 0 && disconnectedSignalProcessVtx->particles_out_size() == 0 ) {
414 //This is a disconnected vertex, add it manually
415 vertices.push_back (std::move(disconnectedSignalProcessVtx));
416 }
417 } else {
418 ATH_MSG_WARNING("Signal process vertex pointer not valid in HepMC Collection for GenEvent #" << cntr << " / " << mcColl->size());
419 }
420
421 // Get the beam particles
422 pair<HepMC::ConstGenParticlePtr,HepMC::ConstGenParticlePtr> beamParticles;
423 bool genEvt_valid_beam_particles=false;
424#ifdef HEPMC3
425 auto beamParticles_vec = genEvt->beams();
426 genEvt_valid_beam_particles=(beamParticles_vec.size()>1);
427 if (genEvt_valid_beam_particles){beamParticles.first=beamParticles_vec[0]; beamParticles.second=beamParticles_vec[1]; }
428 // We want to process particles in barcode order.
429 auto bcmapatt = genEvt->attribute<HepMC::GenEventBarcodes>("barcodes"); // FIXME barcode-based
430 if (!bcmapatt) {
431 ATH_MSG_ERROR("TruthParticleCnvTool.cxx: Event does not contain barcodes attribute");
432 return StatusCode::FAILURE;
433 }
434 std::map<int, HepMC3::ConstGenParticlePtr> bcmap = bcmapatt->barcode_to_particle_map();
435 xTruthParticleContainer->reserve(bcmap.size());
436 for (const auto &[genPartBarcode,part]: bcmap) {
437#else
438 genEvt_valid_beam_particles=genEvt->valid_beam_particles();
439 if ( genEvt_valid_beam_particles ) beamParticles = genEvt->beam_particles();
440 xTruthParticleContainer->reserve(genEvt->particles_size());
441 for (auto part: *genEvt) {
442#endif
443 int genPartUniqueID = HepMC::uniqueID(part);
444 // (a) create TruthParticle
445 xAOD::TruthParticle* xTruthParticle = new xAOD::TruthParticle();
446 // (b) Put particle into container;
447 xTruthParticleContainer->push_back( xTruthParticle );
448 fillParticle(xTruthParticle, part); // (c) Copy HepMC info into the new particle
449 // (d) Build Event<->Particle element link
450 const ElementLink<xAOD::TruthParticleContainer> eltp(*xTruthParticleContainer, xTruthParticleContainer->size()-1);
451 if (isSignalProcess) xTruthEvent->addTruthParticleLink(eltp);
452 if (!isSignalProcess) xTruthPileupEvent->addTruthParticleLink(eltp);
453
454 // Create link between HepMC and xAOD truth
455 if (isSignalProcess) truthLinkVec->push_back(new xAODTruthParticleLink(HepMcParticleLink(genPartUniqueID,0,HepMcParticleLink::IS_POSITION, HepMcParticleLink::IS_ID), eltp));
456 if (!isSignalProcess) truthLinkVec->push_back(new xAODTruthParticleLink(HepMcParticleLink(genPartUniqueID,genEvt->event_number(), HepMcParticleLink::IS_EVENTNUM, HepMcParticleLink::IS_ID), eltp));
457
458 // Is this one of the beam particles?
459 if (genEvt_valid_beam_particles) {
460 if (isSignalProcess) {
461 if (part == beamParticles.first) xTruthEvent->setBeamParticle1Link(eltp);
462 if (part == beamParticles.second) xTruthEvent->setBeamParticle2Link(eltp);
463 }
464 }
465 // (e) Particle's production vertex
466 auto productionVertex = part->production_vertex();
467 // Skip the dummy vertex that HepMC3 adds
468 // Can distinguish it from real vertices because it has
469 // a null event pointer.
470 if (productionVertex && productionVertex->parent_event() != nullptr) {
471 VertexParticles& parts = vertexMap[productionVertex];
472 if (parts.incoming.empty() && parts.outgoing.empty())
473 vertices.push_back (std::move(productionVertex));
474 parts.outgoingEL.push_back(eltp);
475 parts.outgoing.push_back(xTruthParticle);
476 }
477 //
478 // else maybe want to keep track that this is the production vertex
479 //
480 // (f) Particle's decay vertex
481 auto decayVertex = part->end_vertex();
482 if (decayVertex) {
483 VertexParticles& parts = vertexMap[decayVertex];
484 if (parts.incoming.empty() && parts.outgoing.empty())
485 vertices.push_back (std::move(decayVertex));
486 parts.incomingEL.push_back(eltp);
487 parts.incoming.push_back(xTruthParticle);
488 }
489
490 } // end of loop over particles
491
492 // (3) Loop over the map
493 auto signalProcessVtx = HepMC::signal_process_vertex(genEvt); // Get the signal process vertex
494 xTruthVertexContainer->reserve(vertices.size());
495 for (const auto& vertex : vertices) {
496 const auto& parts = vertexMap[vertex];
497 // (a) create TruthVertex
498 xAOD::TruthVertex* xTruthVertex = new xAOD::TruthVertex();
499 // (b) Put particle into container (so has store)
500 xTruthVertexContainer->push_back( xTruthVertex );
501 fillVertex(xTruthVertex, vertex); // (c) Copy HepMC info into the new vertex
502 // (d) Build Event<->Vertex element link
503 ElementLink<xAOD::TruthVertexContainer> eltv(*xTruthVertexContainer, xTruthVertexContainer->size()-1);
504 // Mark if this is the signal process vertex
505 if (vertex == signalProcessVtx && isSignalProcess) xTruthEvent->setSignalProcessVertexLink(eltv);
506 if (isSignalProcess) xTruthEvent->addTruthVertexLink(eltv);
507 if (!isSignalProcess) xTruthPileupEvent->addTruthVertexLink(eltv);
508 // (e) Assign incoming particles to the vertex, from the map
509 xTruthVertex->setIncomingParticleLinks( parts.incomingEL );
510 // (f) Assign outgoing particles to the vertex, from the map
511 xTruthVertex->setOutgoingParticleLinks( parts.outgoingEL );
512 // (g) Set Particle<->Vertex links for incoming particles
513 for (xAOD::TruthParticle* p : parts.incoming) p->setDecayVtxLink(eltv);
514 // (h) Set Particle<->Vertex links for incoming particles
515 for (xAOD::TruthParticle* p : parts.outgoing) p->setProdVtxLink(eltv);
516 } //end of loop over vertices
517
518 } // end of loop over McEventCollection
519
520
521 std::stable_sort(truthLinkVec->begin(), truthLinkVec->end(), SortTruthParticleLink());
522 ATH_MSG_VERBOSE("Summarizing truth link size: " << truthLinkVec->size() );
523
524 return StatusCode::SUCCESS;
525 }
526
527 void xAODTruthCnvAlg::handle(const Incident& incident) {
528 if (m_firstBeginRun && incident.type()==IncidentType::BeginRun) {
529 m_firstBeginRun = false;
530 ServiceHandle<StoreGateSvc> inputStore("StoreGateSvc/InputMetaDataStore", name());
531 if(inputStore.retrieve().isFailure()) {
532 ATH_MSG_ERROR("Failed to retrieve Input Metadata Store");
533 return;
534 }
535 const IOVMetaDataContainer* tagInfo{nullptr};
536 if(inputStore->retrieve(tagInfo,"/TagInfo").isFailure()) {
537 ATH_MSG_WARNING("Failed to retrieve /TagInfo metadata from the input store");
538 return;
539 }
540 if(tagInfo->payloadContainer()->size()>0) {
541 CondAttrListCollection* tagInfoPayload = tagInfo->payloadContainer()->at(0);
542 if(tagInfoPayload->size()>0) {
543 const CondAttrListCollection::AttributeList& al = tagInfoPayload->attributeList(0);
544 if (al.exists("lhefGenerator")){
545 m_metaFields.lhefGenerator = al["lhefGenerator"].data<std::string>();
546 }
547
548 if (al.exists("generators")){
549 m_metaFields.generators = al["generators"].data<std::string>();
550 }
551
552 if (al.exists("evgenProcess")){
553 m_metaFields.evgenProcess = al["evgenProcess"].data<std::string>();
554 }
555
556 if (al.exists("evgenTune")){
557 m_metaFields.evgenTune = al["evgenTune"].data<std::string>();
558 }
559
560 if (al.exists("hardPDF")){
561 m_metaFields.hardPDF = al["hardPDF"].data<std::string>();
562 }
563
564 if (al.exists("softPDF")){
565 m_metaFields.softPDF = al["softPDF"].data<std::string>();
566 }
567 }
568 }
569 else {
570 ATH_MSG_WARNING("Empty Tag Info metadata!");
571 }
572 }
573 }
574
575
576 // A helper to set up a TruthVertex (without filling the ELs)
578 // id was renamed to status in HepMC3.
579 tv->setStatus(HepMC::status(gv)); // For now convert the status back to the old scheme
580 tv->setUid(HepMC::uniqueID(gv)); // FIXME set xAOD::TruthVertex function name appropriately
581 tv->setX(gv->position().x());
582 tv->setY(gv->position().y());
583 tv->setZ(gv->position().z());
584 tv->setT(gv->position().t());
585 }
586
587
588 // A helper to set up a TruthParticle (without filling the ELs)
590 tp->setPdgId(gp->pdg_id());
591 tp->setUid(HepMC::uniqueID(gp)); // FIXME set xAOD::TruthParticle function name appropriately
592 tp->setStatus(HepMC::status(gp)); // For now convert the status back to the old scheme
593
594 auto pol = HepMC::polarization(gp);
595 if (pol.is_defined()) {
596 tp->setPolarizationParameter(pol.theta(), xAOD::TruthParticle::polarizationTheta);
597 tp->setPolarizationParameter(pol.phi(), xAOD::TruthParticle::polarizationPhi);
598 }
599
600 tp->setM(gp->generated_mass());
601 tp->setPx(gp->momentum().px());
602 tp->setPy(gp->momentum().py());
603 tp->setPz(gp->momentum().pz());
604 tp->setE(gp->momentum().e());
605 }
606
607
608 StatusCode
610 const std::string& metaName)
611 {
612 ATH_CHECK( metaStore.retrieve() );
613
614 auto md = std::make_unique<xAOD::TruthMetaDataContainer>();
615 m_tmd = md.get();
616
617 auto aux = std::make_unique<xAOD::TruthMetaDataAuxContainer>();
618 md->setStore( aux.get() );
619
620 // Record the trigger configuration metadata into it:
621 CHECK( metaStore->record( std::move (aux), metaName + "Aux." ) );
622 CHECK( metaStore->record( std::move (md), metaName ) );
623 return StatusCode::SUCCESS;
624 }
625
626
627 StatusCode
629 const HepMC::GenEvent& genEvt,
630 const MetadataFields& metaFields)
631 {
632 // This bit needs to be serialized.
633 lock_t lock (m_mutex);
634
635 //Inserting in a (unordered_)set returns an <iterator, boolean> pair, where the boolean
636 //is used to check if the key already exists (returns false in the case it exists)
637 if( m_existingMetaDataChan.insert(mcChannelNumber).second ) {
638 m_tmd->push_back (std::make_unique <xAOD::TruthMetaData>());
639 xAOD::TruthMetaData* md = m_tmd->back();
640
641#ifdef HEPMC3
643 md->setMcChannelNumber(mcChannelNumber);
644 std::vector<std::string> orderedWeightNameVec;
645 if (!genEvt.run_info()) {
646 for (size_t i=0; i<genEvt.weights().size();i++) orderedWeightNameVec.push_back(std::to_string(i));
647 } else {
648 if (!genEvt.run_info()->weight_names().empty()) {
649 orderedWeightNameVec = genEvt.weight_names();
650 } else {
651 //AV This to be decided. It is always a good idea to have a default weight 1.0.
652 //orderedWeightNameVec.push_back("0");
653 }
654 }
655 md->setWeightNames(orderedWeightNameVec);
656#else
657 // FIXME: class member protection violation here.
658 // This appears to be because WeightContainer has no public methods
659 // to get information about the weight names.
660 const auto& weightNameMap = genEvt.weights().m_names;
661 std::vector<std::string> orderedWeightNameVec;
662 orderedWeightNameVec.reserve( weightNameMap.size() );
663 for (const auto& entry: weightNameMap) {
664 orderedWeightNameVec.push_back(entry.first);
665 }
666
667 //The map from the HepMC record pairs the weight names with a corresponding index,
668 //it is not guaranteed that the indices are ascending when iterating over the map
669 std::sort(orderedWeightNameVec.begin(), orderedWeightNameVec.end(),
670 [&](const std::string& i, const std::string& j){return weightNameMap.at(i) < weightNameMap.at(j);});
671
672 md->setMcChannelNumber(mcChannelNumber);
673 md->setWeightNames( orderedWeightNameVec );
674#endif
675
676 if(!metaFields.lhefGenerator.empty()) {
677 md->setLhefGenerator(metaFields.lhefGenerator);
678 }
679 if(!metaFields.generators.empty()) {
680 md->setGenerators(metaFields.generators);
681 }
682 if(!metaFields.evgenProcess.empty()) {
683 md->setEvgenProcess(metaFields.evgenProcess);
684 }
685 if(!metaFields.evgenTune.empty()) {
686 md->setEvgenTune(metaFields.evgenTune);
687 }
688 if(!metaFields.hardPDF.empty()) {
689 md->setHardPDF(metaFields.hardPDF);
690 }
691 if(!metaFields.softPDF.empty()) {
692 md->setSoftPDF(metaFields.softPDF);
693 }
694 }
695
696 return StatusCode::SUCCESS;
697 }
698
699
700
701} // namespace xAODMaker
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
This file defines the class for a collection of AttributeLists where each one is associated with a ch...
Helpers for checking error return status codes and reporting errors.
#define CHECK(...)
Evaluate an expression and check for errors.
This class is a container for conditions data.
static Double_t al
Handle class for reading from StoreGate.
Handle class for recording to StoreGate.
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
An algorithm that can be simultaneously executed in multiple threads.
This class is a collection of AttributeLists where each one is associated with a channel number.
const AttributeList & attributeList(ChanNum chanNum) const
attribute list for a given channel number
size_type size() const
number of Chan/AttributeList pairs
coral::AttributeList AttributeList
This class is a container for conditions data.
const IOVPayloadContainer * payloadContainer() const
Access to payload container.
size_type size() const
size of payload vector
CondAttrListCollection * at(unsigned int i) const
Element access.
virtual bool isValid() override final
Can the handle be successfully dereferenced?
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
std::mutex m_mutex
Mutex to control access to meta data writing.
StatusCode maybeWrite(uint32_t mcChannelNumber, const HepMC::GenEvent &genEvt, const MetadataFields &metaFields)
std::unordered_set< uint32_t > m_existingMetaDataChan
Set for tracking the mc channels for which we already added meta data.
StatusCode initialize(ServiceHandle< StoreGateSvc > &metaStore, const std::string &metaName)
xAOD::TruthMetaDataContainer * m_tmd
The meta data container to be written out.
SG::WriteHandleKey< xAOD::TruthVertexContainer > m_xaodTruthVertexContainerKey
Gaudi::Property< bool > m_doAllPileUp
Pile-up options.
SG::WriteHandleKey< xAOD::TruthParticleContainer > m_lheTruthParticleContainerKey
SG::WriteHandleKey< xAODTruthParticleLinkVector > m_truthLinkContainerKey
SG::ReadHandleKey< McEventCollection > m_aodContainerKey
The key of the input AOD truth container.
static void fillVertex(xAOD::TruthVertex *tv, const HepMC::ConstGenVertexPtr &gv)
These functions do not set up ELs, just the other variables.
xAODTruthCnvAlg(const std::string &name, ISvcLocator *svcLoc)
Regular algorithm constructor.
SG::WriteHandleKey< xAOD::TruthParticleContainer > m_xaodTruthParticleContainerKey
Gaudi::Property< bool > m_writeMetaData
option to disable writing of metadata (e.g. if running a filter on xAOD in generators)
virtual StatusCode execute(const EventContext &ctx) const override
Function executing the algorithm.
std::map< HepMC::ConstGenVertexPtr, VertexParticles > VertexMap
Convenience handle for a map of vtx ptrs -> connected particles.
ServiceHandle< StoreGateSvc > m_metaStore
Connection to the metadata store.
SG::ReadHandleKey< xAOD::EventInfo > m_evtInfo
Event Info.
Gaudi::Property< bool > m_doInTimePileUp
virtual void handle(const Incident &incident) override
Incident handler.
static void fillParticle(xAOD::TruthParticle *tp, const HepMC::ConstGenParticlePtr &gp)
std::string m_metaName
SG key and name for meta data.
SG::WriteHandleKey< xAOD::TruthEventContainer > m_xaodTruthEventContainerKey
The key for the output xAOD truth containers.
SG::WriteHandleKey< xAOD::TruthPileupEventContainer > m_xaodTruthPUEventContainerKey
virtual StatusCode initialize() override
Function initialising the algorithm.
void addTruthVertexLink(const TruthVertexLink_t &vlink)
Add one truth vertex.
void addTruthParticleLink(const TruthParticleLink_t &plink)
Add one truth particle.
void setWeights(const std::vector< float > &weights)
Set the event weights.
void setSignalProcessVertexLink(const TruthVertexLink_t &link)
Set pointer to a vertex representing the primary beam interaction point.
void setCrossSectionError(float value)
Set the cross-section error.
void setCrossSection(float value)
Set the cross-section.
bool setHeavyIonParameter(int value, HIParam parameter)
Set an integer HI parameter.
bool setPdfInfoParameter(int value, PdfParam parameter)
Set an integer PDF info parameter.
@ NWOUNDEDNWOUNDEDCOLLISIONS
[int]
void setBeamParticle1Link(const TruthParticleLink_t &pcl1)
Set one incoming beam particle.
void setBeamParticle2Link(const TruthParticleLink_t &pcl2)
Set one incoming beam particle.
void setGenerators(const std::string &value)
void setLhefGenerator(const std::string &value)
void setEvgenTune(const std::string &value)
void setHardPDF(const std::string &value)
void setMcChannelNumber(uint32_t value)
void setSoftPDF(const std::string &value)
void setWeightNames(const std::vector< std::string > &value)
void setEvgenProcess(const std::string &value)
void setPy(float value)
Set the y component of the particle's momentum.
void setUid(int value)
Set unique ID.
void setM(float value)
Also store the mass.
void setE(float value)
Set the energy of the particle.
void setPz(float value)
Set the z component of the particle's momentum.
void setStatus(int value)
Set status code.
void setPdgId(int pid)
Set PDG ID code.
void setPx(float value)
Set the x component of the particle's momentum.
@ polarizationPhi
Polarization in ( )
@ polarizationTheta
Polarization in ( )
void setStatus(int value)
Set the vertex status.
void setZ(float value)
Set the vertex's longitudinal distance from the origin.
void setUid(int value)
Set the vertex unique ID.
void setOutgoingParticleLinks(const TPLinks_t &links)
Set all the outgoing particles.
void setIncomingParticleLinks(const TPLinks_t &links)
Set all the incoming particles.
void setT(float value)
Set the vertex time.
void setX(float value)
Set the x displacement of the vertex.
void setY(float value)
Set the y displacement of the vertex.
int signal_process_id(const GenEvent &e)
Definition GenEvent.h:637
int uniqueID(const T &p)
int status(const T &p)
Polarization polarization(const T &a)
const GenParticle * ConstGenParticlePtr
Definition GenParticle.h:38
GenVertex * signal_process_vertex(const GenEvent *e)
Definition GenEvent.h:627
const HepMC::GenVertex * ConstGenVertexPtr
Definition GenVertex.h:60
STL namespace.
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
void stable_sort(DataModel_detail::iterator< DVL > beg, DataModel_detail::iterator< DVL > end)
Specialization of stable_sort for DataVector/List.
TruthVertex_v1 TruthVertex
Typedef to implementation.
Definition TruthVertex.h:15
TruthEvent_v1 TruthEvent
Typedef to implementation.
Definition TruthEvent.h:17
TruthParticle_v1 TruthParticle
Typedef to implementation.
TruthPileupEvent_v1 TruthPileupEvent
Typedef to implementation.
TruthMetaData_v1 TruthMetaData
Typedef to implementation.
Type for tracking particles connected to a single vertex.
bool isSeparatorGenEvent(const HepMC::GenEvent *genEvt)