source: trunk/examples/advanced/microbeam/src/MicrobeamSteppingAction.cc@ 1354

Last change on this file since 1354 was 1342, checked in by garnier, 15 years ago

update ti head

  • Property svn:executable set to *
File size: 6.8 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// -------------------------------------------------------------------
27// $Id: MicrobeamSteppingAction.cc,v 1.10 2010/10/07 14:03:11 sincerti Exp $
28// -------------------------------------------------------------------
29
30#include "G4SteppingManager.hh"
31
32#include "MicrobeamSteppingAction.hh"
33#include "MicrobeamRunAction.hh"
34#include "MicrobeamDetectorConstruction.hh"
35
36//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
37
38MicrobeamSteppingAction::MicrobeamSteppingAction(MicrobeamRunAction* run,MicrobeamDetectorConstruction* det,
39MicrobeamHistoManager* his)
40:Run(run),Detector(det),Histo(his)
41{ }
42
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
44
45MicrobeamSteppingAction::~MicrobeamSteppingAction()
46{ }
47
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49
50void MicrobeamSteppingAction::UserSteppingAction(const G4Step* aStep)
51
52{
53
54// COUNT GAS DETECTOR HITS
55
56if ( ((aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "_CollDet_yoke_")
57 && (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "Isobutane")
58 && (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha"))
59
60 ||
61 ((aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "_CollDet_gap4_")
62 && (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "Isobutane")
63 && (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha"))
64
65 ||
66
67 ((aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "_CollDet_gap5_")
68 && (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "Isobutane")
69 && (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha"))
70
71 )
72 {
73 Run->AddNbOfHitsGas();
74 }
75
76// STOPPING POWER AND BEAM SPOT SIZE AT CELL ENTRANCE
77
78if ( ((aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "Polyprop")
79 && (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "KGM")
80 && (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha"))
81
82 ||
83 ((aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "Polyprop")
84 && (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "physicalCytoplasm")
85 && (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha"))
86 )
87 {
88
89 if( (aStep->GetPreStepPoint()->GetKineticEnergy() - aStep->GetPostStepPoint()->GetKineticEnergy() ) >0)
90 {
91 Histo->FillNtuple(0,0,aStep->GetPreStepPoint()->GetKineticEnergy()/keV);
92 Histo->FillNtuple(0,1,
93 (aStep->GetPreStepPoint()->GetKineticEnergy() -
94 aStep->GetPostStepPoint()->GetKineticEnergy())/ keV/(aStep->GetStepLength()/micrometer));
95 Histo->AddRowNtuple(0);
96 }
97
98 // Average dE over step suggested by Michel Maire
99
100 G4StepPoint* p1 = aStep->GetPreStepPoint();
101 G4ThreeVector coord1 = p1->GetPosition();
102 const G4AffineTransform transformation1 = p1->GetTouchable()->GetHistory()->GetTopTransform();
103 G4ThreeVector localPosition1 = transformation1.TransformPoint(coord1);
104
105 G4StepPoint* p2 = aStep->GetPostStepPoint();
106 G4ThreeVector coord2 = p2->GetPosition();
107 const G4AffineTransform transformation2 = p2->GetTouchable()->GetHistory()->GetTopTransform();
108 G4ThreeVector localPosition2 = transformation2.TransformPoint(coord2);
109
110 G4ThreeVector localPosition = localPosition1 + G4UniformRand()*(localPosition2-localPosition1);
111
112 // end
113
114 Histo->FillNtuple(1,0,localPosition.x()/micrometer);
115 Histo->FillNtuple(1,1,localPosition.y()/micrometer);
116 Histo->AddRowNtuple(1);
117
118 }
119
120// ALPHA RANGE
121
122if (
123
124 (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha")
125
126 &&
127
128 (aStep->GetTrack()->GetKineticEnergy()<1e-6)
129
130 &&
131
132 ( (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "physicalCytoplasm")
133 || (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "KGM")
134 || (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "physicalNucleus") )
135
136 )
137
138 {
139 Histo->FillNtuple(2,0,aStep->GetPostStepPoint()->GetPosition().x()/micrometer);
140 Histo->FillNtuple(2,1,aStep->GetPostStepPoint()->GetPosition().y()/micrometer);
141 Histo->FillNtuple(2,2,aStep->GetPostStepPoint()->GetPosition().z()/micrometer);
142 Histo->AddRowNtuple(2);
143 }
144
145// TOTAL DOSE DEPOSIT AND DOSE DEPOSIT WITHIN A PHANTOM VOXEL
146// FOR ALL PARTICLES
147
148if (aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physicalNucleus")
149
150 {
151 G4double dose = (aStep->GetTotalEnergyDeposit()/joule)/(Run->GetMassNucleus()/kg);
152 Run->AddDoseN(dose);
153
154 G4ThreeVector v;
155 Run->AddDoseBox(aStep->GetPreStepPoint()->GetTouchableHandle()->GetReplicaNumber(),
156 aStep->GetTotalEnergyDeposit()/eV);
157 }
158
159
160if (aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physicalCytoplasm")
161
162 {
163 G4double dose = (aStep->GetTotalEnergyDeposit()/joule)/(Run->GetMassCytoplasm()/kg);
164 Run->AddDoseC(dose);
165
166 G4ThreeVector v;
167 Run->AddDoseBox(aStep->GetPreStepPoint()->GetTouchableHandle()->GetReplicaNumber(),
168 aStep->GetTotalEnergyDeposit()/eV);
169 }
170}
Note: See TracBrowser for help on using the repository browser.