Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 
36 
37 #include "G4SteppingManager.hh"
38 #include "G4TrackVector.hh"
40 #include "G4ios.hh"
41 #include "G4SteppingManager.hh"
42 #include "G4Track.hh"
43 #include "G4Step.hh"
44 #include "G4StepPoint.hh"
45 #include "G4TrackStatus.hh"
46 #include "G4TrackVector.hh"
47 #include "G4ParticleDefinition.hh"
48 #include "G4ParticleTypes.hh"
49 #include "G4UserEventAction.hh"
50 
52 
54 
57 {
58  runAction = run;
59 }
60 
63 {
64 }
65 
68 {
69  /*
70  // USEFULL METHODS TO RETRIEVE INFORMATION DURING THE STEPS
71  if( (aStep->GetTrack()->GetVolume()->GetName() == "DetectorPhys")
72  && aStep->GetTrack()->GetDefinition()->GetParticleName() == "proton")
73  //G4int evtNb = G4RunManager::GetRunManager()->GetCurrentEvent() -> GetEventID();
74  {
75  G4cout << "ENERGIA: " << aStep->GetTrack()->GetKineticEnergy()
76  << " VOLUME " << aStep->GetTrack()->GetVolume()->GetName()
77  << " MATERIALE " << aStep -> GetTrack() -> GetMaterial() -> GetName()
78  << " EVENTO " << G4RunManager::GetRunManager()->GetCurrentEvent() -> GetEventID()
79  << " POS " << aStep->GetTrack()->GetPosition().x()
80  << G4endl;
81  }
82  */
83 
84  if( aStep->GetTrack()->GetVolume()->GetName() == "NewDetectorPhys"){
85 #ifdef G4ANALYSIS_USE_ROOT
86  G4ParticleDefinition *def = aStep->GetTrack()->GetDefinition();
87  G4double secondaryParticleKineticEnergy = aStep->GetTrack()->GetKineticEnergy();
88  G4String particleType = def->GetParticleType(); // particle type = nucleus for d, t, He3, alpha, and heavier nuclei
89  G4String particleName = def->GetParticleName(); // e.g. for alpha: the name = "alpha" and type = "nucleus"
90  if(particleType == "nucleus") {
91  G4int A = def->GetBaryonNumber();
92  G4double Z = def->GetPDGCharge();
93  G4double posX = aStep->GetTrack()->GetPosition().x() / CLHEP::cm;
94  G4double posY = aStep->GetTrack()->GetPosition().y() / CLHEP::cm;
95  G4double posZ = aStep->GetTrack()->GetPosition().z() / CLHEP::cm;
96  G4double energy = secondaryParticleKineticEnergy / A / CLHEP::MeV;
97 
99  analysisMgr->FillFragmentTuple(A, Z, energy, posX, posY, posZ);
100  } else if(particleName == "proton") { // proton (hydrogen-1) is a special case
101  G4double posX = aStep->GetTrack()->GetPosition().x() / CLHEP::cm ;
102  G4double posY = aStep->GetTrack()->GetPosition().y() / CLHEP::cm ;
103  G4double posZ = aStep->GetTrack()->GetPosition().z() / CLHEP::cm ;
104  G4double energy = secondaryParticleKineticEnergy * CLHEP::MeV; // Hydrogen-1: A = 1, Z = 1
105  HadrontherapyAnalysisManager::GetInstance()->FillFragmentTuple(1, 1.0, energy, posX, posY, posZ);
106  }
107 
108  G4String secondaryParticleName = def -> GetParticleName();
109  //G4cout <<"Particle: " << secondaryParticleName << G4endl;
110  //G4cout <<"Energy: " << secondaryParticleKineticEnergy << G4endl;
112  //There is a bunch of stuff recorded with the energy 0, something should perhaps be done about this.
113  if(secondaryParticleName == "proton") {
114  analysis->hydrogenEnergy(secondaryParticleKineticEnergy / CLHEP::MeV);
115  }
116  if(secondaryParticleName == "deuteron") {
117  analysis->hydrogenEnergy((secondaryParticleKineticEnergy/2) / CLHEP::MeV);
118  }
119  if(secondaryParticleName == "triton") {
120  analysis->hydrogenEnergy((secondaryParticleKineticEnergy/3) / CLHEP::MeV);
121  }
122  if(secondaryParticleName == "alpha") {
123  analysis->heliumEnergy((secondaryParticleKineticEnergy/4) / CLHEP::MeV);
124  }
125  if(secondaryParticleName == "He3"){
126  analysis->heliumEnergy((secondaryParticleKineticEnergy/3) / CLHEP::MeV);
127  }
128 #endif
129 
131  }
132 
133  // Electromagnetic and hadronic processes of primary particles in the phantom
134  //setting phantomPhys correctly will break something here fixme
135  if ((aStep -> GetTrack() -> GetTrackID() == 1) &&
136  (aStep -> GetTrack() -> GetVolume() -> GetName() == "PhantomPhys") &&
137  (aStep -> GetPostStepPoint() -> GetProcessDefinedStep() != NULL))
138  {
139  G4String process = aStep -> GetPostStepPoint() ->
140  GetProcessDefinedStep() -> GetProcessName();
141 
142  if ((process == "Transportation") || (process == "StepLimiter")) {;}
143  else
144  {
145  if ((process == "msc") || (process == "hLowEIoni") || (process == "hIoni"))
146  {
147  runAction -> AddEMProcess();
148  }
149  else
150  {
151  runAction -> AddHadronicProcess();
152 
153  if ( (process != "LElastic") && (process != "ProtonInelastic") && (process != "hElastic") )
154  G4cout << "Warning! Unknown proton process: "<< process << G4endl;
155  }
156  }
157  }
158 
159  // Retrieve information about the secondary particles originated in the phantom
160 
161  G4SteppingManager* steppingManager = fpSteppingManager;
162 
163  // check if it is alive
164  //if(theTrack-> GetTrackStatus() == fAlive) { return; }
165 
166  // Retrieve the secondary particles
167  G4TrackVector* fSecondary = steppingManager -> GetfSecondary();
168 
169  for(size_t lp1=0;lp1<(*fSecondary).size(); lp1++)
170  {
171  G4String volumeName = (*fSecondary)[lp1] -> GetVolume() -> GetName();
172 
173  if (volumeName == "phantomPhys")
174  {
175 #ifdef G4ANALYSIS_USE_ROOT
176  G4String secondaryParticleName = (*fSecondary)[lp1]->GetDefinition() -> GetParticleName();
177  G4double secondaryParticleKineticEnergy = (*fSecondary)[lp1] -> GetKineticEnergy();
178 
180 
181  if (secondaryParticleName == "e-")
182  analysis -> electronEnergyDistribution(secondaryParticleKineticEnergy/CLHEP::MeV);
183 
184  if (secondaryParticleName == "gamma")
185  analysis -> gammaEnergyDistribution(secondaryParticleKineticEnergy/CLHEP::MeV);
186 
187  if (secondaryParticleName == "deuteron")
188  analysis -> deuteronEnergyDistribution(secondaryParticleKineticEnergy/CLHEP::MeV);
189 
190  if (secondaryParticleName == "triton")
191  analysis -> tritonEnergyDistribution(secondaryParticleKineticEnergy/CLHEP::MeV);
192 
193  if (secondaryParticleName == "alpha")
194  analysis -> alphaEnergyDistribution(secondaryParticleKineticEnergy/CLHEP::MeV);
195 
196  G4double z = (*fSecondary)[lp1]-> GetDynamicParticle() -> GetDefinition() -> GetPDGCharge();
197  if (z > 0.)
198  {
199  G4int a = (*fSecondary)[lp1]-> GetDynamicParticle() -> GetDefinition() -> GetBaryonNumber();
200  G4int electronOccupancy = (*fSecondary)[lp1] -> GetDynamicParticle() -> GetTotalOccupancy();
201 
202  // If a generic ion is originated in the detector, its baryonic number, PDG charge,
203  // total number of electrons in the orbitals are stored in a ntuple
204  analysis -> genericIonInformation(a, z, electronOccupancy, secondaryParticleKineticEnergy/CLHEP::MeV);
205  }
206 #endif
207  }
208  }
209 }
210 
211 
212 
213 
214