ATLAS Offline Software
Loading...
Searching...
No Matches
Geant4TruthIncident.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef ISF_GEANT4TOOLS_Geant4TruthIncident_H
6#define ISF_GEANT4TOOLS_Geant4TruthIncident_H
7
8// ISF includes
9#include "ISF_Event/ITruthIncident.h" //inheritance
10
11//G4 includes
12#include "G4ThreeVector.hh" //typedef
13
14// HepMC includes
15#include "AtlasHepMC/SimpleVector.h" //typedef FourVector
16#include "AtlasHepMC/GenParticle_fwd.h" //typedef GenParticlePtr
17
18#include "AtlasDetDescr/AtlasRegion.h" //enum
19
20#include "CxxUtils/checker_macros.h" //ATLAS_THREAD_SAFE
21// std
22#include <vector>
23
24
25// forward declarations
26class G4Step;
27class G4Track;
29
30namespace ISF {
31 class ISFParticle;
32}
33
34
35namespace iGeant4 {
36
43
45 public:
46 Geant4TruthIncident( const G4Step*,
47 const ISF::ISFParticle& baseISP,
49 AtlasG4EventUserInfo* atlasG4EvtUserInfo);
51
53 const HepMC::FourVector& position() const override final;
54
56 int physicsProcessCategory() const override final;
58 int physicsProcessCode() const override final;
59
61 double parentP2() const override final;
63 double parentPt2() const override final;
65 double parentEkin() const override final;
67 int parentPdgCode() const override final;
69 int parentBarcode() override final; // TODO Remove this method
71 int parentUniqueID() override final;
73 int parentStatus() override final;
75 bool parentSurvivesIncident() const override final;
78 HepMC::GenParticlePtr parentParticleAfterIncident(int newBC) override final;
79
81 const G4ThreeVector childP(unsigned short index) const;
83 double childP2(unsigned short index) const override final;
85 double childPt2(unsigned short index) const override final;
87 double childEkin(unsigned short index) const override final;
89 int childPdgCode(unsigned short index) const override final;
91 int childBarcode(unsigned short index) const override final; // TODO Remove - only used in one place in TruthSvc
92
100 ISF::InteractionClass_t interactionClassification() const override final;
101
102 // only called once accepted
103
105 HepMC::GenParticlePtr parentParticle() override final;
108 HepMC::GenParticlePtr childParticle(unsigned short index,
109 int bc) override final;
110 private:
113 inline void prepareChildren();
114
116 inline bool particleAlive(const G4Track *track) const;
117
118 HepMC::GenParticlePtr convert(const G4Track *particle, const int barcode, const bool secondary) const;
119
121 HepMC::FourVector m_position;
122 const G4Step* m_step{};
124
126 std::vector<const G4Track*> m_children;
127
129 };
130
131}
132
133#endif // ISF_GEANT4TOOLS_Geant4TruthIncident_H
Define macros for attributes used to control the static checker.
This class is attached to G4Event objects as UserInformation.
The generic ISF particle definition,.
Definition ISFParticle.h:42
ISF interface class for TruthIncidents.
AtlasDetDescr::AtlasRegion geoID()
Return the SimGeoID corresponding to the vertex of the truth incident.
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.
AtlasG4EventUserInfo *m_atlasG4EvtUserInfo ATLAS_THREAD_SAFE
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.
Geant4TruthIncident(const G4Step *, const ISF::ISFParticle &baseISP, AtlasDetDescr::AtlasRegion geoID, AtlasG4EventUserInfo *atlasG4EvtUserInfo)
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
AtlasRegion
A simple enum of ATLAS regions and sub-detectors.
Definition AtlasRegion.h:21
GenParticle * GenParticlePtr
Definition GenParticle.h:37
ISFParticleOrderedQueue.
Definition index.py:1
#define private