11 #include "montecarlo/dataobjs/GGSTHadrIntInfo.h"
12 #include "montecarlo/dataobjs/GGSTParticle.h"
16 #include "G4HadronicProcess.hh"
17 #include "G4SystemOfUnits.hh"
26 : _energyFraction(0.9), _outputProducts(false), _outBase(
""), _outTreeName(
""), _outRootFile(nullptr),
30 _quasiElasticInteractions =
new TClonesArray(
"GGSTHadrIntInfo", 10);
31 _inelasticInteractions =
new TClonesArray(
"GGSTHadrIntInfo", 1);
39 delete _quasiElasticInteractions;
40 delete _inelasticInteractions;
48 G4PrimaryParticle *primary = track->GetDynamicParticle()->GetPrimaryParticle();
49 if (primary !=
nullptr)
51 if (primary->GetParticleDefinition()->GetBaryonNumber() > 0)
52 _tracks[track] = track->GetTrackID();
59 if (_tracks.size() > 0) {
61 std::map<const G4Track *, int>::iterator trackIter = _tracks.find(track);
62 if (trackIter != _tracks.end()) {
64 const G4Step *step = track->GetStep();
65 int pType = step->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessType();
66 if (pType == fHadronic) {
67 int primaryPDG = track->GetDefinition()->GetPDGEncoding();
70 float maxSecondaryEnergy = step->GetPreStepPoint()->GetKineticEnergy() * _energyFraction;
71 auto *secondaries = step->GetSecondaryInCurrentStep();
72 unsigned int nSecondaries = secondaries->size();
73 bool inelastic =
true;
75 for (
unsigned int is = 0; is < nSecondaries; is++) {
76 auto *currentSecondaryTrack = (*secondaries)[is];
77 if (currentSecondaryTrack->GetDefinition()->GetPDGEncoding() == primaryPDG &&
78 currentSecondaryTrack->GetKineticEnergy() > maxSecondaryEnergy) {
82 (
GGSTHadrIntInfo *)(_quasiElasticInteractions->ConstructedAt(_quasiElasticInteractions->GetEntries()));
86 _tracks[currentSecondaryTrack] = trackIter->second;
96 newElement = (
GGSTHadrIntInfo *)(_inelasticInteractions->ConstructedAt(_inelasticInteractions->GetEntries()));
103 for (
int i = 0; i < 3; i++) {
104 newElement->
primary.
pos[i] = track->GetPosition()[i] / cm;
105 newElement->
primary.
mom[i] = track->GetMomentum()[i] / GeV;
107 newElement->
primary.
time = track->GetGlobalTime() / ns;
108 const G4Isotope *target =
109 ((G4HadronicProcess *)(step->GetPostStepPoint()->GetProcessDefinedStep()))->GetTargetIsotope();
110 newElement->
targetPDGCode = 1000000000 + target->GetZ() * 10000 + target->GetN() * 10;
112 newElement->
processName = step->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName().data();
113 newElement->
volumeName = step->GetPostStepPoint()->GetPhysicalVolume()->GetName();
114 if (_outputProducts) {
118 newElement->
products =
new TClonesArray(
"GGSTParticle");
120 for (
unsigned int is = 0; is < nSecondaries; is++) {
122 auto *currentSecondaryTrack = (*secondaries)[is];
123 sec->
PDGCode = currentSecondaryTrack->GetDefinition()->GetPDGEncoding();
124 sec->
trackID = currentSecondaryTrack->GetTrackID();
125 for (
int i = 0; i < 3; i++) {
126 sec->
pos[i] = currentSecondaryTrack->GetPosition()[i] / cm;
127 sec->
mom[i] = currentSecondaryTrack->GetMomentum()[i] / GeV;
129 sec->
time = currentSecondaryTrack->GetGlobalTime() / ns;
139 _tracks.erase(trackIter);
148 _inelasticInteractions->Clear(
"C");
149 _quasiElasticInteractions->Clear(
"C");
156 if (_outTreeName !=
"") {
168 if (_outTreeName ==
"") {
170 _outTree->SetTitle(TString(_outTree->GetTitle()) +
"HadrInt ");
172 _outTree =
new TTree(_outTreeName.data(),
"Hadronic interactions of primary particles");
175 _outTree->Branch(
"inelastic",
"TClonesArray", &_inelasticInteractions);
178 _outTree->Branch(
"quasiElastic",
"TClonesArray", &_quasiElasticInteractions);
185 if (_outTreeName !=
"") {
191 _outRootFile =
nullptr;
Float_t mom[3]
Momentum at generation [GeV].
Int_t originalTrackID
Track ID of the original primary particle.
GGSTParticle primary
The particle originating the interaction.
void BeginOfEventAction(const G4Event *event)
Clears the interaction arrays.
GGSHadrIntAction()
Constructor.
TFile * GetFileForThisRun(const std::filesystem::path &baseName, const G4Run *run)
Opens a file for a given run and returns a pointer to it.
void BeginOfRunAction(const G4Run *run)
Opens the output file for the current run and prepares the output tree.
Int_t PDGCode
PDG code of particle (see http://www3.nd.edu/~avillano/geant4/geant4_pid.html).
void CloseFileForThisRun(const std::filesystem::path &baseName, const G4Run *run)
Closes the ROOT output file.
void PostUserTrackingAction(const G4Track *track)
Check the track end for a hadronic interaction.
Float_t time
Time of generation [ns].
TClonesArray * products
Particles produced in the interaction.
TString volumeName
The name of the physical volume where the interaction took place.
Action which finds the hadronic interaction points for each primary particle.
Float_t pos[3]
Point of generation [cm].
void EndOfRunAction(const G4Run *run)
Closes the output file for the current run.
A simple class to carry informations about hadronic interactions.
Class to store G4 particles.
A control messenger for GGSHadrIntAction.
TTree * GetDefaultTree(TFile *file)
Gets the default tree for this file.
void PreUserTrackingAction(const G4Track *track)
Initialization of primary track.
static GGSRootFileService & GetInstance()
Get reference to GGSRootFileService unique instance.
Int_t targetPDGCode
PDG code of target nucleus.
void EndOfEventAction(const G4Event *event)
Fills the output tree.
#define RegisterUA(uaClassName)
Macro for registration of user actions classes.
void Clear(Option_t *="")
Resets all the members.
TString processName
Hadronic process name.
virtual ~GGSHadrIntAction()
Destructor.