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

Last change on this file since 810 was 807, checked in by garnier, 16 years ago

update

  • Property svn:executable set to *
File size: 6.6 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.8 2007/08/22 13:58:33 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)
39:Run(run),Detector(det)
40{ }
41
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
43
44MicrobeamSteppingAction::~MicrobeamSteppingAction()
45{ }
46
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48
49void MicrobeamSteppingAction::UserSteppingAction(const G4Step* aStep)
50 
51{ 
52
53// COUNT GAS DETECTOR HITS
54
55if (       ((aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "_CollDet_yoke_")
56        &&  (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "Isobutane")
57        &&  (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha"))
58       
59        || 
60           ((aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "_CollDet_gap4_")
61        &&  (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "Isobutane")
62        &&  (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha"))
63
64        ||
65
66           ((aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "_CollDet_gap5_")
67        &&  (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "Isobutane")
68        &&  (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha"))
69
70    )
71        {
72          Run->AddNbOfHitsGas();       
73        }
74       
75// STOPPING POWER AND BEAM SPOT SIZE AT CELL ENTRANCE
76
77if (       ((aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "Polyprop")
78        &&  (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "KGM")
79        &&  (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha"))
80       
81        || 
82           ((aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "Polyprop")
83        &&  (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "physicalCytoplasm")
84        &&  (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha"))
85    )
86        {
87       
88        if(aStep->GetTotalEnergyDeposit()>0) 
89         {
90          FILE *myFile;
91          myFile=fopen("stoppingPower.txt","a");       
92          fprintf(myFile,"%e \n",
93          (aStep->GetTotalEnergyDeposit()/keV)/(aStep->GetStepLength()/micrometer));   
94          fclose (myFile);
95         }
96
97         // Average dE over step syggested by Michel Maire
98
99         G4StepPoint* p1 = aStep->GetPreStepPoint();
100         G4ThreeVector coord1 = p1->GetPosition();
101         const G4AffineTransform transformation1 = p1->GetTouchable()->GetHistory()->GetTopTransform();
102         G4ThreeVector localPosition1 = transformation1.TransformPoint(coord1);
103
104         G4StepPoint* p2 = aStep->GetPostStepPoint();
105         G4ThreeVector coord2 = p2->GetPosition();
106         const G4AffineTransform transformation2 = p2->GetTouchable()->GetHistory()->GetTopTransform();
107         G4ThreeVector localPosition2 = transformation2.TransformPoint(coord2);
108
109         G4ThreeVector localPosition = localPosition1 + G4UniformRand()*(localPosition2-localPosition1);
110         
111         // end
112
113         FILE *myFile;
114         myFile=fopen("beamPosition.txt","a");
115         fprintf
116         ( 
117          myFile,"%e %e \n",
118          localPosition.x()/micrometer,
119          localPosition.y()/micrometer
120         );
121         fclose (myFile);                               
122
123        }
124
125// ALPHA RANGE
126
127if (
128
129        (aStep->GetTrack()->GetDynamicParticle()->GetDefinition() ->GetParticleName() == "alpha")
130       
131        &&
132       
133        (aStep->GetTrack()->GetKineticEnergy()<1e-6)
134       
135        &&
136                       
137          ( (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "physicalCytoplasm")
138        ||  (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "KGM")       
139        ||  (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() == "physicalNucleus") ) 
140               
141   )
142       
143        {               
144         FILE *myFile;
145         myFile=fopen("range.txt","a");
146         fprintf
147         ( myFile,"%e %e %e\n",
148          (aStep->GetTrack()->GetPosition().x())/micrometer,
149          (aStep->GetTrack()->GetPosition().y())/micrometer,
150          (aStep->GetTrack()->GetPosition().z())/micrometer
151         );
152         fclose (myFile);                               
153        }
154
155// TOTAL DOSE DEPOSIT AND DOSE DEPOSIT WITHIN A PHANTOM VOXEL
156
157if (aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName()  == "physicalNucleus")
158
159        { 
160         G4double dose = (e_SI*(aStep->GetTotalEnergyDeposit()/eV))/(Run->GetMassNucleus());
161         Run->AddDoseN(dose);
162
163         G4ThreeVector v;
164         Run->AddDoseBox(aStep->GetPreStepPoint()->GetTouchableHandle()->GetReplicaNumber(),
165          aStep->GetTotalEnergyDeposit()/eV);
166        }
167
168
169if (aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName()  == "physicalCytoplasm")
170
171        { 
172         G4double dose = (e_SI*(aStep->GetTotalEnergyDeposit()/eV))/(Run->GetMassCytoplasm());
173         Run->AddDoseC(dose);
174
175         G4ThreeVector v;
176         Run->AddDoseBox(aStep->GetPreStepPoint()->GetTouchableHandle()->GetReplicaNumber(),
177          aStep->GetTotalEnergyDeposit()/eV);
178        }
179}
Note: See TracBrowser for help on using the repository browser.