Loading...
Searching...
No Matches
LXeSteppingAction.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//
27/// \file optical/LXe/src/LXeSteppingAction.cc
28/// \brief Implementation of the LXeSteppingAction class
29//
30//
31#include "LXeSteppingAction.hh"
32
33#include "LXeEventAction.hh"
34#include "LXePMTSD.hh"
36#include "LXeTrajectory.hh"
38
39#include "G4OpticalPhoton.hh"
40#include "G4ProcessManager.hh"
41#include "G4SDManager.hh"
42#include "G4Step.hh"
43#include "G4SteppingManager.hh"
44#include "G4StepPoint.hh"
45#include "G4Track.hh"
46#include "G4TrackStatus.hh"
47#include "G4VPhysicalVolume.hh"
48
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
50
56
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62
64{
65 G4Track* theTrack = theStep->GetTrack();
66 const G4ParticleDefinition* part = theTrack->GetDefinition();
67 G4int pdg = part->GetPDGEncoding();
68
69 if(theTrack->GetCurrentStepNumber() == 1)
70 fExpectedNextStatus = Undefined;
71
72 auto trackInformation =
73 static_cast<LXeUserTrackInformation*>(theTrack->GetUserInformation());
74
75 G4StepPoint* thePrePoint = theStep->GetPreStepPoint();
76 G4VPhysicalVolume* thePrePV = thePrePoint->GetPhysicalVolume();
77
78 G4StepPoint* thePostPoint = theStep->GetPostStepPoint();
79 G4VPhysicalVolume* thePostPV = thePostPoint->GetPhysicalVolume();
80
81 G4OpBoundaryProcessStatus boundaryStatus = Undefined;
82
83 // find the boundary process only once
84 if(nullptr == fBoundary && pdg == -22)
85 {
86 G4ProcessManager* pm = part->GetProcessManager();
87 G4int nprocesses = pm->GetProcessListLength();
88 G4ProcessVector* pv = pm->GetProcessList();
89 for(G4int i = 0; i < nprocesses; ++i)
90 {
91 if(nullptr != (*pv)[i] && (*pv)[i]->GetProcessName() == "OpBoundary")
92 {
93 fBoundary = dynamic_cast<G4OpBoundaryProcess*>((*pv)[i]);
94 break;
95 }
96 }
97 }
98
99 if(theTrack->GetParentID() == 0)
100 {
101 // This is a primary track
102 auto secondaries = theStep->GetSecondaryInCurrentStep();
103 // If we haven't already found the conversion position and there were
104 // secondaries generated, then search for it
105 // since this is happening before the secondary is being tracked,
106 // the vertex position has not been set yet (set in initial step)
107 if(nullptr != secondaries && !fEventAction->IsConvPosSet())
108 {
109 if(!secondaries->empty())
110 {
111 for(auto & tr : *secondaries)
112 {
113 const G4VProcess* creator = tr->GetCreatorProcess();
114 if(nullptr != creator)
115 {
116 G4int type = creator->GetProcessSubType();
117 // 12 - photoeffect
118 // 13 - Compton scattering
119 // 14 - gamma conversion
120 if(type >= 12 && type <= 14)
121 {
122 fEventAction->SetConvPos(tr->GetPosition());
123 }
124 }
125 }
126 }
127 }
128 if(fOneStepPrimaries && thePrePV->GetName() == "scintillator")
129 theTrack->SetTrackStatus(fStopAndKill);
130 }
131
132 if(nullptr == thePostPV)
133 { // out of world
134 fExpectedNextStatus = Undefined;
135 return;
136 }
137
138 // Optical photon only
139 if(pdg == -22)
140 {
141 if(thePrePV->GetName() == "Slab")
142 // force drawing of photons in WLS slab
143 trackInformation->SetForceDrawTrajectory(true);
144 else if(thePostPV->GetName() == "expHall")
145 // Kill photons entering expHall from something other than Slab
146 theTrack->SetTrackStatus(fStopAndKill);
147
148 // Was the photon absorbed by the absorption process
149 auto proc = thePostPoint->GetProcessDefinedStep();
150 if(nullptr != proc && proc->GetProcessName() == "OpAbsorption")
151 {
153 trackInformation->AddTrackStatusFlag(absorbed);
154 }
155 if(nullptr != fBoundary)
156 boundaryStatus = fBoundary->GetStatus();
157
158 if(thePostPoint->GetStepStatus() == fGeomBoundary)
159 {
160 // Check to see if the particle was actually at a boundary
161 // Otherwise the boundary status may not be valid
162 if(fExpectedNextStatus == StepTooSmall)
163 {
164 if(boundaryStatus != StepTooSmall)
165 {
166 G4cout << "LXeSteppingAction::UserSteppingAction(): "
167 << "trackID=" << theTrack->GetTrackID()
168 << " parentID=" << theTrack->GetParentID()
169 << " " << part->GetParticleName()
170 << " E(MeV)=" << theTrack->GetKineticEnergy()
171 << "n/ at " << theTrack->GetPosition()
172 << " prePV: " << thePrePV->GetName()
173 << " postPV: " << thePostPV->GetName()
174 << G4endl;
175 G4ExceptionDescription ed;
176 ed << "LXeSteppingAction: "
177 << "No reallocation step after reflection!"
178 << "Something is wrong with the surface normal or geometry";
179 G4Exception("LXeSteppingAction:", "LXeExpl01", JustWarning, ed, "");
180 return;
181 }
182 }
183 fExpectedNextStatus = Undefined;
184 switch(boundaryStatus)
185 {
186 case Absorption:
187 trackInformation->AddTrackStatusFlag(boundaryAbsorbed);
189 break;
190 case Detection: // Note, this assumes that the volume causing detection
191 // is the photocathode because it is the only one with
192 // non-zero efficiency
193 {
194 // Trigger sensitive detector manually since photon is
195 // absorbed but status was Detection
196 G4SDManager* SDman = G4SDManager::GetSDMpointer();
197 G4String sdName = "/LXeDet/pmtSD";
198 LXePMTSD* pmtSD = (LXePMTSD*) SDman->FindSensitiveDetector(sdName);
199 if(pmtSD)
200 pmtSD->ProcessHits_boundary(theStep, nullptr);
201 trackInformation->AddTrackStatusFlag(hitPMT);
202 break;
203 }
204 case FresnelReflection:
205 case TotalInternalReflection:
206 case LambertianReflection:
207 case LobeReflection:
208 case SpikeReflection:
209 case BackScattering:
210 trackInformation->IncReflections();
211 fExpectedNextStatus = StepTooSmall;
212 break;
213 default:
214 break;
215 }
216 if(thePostPV->GetName() == "sphere")
217 trackInformation->AddTrackStatusFlag(hitSphere);
218 }
219 }
220}
Definition of the LXeEventAction class.
Definition of the LXePMTSD class.
Definition of the LXeSteppingAction class.
Definition of the LXeSteppingMessenger class.
Definition of the LXeTrajectory class.
Definition of the LXeUserTrackInformation class.
G4double IsConvPosSet()
void SetConvPos(const G4ThreeVector &p)
void IncBoundaryAbsorption()
G4bool ProcessHits_boundary(const G4Step *, G4TouchableHistory *)
Definition LXePMTSD.cc:104
~LXeSteppingAction() override
LXeEventAction * fEventAction
void UserSteppingAction(const G4Step *) override
G4OpBoundaryProcess * fBoundary
LXeSteppingAction(LXeEventAction *)
LXeSteppingMessenger * fSteppingMessenger
G4OpBoundaryProcessStatus fExpectedNextStatus

Applications | User Support | Publications | Collaboration