source: trunk/examples/advanced/hadrontherapy/src/HadrontherapySteppingAction.cc @ 1272

Last change on this file since 1272 was 1230, checked in by garnier, 15 years ago

update to geant4.9.3

File size: 8.1 KB
Line 
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// HadrontherapyProtonSteppingAction.cc;
27// See more at: http://g4advancedexamples.lngs.infn.it/Examples/hadrontherapy
28
29#include "G4SteppingManager.hh"
30#include "G4TrackVector.hh"
31#include "HadrontherapySteppingAction.hh"
32#include "G4ios.hh"
33#include "G4SteppingManager.hh"
34#include "G4Track.hh"
35#include "G4Step.hh"
36#include "G4StepPoint.hh"
37#include "G4TrackStatus.hh"
38#include "G4TrackVector.hh"
39#include "G4ParticleDefinition.hh"
40#include "G4ParticleTypes.hh"
41
42#include "HadrontherapyAnalysisManager.hh"
43
44#include "HadrontherapyRunAction.hh"
45
46/////////////////////////////////////////////////////////////////////////////
47HadrontherapySteppingAction::HadrontherapySteppingAction( HadrontherapyRunAction *run)
48{
49  runAction = run;
50}
51
52/////////////////////////////////////////////////////////////////////////////
53HadrontherapySteppingAction::~HadrontherapySteppingAction()
54{
55}
56
57/////////////////////////////////////////////////////////////////////////////
58void HadrontherapySteppingAction::UserSteppingAction(const G4Step* aStep)
59{ 
60    if( aStep->GetTrack()->GetVolume()->GetName() == "NewDetectorPhys"){
61#ifdef ANALYSIS_USE
62      G4ParticleDefinition *def = aStep->GetTrack()->GetDefinition();
63      G4double secondaryParticleKineticEnergy =  aStep->GetTrack()->GetKineticEnergy();     
64      G4String particleType = def->GetParticleType(); // particle type = nucleus for d, t, He3, alpha, and heavier nuclei
65      G4String particleName = def->GetParticleName(); // e.g. for alpha: the name = "alpha" and type = "nucleus"
66      if(particleType == "nucleus") {
67        G4int A = def->GetBaryonNumber();
68        G4double Z = def->GetPDGCharge();
69        G4double posX = aStep->GetTrack()->GetPosition().x() / cm;
70        G4double posY = aStep->GetTrack()->GetPosition().y() / cm;
71        G4double posZ = aStep->GetTrack()->GetPosition().z() / cm;
72        G4double energy = secondaryParticleKineticEnergy / A / MeV;
73
74        HadrontherapyAnalysisManager* analysisMgr =  HadrontherapyAnalysisManager::getInstance();   
75        analysisMgr->fillFragmentTuple(A, Z, energy, posX, posY, posZ);
76      } else if(particleName == "proton") {   // proton (hydrogen-1) is a special case
77        G4double posX = aStep->GetTrack()->GetPosition().x() / cm ;
78        G4double posY = aStep->GetTrack()->GetPosition().y() / cm ;
79        G4double posZ = aStep->GetTrack()->GetPosition().z() / cm ;
80        G4double energy = secondaryParticleKineticEnergy * MeV;    // Hydrogen-1: A = 1, Z = 1
81        HadrontherapyAnalysisManager::getInstance()->fillFragmentTuple(1, 1.0, energy, posX, posY, posZ);
82      }
83
84      G4String secondaryParticleName =  def -> GetParticleName(); 
85      //G4cout <<"Particle: " << secondaryParticleName << G4endl;
86      //G4cout <<"Energy: " << secondaryParticleKineticEnergy << G4endl;
87        HadrontherapyAnalysisManager* analysis =  HadrontherapyAnalysisManager::getInstance();   
88        //There is a bunch of stuff recorded with the energy 0, something should perhaps be done about this.
89        if(secondaryParticleName == "proton") {
90          analysis->hydrogenEnergy(secondaryParticleKineticEnergy / MeV);
91        }
92        if(secondaryParticleName == "deuteron") {
93          analysis->hydrogenEnergy((secondaryParticleKineticEnergy/2) / MeV);
94        }
95        if(secondaryParticleName == "triton") {
96          analysis->hydrogenEnergy((secondaryParticleKineticEnergy/3) / MeV);
97        }
98        if(secondaryParticleName == "alpha") {
99          analysis->heliumEnergy((secondaryParticleKineticEnergy/4) / MeV);
100        }
101        if(secondaryParticleName == "He3"){
102          analysis->heliumEnergy((secondaryParticleKineticEnergy/3) / MeV);             
103        }
104#endif
105
106        aStep->GetTrack()->SetTrackStatus(fKillTrackAndSecondaries);
107    }
108
109  // Electromagnetic and hadronic processes of primary particles in the phantom
110  //setting phantomPhys correctly will break something here fixme
111  if ((aStep -> GetTrack() -> GetTrackID() == 1) &&
112    (aStep -> GetTrack() -> GetVolume() -> GetName() == "PhantomPhys") &&
113    (aStep -> GetPostStepPoint() -> GetProcessDefinedStep() != NULL))
114            {
115              G4String process = aStep -> GetPostStepPoint() -> 
116                GetProcessDefinedStep() -> GetProcessName();
117 
118              if ((process == "Transportation") || (process == "StepLimiter")) {;}
119              else {
120                if ((process == "msc") || (process == "hLowEIoni") || (process == "hIoni")) 
121                  { 
122                    runAction -> AddEMProcess();
123                  } 
124                else 
125                 {
126                   runAction -> AddHadronicProcess();
127
128                   if ( (process != "LElastic") && (process != "ProtonInelastic") && (process != "hElastic") )
129                     G4cout << "Warning! Unknown proton process: "<< process << G4endl;
130                 }
131              }         
132            }
133
134 // Retrieve information about the secondaries originated in the phantom
135
136 G4SteppingManager*  steppingManager = fpSteppingManager;
137 
138  // check if it is alive
139  //if(theTrack-> GetTrackStatus() == fAlive) { return; }
140
141  // Retrieve the secondary particles
142  G4TrackVector* fSecondary = steppingManager -> GetfSecondary();
143     
144  for(size_t lp1=0;lp1<(*fSecondary).size(); lp1++)
145    { 
146      G4String volumeName = (*fSecondary)[lp1] -> GetVolume() -> GetName(); 
147 
148      if (volumeName == "phantomPhys")
149        {
150#ifdef ANALYSIS_USE   
151          G4String secondaryParticleName =  (*fSecondary)[lp1]->GetDefinition() -> GetParticleName(); 
152          G4double secondaryParticleKineticEnergy =  (*fSecondary)[lp1] -> GetKineticEnergy();     
153
154          HadrontherapyAnalysisManager* analysis =  HadrontherapyAnalysisManager::getInstance();   
155       
156          if (secondaryParticleName == "e-")
157            analysis -> electronEnergyDistribution(secondaryParticleKineticEnergy/MeV);
158                   
159          if (secondaryParticleName == "gamma")
160            analysis -> gammaEnergyDistribution(secondaryParticleKineticEnergy/MeV);
161
162          if (secondaryParticleName == "deuteron")
163            analysis -> deuteronEnergyDistribution(secondaryParticleKineticEnergy/MeV);
164       
165          if (secondaryParticleName == "triton")
166            analysis -> tritonEnergyDistribution(secondaryParticleKineticEnergy/MeV);
167           
168          if (secondaryParticleName == "alpha")
169            analysis -> alphaEnergyDistribution(secondaryParticleKineticEnergy/MeV);
170       
171          G4double z = (*fSecondary)[lp1]-> GetDynamicParticle() -> GetDefinition() -> GetPDGCharge();
172          if (z > 0.)
173            {     
174              G4int a = (*fSecondary)[lp1]-> GetDynamicParticle() -> GetDefinition() -> GetBaryonNumber();
175              G4int electronOccupancy = (*fSecondary)[lp1] ->  GetDynamicParticle() -> GetTotalOccupancy(); 
176             
177              // If a generic ion is originated in the detector, its baryonic number, PDG charge,
178              // total number of electrons in the orbitals are stored in a ntuple
179              analysis -> genericIonInformation(a, z, electronOccupancy, secondaryParticleKineticEnergy/MeV);
180            }
181#endif
182        }
183    }
184}
185
186
187
188
189
Note: See TracBrowser for help on using the repository browser.