1 | #include "ELYSE/SteppingAction.hh" |
---|
2 | |
---|
3 | //GEANT 4 |
---|
4 | #include "G4SteppingManager.hh" |
---|
5 | #include "G4SDManager.hh" |
---|
6 | #include "G4EventManager.hh" |
---|
7 | #include "G4ProcessManager.hh" |
---|
8 | #include "G4Track.hh" |
---|
9 | #include "G4Step.hh" |
---|
10 | #include "G4Event.hh" |
---|
11 | #include "G4StepPoint.hh" |
---|
12 | #include "G4TrackStatus.hh" |
---|
13 | #include "G4VPhysicalVolume.hh" |
---|
14 | #include "G4ParticleDefinition.hh" |
---|
15 | #include "G4ParticleTypes.hh" |
---|
16 | #include "G4OpBoundaryProcess.hh" |
---|
17 | |
---|
18 | //ELYSE |
---|
19 | #include "ELYSE/EventAction.hh" |
---|
20 | #include "ELYSE/TrackingAction.hh" |
---|
21 | #include "ELYSE/Trajectory.hh" |
---|
22 | #include "ELYSE/TrappingVolume.hh" |
---|
23 | #include "ELYSE/TrackInformation.hh" |
---|
24 | |
---|
25 | |
---|
26 | ELYSE::SteppingAction::SteppingAction() { |
---|
27 | }//Ctor |
---|
28 | |
---|
29 | //--------------------------------------------------------------------------------------------------- |
---|
30 | |
---|
31 | ELYSE::SteppingAction::~SteppingAction() { |
---|
32 | }//Dtor |
---|
33 | |
---|
34 | //--------------------------------------------------------------------------------------------------- |
---|
35 | void ELYSE::SteppingAction::UserSteppingAction(const G4Step * theStep){ |
---|
36 | |
---|
37 | G4Track* theTrack = theStep->GetTrack(); |
---|
38 | |
---|
39 | G4StepPoint* thePostPoint = theStep->GetPostStepPoint(); |
---|
40 | G4VPhysicalVolume* thePostPV = thePostPoint->GetPhysicalVolume(); |
---|
41 | |
---|
42 | G4OpBoundaryProcessStatus boundaryStatus=Undefined; |
---|
43 | static G4OpBoundaryProcess* boundary=NULL; |
---|
44 | |
---|
45 | if(!thePostPV){//out of world |
---|
46 | return;} |
---|
47 | |
---|
48 | G4ParticleDefinition* particleType = theTrack->GetDefinition(); |
---|
49 | G4String particleName = particleType->GetParticleName(); |
---|
50 | |
---|
51 | //-----------only optical photons |
---|
52 | |
---|
53 | if(particleName == "opticalphoton") { |
---|
54 | |
---|
55 | //------------find the boundary process only once |
---|
56 | if( !boundary ) { |
---|
57 | G4ProcessManager* pm = theStep->GetTrack()->GetDefinition()->GetProcessManager(); |
---|
58 | G4int nprocesses = pm->GetProcessListLength(); |
---|
59 | G4ProcessVector* pv = pm->GetProcessList(); |
---|
60 | G4int i; |
---|
61 | for( i=0;i<nprocesses;i++){ |
---|
62 | if((*pv)[i]->GetProcessName()=="OpBoundary"){ |
---|
63 | boundary = (G4OpBoundaryProcess*)(*pv)[i]; |
---|
64 | break; |
---|
65 | } |
---|
66 | } |
---|
67 | }//set boundary once |
---|
68 | |
---|
69 | //if the OpticalBoundary process is not registered for Optical Photon => FATAL error |
---|
70 | if(!boundary) { |
---|
71 | G4cout << "(JEC) UserSteppingAction: FATAL: NO OpBoundary FOUND" << G4endl; |
---|
72 | exit(0); |
---|
73 | } |
---|
74 | |
---|
75 | //------------Kill photons entering expHall |
---|
76 | if(thePostPV->GetName()=="expHall") { |
---|
77 | theTrack->SetTrackStatus(fStopAndKill);} |
---|
78 | |
---|
79 | //----------- |
---|
80 | |
---|
81 | boundaryStatus = boundary->GetStatus(); |
---|
82 | |
---|
83 | // Check to see if the particle was actually at a boundary |
---|
84 | // Otherwise the boundary status may not be valid |
---|
85 | // Prior to Geant4.6.0-p1 this would not have been enough to check |
---|
86 | |
---|
87 | if(thePostPoint->GetStepStatus()==fGeomBoundary) { |
---|
88 | |
---|
89 | //possible G4OpBoundaryProcessStatus { Undefined, |
---|
90 | // FresnelRefraction, FresnelReflection, |
---|
91 | // TotalInternalReflection, |
---|
92 | // LambertianReflection, LobeReflection, |
---|
93 | // SpikeReflection, BackScattering, |
---|
94 | // Absorption, Detection, NotAtBoundary, |
---|
95 | // SameMaterial, StepTooSmall, NoRINDEX} |
---|
96 | switch(boundaryStatus){ |
---|
97 | case FresnelRefraction: |
---|
98 | G4cout << "(JEC) Status of PostPoint: FresnelRefraction" << G4endl; break; |
---|
99 | case FresnelReflection: |
---|
100 | G4cout << "(JEC) Status of PostPoint: FresnelReflection" << G4endl; break; |
---|
101 | case TotalInternalReflection: |
---|
102 | G4cout << "(JEC) Status of PostPoint: TotalInternalReflection" << G4endl; break; |
---|
103 | case LambertianReflection: |
---|
104 | G4cout << "(JEC) Status of PostPoint: LambertianReflection" << G4endl; break; |
---|
105 | case LobeReflection: |
---|
106 | G4cout << "(JEC) Status of PostPoint: LobeReflection" << G4endl; break; |
---|
107 | case SpikeReflection: |
---|
108 | G4cout << "(JEC) Status of PostPoint: SpikeReflection" << G4endl; break; |
---|
109 | case Absorption: |
---|
110 | G4cout << "(JEC) Status of PostPoint: Absorption" << G4endl; break; |
---|
111 | |
---|
112 | case Detection: // Note, this assumes that the volume causing detection |
---|
113 | // is the interface between Tyvek and Air because it is the only one with |
---|
114 | // 100% efficiency |
---|
115 | { |
---|
116 | // Triger sensitive detector manually since photon is |
---|
117 | // absorbed but status was Detection |
---|
118 | G4SDManager* SDman = G4SDManager::GetSDMpointer(); |
---|
119 | //see DetectorConstruction::Construct for the PATH of the Sensitive Volume |
---|
120 | TrappingVolume* trapVol = 0; |
---|
121 | if(SDman) trapVol = (TrappingVolume*)SDman->FindSensitiveDetector("/ELYSE/TrappingVolume"); |
---|
122 | |
---|
123 | if(trapVol) { |
---|
124 | trapVol->ProcessHits_constStep(theStep,0); |
---|
125 | |
---|
126 | //------------Save the hit into the track info |
---|
127 | |
---|
128 | TrackInformation* aInfo = (TrackInformation*)theTrack->GetUserInformation(); |
---|
129 | G4cout << "(JEC) SteppingAction::opticalphoton detection = " |
---|
130 | << " " << theTrack << " " << aInfo << G4endl; |
---|
131 | if(aInfo)aInfo->AddTrackStatusFlag(hitProduced); |
---|
132 | |
---|
133 | } else { |
---|
134 | G4cout << "(JEC) UserSteppingAction: FATAL: NO Sensitive Det. FOUND" << G4endl; |
---|
135 | exit(0); |
---|
136 | } |
---|
137 | break; |
---|
138 | }//case Detection |
---|
139 | |
---|
140 | case NotAtBoundary: |
---|
141 | G4cout << "(JEC) Status of PostPoint: NotAtBoundary" << G4endl; break; |
---|
142 | case SameMaterial: |
---|
143 | G4cout << "(JEC) Status of PostPoint: SameMaterial" << G4endl; break; |
---|
144 | case StepTooSmall: |
---|
145 | G4cout << "(JEC) Status of PostPoint: StepTooSmall" << G4endl; break; |
---|
146 | case NoRINDEX: |
---|
147 | G4cout << "(JEC) Status of PostPoint: NoRINDEX" << G4endl; break; |
---|
148 | |
---|
149 | default: |
---|
150 | break; |
---|
151 | |
---|
152 | }//eo switch on boundary values |
---|
153 | |
---|
154 | }//test on GeomBoundary |
---|
155 | |
---|
156 | }//Optical photon |
---|
157 | |
---|
158 | |
---|
159 | }//UserSteppingAction |
---|
160 | |
---|
161 | |
---|
162 | |
---|
163 | |
---|
164 | |
---|
165 | |
---|
166 | |
---|
167 | |
---|
168 | |
---|
169 | |
---|
170 | |
---|
171 | |
---|
172 | |
---|
173 | |
---|
174 | |
---|
175 | |
---|
176 | |
---|
177 | |
---|