Geant4  10.02.p01
SteppingAction.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 example is provided by the Geant4-DNA collaboration
27 // Any report or published results obtained using the Geant4-DNA software
28 // shall cite the following Geant4-DNA collaboration publication:
29 // Med. Phys. 37 (2010) 4692-4708
30 // The Geant4-DNA web site is available at http://geant4-dna.org
31 //
32 // If you use this example, please cite the following publication:
33 // Rad. Prot. Dos. 133 (2009) 2-11
34 
35 #include "G4SystemOfUnits.hh"
36 #include "G4SteppingManager.hh"
37 
38 #include "SteppingAction.hh"
39 #include "RunAction.hh"
40 #include "DetectorConstruction.hh"
41 #include "Analysis.hh"
42 
43 #include "G4Alpha.hh"
44 #include "G4Electron.hh"
45 
46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47 
49 :fRun(run),fDetector(det)
50 { }
51 
52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
53 
55 { }
56 
57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
58 
60 
61 {
62  // Analysis manager
63 
64  G4AnalysisManager* man = G4AnalysisManager::Instance();
65 
66  // Read phantom - Singleton
67 
69 
70  //
71 
72  // Material : 1 is cytoplasm, 2 is nucleus
73 
74  G4int matVoxelPRE = -1;
75  G4int matVoxelPOST = -1;
76  G4int tmp=-1;
77 
79 
80  if (tmp>0)
81  {
82  matVoxelPRE = fMyCellParameterisation->GetTissueType(tmp);
83  }
84 
86 
87  if (tmp>0)
88  {
89  matVoxelPOST = fMyCellParameterisation->GetTissueType(tmp);
90  }
91 
92 // COUNT GAS DETECTOR HITS
93 
97 
98  ||
102 
103  ||
104 
108 
109  )
110  {
111  fRun->AddNbOfHitsGas();
112  }
113 
114 // STOPPING POWER AND BEAM SPOT SIZE AT CELL ENTRANCE
115 
119 
120  ||
122  && (matVoxelPOST == 1)
124  )
125  {
126 
127  if( (aStep->GetPreStepPoint()->GetKineticEnergy() - aStep->GetPostStepPoint()->GetKineticEnergy() ) >0)
128  {
129  //Fill ntupleid=1
130  man->FillNtupleDColumn(1,0,aStep->GetPreStepPoint()->GetKineticEnergy()/keV);
131  man->FillNtupleDColumn(1,1,
132  (aStep->GetPreStepPoint()->GetKineticEnergy() -
133  aStep->GetPostStepPoint()->GetKineticEnergy())/
134  keV/(aStep->GetStepLength()/micrometer));
135  man->AddNtupleRow(1);
136  }
137 
138  // Average dE over step suggested by Michel Maire
139 
140  G4StepPoint* p1 = aStep->GetPreStepPoint();
141  G4ThreeVector coord1 = p1->GetPosition();
142  const G4AffineTransform transformation1 = p1->GetTouchable()->GetHistory()->GetTopTransform();
143  G4ThreeVector localPosition1 = transformation1.TransformPoint(coord1);
144 
145  G4StepPoint* p2 = aStep->GetPostStepPoint();
146  G4ThreeVector coord2 = p2->GetPosition();
147  const G4AffineTransform transformation2 = p2->GetTouchable()->GetHistory()->GetTopTransform();
148  G4ThreeVector localPosition2 = transformation2.TransformPoint(coord2);
149 
150  G4ThreeVector localPosition = localPosition1 + G4UniformRand()*(localPosition2-localPosition1);
151 
152  // end
153 
154  //Fill ntupleid=2
155  man->FillNtupleDColumn(2,0,localPosition.x()/micrometer);
156  man->FillNtupleDColumn(2,1,localPosition.y()/micrometer);
157  man->AddNtupleRow(2);
158  }
159 
160 // ALPHA RANGE
161 
162 if (
163 
165 
166  &&
167 
168  (aStep->GetTrack()->GetKineticEnergy()<1e-6)
169 
170  &&
171 
172  ( (matVoxelPOST==1)
174  || (matVoxelPOST==2) )
175 
176  )
177 
178  {
179  //Fill ntupleid=3
180  man->FillNtupleDColumn(3,0,
181  aStep->GetPostStepPoint()->GetPosition().x()/micrometer);
182  man->FillNtupleDColumn(3,1,
183  aStep->GetPostStepPoint()->GetPosition().y()/micrometer);
184  man->FillNtupleDColumn(3,2,
185  aStep->GetPostStepPoint()->GetPosition().z()/micrometer);
186  man->AddNtupleRow(3);
187  }
188 
189 // TOTAL DOSE DEPOSIT AND DOSE DEPOSIT WITHIN A PHANTOM VOXEL
190 // FOR ALL PARTICLES
191 
192 if (matVoxelPRE == 2)
193 
194  {
195  G4double dose = (aStep->GetTotalEnergyDeposit()/joule)/(fRun->GetMassNucleus()/kg);
196  fRun->AddDoseN(dose);
197 
198  G4ThreeVector v;
200  aStep->GetTotalEnergyDeposit()/eV);
201  }
202 
203 
204 if (matVoxelPRE == 1)
205 
206  {
208  fRun->AddDoseC(dose);
209 
210  G4ThreeVector v;
212  aStep->GetTotalEnergyDeposit()/eV);
213  }
214 
215 // PROTECTION AGAINST MSC LOOPS FOR e-
216 
217 if ( aStep->GetTotalEnergyDeposit()/MeV<1e-25
219 {
221  /*
222  G4cout << "*** Warning *** : msc loop for "
223  << aStep->GetTrack()->GetDefinition()->GetParticleName()
224  << " in " <<
225  aStep->GetPostStepPoint()->GetTouchableHandle()->GetVolume()->GetName() << G4endl;
226  */
227 }
228 
229 }
void SetTrackStatus(const G4TrackStatus aTrackStatus)
G4ParticleDefinition * GetDefinition() const
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static const double MeV
Definition: G4SIunits.hh:211
void AddDoseN(G4double dose)
Definition: RunAction.hh:60
CLHEP::Hep3Vector G4ThreeVector
G4double GetStepLength() const
const G4DynamicParticle * GetDynamicParticle() const
G4LogicalVolume * GetLogicalCollDetYoke()
virtual const G4NavigationHistory * GetHistory() const
Definition: G4VTouchable.cc:86
static const double joule
Definition: G4SIunits.hh:201
G4ParticleDefinition * GetDefinition() const
const G4VTouchable * GetTouchable() const
void UserSteppingAction(const G4Step *)
int G4int
Definition: G4Types.hh:78
G4StepPoint * GetPreStepPoint() const
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4LogicalVolume * GetLogicalIsobutane()
const G4ThreeVector & GetPosition() const
static const double kg
Definition: G4SIunits.hh:179
ExG4HbookAnalysisManager G4AnalysisManager
Definition: g4hbook_defs.hh:66
Definition: G4Step.hh:76
void AddNbOfHitsGas()
Definition: RunAction.hh:72
DetectorConstruction * fDetector
G4double GetTotalEnergyDeposit() const
static const double micrometer
Definition: G4SIunits.hh:99
static const double eV
Definition: G4SIunits.hh:212
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4LogicalVolume * GetLogicalVolume() const
CellParameterisation * fMyCellParameterisation
G4double GetMassNucleus()
Definition: RunAction.hh:75
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
void AddDoseC(G4double dose)
Definition: RunAction.hh:64
G4StepPoint * GetPostStepPoint() const
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:58
G4LogicalVolume * GetLogicalKgm()
Detector construction class to demonstrate various ways of placement.
const G4AffineTransform & GetTopTransform() const
static const double keV
Definition: G4SIunits.hh:213
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
G4LogicalVolume * GetLogicalCollDetGap4()
void AddDoseBox(G4int i, G4double x)
Definition: RunAction.hh:80
Run action class.
Definition: RunAction.hh:45
G4LogicalVolume * GetLogicalPolyprop()
G4Track * GetTrack() const
G4double GetMassCytoplasm()
Definition: RunAction.hh:78
static CellParameterisation * Instance()
RunAction * fRun
static G4Alpha * AlphaDefinition()
Definition: G4Alpha.cc:84
const G4TouchableHandle & GetTouchableHandle() const