GGS(GenericGEANT4Simulation)Software  2.99.0
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Macros
GGSHadrIntAction.cpp
Go to the documentation of this file.
1 /*
2  * GGSHadrIntAction.cpp
3  *
4  * Created on: 31 May 2011
5  * Author: Nicola mori
6  */
7 
11 #include "montecarlo/dataobjs/GGSTHadrIntInfo.h"
12 #include "montecarlo/dataobjs/GGSTParticle.h"
15 
16 #include "G4HadronicProcess.hh"
17 #include "G4SystemOfUnits.hh"
18 
19 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
20 
22 
23 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
24 
26  : _energyFraction(0.9), _outputProducts(false), _outBase(""), _outTreeName(""), _outRootFile(nullptr),
27  _outTree(nullptr) {
28 
29  _messenger = new GGSHadrIntMessenger(this);
30  _quasiElasticInteractions = new TClonesArray("GGSTHadrIntInfo", 10);
31  _inelasticInteractions = new TClonesArray("GGSTHadrIntInfo", 1);
32 }
33 
34 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
35 
37 
38  delete _messenger;
39  delete _quasiElasticInteractions;
40  delete _inelasticInteractions;
41 }
42 
43 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
44 
45 void GGSHadrIntAction::PreUserTrackingAction(const G4Track *track) {
46 
47  // Check whether it is a primary particle
48  G4PrimaryParticle *primary = track->GetDynamicParticle()->GetPrimaryParticle();
49  if (primary != nullptr)
50  // Check whether it is a baryon
51  if (primary->GetParticleDefinition()->GetBaryonNumber() > 0)
52  _tracks[track] = track->GetTrackID();
53 }
54 
55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56 
57 void GGSHadrIntAction::PostUserTrackingAction(const G4Track *track) {
58 
59  if (_tracks.size() > 0) {
60  // See if the track is in the tracking list
61  std::map<const G4Track *, int>::iterator trackIter = _tracks.find(track);
62  if (trackIter != _tracks.end()) {
63  // Check process type
64  const G4Step *step = track->GetStep();
65  int pType = step->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessType();
66  if (pType == fHadronic) {
67  int primaryPDG = track->GetDefinition()->GetPDGEncoding();
68 
69  // Check if interaction is quasi-elastic
70  float maxSecondaryEnergy = step->GetPreStepPoint()->GetKineticEnergy() * _energyFraction;
71  auto *secondaries = step->GetSecondaryInCurrentStep();
72  unsigned int nSecondaries = secondaries->size();
73  bool inelastic = true;
74  GGSTHadrIntInfo *newElement = nullptr;
75  for (unsigned int is = 0; is < nSecondaries; is++) {
76  auto *currentSecondaryTrack = (*secondaries)[is];
77  if (currentSecondaryTrack->GetDefinition()->GetPDGEncoding() == primaryPDG &&
78  currentSecondaryTrack->GetKineticEnergy() > maxSecondaryEnergy) {
79 
80  // Quasi-elastic interaction
81  newElement =
82  (GGSTHadrIntInfo *)(_quasiElasticInteractions->ConstructedAt(_quasiElasticInteractions->GetEntries()));
83  newElement->Clear();
84 
85  // Add an entry in the track list for the new "primary"
86  _tracks[currentSecondaryTrack] = trackIter->second;
87  // Set flag to avoid treating it as inelastic
88  inelastic = false;
89  // End loop
90  break;
91  }
92  }
93 
94  if (inelastic) {
95  // Inelastic interaction
96  newElement = (GGSTHadrIntInfo *)(_inelasticInteractions->ConstructedAt(_inelasticInteractions->GetEntries()));
97  newElement->Clear();
98  }
99 
100  // Fill the new element
101  newElement->primary.PDGCode = primaryPDG;
102  newElement->primary.trackID = step->GetTrack()->GetTrackID();
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;
106  }
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;
111  newElement->originalTrackID = trackIter->second;
112  newElement->processName = step->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName().data();
113  newElement->volumeName = step->GetPostStepPoint()->GetPhysicalVolume()->GetName();
114  if (_outputProducts) {
115  if (newElement->products) {
116  newElement->products->Clear("C");
117  } else {
118  newElement->products = new TClonesArray("GGSTParticle");
119  }
120  for (unsigned int is = 0; is < nSecondaries; is++) {
121  GGSTParticle *sec = new ((*(newElement->products))[is]) GGSTParticle;
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;
128  }
129  sec->time = currentSecondaryTrack->GetGlobalTime() / ns;
130  }
131 
132  } else {
133  delete newElement->products;
134  newElement->products = nullptr;
135  } // if(_outputProducts)
136  } // if (pType == fHadronic)
137 
138  // Delete the currently followed particle from list
139  _tracks.erase(trackIter);
140  } // if (trackIter != _tracks.end())
141  }
142 }
143 
144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145 
147 
148  _inelasticInteractions->Clear("C");
149  _quasiElasticInteractions->Clear("C");
150 }
151 
152 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
153 
154 void GGSHadrIntAction::EndOfEventAction(const G4Event *) {
155 
156  if (_outTreeName != "") {
157  // Call fill only if not using the default tree
158  _outRootFile->cd();
159  _outTree->Fill();
160  }
161 }
162 
163 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
164 
165 void GGSHadrIntAction::BeginOfRunAction(const G4Run *run) {
166 
167  _outRootFile = GGSRootFileService::GetInstance().GetFileForThisRun(_outBase, run);
168  if (_outTreeName == "") {
169  _outTree = GGSRootFileService::GetInstance().GetDefaultTree(_outRootFile);
170  _outTree->SetTitle(TString(_outTree->GetTitle()) + "HadrInt ");
171  } else
172  _outTree = new TTree(_outTreeName.data(), "Hadronic interactions of primary particles");
173 
174  // Create inelastic interactions branch
175  _outTree->Branch("inelastic", "TClonesArray", &_inelasticInteractions);
176 
177  // Create quasi-elastic interactions branch
178  _outTree->Branch("quasiElastic", "TClonesArray", &_quasiElasticInteractions);
179 }
180 
181 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
182 
183 void GGSHadrIntAction::EndOfRunAction(const G4Run *run) {
184 
185  if (_outTreeName != "") {
186  // Write the tree if we are not using the default one
187  _outRootFile->cd();
188  _outTree->Write();
189  }
191  _outRootFile = nullptr;
192  _outTree = nullptr;
193 }
Float_t mom[3]
Momentum at generation [GeV].
Definition: GGSTParticle.h:24
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).
Definition: GGSTParticle.h:21
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].
Definition: GGSTParticle.h:25
Int_t trackID
Track ID.
Definition: GGSTParticle.h:22
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].
Definition: GGSTParticle.h:23
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.
Definition: GGSTParticle.h:19
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.