GGS(GenericGEANT4Simulation)Software  2.6.0
 All Data Structures Namespaces Files Functions Variables Typedefs Macros
Public Member Functions
GGSHadrIntAction Class Reference

Action which finds the hadronic interaction points for each primary particle. More...

#include <GGSHadrIntAction.h>

Inheritance diagram for GGSHadrIntAction:
Inheritance graph
[legend]
Collaboration diagram for GGSHadrIntAction:
Collaboration graph
[legend]

Public Member Functions

 GGSHadrIntAction ()
 Constructor.
 
virtual ~GGSHadrIntAction ()
 Destructor.
 
void PreUserTrackingAction (const G4Track *track)
 Initialization of primary track. More...
 
void PostUserTrackingAction (const G4Track *track)
 Check the track end for a hadronic interaction. More...
 
void BeginOfEventAction (const G4Event *event)
 Clears the interaction arrays. More...
 
void EndOfEventAction (const G4Event *event)
 Fills the output tree. More...
 
void BeginOfRunAction (const G4Run *run)
 Opens the output file for the current run and prepares the output tree. More...
 
void EndOfRunAction (const G4Run *run)
 Closes the output file for the current run. More...
 
void SetOutputFileBase (const std::string &outFileBase)
 Sets the output file base name. More...
 
void SetOutputTreeName (const std::string &outTreeName)
 Sets the output tree name. More...
 
void SetEnergyFraction (double fraction)
 Sets the maximum secondary energy fraction for hadronic interactions. More...
 
void SetOutputOfProducts (bool output)
 Toggles the output of hadronic interaction products. More...
 
- Public Member Functions inherited from GGSUserAction
 GGSUserAction ()
 Constructor. More...
 
virtual ~GGSUserAction ()
 Destructor.
 
G4ClassificationOfNewTrack ClassifyNewTrack (const G4Track *)
 Override of the ClassifyNewTrack method. More...
 

Detailed Description

Action which finds the hadronic interaction points for each primary particle.

This class distinguish quasi-elastic interactions from inelastic ones, providing a list of information for both interactions categories for each primary particle.

Definition at line 39 of file GGSHadrIntAction.h.

Member Function Documentation

void GGSHadrIntAction::BeginOfEventAction ( const G4Event *  event)

Clears the interaction arrays.

Parameters
eventThe current event.

Definition at line 148 of file GGSHadrIntAction.cpp.

148  {
149 
150  _inelasticInteractions->Clear("C");
151  _quasiElasticInteractions->Clear("C");
152 }
void GGSHadrIntAction::BeginOfRunAction ( const G4Run *  run)

Opens the output file for the current run and prepares the output tree.

Parameters
runThe current run.

Definition at line 168 of file GGSHadrIntAction.cpp.

168  {
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 }
TFile * GetFileForThisRun(const path &baseName, const G4Run *run)
Opens a file for a given run and returns a pointer to it.
TTree * GetDefaultTree(TFile *file)
Gets the default tree for this file.
static GGSRootFileService & GetInstance()
Get reference to GGSRootFileService unique instance.
void GGSHadrIntAction::EndOfEventAction ( const G4Event *  event)

Fills the output tree.

Parameters
eventThe current event.

Definition at line 156 of file GGSHadrIntAction.cpp.

156  {
157 
158  if (_outTreeName != "") {
159  // Call fill only if not using the default tree
160  _outRootFile->cd();
161  _outTree->Fill();
162  }
163 
164 }
void GGSHadrIntAction::EndOfRunAction ( const G4Run *  run)

Closes the output file for the current run.

Parameters
runThe current run.

Definition at line 187 of file GGSHadrIntAction.cpp.

187  {
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 }
void CloseFileForThisRun(const path &baseName)
Closes the ROOT output file.
static GGSRootFileService & GetInstance()
Get reference to GGSRootFileService unique instance.
void GGSHadrIntAction::PostUserTrackingAction ( const G4Track *  track)

Check the track end for a hadronic interaction.

For each followed particle, this routine will check if a hadronic interaction has happened. If no particle identical to the primary and carrying more than a specified fraction of the energy of the primary is produced, the interaction point is recorded as "inelastic". Otherwise, it is recorded as "quasi-elastic", the secondary carrying most of the energy is tagged as the new primary and followed until an inelastic interaction or the end of the event.

Parameters
trackThe current track.

Definition at line 57 of file GGSHadrIntAction.cpp.

57  {
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 }
Float_t mom[3]
Momentum at generation [GeV].
Definition: GGSTParticle.h:25
Int_t originalTrackID
Track ID of the original primary particle.
GGSTParticle primary
The particle originating the interaction.
Int_t PDGCode
PDG code of particle (see http://www3.nd.edu/~avillano/geant4/geant4_pid.html).
Definition: GGSTParticle.h:22
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.
Float_t pos[3]
Point of generation [cm].
Definition: GGSTParticle.h:24
A simple class to carry informations about hadronic interactions.
Class to store G4 particles.
Definition: GGSTParticle.h:19
Int_t targetPDGCode
PDG code of target nucleus.
void Clear(Option_t *="")
Resets all the members.
TString processName
Hadronic process name.
void GGSHadrIntAction::PreUserTrackingAction ( const G4Track *  track)

Initialization of primary track.

This routine checks if the current track belongs to a primary particle, and in affirmative case it adds it to the followed particles vector

Parameters
trackThe current track.

Definition at line 44 of file GGSHadrIntAction.cpp.

44  {
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 }
void GGSHadrIntAction::SetEnergyFraction ( double  fraction)
inline

Sets the maximum secondary energy fraction for hadronic interactions.

To to tag "inelastic" interactions, the produced secondaries which are identical to the primary are requested to carry less than a fraction of the primary's energy. This routine sets this fraction, which default value is set in constructor.

Parameters
fractionThe maximum secondary energy fraction.

Definition at line 126 of file GGSHadrIntAction.h.

126  {
127  _energyFraction = fraction;
128  }
void GGSHadrIntAction::SetOutputFileBase ( const std::string &  outFileBase)
inline

Sets the output file base name.

The file base name can be with or without extension (.root will be automatically used as extension). For each run, the run number will be appended to the base name before the .root extension. If no value is provided the file base name will fallback to the default value set in GSRootFileservice. *

Parameters
outFileBaseThe output file base name.
See Also
GGSRootFileservice

Definition at line 104 of file GGSHadrIntAction.h.

104  {
105  _outBase = outFileBase;
106  }
void GGSHadrIntAction::SetOutputOfProducts ( bool  output)
inline

Toggles the output of hadronic interaction products.

Parameters
outputIf true the products of the interaction will be saved in the file.

Definition at line 134 of file GGSHadrIntAction.h.

134  {
135  _outputProducts = output;
136  }
void GGSHadrIntAction::SetOutputTreeName ( const std::string &  outTreeName)
inline

Sets the output tree name.

This name will be used for the TTree object where the hits for each event will be stored. By default, this value is set to "HadrInt" in constructor.

Parameters
outTreeNameThe output tree name.

Definition at line 114 of file GGSHadrIntAction.h.

114  {
115  _outTreeName = outTreeName;
116  }

The documentation for this class was generated from the following files: