Geant4  10.00.p02
HadrontherapySteppingAction.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // This is the *BASIC* version of Hadrontherapy, a Geant4-based application
27 // See more at: http://g4advancedexamples.lngs.infn.it/Examples/hadrontherapy
28 //
29 // Visit the Hadrontherapy web site (http://www.lns.infn.it/link/Hadrontherapy) to request
30 // the *COMPLETE* version of this program, together with its documentation;
31 // Hadrontherapy (both basic and full version) are supported by the Italian INFN
32 // Institute in the framework of the MC-INFN Group
33 //
34 
35 #include "G4SteppingManager.hh"
36 #include "G4TrackVector.hh"
38 #include "G4ios.hh"
39 #include "G4SteppingManager.hh"
40 #include "G4Track.hh"
41 #include "G4Step.hh"
42 #include "G4StepPoint.hh"
43 #include "G4TrackStatus.hh"
44 #include "G4TrackVector.hh"
45 #include "G4ParticleDefinition.hh"
46 #include "G4ParticleTypes.hh"
47 #include "G4UserEventAction.hh"
49 #include "G4VSensitiveDetector.hh"
52 #include "G4SystemOfUnits.hh"
53 
56 {
57  runAction = run;
58 }
59 
62 {
63 }
64 
67 {
68  // G4TransportationManager* tManager = G4TransportationManager::GetTransportationManager();
69  //G4VPhysicalVolume* pW = tManager->GetParallelWorld ("DetectorROGeometry");
70 
71  //G4Navigator* gNav = tManager->GetNavigator(pW);
72 
73  // G4VPhysicalVolume* currentVol = gNav->LocateGlobalPointAndSetup(aStep->GetTrack()->GetPosition());
74 
75  // G4cout << "Step: " << currentVol->GetName() << " " << aStep->GetTrack()->GetPosition() << G4endl;
76  //G4cout << "G4LogicalVolume: = " << currentVol->GetLogicalVolume()->GetName() << G4endl;
77  //if (currentVol->GetLogicalVolume()->GetSensitiveDetector())
78  // G4cout << "Sensitive Detector: " << currentVol->GetLogicalVolume()->GetSensitiveDetector()->GetName() << G4endl;
79 
80 
81  /*
82  // USEFULL METHODS TO RETRIEVE INFORMATION DURING THE STEPS
83  if( (aStep->GetTLocateGlobalPointAndSeturack()->GetVolume()->GetName() == "DetectorPhys")
84  && aStep->GetTrack()->GetDefinition()->GetParticleName() == "proton")
85  //G4int evtNb = G4RunManager::GetRunManager()->GetCurrentEvent() -> GetEventID();
86  {
87  G4cout << "ENERGIA: " << aStep->GetTrack()->GetKineticEnergy()
88  << " VOLUME " << aStep->GetTrack()->GetVolume()->GetName()
89  << " MATERIALE " << aStep -> GetTrack() -> GetMaterial() -> GetName()
90  << " EVENTO " << G4RunManager::GetRunManager()->GetCurrentEvent() -> GetEventID()
91  << " POS " << aStep->GetTrack()->GetPosition().x()
92  << G4endl;
93  }
94  */
95 
96  if( aStep->GetTrack()->GetVolume()->GetName() == "NewDetectorPhys"){
97 #ifdef G4ANALYSIS_USE_ROOT
98  G4ParticleDefinition *def = aStep->GetTrack()->GetDefinition();
99  G4double secondaryParticleKineticEnergy = aStep->GetTrack()->GetKineticEnergy();
100  G4String particleType = def->GetParticleType(); // particle type = nucleus for d, t, He3, alpha, and heavier nuclei
101  G4String particleName = def->GetParticleName(); // e.g. for alpha: the name = "alpha" and type = "nucleus"
102  if(particleType == "nucleus") {
103  G4int A = def->GetBaryonNumber();
104  G4double Z = def->GetPDGCharge();
105  G4double posX = aStep->GetTrack()->GetPosition().x() / cm;
106  G4double posY = aStep->GetTrack()->GetPosition().y() / cm;
107  G4double posZ = aStep->GetTrack()->GetPosition().z() / cm;
108  G4double energy = secondaryParticleKineticEnergy / A / MeV;
109 
111  analysisMgr->FillFragmentTuple(A, Z, energy, posX, posY, posZ);
112  } else if(particleName == "proton") { // proton (hydrogen-1) is a special case
113  G4double posX = aStep->GetTrack()->GetPosition().x() / cm ;
114  G4double posY = aStep->GetTrack()->GetPosition().y() / cm ;
115  G4double posZ = aStep->GetTrack()->GetPosition().z() / cm ;
116  G4double energy = secondaryParticleKineticEnergy * MeV; // Hydrogen-1: A = 1, Z = 1
117  HadrontherapyAnalysisManager::GetInstance()->FillFragmentTuple(1, 1.0, energy, posX, posY, posZ);
118  }
119 
120  G4String secondaryParticleName = def -> GetParticleName();
121  //G4cout <<"Particle: " << secondaryParticleName << G4endl;
122  //G4cout <<"Energy: " << secondaryParticleKineticEnergy << G4endl;
124  //There is a bunch of stuff recorded with the energy 0, something should perhaps be done about this.
125  if(secondaryParticleName == "proton") {
126  analysis->hydrogenEnergy(secondaryParticleKineticEnergy / MeV);
127  }
128  if(secondaryParticleName == "deuteron") {
129  analysis->hydrogenEnergy((secondaryParticleKineticEnergy/2) / MeV);
130  }
131  if(secondaryParticleName == "triton") {
132  analysis->hydrogenEnergy((secondaryParticleKineticEnergy/3) / MeV);
133  }
134  if(secondaryParticleName == "alpha") {
135  analysis->heliumEnergy((secondaryParticleKineticEnergy/4) / MeV);
136  }
137  if(secondaryParticleName == "He3"){
138  analysis->heliumEnergy((secondaryParticleKineticEnergy/3) / MeV);
139  }
140 #endif
141 
143  }
144 
145  // Electromagnetic and hadronic processes of primary particles in the phantom
146  //setting phantomPhys correctly will break something here fixme
147  if ((aStep -> GetTrack() -> GetTrackID() == 1) &&
148  (aStep -> GetTrack() -> GetVolume() -> GetName() == "PhantomPhys") &&
149  (aStep -> GetPostStepPoint() -> GetProcessDefinedStep() != NULL))
150  {
151  G4String process = aStep -> GetPostStepPoint() ->
152  GetProcessDefinedStep() -> GetProcessName();
153 
154  if ((process == "Transportation") || (process == "StepLimiter")) {;}
155  else
156  {
157  if ((process == "msc") || (process == "hLowEIoni") || (process == "hIoni"))
158  {
159  runAction -> AddEMProcess();
160  }
161  else
162  {
163  runAction -> AddHadronicProcess();
164 
165  if ( (process != "LElastic") && (process != "ProtonInelastic") && (process != "hElastic") )
166  G4cout << "Warning! Unknown proton process: "<< process << G4endl;
167  }
168  }
169  }
170 
171  // Retrieve information about the secondary particles originated in the phantom
172 
173  G4SteppingManager* steppingManager = fpSteppingManager;
174 
175  // check if it is alive
176  //if(theTrack-> GetTrackStatus() == fAlive) { return; }
177 
178  // Retrieve the secondary particles
179  G4TrackVector* fSecondary = steppingManager -> GetfSecondary();
180 
181  for(size_t lp1=0;lp1<(*fSecondary).size(); lp1++)
182  {
183  G4String volumeName = (*fSecondary)[lp1] -> GetVolume() -> GetName();
184 
185  if (volumeName == "phantomPhys")
186  {
187 #ifdef G4ANALYSIS_USE_ROOT
188  G4String secondaryParticleName = (*fSecondary)[lp1]->GetDefinition() -> GetParticleName();
189  G4double secondaryParticleKineticEnergy = (*fSecondary)[lp1] -> GetKineticEnergy();
190 
192 
193  if (secondaryParticleName == "e-")
194  analysis -> electronEnergyDistribution(secondaryParticleKineticEnergy/MeV);
195 
196  if (secondaryParticleName == "gamma")
197  analysis -> gammaEnergyDistribution(secondaryParticleKineticEnergy/MeV);
198 
199  if (secondaryParticleName == "deuteron")
200  analysis -> deuteronEnergyDistribution(secondaryParticleKineticEnergy/MeV);
201 
202  if (secondaryParticleName == "triton")
203  analysis -> tritonEnergyDistribution(secondaryParticleKineticEnergy/MeV);
204 
205  if (secondaryParticleName == "alpha")
206  analysis -> alphaEnergyDistribution(secondaryParticleKineticEnergy/MeV);
207 
208  G4double z = (*fSecondary)[lp1]-> GetDynamicParticle() -> GetDefinition() -> GetPDGCharge();
209  if (z > 0.)
210  {
211  G4int a = (*fSecondary)[lp1]-> GetDynamicParticle() -> GetDefinition() -> GetBaryonNumber();
212  G4int electronOccupancy = (*fSecondary)[lp1] -> GetDynamicParticle() -> GetTotalOccupancy();
213 
214  // If a generic ion is originated in the detector, its baryonic number, PDG charge,
215  // total number of electrons in the orbitals are stored in a ntuple
216  analysis -> genericIonInformation(a, z, electronOccupancy, secondaryParticleKineticEnergy/MeV);
217  }
218 #endif
219  }
220  }
221 }
222 
223 
224 
225 
226 
static const double cm
Definition: G4SIunits.hh:106
void SetTrackStatus(const G4TrackStatus aTrackStatus)
G4ParticleDefinition * GetDefinition() const
static HadrontherapyAnalysisManager * GetInstance()
Get the pointer to the analysis manager.
static const double MeV
Definition: G4SIunits.hh:193
HadrontherapySteppingAction(HadrontherapyRunAction *)
G4double z
Definition: TRTMaterials.hh:39
const G4ThreeVector & GetPosition() const
G4double a
Definition: TRTMaterials.hh:39
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
G4SteppingManager * fpSteppingManager
const G4String & GetParticleType() const
Definition: G4Step.hh:76
static const G4double A[nN]
std::vector< G4Track * > G4TrackVector
G4double energy(const ThreeVector &p, const G4double m)
G4VPhysicalVolume * GetVolume() const
A class for connecting the simulation to an analysis package.
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4Track * GetTrack() const
G4double GetPDGCharge() const