GGS(GenericGEANT4Simulation)Software  2.6.0
 All Data Structures Namespaces Files Functions Variables Typedefs 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(NULL), _outTree(NULL) {
27 
28  _messenger = new GGSHadrIntMessenger(this);
29  _quasiElasticInteractions = new TClonesArray("GGSTHadrIntInfo", 10);
30  _inelasticInteractions = new TClonesArray("GGSTHadrIntInfo", 1);
31 }
32 
33 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
34 
36 
37  delete _messenger;
38  delete _quasiElasticInteractions;
39  delete _inelasticInteractions;
40 }
41 
42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
43 
44 void GGSHadrIntAction::PreUserTrackingAction(const G4Track *track) {
45 
46  // Check whether it is a primary particle
47  G4PrimaryParticle *primary = track->GetDynamicParticle()->GetPrimaryParticle();
48  if (primary != NULL)
49  // Check whether it is a baryon
50  if (primary->GetParticleDefinition()->GetBaryonNumber() > 0)
51  _tracks[track] = track->GetTrackID();
52 
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  const G4TrackVector *secondaries = step->GetSecondary();
72  unsigned int nSecondaries = secondaries->size();
73  bool inelastic = true;
74  GGSTHadrIntInfo *newElement = NULL;
75  for (unsigned int is = 0; is < nSecondaries; is++) {
76  G4Track *currentSecondaryTrack = (*secondaries)[is];
77  if (currentSecondaryTrack->GetDefinition()->GetPDGEncoding() == primaryPDG
78  && currentSecondaryTrack->GetKineticEnergy() > maxSecondaryEnergy) {
79 
80  // Quasi-elastic interaction
81  newElement = (GGSTHadrIntInfo*) (_quasiElasticInteractions->ConstructedAt(
82  _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  if (_outputProducts) {
114  if (newElement->products) {
115  newElement->products->Clear("C");
116  }
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  G4Track *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  }
133  else {
134  delete newElement->products;
135  newElement->products = NULL;
136  } //if(_outputProducts)
137  } //if (pType == fHadronic)
138 
139  // Delete the currently followed particle from list
140  _tracks.erase(trackIter);
141  } //if (trackIter != _tracks.end())
142 
143  }
144 }
145 
146 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
147 
149 
150  _inelasticInteractions->Clear("C");
151  _quasiElasticInteractions->Clear("C");
152 }
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
155 
156 void GGSHadrIntAction::EndOfEventAction(const G4Event *) {
157 
158  if (_outTreeName != "") {
159  // Call fill only if not using the default tree
160  _outRootFile->cd();
161  _outTree->Fill();
162  }
163 
164 }
165 
166 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
167 
168 void GGSHadrIntAction::BeginOfRunAction(const G4Run *run) {
169 
170  _outRootFile = GGSRootFileService::GetInstance().GetFileForThisRun(_outBase, run);
171  if (_outTreeName == "") {
172  _outTree = GGSRootFileService::GetInstance().GetDefaultTree(_outRootFile);
173  _outTree->SetTitle(TString(_outTree->GetTitle()) + "HadrInt ");
174  }
175  else
176  _outTree = new TTree(_outTreeName.data(), "Hadronic interactions of primary particles");
177 
178  // Create inelastic interactions branch
179  _outTree->Branch("inelastic", "TClonesArray", &_inelasticInteractions);
180 
181  // Create quasi-elastic interactions branch
182  _outTree->Branch("quasiElastic", "TClonesArray", &_quasiElasticInteractions);
183 }
184 
185 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
186 
187 void GGSHadrIntAction::EndOfRunAction(const G4Run *) {
188 
189  if (_outTreeName != "") {
190  // Write the tree if we are not using the default one
191  _outRootFile->cd();
192  _outTree->Write();
193  }
195  _outRootFile = NULL;
196  _outTree = NULL;
197 }
Float_t mom[3]
Momentum at generation [GeV].
Definition: GGSTParticle.h:25
Int_t originalTrackID
Track ID of the original primary particle.
TFile * GetFileForThisRun(const path &baseName, const G4Run *run)
Opens a file for a given run and returns a pointer to it.
GGSTParticle primary
The particle originating the interaction.
void BeginOfEventAction(const G4Event *event)
Clears the interaction arrays.
GGSHadrIntAction()
Constructor.
void BeginOfRunAction(const G4Run *run)
Opens the output file for the current run and prepares the output tree.
void CloseFileForThisRun(const path &baseName)
Closes the ROOT output file.
Int_t PDGCode
PDG code of particle (see http://www3.nd.edu/~avillano/geant4/geant4_pid.html).
Definition: GGSTParticle.h:22
void PostUserTrackingAction(const G4Track *track)
Check the track end for a hadronic interaction.
Float_t time
Time of generation [ns].
Definition: GGSTParticle.h:26
Int_t trackID
Track ID.
Definition: GGSTParticle.h:23
TClonesArray * products
Particles produced in the interaction.
Action which finds the hadronic interaction points for each primary particle.
Float_t pos[3]
Point of generation [cm].
Definition: GGSTParticle.h:24
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.