Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
advanced/microbeam/src/SteppingAction.cc
이 파일의 문서화 페이지로 가기
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  // Material : 1 is cytoplasm, 2 is nucleus
71 
72  G4int matVoxelPRE = -1;
73  G4int matVoxelPOST = -1;
74 
75  const G4StepPoint* preStep = aStep->GetPreStepPoint();
76  const G4StepPoint* postStep = aStep->GetPostStepPoint();
77  const G4Track* track = aStep->GetTrack();
78 
79  const G4LogicalVolume* preVolume =
80  preStep->GetPhysicalVolume()->GetLogicalVolume();
81 
82  const G4LogicalVolume* postVolume = nullptr;
83  if(postStep->GetPhysicalVolume())
84  {
85  postVolume = postStep->GetPhysicalVolume()->GetLogicalVolume();
86  }
87  const G4ParticleDefinition* particle =
89 
90  G4int preReplicaNumber = preStep->GetTouchableHandle()->GetReplicaNumber();
92 
93  if (preReplicaNumber>0)
94  {
95  matVoxelPRE = fMyCellParameterisation->GetTissueType(preReplicaNumber);
96  }
97 
98  if(postVolume)
99  {
100  G4int postReplicaNumber = postStep->GetTouchableHandle()->GetReplicaNumber();
101  if (postReplicaNumber>0)
102  {
103  matVoxelPOST = fMyCellParameterisation->GetTissueType(postReplicaNumber);
104  }
105  }
106 
107  // COUNT GAS DETECTOR HITS
108 
109  if (particle == G4Alpha::AlphaDefinition())
110  {
111  if(postVolume == fDetector->GetLogicalIsobutane() &&
112  ((preVolume == fDetector->GetLogicalCollDetYoke())
113  ||
114  (preVolume == fDetector->GetLogicalCollDetGap4())
115  ||
116  (preVolume == fDetector->GetLogicalCollDetGap4())))
117  {
118  fRun->AddNbOfHitsGas();
119  }
120 
121  // STOPPING POWER AND BEAM SPOT SIZE AT CELL ENTRANCE
122  if(preVolume == fDetector->GetLogicalPolyprop() &&
123  ( (postVolume == fDetector->GetLogicalKgm()) ||
124  (matVoxelPOST == 1)) )
125  {
126  G4double deltaE = preStep->GetKineticEnergy()
127  - postStep->GetKineticEnergy();
128  if(deltaE > 0.0)
129  {
130  //Fill ntupleid=1
131  man->FillNtupleDColumn(1,0,preStep->GetKineticEnergy()/keV);
132  man->FillNtupleDColumn(1,1,deltaE*micrometer/(keV*aStep->GetStepLength()));
133  man->AddNtupleRow(1);
134  }
135 
136  // Average dE over step suggested by Michel Maire
137  G4ThreeVector coord1 = preStep->GetPosition();
138  const G4AffineTransform transformation1 =
139  preStep->GetTouchable()->GetHistory()->GetTopTransform();
140  G4ThreeVector localPosition1 = transformation1.TransformPoint(coord1);
141 
142  G4ThreeVector coord2 = postStep->GetPosition();
143  const G4AffineTransform transformation2 =
144  postStep->GetTouchable()->GetHistory()->GetTopTransform();
145  G4ThreeVector localPosition2 = transformation2.TransformPoint(coord2);
146 
147  G4ThreeVector localPosition =
148  localPosition1 + G4UniformRand()*(localPosition2-localPosition1);
149 
150  //Fill ntupleid=2
151  man->FillNtupleDColumn(2,0,localPosition.x()/micrometer);
152  man->FillNtupleDColumn(2,1,localPosition.y()/micrometer);
153  man->AddNtupleRow(2);
154  }
155 
156  // ALPHA RANGE
157  if (postStep->GetKineticEnergy() < eV &&
158  ( (matVoxelPOST==1) ||
159  (postVolume == fDetector->GetLogicalKgm()) ||
160  (matVoxelPOST==2) ) )
161  {
162  //Fill ntupleid=3
163  man->FillNtupleDColumn(3,0,postStep->GetPosition().x()/micrometer);
164  man->FillNtupleDColumn(3,1,postStep->GetPosition().y()/micrometer);
165  man->FillNtupleDColumn(3,2,postStep->GetPosition().z()/micrometer);
166  man->AddNtupleRow(3);
167  }
168 
169  // TOTAL DOSE DEPOSIT AND DOSE DEPOSIT WITHIN A PHANTOM VOXEL
170  // FOR ALL PARTICLES
171  }
172 
173  if (matVoxelPRE == 2)
174  {
175  G4double dose = (edep/joule)/(fRun->GetMassNucleus()/kg);
176  fRun->AddDoseN(dose);
177  fRun->AddDoseBox(preReplicaNumber, edep/eV);
178  }
179  else if (matVoxelPRE == 1)
180  {
181  G4double dose = (edep/joule)/(fRun->GetMassCytoplasm()/kg);
182  fRun->AddDoseC(dose);
183  fRun->AddDoseBox(preReplicaNumber, edep/eV);
184  }
185 
186  // PROTECTION AGAINST POSSIBLE MSC LOOPS FOR e-
187 
188  // if ( edep/MeV<1e-25 && particle == G4Electron::Electron())
189  // {
190  //aStep->GetTrack()->SetTrackStatus(fStopAndKill);
191  /*
192  G4cout << "*** Warning *** : msc loop for "
193  << track->GetDefinition()->GetParticleName()
194  << " in " <<
195  postPoint->GetTouchableHandle()->GetVolume()->GetName() << G4endl;
196  */
197  // }
198 }
const DetectorConstruction * fDetector
static constexpr double micrometer
Definition: G4SIunits.hh:100
const G4VTouchable * GetTouchable() const
G4LogicalVolume * GetLogicalVolume() const
G4StepPoint * GetPreStepPoint() const
Float_t dose
static constexpr double keV
Definition: G4SIunits.hh:216
void AddDoseBox(G4int i, G4double x)
double z() const
G4VPhysicalVolume * GetPhysicalVolume() const
static constexpr double kg
Definition: G4SIunits.hh:182
const G4TouchableHandle & GetTouchableHandle() const
double G4double
Definition: G4Types.hh:76
G4ParticleDefinition * GetDefinition() const
G4double GetStepLength() const
const G4ThreeVector & GetPosition() const
static CellParameterisation * Instance()
#define G4UniformRand()
Definition: Randomize.hh:53
const G4AffineTransform & GetTopTransform() const
G4Track * GetTrack() const
G4double GetTotalEnergyDeposit() const
static constexpr double eV
Definition: G4SIunits.hh:215
G4CsvAnalysisManager G4AnalysisManager
Definition: g4csv_defs.hh:77
Definition: G4Step.hh:76
G4StepPoint * GetPostStepPoint() const
void UserSteppingAction(const G4Step *)
G4double GetKineticEnergy() const
Double_t edep
const CellParameterisation * fMyCellParameterisation
static G4Alpha * AlphaDefinition()
Definition: G4Alpha.cc:84
int G4int
Definition: G4Types.hh:78
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:58
G4int GetTissueType(G4int i) const
static constexpr double joule
Definition: G4SIunits.hh:204
double x() const
double y() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
virtual const G4NavigationHistory * GetHistory() const
Definition: G4VTouchable.cc:86
const G4DynamicParticle * GetDynamicParticle() const
Simple detector construction with a box volume placed in a world.