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

Last change on this file since 1332 was 1313, checked in by garnier, 14 years ago

geant4.9.4 beta rc0

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