Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
IORTSteppingAction Class Reference

#include <IORTSteppingAction.hh>

Inheritance diagram for IORTSteppingAction:
Collaboration diagram for IORTSteppingAction:

Public Member Functions

 IORTSteppingAction (IORTRunAction *)
 
 ~IORTSteppingAction ()
 
void UserSteppingAction (const G4Step *)
 
- Public Member Functions inherited from G4UserSteppingAction
 G4UserSteppingAction ()
 
virtual ~G4UserSteppingAction ()
 
virtual void SetSteppingManagerPointer (G4SteppingManager *pValue)
 

Additional Inherited Members

- Protected Attributes inherited from G4UserSteppingAction
G4SteppingManagerfpSteppingManager
 

Detailed Description

Definition at line 58 of file IORTSteppingAction.hh.

Constructor & Destructor Documentation

IORTSteppingAction::IORTSteppingAction ( IORTRunAction run)

Definition at line 63 of file IORTSteppingAction.cc.

64 {
65  runAction = run;
66 }
IORTSteppingAction::~IORTSteppingAction ( )

Definition at line 69 of file IORTSteppingAction.cc.

70 {
71 }

Member Function Documentation

void IORTSteppingAction::UserSteppingAction ( const G4Step aStep)
virtual

Reimplemented from G4UserSteppingAction.

Definition at line 74 of file IORTSteppingAction.cc.

75 {
76  /*
77  // USEFULL METHODS TO RETRIEVE INFORMATION DURING THE STEPS
78  if( (aStep->GetTrack()->GetVolume()->GetName() == "DetectorPhys")
79  && aStep->GetTrack()->GetDefinition()->GetParticleName() == "proton")
80  //G4int evtNb = G4RunManager::GetRunManager()->GetCurrentEvent() -> GetEventID();
81  {
82  G4cout << "ENERGIA: " << aStep->GetTrack()->GetKineticEnergy()
83  << " VOLUME " << aStep->GetTrack()->GetVolume()->GetName()
84  << " MATERIALE " << aStep -> GetTrack() -> GetMaterial() -> GetName()
85  << " EVENTO " << G4RunManager::GetRunManager()->GetCurrentEvent() -> GetEventID()
86  << " POS " << aStep->GetTrack()->GetPosition().x()
87  << G4endl;
88  }
89  */
90 
91  if( aStep->GetTrack()->GetVolume()->GetName() == "NewDetectorPhys"){
92  G4ParticleDefinition *def = aStep->GetTrack()->GetDefinition();
93  G4double secondaryParticleKineticEnergy = aStep->GetTrack()->
94  GetKineticEnergy();
95  G4String particleType = def->GetParticleType(); // particle type = nucleus for d, t, He3, alpha, and heavier nuclei
96  G4String particleName = def->GetParticleName(); // e.g. for alpha: the name = "alpha" and type = "nucleus"
97  if(particleType == "nucleus") {
98  G4int A = def->GetBaryonNumber();
99  G4double Z = def->GetPDGCharge();
100  G4double posX = aStep->GetTrack()->GetPosition().x() /cm;
101  G4double posY = aStep->GetTrack()->GetPosition().y() /cm;
102  G4double posZ = aStep->GetTrack()->GetPosition().z() /cm;
103  G4double energy = secondaryParticleKineticEnergy / A /MeV;
104 
106  analysisMgr->FillFragmentTuple(A, Z, energy, posX, posY, posZ);
107  } else if(particleName == "proton") { // proton (hydrogen-1) is a special case
108  G4double posX = aStep->GetTrack()->GetPosition().x() /cm ;
109  G4double posY = aStep->GetTrack()->GetPosition().y() /cm ;
110  G4double posZ = aStep->GetTrack()->GetPosition().z() / cm ;
111  G4double energy = secondaryParticleKineticEnergy / MeV; // Hydrogen-1: A = 1, Z = 1
112  IORTAnalysisManager::GetInstance()->FillFragmentTuple(1, 1.0, energy, posX, posY, posZ);
113  }
114 
115  G4String secondaryParticleName = def -> GetParticleName();
116  //G4cout <<"Particle: " << secondaryParticleName << G4endl;
117  //G4cout <<"Energy: " << secondaryParticleKineticEnergy << G4endl;
119  //There is a bunch of stuff recorded with the energy 0, something should perhaps be done about this.
120  if(secondaryParticleName == "proton") {
121  analysis->hydrogenEnergy(secondaryParticleKineticEnergy /MeV);
122  }
123  if(secondaryParticleName == "deuteron") {
124  analysis->hydrogenEnergy((secondaryParticleKineticEnergy/2) / MeV);
125  }
126  if(secondaryParticleName == "triton") {
127  analysis->hydrogenEnergy((secondaryParticleKineticEnergy/3) / MeV);
128  }
129  if(secondaryParticleName == "alpha") {
130  analysis->heliumEnergy((secondaryParticleKineticEnergy/4) / MeV);
131  }
132  if(secondaryParticleName == "He3"){
133  analysis->heliumEnergy((secondaryParticleKineticEnergy/3) / MeV);
134  }
135 
137  }
138 
139  // Electromagnetic and hadronic processes of primary particles in the phantom
140  //setting phantomPhys correctly will break something here fixme
141  if ((aStep -> GetTrack() -> GetTrackID() == 1) &&
142  (aStep -> GetTrack() -> GetVolume() -> GetName() == "PhantomPhys") &&
143  (aStep -> GetPostStepPoint() -> GetProcessDefinedStep() != NULL))
144  {
145  G4String process = aStep -> GetPostStepPoint() ->
146  GetProcessDefinedStep() -> GetProcessName();
147 
148  if ((process == "Transportation") || (process == "StepLimiter")) {;}
149  else
150  {
151  if ((process == "msc") || (process == "hLowEIoni") || (process == "hIoni"))
152  {
153  runAction -> AddEMProcess();
154  }
155  else
156  {
157  runAction -> AddHadronicProcess();
158 
159  if ( (process != "LElastic") && (process != "ProtonInelastic") && (process != "hElastic") )
160  G4cout << "Warning! Unknown proton process: "<< process << G4endl;
161  }
162  }
163  }
164 
165  // Retrieve information about the secondary particles originated in the phantom
166 
167  G4SteppingManager* steppingManager = fpSteppingManager;
168 
169  // check if it is alive
170  //if(theTrack-> GetTrackStatus() == fAlive) { return; }
171 
172  // Retrieve the secondary particles
173  G4TrackVector* fSecondary = steppingManager -> GetfSecondary();
174 
175  for(size_t lp1=0;lp1<(*fSecondary).size(); lp1++)
176  {
177  G4String volumeName = (*fSecondary)[lp1] -> GetVolume() -> GetName();
178 
179  if (volumeName == "phantomPhys")
180  {
181  G4String secondaryParticleName = (*fSecondary)[lp1]->GetDefinition() -> GetParticleName();
182  G4double secondaryParticleKineticEnergy = (*fSecondary)[lp1] -> GetKineticEnergy();
183 
185 
186  if (secondaryParticleName == "e-")
187  analysis -> electronEnergyDistribution(secondaryParticleKineticEnergy/MeV);
188 
189  if (secondaryParticleName == "gamma")
190  analysis -> gammaEnergyDistribution(secondaryParticleKineticEnergy/MeV);
191 
192  if (secondaryParticleName == "deuteron")
193  analysis -> deuteronEnergyDistribution(secondaryParticleKineticEnergy/MeV);
194 
195  if (secondaryParticleName == "triton")
196  analysis -> tritonEnergyDistribution(secondaryParticleKineticEnergy/MeV);
197 
198  if (secondaryParticleName == "alpha")
199  analysis -> alphaEnergyDistribution(secondaryParticleKineticEnergy/MeV);
200 
201  G4double z = (*fSecondary)[lp1]-> GetDynamicParticle() -> GetDefinition() -> GetPDGCharge();
202  if (z > 0.)
203  {
204  G4int a = (*fSecondary)[lp1]-> GetDynamicParticle() -> GetDefinition() -> GetBaryonNumber();
205  G4int electronOccupancy = (*fSecondary)[lp1] -> GetDynamicParticle() -> GetTotalOccupancy();
206 
207  // If a generic ion is originated in the detector, its baryonic number, PDG charge,
208  // total number of electrons in the orbitals are stored in a ntuple
209  analysis -> genericIonInformation(a, z, electronOccupancy, secondaryParticleKineticEnergy/MeV);
210  }
211  }
212  }
213 }
void SetTrackStatus(const G4TrackStatus aTrackStatus)
G4ParticleDefinition * GetDefinition() const
double x() const
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
const G4ThreeVector & GetPosition() const
void hydrogenEnergy(G4double secondaryParticleKineticEnergy)
Energy distribution of the hydrogen (proton, d, t) particles after the phantom.
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
double z() const
void FillFragmentTuple(G4int A, G4double Z, G4double energy, G4double posX, G4double posY, G4double posZ)
Energy ntuple.
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
const G4String & GetName() const
G4SteppingManager * fpSteppingManager
static constexpr double cm
Definition: G4SIunits.hh:119
const G4String & GetParticleType() const
static IORTAnalysisManager * GetInstance()
std::vector< G4Track * > G4TrackVector
G4double energy(const ThreeVector &p, const G4double m)
double y() const
G4VPhysicalVolume * GetVolume() const
tuple z
Definition: test.py:28
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
G4Track * GetTrack() const
G4double GetPDGCharge() const
void heliumEnergy(G4double secondaryParticleKineticEnergy)
Energy distribution of the helium (He3 and alpha) particles after the phantom.

Here is the call graph for this function:


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