ATLAS Offline Software
Geant4TruthIncident.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2024 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
21 #include "AtlasHepMC/GenParticle.h"
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
39 #include "ISF_Event/ISFParticle.h"
40 
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  - retrives 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 
82 const 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 
116 int iGeant4::Geant4TruthIncident::parentBarcode() { // TODO Remove this method
117  auto parent = parentParticle();
118 
120 }
121 
123  auto parent = parentParticle();
124 
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)
159  if ( !m_parentParticleAfterIncident ) {
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
185  return m_parentParticleAfterIncident;
186 }
187 
188 double iGeant4::Geant4TruthIncident::childP2(unsigned short i) const {
189  const G4ThreeVector & mom= m_children[i]->GetMomentum();
190  return mom.mag2();
191 }
192 
193 const G4ThreeVector iGeant4::Geant4TruthIncident::childP(unsigned short i) const {
194  const G4ThreeVector & mom= m_children[i]->GetMomentum();
195  return mom;
196 }
197 
198 double iGeant4::Geant4TruthIncident::childPt2(unsigned short i) const {
199  const G4ThreeVector & mom= m_children[i]->GetMomentum();
200  return mom.perp2();
201 }
202 
203 double 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 
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
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
274 HepMC::GenParticlePtr iGeant4::Geant4TruthIncident::convert(const G4Track *track, const int, const bool secondary) const {
275 #else
276 HepMC::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
287  if (this->interactionClassification() == ISF::QS_PREDEF_VTX) {
288  const G4DynamicParticle* dynPart = track->GetDynamicParticle();
289  bool hasPredefinedDecay = (dynPart && (nullptr!=(dynPart->GetPreAssignedDecayProducts())));
290  status = (hasPredefinedDecay)? 2 : 1;
291  }
292  HepMC::GenParticlePtr newParticle = HepMC::newGenParticlePtr(fourMomentum, pdgCode, status);
293 
294 #ifndef HEPMC3
295  // This should be a *secondary* track. If it has a primary, it was a decay and
296  // we are running with quasi-stable particle simulation. Note that if the primary
297  // track is passed in as a secondary that survived the interaction, then this was
298  // *not* a decay and we should not treat it in this way
299  if (secondary &&
300  track->GetDynamicParticle() &&
301  track->GetDynamicParticle()->GetPrimaryParticle() &&
302  track->GetDynamicParticle()->GetPrimaryParticle()->GetUserInformation()){
303  // Then the new particle should use the same barcode as the old one!!
304  PrimaryParticleInformation* primaryPartInfo = dynamic_cast<PrimaryParticleInformation*>( track->GetDynamicParticle()->GetPrimaryParticle()->GetUserInformation() );
305  HepMC::suggest_barcode( newParticle, primaryPartInfo->GetParticleBarcode() );
306  } else {
307  HepMC::suggest_barcode( newParticle, barcode );
308  }
309 #endif
310 
311  return newParticle;
312 }
313 
314 
316 
317  const std::vector<const G4Track*>* tracks = m_step->GetSecondaryInCurrentStep();
318  const int iSize = tracks->size();
319  const int iLast = iSize - numberOfChildren() - 1; //NB can be -1.
320  for(int i=iSize-1;i>iLast;i--) {
321  m_children.push_back((*tracks)[i]);
322  }
323 }
324 
333  G4Track* track=m_step->GetTrack();
334  const G4DynamicParticle* dynPart = track->GetDynamicParticle();
335  bool parentIsQuasiStable = (nullptr!=(dynPart->GetPreAssignedDecayProducts()));
336  const G4VProcess *process = m_step->GetPostStepPoint()->GetProcessDefinedStep();
337  const int processType = process->GetProcessType();
338  const int processSubType = process->GetProcessSubType();
339  ISF::InteractionClass_t classification(ISF::STD_VTX);
340  if(parentIsQuasiStable) {
341  if(this->parentSurvivesIncident()) {
342  classification = ISF::QS_SURV_VTX;
343  }
344  else if(processType==6 && processSubType==201) {
345  classification = ISF::QS_PREDEF_VTX;
346  }
347  else {
348  classification = ISF::QS_DEST_VTX;
349  }
350  }
351  return classification;
352  // end of code to determine interaction type
353 }
iGeant4::Geant4TruthIncident::parentP2
double parentP2() const override final
Return p^2 of the parent particle.
Definition: Geant4TruthIncident.cxx:98
iGeant4::Geant4TruthIncident::parentBarcode
int parentBarcode() override final
Return the barcode of the parent particle.
Definition: Geant4TruthIncident.cxx:116
iGeant4::Geant4TruthIncident::parentParticleAfterIncident
HepMC::GenParticlePtr parentParticleAfterIncident(int newBC) override final
Return the parent particle after the TruthIncident vertex (and give it a new barcode)
Definition: Geant4TruthIncident.cxx:149
iGeant4::Geant4TruthIncident::parentParticle
HepMC::GenParticlePtr parentParticle() override final
Return the parent particle as a HepMC particle type.
Definition: Geant4TruthIncident.cxx:134
HepMC::suggest_barcode
bool suggest_barcode(T &p, int i)
Definition: GenEvent.h:548
AtlasG4EventUserInfo
This class is attached to G4Event objects as UserInformation. It holds a pointer to the HepMC::GenEve...
Definition: AtlasG4EventUserInfo.h:21
iGeant4::Geant4TruthIncident::parentStatus
int parentStatus() override final
Return the status of the parent particle.
Definition: Geant4TruthIncident.cxx:128
iGeant4::Geant4TruthIncident::parentUniqueID
int parentUniqueID() override final
Return the unique ID of the parent particle.
Definition: Geant4TruthIncident.cxx:122
ISF::STD_VTX
@ STD_VTX
Definition: ITruthIncident.h:28
index
Definition: index.py:1
HepMC::SIM_STATUS_INCREMENT
constexpr int SIM_STATUS_INCREMENT
Constant defining the barcode threshold for regenerated particles, i.e. particles surviving an intera...
Definition: MagicNumbers.h:35
TrackHelper.h
AtlasDetDescr::AtlasRegion
AtlasRegion
Definition: AtlasRegion.h:27
ISF::QS_PREDEF_VTX
@ QS_PREDEF_VTX
Definition: ITruthIncident.h:31
HepMC::GenParticlePtr
GenParticle * GenParticlePtr
Definition: GenParticle.h:37
python.SystemOfUnits.MeV
int MeV
Definition: SystemOfUnits.py:154
ISF::ISFParticle
Definition: ISFParticle.h:42
iGeant4::Geant4TruthIncident::parentEkin
double parentEkin() const override final
Return Ekin of the parent particle.
Definition: Geant4TruthIncident.cxx:108
iGeant4::Geant4TruthIncident::convert
HepMC::GenParticlePtr convert(const G4Track *particle, const int barcode, const bool secondary) const
Definition: Geant4TruthIncident.cxx:276
GenParticle.h
SUSY_SimplifiedModel_PostInclude.process
string process
Definition: SUSY_SimplifiedModel_PostInclude.py:42
iGeant4::Geant4TruthIncident::m_step
const G4Step * m_step
Definition: Geant4TruthIncident.h:122
iGeant4::Geant4TruthIncident::m_position
HepMC::FourVector m_position
Definition: Geant4TruthIncident.h:121
PrimaryParticleInformation::GetParticleBarcode
int GetParticleBarcode() const
Definition: PrimaryParticleInformation.cxx:18
iGeant4::Geant4TruthIncident::childPt2
double childPt2(unsigned short index) const override final
Return pT^2 of the i-th child particle.
Definition: Geant4TruthIncident.cxx:198
iGeant4::Geant4TruthIncident::physicsProcessCode
int physicsProcessCode() const override final
Return specific physics process code of the truth incident (eg ionisation, bremsstrahlung,...
Definition: Geant4TruthIncident.cxx:92
python.setupRTTAlg.size
int size
Definition: setupRTTAlg.py:39
TrackHelper
Definition: TrackHelper.h:14
ISFParticle.h
ParticleGun_EoverP_Config.mom
mom
Definition: ParticleGun_EoverP_Config.py:63
ParticleGun_FastCalo_ChargeFlip_Config.energy
energy
Definition: ParticleGun_FastCalo_ChargeFlip_Config.py:78
iGeant4::Geant4TruthIncident::position
const HepMC::FourVector & position() const override final
Return HepMC position of the truth vertex.
Definition: Geant4TruthIncident.cxx:82
lumiFormat.i
int i
Definition: lumiFormat.py:92
iGeant4::Geant4TruthIncident::physicsProcessCategory
int physicsProcessCategory() const override final
Return category of the physics process represented by the truth incident (eg hadronic,...
Definition: Geant4TruthIncident.cxx:86
iGeant4::Geant4TruthIncident::childEkin
double childEkin(unsigned short index) const override final
Return Ekin of the i-th child particle.
Definition: Geant4TruthIncident.cxx:203
HepMC::barcode
int barcode(const T *p)
Definition: Barcode.h:16
ISF::QS_SURV_VTX
@ QS_SURV_VTX
Definition: ITruthIncident.h:29
VTrackInformation::SetClassification
void SetClassification(TrackClassification tc)
update the classification of the currently tracked particle, usually called when a new G4Track is cre...
Definition: VTrackInformation.h:45
ISF::QS_DEST_VTX
@ QS_DEST_VTX
Definition: ITruthIncident.h:30
HepMC::uniqueID
int uniqueID(const T &p)
Definition: MagicNumbers.h:113
test_pyathena.parent
parent
Definition: test_pyathena.py:15
TrackInformation::SetRegenerationNr
void SetRegenerationNr(int i)
update the number of times the particle represented by the G4Track has undergone a non-destructive in...
Definition: TrackInformation.h:96
iGeant4::Geant4TruthIncident::childP2
double childP2(unsigned short index) const override final
Return p^2 of the i-th child particle.
Definition: Geant4TruthIncident.cxx:188
trackInfo
Definition: TrigInDetUtils.h:13
TrackInformation
Implementation of VTrackInformation. Instances of this class are attached as UserInformation to G4Tra...
Definition: TrackInformation.h:41
HepMC::UNDEFINED_ID
constexpr int UNDEFINED_ID
Definition: MagicNumbers.h:55
TrackInformation::GetRegenerationNr
int GetRegenerationNr() const
return the number of times the particle represented by the G4Track has undergone a non-destructive in...
Definition: TrackInformation.h:90
iGeant4::Geant4TruthIncident::childP
const G4ThreeVector childP(unsigned short index) const
Return p of the i-th child particle.
Definition: Geant4TruthIncident.cxx:193
TrackHelper::GetTrackInformation
TrackInformation * GetTrackInformation()
Definition: TrackHelper.h:24
iGeant4::Geant4TruthIncident::prepareChildren
void prepareChildren()
prepare the child particles
Definition: Geant4TruthIncident.cxx:315
runSelector.processType
processType
Definition: runSelector.py:225
HepMC::SIM_STATUS_THRESHOLD
constexpr int SIM_STATUS_THRESHOLD
Constant definiting the status threshold for simulated particles, eg. can be used to separate generat...
Definition: MagicNumbers.h:38
MagicNumbers.h
TrackInformation.h
ISF::InteractionClass_t
InteractionClass_t
The interaction classifications are described as follows: STD_VTX: interaction of a particle without ...
Definition: ITruthIncident.h:27
python.PhysicalConstants.c_light
float c_light
Definition: PhysicalConstants.py:63
TMVAToMVAUtils::convert
std::unique_ptr< MVAUtils::BDT > convert(TMVA::MethodBDT *bdt, bool isRegression=true, bool useYesNoLeaf=false)
Definition: TMVAToMVAUtils.h:114
iGeant4::Geant4TruthIncident::parentSurvivesIncident
bool parentSurvivesIncident() const override final
Return a boolean whether or not the parent particle survives the incident.
Definition: Geant4TruthIncident.cxx:138
python.LumiBlobConversion.pos
pos
Definition: LumiBlobConversion.py:18
PrimaryParticleInformation
This class is attached to G4PrimaryParticle objects as UserInformation. The member variable m_thePart...
Definition: PrimaryParticleInformation.h:39
DeMoScan.index
string index
Definition: DeMoScan.py:362
iGeant4::Geant4TruthIncident::particleAlive
bool particleAlive(const G4Track *track) const
check if the given G4Track represents a particle that is alive in ISF or ISF-G4
Definition: Geant4TruthIncident.cxx:255
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
AtlasG4EventUserInfo.h
CaloSwCorrections.time
def time(flags, cells_name, *args, **kw)
Definition: CaloSwCorrections.py:242
iGeant4::Geant4TruthIncident::childParticle
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...
Definition: Geant4TruthIncident.cxx:228
TrackHelper::IsPrimary
bool IsPrimary() const
Definition: TrackHelper.cxx:15
PrimaryParticleInformation.h
Geant4TruthIncident.h
LArCellBinning.step
step
Definition: LArCellBinning.py:158
iGeant4::Geant4TruthIncident::childPdgCode
int childPdgCode(unsigned short index) const override final
Return the PDG Code of the i-th child particle.
Definition: Geant4TruthIncident.cxx:207
iGeant4::Geant4TruthIncident::parentPdgCode
int parentPdgCode() const override final
Return the PDG Code of the parent particle.
Definition: Geant4TruthIncident.cxx:112
merge.status
status
Definition: merge.py:17
xAOD::track
@ track
Definition: TrackingPrimitives.h:512
iGeant4::Geant4TruthIncident::interactionClassification
ISF::InteractionClass_t interactionClassification() const override final
The interaction classifications are described as follows: STD_VTX: interaction of a particle without ...
Definition: Geant4TruthIncident.cxx:332
ISFG4Helper.h
iGeant4::ISFG4Helper::getISFTrackInfo
static VTrackInformation * getISFTrackInfo(const G4Track &aTrack)
return a valid UserInformation object of the G4Track for use within the ISF
Definition: ISFG4Helper.cxx:69
VTrackInformation::RegisteredSecondary
@ RegisteredSecondary
Definition: VTrackInformation.h:32
iGeant4::Geant4TruthIncident::childBarcode
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...
Definition: Geant4TruthIncident.cxx:211
iGeant4::Geant4TruthIncident::Geant4TruthIncident
Geant4TruthIncident()
iGeant4::Geant4TruthIncident::parentPt2
double parentPt2() const override final
Return pT^2 of the parent particle.
Definition: Geant4TruthIncident.cxx:103
VTrackInformation::RegeneratedPrimary
@ RegeneratedPrimary
Definition: VTrackInformation.h:32