ATLAS Offline Software
Loading...
Searching...
No Matches
Geant4TruthIncident.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// class header
7
8// package includes
10
11// Atlas G4 Helpers
13#include "MCTruth/TrackHelper.h"
16
17// Units
18#include "GaudiKernel/PhysicalConstants.h"
19
20// HepMC includes
23
24// Geant4 includes
25#include "G4Step.hh"
26#include "G4Track.hh"
27#include "G4VProcess.hh"
28
29#include "G4TrackStatus.hh"
30#include "G4ProcessType.hh"
31#include "G4EmProcessSubType.hh"
32#include "G4DynamicParticle.hh"
33#include "G4PrimaryParticle.hh"
34
35#include "G4EventManager.hh"
36#include "G4Event.hh"
37
38// ISF includes
40#include <limits>
41/*
42 Comments:
43 what about parent particle surviving (e.g. bremstrahlung)
44
45
46 Existing code:
47 Simulation/G4Sim/SimHelpers/src/StepHelper.cxx
48 - provids convenient information of a G4Step
49 Simulation/G4Utilities/G4TruthStrategies/src/BremsstrahlungStrategy.cxx
50 - retrieves information from a G4Step (via StepHelper)
51 Simulation/G4Sim/MCTruth/MCTruth/TruthStrategy.h
52 - common base for different truth strategies
53 Simulation/G4Sim/MCTruth/src/AtlasG4EventUserInfo.cxx
54 - stores HepMCevent in G4
55 Simulation/G4Sim/MCTruth/src/TrackInformation.cxx
56 Simulation/G4Sim/MCTruth/src/TrackHelper.cxx
57 - store/manage barcode
58*/
59
60
62 const ISF::ISFParticle& baseISP,
64 AtlasG4EventUserInfo *atlasG4EvtUserInfo) :
65 ITruthIncident(geoID, step->GetSecondaryInCurrentStep()->size()), // switch to G4Step::GetNumberOfSecondariesInCurrentStep() once we're using G4 10.2 or later
66 m_positionSet(false),
67 m_position(),
68 m_step(step),
69 m_baseISP(baseISP),
70 m_atlasG4EvtUserInfo(atlasG4EvtUserInfo)
71{
72 // prepare children:
74
75 // calculate position:
76 const G4StepPoint *postStepPoint = m_step->GetPostStepPoint();
77 const G4ThreeVector &pos = postStepPoint->GetPosition();
78 const G4double time = postStepPoint->GetGlobalTime()*Gaudi::Units::c_light;
79 m_position.set( pos.x(), pos.y(), pos.z(), time );
80}
81
82const HepMC::FourVector& iGeant4::Geant4TruthIncident::position() const {
83 return m_position;
84}
85
87 const G4VProcess *process = m_step->GetPostStepPoint()->GetProcessDefinedStep();
88 // TODO: need to check that G4ProcessSubTypes Match int
89 return process->GetProcessType();
90}
91
93 const G4VProcess *process = m_step->GetPostStepPoint()->GetProcessDefinedStep();
94 // TODO: need to check that G4ProcessSubTypes Match int
95 return process->GetProcessSubType();
96}
97
99 const G4ThreeVector & mom= m_step->GetPreStepPoint()->GetMomentum();
100 return mom.mag2();
101}
102
104 const G4ThreeVector & mom= m_step->GetPreStepPoint()->GetMomentum();
105 return mom.perp2();
106}
107
109 return (m_step->GetPreStepPoint()->GetKineticEnergy()/CLHEP::MeV);
110}
111
113 return m_step->GetTrack()->GetDefinition()->GetPDGEncoding();
114}
115
116int iGeant4::Geant4TruthIncident::parentBarcode() { // TODO Remove this method
117 auto parent = parentParticle();
118
119 return (parent) ? HepMC::barcode(parent) : HepMC::UNDEFINED_ID;
120}
121
123 auto parent = parentParticle();
124
125 return (parent) ? HepMC::uniqueID(parent) : HepMC::UNDEFINED_ID;
126}
127
129 auto parent = parentParticle();
130
131 return (parent) ? parent->status() : 0;
132}
133
135 return m_atlasG4EvtUserInfo->GetCurrentGenParticle();
136}
137
139 const G4Track *track = m_step->GetTrack();
140
141 // check if particle is a alive in G4 or in ISF
142 if ( particleAlive(track) ) {
143 return true;
144 } else {
145 return false;
146 }
147}
148
150 const G4Track *track = m_step->GetTrack();
151
152 // check if particle is a alive in G4 or in ISF
153 if ( !parentSurvivesIncident() ) {
154 return nullptr;
155 }
156
157 // internally cache the HepMC representation of the parent particle
158 // (to allow multiple calls to this method on the same instance)
160 // create new HepMC particle, using momentum and energy
161 // from G4DynamicParticle (which should be equivalent to postStep)
162 m_parentParticleAfterIncident = convert(track, newBarcode, false);
163
164 m_atlasG4EvtUserInfo->SetCurrentGenParticle( m_parentParticleAfterIncident );
165
166 // store (new) hepmc particle in track's UserInformation
167 TrackHelper tHelper(track);
168 TrackInformation *tInfo = tHelper.GetTrackInformation();
169 if (tInfo) {
170 // do NOT update the TrackInformation for regenerated particles!
171 // (most recent truth info is kept in AtlasG4EventUserInfo)
172 //tInfo->SetCurrentGenParticle( m_parentParticleAfterIncident );
173 int regenerationNr = tInfo->GetRegenerationNr();
174 regenerationNr++;
175 tInfo->SetRegenerationNr(regenerationNr);
176 if ( tHelper.IsPrimary() ) {
178 }
179 }
180
181 } // <-- end if m_parentParticleAfterIncident is not filled
182
183
184 // return the HepMC particle representation of the parent
186}
187
188double iGeant4::Geant4TruthIncident::childP2(unsigned short i) const {
189 const G4ThreeVector & mom= m_children[i]->GetMomentum();
190 return mom.mag2();
191}
192
193const G4ThreeVector iGeant4::Geant4TruthIncident::childP(unsigned short i) const {
194 const G4ThreeVector & mom= m_children[i]->GetMomentum();
195 return mom;
196}
197
198double iGeant4::Geant4TruthIncident::childPt2(unsigned short i) const {
199 const G4ThreeVector & mom= m_children[i]->GetMomentum();
200 return mom.perp2();
201}
202
203double iGeant4::Geant4TruthIncident::childEkin(unsigned short i) const {
204 return (m_children[i]->GetKineticEnergy()/CLHEP::MeV);
205}
206
208 return m_children[i]->GetDefinition()->GetPDGEncoding();
209}
210
212 // the G4Track instance for the current child particle
213 const G4Track* track = m_children[index];
214 // This should be a *secondary* track. If it has a primary, it was a decay and
215 // we are running with quasi-stable particle simulation. Note that if the primary
216 // track is passed in as a secondary that survived the interaction, then this was
217 // *not* a decay and we should not treat it in this way
218 if (track->GetDynamicParticle() &&
219 track->GetDynamicParticle()->GetPrimaryParticle() &&
220 track->GetDynamicParticle()->GetPrimaryParticle()->GetUserInformation()){
221 // Then the new particle should use the same barcode as the old one!!
222 PrimaryParticleInformation* primaryPartInfo = dynamic_cast<PrimaryParticleInformation*>( track->GetDynamicParticle()->GetPrimaryParticle()->GetUserInformation() );
223 return primaryPartInfo->GetParticleBarcode();
224 }
225 return 0;
226}
227
229 int newBarcode) {
230 // the G4Track instance for the current child particle
231 const G4Track* thisChildTrack = m_children[i];
232
233 // NB: NOT checking if secondary is actually alive. Even with zero momentum,
234 // secondary could decay right away and create further particles which pass the
235 // truth strategies.
236
237 HepMC::GenParticlePtr hepParticle = convert( thisChildTrack, newBarcode, true );
238 TrackHelper tHelper(thisChildTrack);
240
241 // needed to make AtlasG4 work with ISF TruthService
242 if (trackInfo==nullptr) {
243 trackInfo = new TrackInformation( hepParticle );
244 thisChildTrack->SetUserInformation( trackInfo );
245 }
246
247 trackInfo->SetCurrentGenParticle(hepParticle);
249 trackInfo->SetRegenerationNr(0);
250
251 return hepParticle;
252}
253
254
255bool iGeant4::Geant4TruthIncident::particleAlive(const G4Track *track) const {
256 G4TrackStatus trackStatus = track->GetTrackStatus();
257
258 if ( trackStatus != fAlive && trackStatus != fStopButAlive ) {
259 // parent does not exist in G4 anymore after this step
260
261 // check whether the particle was returned to ISF
262 auto* trackInfo = ISFG4Helper::getISFTrackInfo( *track );
263 bool returnedToISF = trackInfo ? trackInfo->GetReturnedToISF() : false;
264 if ( !returnedToISF ) {
265 // particle was not sent to ISF either
266 return false;
267 }
268 }
269
270 return true;
271}
272
273#ifdef HEPMC3
274HepMC::GenParticlePtr iGeant4::Geant4TruthIncident::convert(const G4Track *track, const int, const bool secondary) const {
275#else
276HepMC::GenParticlePtr iGeant4::Geant4TruthIncident::convert(const G4Track *track, const int barcode, const bool secondary) const {
277#endif
278
279 const G4ThreeVector & mom = track->GetMomentum();
280 const double energy = track->GetTotalEnergy();
281 const int pdgCode = track->GetDefinition()->GetPDGEncoding();
282 const HepMC::FourVector fourMomentum( mom.x(), mom.y(), mom.z(), energy);
283
284 const HepMC::GenParticlePtr parent = m_atlasG4EvtUserInfo->GetCurrentGenParticle();
285 int status = (secondary) ? 1 + HepMC::SIM_STATUS_THRESHOLD : parent->status() + HepMC::SIM_STATUS_INCREMENT;
286 // Treat child particles of pre-defined decays differently
288 const G4DynamicParticle* dynPart = track->GetDynamicParticle();
289 bool hasPredefinedDecay = (dynPart && (nullptr!=(dynPart->GetPreAssignedDecayProducts())));
290 status = (hasPredefinedDecay)? 2 : 1;
291 }
292 // Treat particles with a pre-defined decay which have survived an
293 // interaction differently for now.
295 status = status%HepMC::SIM_STATUS_INCREMENT;
296 // Such particles were previously assigned barcodes below
297 // HepMC::SIM_BARCODE_OFFSET and therefore
298 // HepMC::is_simulation_particle(barcode) would return false. Not
299 // applying HepMC::SIM_STATUS_INCREMENT means that
300 // HepMC::is_simulation_particle(status) will also return false.
301 // TODO This should be revisited in the future should we decide
302 // that we want to be able to easily identify such particles.
303 }
304 HepMC::GenParticlePtr newParticle = HepMC::newGenParticlePtr(fourMomentum, pdgCode, status);
305
306#ifndef HEPMC3
307 // This should be a *secondary* track. If it has a primary, it was a decay and
308 // we are running with quasi-stable particle simulation. Note that if the primary
309 // track is passed in as a secondary that survived the interaction, then this was
310 // *not* a decay and we should not treat it in this way
311 if (secondary &&
312 track->GetDynamicParticle() &&
313 track->GetDynamicParticle()->GetPrimaryParticle() &&
314 track->GetDynamicParticle()->GetPrimaryParticle()->GetUserInformation()){
315 // Then the new particle should use the same barcode as the old one!!
316 PrimaryParticleInformation* primaryPartInfo = dynamic_cast<PrimaryParticleInformation*>( track->GetDynamicParticle()->GetPrimaryParticle()->GetUserInformation() );
317 HepMC::suggest_barcode( newParticle, primaryPartInfo->GetParticleBarcode() );
318 } else {
319 HepMC::suggest_barcode( newParticle, barcode );
320 }
321#endif
322
323 return newParticle;
324}
325
326
328
329 const std::vector<const G4Track*>* tracks = m_step->GetSecondaryInCurrentStep();
330 const int iSize = tracks->size();
331 const int iLast = iSize - numberOfChildren() - 1; //NB can be -1.
332 for(int i=iSize-1;i>iLast;i--) {
333 m_children.push_back((*tracks)[i]);
334 }
335}
336
345 G4Track* track=m_step->GetTrack();
346 const G4DynamicParticle* dynPart = track->GetDynamicParticle();
347 const double preAsgnDecPropTime = dynPart->GetPreAssignedDecayProperTime();
348 const double trackPropTime = track->GetProperTime();
349 bool parentIsQuasiStable = (nullptr!=(dynPart->GetPreAssignedDecayProducts()));
350 const G4VProcess *process = m_step->GetPostStepPoint()->GetProcessDefinedStep();
351 const int processType = process->GetProcessType();
353 if(parentIsQuasiStable) {
354 if(this->parentSurvivesIncident()) {
355 classification = ISF::QS_SURV_VTX;
356 }
357 else if(processType==6 && (trackPropTime-preAsgnDecPropTime) < -std::numeric_limits<double>::epsilon()) {
358 // Particle decayed before its expected pre-defined decay time
359 classification = ISF::QS_DEST_VTX;
360 }
361 else {
362 classification = ISF::QS_PREDEF_VTX;
363 }
364 }
365 return classification;
366 // end of code to determine interaction type
367}
This class is attached to G4Event objects as UserInformation.
The generic ISF particle definition,.
Definition ISFParticle.h:42
unsigned short numberOfChildren() const
Return total number of child particles.
AtlasDetDescr::AtlasRegion geoID()
Return the SimGeoID corresponding to the vertex of the truth incident.
ITruthIncident(AtlasDetDescr::AtlasRegion geoID, unsigned short numChildren)
This class is attached to G4PrimaryParticle objects as UserInformation.
TrackInformation * GetTrackInformation()
Definition TrackHelper.h:28
bool IsPrimary() const
Implementation of VTrackInformation.
int GetRegenerationNr() const
return the number of times the particle represented by the G4Track has undergone a non-destructive in...
void SetRegenerationNr(int i)
update the number of times the particle represented by the G4Track has undergone a non-destructive in...
void SetClassification(TrackClassification tc)
update the classification of the currently tracked particle, usually called when a new G4Track is cre...
int physicsProcessCode() const override final
Return specific physics process code of the truth incident (eg ionisation, bremsstrahlung,...
double childEkin(unsigned short index) const override final
Return Ekin of the i-th child particle.
const ISF::ISFParticle & m_baseISP
int parentBarcode() override final
Return the barcode of the parent particle.
HepMC::GenParticlePtr m_parentParticleAfterIncident
double parentPt2() const override final
Return pT^2 of the parent particle.
bool particleAlive(const G4Track *track) const
check if the given G4Track represents a particle that is alive in ISF or ISF-G4
int childBarcode(unsigned short index) const override final
Return the barcode of the i-th child particle (if defined as part of the TruthIncident) otherwise ret...
double childP2(unsigned short index) const override final
Return p^2 of the i-th child particle.
int parentStatus() override final
Return the status of the parent particle.
HepMC::GenParticlePtr parentParticleAfterIncident(int newBC) override final
Return the parent particle after the TruthIncident vertex (and give it a new barcode)
ISF::InteractionClass_t interactionClassification() const override final
The interaction classifications are described as follows: STD_VTX: interaction of a particle without ...
int physicsProcessCategory() const override final
Return category of the physics process represented by the truth incident (eg hadronic,...
const HepMC::FourVector & position() const override final
Return HepMC position of the truth vertex.
bool parentSurvivesIncident() const override final
Return a boolean whether or not the parent particle survives the incident.
int parentPdgCode() const override final
Return the PDG Code of the parent particle.
HepMC::GenParticlePtr childParticle(unsigned short index, int bc) override final
Return the i-th child as a HepMC particle type and assign the given Barcode to the simulator particle...
int childPdgCode(unsigned short index) const override final
Return the PDG Code of the i-th child particle.
int parentUniqueID() override final
Return the unique ID of the parent particle.
void prepareChildren()
prepare the child particles
double parentP2() const override final
Return p^2 of the parent particle.
std::vector< const G4Track * > m_children
double childPt2(unsigned short index) const override final
Return pT^2 of the i-th child particle.
double parentEkin() const override final
Return Ekin of the parent particle.
HepMC::GenParticlePtr parentParticle() override final
Return the parent particle as a HepMC particle type.
const G4ThreeVector childP(unsigned short index) const
Return p of the i-th child particle.
HepMC::GenParticlePtr convert(const G4Track *particle, const int barcode, const bool secondary) const
static VTrackInformation * getISFTrackInfo(const G4Track &aTrack)
return a valid UserInformation object of the G4Track for use within the ISF
const std::string process
AtlasRegion
A simple enum of ATLAS regions and sub-detectors.
Definition AtlasRegion.h:21
int barcode(const T *p)
Definition Barcode.h:16
int uniqueID(const T &p)
constexpr int UNDEFINED_ID
constexpr int SIM_STATUS_INCREMENT
Constant defining the barcode threshold for regenerated particles, i.e. particles surviving an intera...
constexpr int SIM_STATUS_THRESHOLD
Constant definiting the status threshold for simulated particles, eg. can be used to separate generat...
bool suggest_barcode(T &p, int i)
Definition GenEvent.h:672
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
GenParticle * GenParticlePtr
Definition GenParticle.h:37
InteractionClass_t
The interaction classifications are described as follows: STD_VTX: interaction of a particle without ...
@ QS_SURV_VTX
@ QS_PREDEF_VTX
@ QS_DEST_VTX
Definition index.py:1