Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MicrobeamSteppingAction.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 // -------------------------------------------------------------------
27 // $Id$
28 // -------------------------------------------------------------------
29 
30 #include "G4SystemOfUnits.hh"
31 #include "G4SteppingManager.hh"
32 
34 #include "MicrobeamRunAction.hh"
36 
37 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
38 
41 :Run(run),Detector(det),Histo(his)
42 { }
43 
44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45 
47 { }
48 
49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
50 
52 
53 {
54 
55 // COUNT GAS DETECTOR HITS
56 
57 if ( ((aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "_CollDet_yoke_")
58  && (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "Isobutane")
59  && (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha"))
60 
61  ||
62  ((aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "_CollDet_gap4_")
63  && (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "Isobutane")
64  && (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha"))
65 
66  ||
67 
68  ((aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "_CollDet_gap5_")
69  && (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "Isobutane")
70  && (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha"))
71 
72  )
73  {
74  Run->AddNbOfHitsGas();
75  }
76 
77 // STOPPING POWER AND BEAM SPOT SIZE AT CELL ENTRANCE
78 
79 if ( ((aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "Polyprop")
80  && (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "KGM")
81  && (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha"))
82 
83  ||
84  ((aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "Polyprop")
85  && (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "physicalCytoplasm")
86  && (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha"))
87  )
88  {
89 
90  if( (aStep->GetPreStepPoint()->GetKineticEnergy() - aStep->GetPostStepPoint()->GetKineticEnergy() ) >0)
91  {
92  Histo->FillNtuple(0,0,aStep->GetPreStepPoint()->GetKineticEnergy()/keV);
93  Histo->FillNtuple(0,1,
94  (aStep->GetPreStepPoint()->GetKineticEnergy() -
96  Histo->AddRowNtuple(0);
97  }
98 
99  // Average dE over step suggested by Michel Maire
100 
101  G4StepPoint* p1 = aStep->GetPreStepPoint();
102  G4ThreeVector coord1 = p1->GetPosition();
103  const G4AffineTransform transformation1 = p1->GetTouchable()->GetHistory()->GetTopTransform();
104  G4ThreeVector localPosition1 = transformation1.TransformPoint(coord1);
105 
106  G4StepPoint* p2 = aStep->GetPostStepPoint();
107  G4ThreeVector coord2 = p2->GetPosition();
108  const G4AffineTransform transformation2 = p2->GetTouchable()->GetHistory()->GetTopTransform();
109  G4ThreeVector localPosition2 = transformation2.TransformPoint(coord2);
110 
111  G4ThreeVector localPosition = localPosition1 + G4UniformRand()*(localPosition2-localPosition1);
112 
113  // end
114 
115  Histo->FillNtuple(1,0,localPosition.x()/micrometer);
116  Histo->FillNtuple(1,1,localPosition.y()/micrometer);
117  Histo->AddRowNtuple(1);
118 
119  }
120 
121 // ALPHA RANGE
122 
123 if (
124 
125  (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha")
126 
127  &&
128 
129  (aStep->GetTrack()->GetKineticEnergy()<1e-6)
130 
131  &&
132 
133  ( (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "physicalCytoplasm")
134  || (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "KGM")
135  || (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "physicalNucleus") )
136 
137  )
138 
139  {
140  Histo->FillNtuple(2,0,aStep->GetPostStepPoint()->GetPosition().x()/micrometer);
141  Histo->FillNtuple(2,1,aStep->GetPostStepPoint()->GetPosition().y()/micrometer);
142  Histo->FillNtuple(2,2,aStep->GetPostStepPoint()->GetPosition().z()/micrometer);
143  Histo->AddRowNtuple(2);
144  }
145 
146 // TOTAL DOSE DEPOSIT AND DOSE DEPOSIT WITHIN A PHANTOM VOXEL
147 // FOR ALL PARTICLES
148 
149 if (aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physicalNucleus")
150 
151  {
153  Run->AddDoseN(dose);
154 
157  aStep->GetTotalEnergyDeposit()/eV);
158  }
159 
160 
161 if (aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physicalCytoplasm")
162 
163  {
165  Run->AddDoseC(dose);
166 
169  aStep->GetTotalEnergyDeposit()/eV);
170  }
171 }