GGS(GenericGEANT4Simulation)Software  2.7.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 38 of file GGSHadrIntAction.h.

Member Function Documentation

void GGSHadrIntAction::BeginOfEventAction ( const G4Event *  event)

Clears the interaction arrays.

Parameters
eventThe current event.

Definition at line 145 of file GGSHadrIntAction.cpp.

145  {
146 
147  _inelasticInteractions->Clear("C");
148  _quasiElasticInteractions->Clear("C");
149 }
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 164 of file GGSHadrIntAction.cpp.

164  {
165 
166  _outRootFile = GGSRootFileService::GetInstance().GetFileForThisRun(_outBase, run);
167  if (_outTreeName == "") {
168  _outTree = GGSRootFileService::GetInstance().GetDefaultTree(_outRootFile);
169  _outTree->SetTitle(TString(_outTree->GetTitle()) + "HadrInt ");
170  } else
171  _outTree = new TTree(_outTreeName.data(), "Hadronic interactions of primary particles");
172 
173  // Create inelastic interactions branch
174  _outTree->Branch("inelastic", "TClonesArray", &_inelasticInteractions);
175 
176  // Create quasi-elastic interactions branch
177  _outTree->Branch("quasiElastic", "TClonesArray", &_quasiElasticInteractions);
178 }
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 153 of file GGSHadrIntAction.cpp.

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

Closes the output file for the current run.

Parameters
runThe current run.

Definition at line 182 of file GGSHadrIntAction.cpp.

182  {
183 
184  if (_outTreeName != "") {
185  // Write the tree if we are not using the default one
186  _outRootFile->cd();
187  _outTree->Write();
188  }
190  _outRootFile = NULL;
191  _outTree = NULL;
192 }
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 56 of file GGSHadrIntAction.cpp.

56  {
57 
58  if (_tracks.size() > 0) {
59  // See if the track is in the tracking list
60  std::map<const G4Track *, int>::iterator trackIter = _tracks.find(track);
61  if (trackIter != _tracks.end()) {
62  // Check process type
63  const G4Step *step = track->GetStep();
64  int pType = step->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessType();
65  if (pType == fHadronic) {
66  int primaryPDG = track->GetDefinition()->GetPDGEncoding();
67 
68  // Check if interaction is quasi-elastic
69  float maxSecondaryEnergy = step->GetPreStepPoint()->GetKineticEnergy() * _energyFraction;
70  auto *secondaries = step->GetSecondaryInCurrentStep();
71  unsigned int nSecondaries = secondaries->size();
72  bool inelastic = true;
73  GGSTHadrIntInfo *newElement = NULL;
74  for (unsigned int is = 0; is < nSecondaries; is++) {
75  auto *currentSecondaryTrack = (*secondaries)[is];
76  if (currentSecondaryTrack->GetDefinition()->GetPDGEncoding() == primaryPDG &&
77  currentSecondaryTrack->GetKineticEnergy() > maxSecondaryEnergy) {
78 
79  // Quasi-elastic interaction
80  newElement =
81  (GGSTHadrIntInfo *)(_quasiElasticInteractions->ConstructedAt(_quasiElasticInteractions->GetEntries()));
82  newElement->Clear();
83 
84  // Add an entry in the track list for the new "primary"
85  _tracks[currentSecondaryTrack] = trackIter->second;
86  // Set flag to avoid treating it as inelastic
87  inelastic = false;
88  // End loop
89  break;
90  }
91  }
92 
93  if (inelastic) {
94  // Inelastic interaction
95  newElement = (GGSTHadrIntInfo *)(_inelasticInteractions->ConstructedAt(_inelasticInteractions->GetEntries()));
96  newElement->Clear();
97  }
98 
99  // Fill the new element
100  newElement->primary.PDGCode = primaryPDG;
101  newElement->primary.trackID = step->GetTrack()->GetTrackID();
102  for (int i = 0; i < 3; i++) {
103  newElement->primary.pos[i] = track->GetPosition()[i] / cm;
104  newElement->primary.mom[i] = track->GetMomentum()[i] / GeV;
105  }
106  newElement->primary.time = track->GetGlobalTime() / ns;
107  const G4Isotope *target =
108  ((G4HadronicProcess *)(step->GetPostStepPoint()->GetProcessDefinedStep()))->GetTargetIsotope();
109  newElement->targetPDGCode = 1000000000 + target->GetZ() * 10000 + target->GetN() * 10;
110  newElement->originalTrackID = trackIter->second;
111  newElement->processName = step->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName().data();
112  newElement->volumeName = step->GetPostStepPoint()->GetPhysicalVolume()->GetName();
113  if (_outputProducts) {
114  if (newElement->products) {
115  newElement->products->Clear("C");
116  } else {
117  newElement->products = new TClonesArray("GGSTParticle");
118  }
119  for (unsigned int is = 0; is < nSecondaries; is++) {
120  GGSTParticle *sec = new ((*(newElement->products))[is]) GGSTParticle;
121  auto *currentSecondaryTrack = (*secondaries)[is];
122  sec->PDGCode = currentSecondaryTrack->GetDefinition()->GetPDGEncoding();
123  sec->trackID = currentSecondaryTrack->GetTrackID();
124  for (int i = 0; i < 3; i++) {
125  sec->pos[i] = currentSecondaryTrack->GetPosition()[i] / cm;
126  sec->mom[i] = currentSecondaryTrack->GetMomentum()[i] / GeV;
127  }
128  sec->time = currentSecondaryTrack->GetGlobalTime() / ns;
129  }
130 
131  } else {
132  delete newElement->products;
133  newElement->products = NULL;
134  } // if(_outputProducts)
135  } // if (pType == fHadronic)
136 
137  // Delete the currently followed particle from list
138  _tracks.erase(trackIter);
139  } // if (trackIter != _tracks.end())
140  }
141 }
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.
Int_t PDGCode
PDG code of particle (see http://www3.nd.edu/~avillano/geant4/geant4_pid.html).
Definition: GGSTParticle.h:21
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.
Float_t pos[3]
Point of generation [cm].
Definition: GGSTParticle.h:23
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 }
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 120 of file GGSHadrIntAction.h.

120 { _energyFraction = fraction; }
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 102 of file GGSHadrIntAction.h.

102 { _outBase = outFileBase; }
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 126 of file GGSHadrIntAction.h.

126 { _outputProducts = output; }
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 110 of file GGSHadrIntAction.h.

110 { _outTreeName = outTreeName; }

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