Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
advanced/microbeam/src/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 
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46 
48 :fRun(run),fDetector(det)
49 { }
50 
51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
52 
54 { }
55 
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
57 
59 
60 {
61  G4AnalysisManager* man = G4AnalysisManager::Instance();
62 
63 // COUNT GAS DETECTOR HITS
64 
68 
69  ||
73 
74  ||
75 
79 
80  )
81  {
82  fRun->AddNbOfHitsGas();
83  }
84 
85 // STOPPING POWER AND BEAM SPOT SIZE AT CELL ENTRANCE
86 
87 if ( ((aStep->GetPreStepPoint()->GetTouchableHandle()->GetVolume()->GetLogicalVolume() == fDetector->GetLogicalPolyprop())
88  && (aStep->GetPostStepPoint()->GetTouchableHandle()->GetVolume()->GetLogicalVolume() == fDetector->GetLogicalKgm())
90 
91  ||
93  && (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "physicalCytoplasm")
95  )
96  {
97 
98  if( (aStep->GetPreStepPoint()->GetKineticEnergy() - aStep->GetPostStepPoint()->GetKineticEnergy() ) >0)
99  {
100  //Fill ntupleid=1
101  man->FillNtupleDColumn(1,0,aStep->GetPreStepPoint()->GetKineticEnergy()/keV);
102  man->FillNtupleDColumn(1,1,
103  (aStep->GetPreStepPoint()->GetKineticEnergy() -
104  aStep->GetPostStepPoint()->GetKineticEnergy())/
105  keV/(aStep->GetStepLength()/micrometer));
106  man->AddNtupleRow(1);
107  }
108 
109  // Average dE over step suggested by Michel Maire
110 
111  G4StepPoint* p1 = aStep->GetPreStepPoint();
112  G4ThreeVector coord1 = p1->GetPosition();
113  const G4AffineTransform transformation1 = p1->GetTouchable()->GetHistory()->GetTopTransform();
114  G4ThreeVector localPosition1 = transformation1.TransformPoint(coord1);
115 
116  G4StepPoint* p2 = aStep->GetPostStepPoint();
117  G4ThreeVector coord2 = p2->GetPosition();
118  const G4AffineTransform transformation2 = p2->GetTouchable()->GetHistory()->GetTopTransform();
119  G4ThreeVector localPosition2 = transformation2.TransformPoint(coord2);
120 
121  G4ThreeVector localPosition = localPosition1 + G4UniformRand()*(localPosition2-localPosition1);
122 
123  // end
124  //Fill ntupleid=2
125  man->FillNtupleDColumn(2,0,localPosition.x()/micrometer);
126  man->FillNtupleDColumn(2,1,localPosition.y()/micrometer);
127  man->AddNtupleRow(2);
128  }
129 
130 // ALPHA RANGE
131 
132 if (
133 
135 
136  &&
137 
138  (aStep->GetTrack()->GetKineticEnergy()<1e-6)
139 
140  &&
141 
142  ( (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "physicalCytoplasm")
143  || (aStep->GetPostStepPoint()->GetTouchableHandle()->GetVolume()->GetLogicalVolume() == fDetector->GetLogicalKgm())
144  || (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "physicalNucleus") )
145 
146  )
147 
148  {
149  //Fill ntupleid=3
150  man->FillNtupleDColumn(3,0,
151  aStep->GetPostStepPoint()->GetPosition().x()/micrometer);
152  man->FillNtupleDColumn(3,1,
153  aStep->GetPostStepPoint()->GetPosition().y()/micrometer);
154  man->FillNtupleDColumn(3,2,
155  aStep->GetPostStepPoint()->GetPosition().z()/micrometer);
156  man->AddNtupleRow(3);
157  }
158 
159 // TOTAL DOSE DEPOSIT AND DOSE DEPOSIT WITHIN A PHANTOM VOXEL
160 // FOR ALL PARTICLES
161 
162 if (aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physicalNucleus")
163 
164  {
165  G4double dose = (aStep->GetTotalEnergyDeposit()/joule)/(fRun->GetMassNucleus()/kg);
166  fRun->AddDoseN(dose);
167 
170  aStep->GetTotalEnergyDeposit()/eV);
171  }
172 
173 
174 if (aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physicalCytoplasm")
175 
176  {
177  G4double dose = (aStep->GetTotalEnergyDeposit()/joule)/(fRun->GetMassCytoplasm()/kg);
178  fRun->AddDoseC(dose);
179 
182  aStep->GetTotalEnergyDeposit()/eV);
183  }
184 }
G4double GetStepLength() const
double x() const
const G4DynamicParticle * GetDynamicParticle() const
virtual const G4NavigationHistory * GetHistory() const
Definition: G4VTouchable.cc:86
G4ParticleDefinition * GetDefinition() const
const G4VTouchable * GetTouchable() const
void UserSteppingAction(const G4Step *)
double z() const
G4StepPoint * GetPreStepPoint() const
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4VPhysicalVolume * GetPhysicalVolume() const
const G4String & GetName() const
const G4ThreeVector & GetPosition() const
ExG4HbookAnalysisManager G4AnalysisManager
Definition: g4hbook_defs.hh:46
Definition: G4Step.hh:76
G4double GetTotalEnergyDeposit() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
G4StepPoint * GetPostStepPoint() const
double y() const
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:58
const G4AffineTransform & GetTopTransform() const
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
G4Track * GetTrack() const
int micrometer
Definition: hepunit.py:34
static G4Alpha * AlphaDefinition()
Definition: G4Alpha.cc:84
const G4TouchableHandle & GetTouchableHandle() const