Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 //
28 //
29 //
30 #include "LXeSteppingAction.hh"
31 #include "LXeEventAction.hh"
32 #include "LXeTrackingAction.hh"
33 #include "LXeTrajectory.hh"
34 #include "LXePMTSD.hh"
37 #include "LXeSteppingMessenger.hh"
38 #include "LXeRecorderBase.hh"
39 
40 #include "G4SteppingManager.hh"
41 #include "G4SDManager.hh"
42 #include "G4EventManager.hh"
43 #include "G4ProcessManager.hh"
44 #include "G4Track.hh"
45 #include "G4Step.hh"
46 #include "G4Event.hh"
47 #include "G4StepPoint.hh"
48 #include "G4TrackStatus.hh"
49 #include "G4VPhysicalVolume.hh"
50 #include "G4ParticleDefinition.hh"
51 #include "G4ParticleTypes.hh"
52 
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54 
56  : fRecorder(r),fOneStepPrimaries(false)
57 {
58  fSteppingMessenger = new LXeSteppingMessenger(this);
59 
60  fExpectedNextStatus = Undefined;
61 }
62 
63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64 
66 
67 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68 
70  G4Track* theTrack = theStep->GetTrack();
71 
72  LXeUserTrackInformation* trackInformation
74  LXeUserEventInformation* eventInformation
77 
78  G4StepPoint* thePrePoint = theStep->GetPreStepPoint();
79  G4VPhysicalVolume* thePrePV = thePrePoint->GetPhysicalVolume();
80 
81  G4StepPoint* thePostPoint = theStep->GetPostStepPoint();
82  G4VPhysicalVolume* thePostPV = thePostPoint->GetPhysicalVolume();
83 
84  G4OpBoundaryProcessStatus boundaryStatus=Undefined;
85  static G4OpBoundaryProcess* boundary=NULL;
86 
87  //find the boundary process only once
88  if(!boundary){
90  = theStep->GetTrack()->GetDefinition()->GetProcessManager();
91  G4int nprocesses = pm->GetProcessListLength();
92  G4ProcessVector* pv = pm->GetProcessList();
93  G4int i;
94  for( i=0;i<nprocesses;i++){
95  if((*pv)[i]->GetProcessName()=="OpBoundary"){
96  boundary = (G4OpBoundaryProcess*)(*pv)[i];
97  break;
98  }
99  }
100  }
101 
102  if(theTrack->GetParentID()==0){
103  //This is a primary track
104 
109 
110  //If we havent already found the conversion position and there were
111  //secondaries generated, then search for it
112  if(!eventInformation->IsConvPosSet() && tN2ndariesTot>0 ){
113  for(size_t lp1=(*fSecondary).size()-tN2ndariesTot;
114  lp1<(*fSecondary).size(); lp1++){
115  const G4VProcess* creator=(*fSecondary)[lp1]->GetCreatorProcess();
116  if(creator){
117  G4String creatorName=creator->GetProcessName();
118  if(creatorName=="phot"||creatorName=="compt"||creatorName=="conv"){
119  //since this is happening before the secondary is being tracked
120  //the Vertex position has not been set yet(set in initial step)
121  eventInformation->SetConvPos((*fSecondary)[lp1]->GetPosition());
122  }
123  }
124  }
125  }
126 
127  if(fOneStepPrimaries&&thePrePV->GetName()=="scintillator")
128  theTrack->SetTrackStatus(fStopAndKill);
129  }
130 
131  if(!thePostPV){//out of world
132  return;
133  }
134 
135  G4ParticleDefinition* particleType = theTrack->GetDefinition();
136  if(particleType==G4OpticalPhoton::OpticalPhotonDefinition()){
137  //Optical photon only
138 
139  if(thePrePV->GetName()=="Slab")
140  //force drawing of photons in WLS slab
141  trackInformation->SetForceDrawTrajectory(true);
142  else if(thePostPV->GetName()=="expHall")
143  //Kill photons entering expHall from something other than Slab
144  theTrack->SetTrackStatus(fStopAndKill);
145 
146  //Was the photon absorbed by the absorption process
147  if(thePostPoint->GetProcessDefinedStep()->GetProcessName()
148  =="OpAbsorption"){
149  eventInformation->IncAbsorption();
150  trackInformation->AddTrackStatusFlag(absorbed);
151  }
152 
153  boundaryStatus=boundary->GetStatus();
154 
155  //Check to see if the partcile was actually at a boundary
156  //Otherwise the boundary status may not be valid
157  //Prior to Geant4.6.0-p1 this would not have been enough to check
158  if(thePostPoint->GetStepStatus()==fGeomBoundary){
159  if(fExpectedNextStatus==StepTooSmall){
160  if(boundaryStatus!=StepTooSmall){
162  ed << "LXeSteppingAction::UserSteppingAction(): "
163  << "No reallocation step after reflection!"
164  << G4endl;
165  G4Exception("LXeSteppingAction::UserSteppingAction()", "LXeExpl01",
166  FatalException,ed,
167  "Something is wrong with the surface normal or geometry");
168  }
169  }
170  fExpectedNextStatus=Undefined;
171  switch(boundaryStatus){
172  case Absorption:
173  trackInformation->AddTrackStatusFlag(boundaryAbsorbed);
174  eventInformation->IncBoundaryAbsorption();
175  break;
176  case Detection: //Note, this assumes that the volume causing detection
177  //is the photocathode because it is the only one with
178  //non-zero efficiency
179  {
180  //Triger sensitive detector manually since photon is
181  //absorbed but status was Detection
183  G4String sdName="/LXeDet/pmtSD";
184  LXePMTSD* pmtSD = (LXePMTSD*)SDman->FindSensitiveDetector(sdName);
185  if(pmtSD)pmtSD->ProcessHits_constStep(theStep,NULL);
186  trackInformation->AddTrackStatusFlag(hitPMT);
187  break;
188  }
189  case FresnelReflection:
192  case LobeReflection:
193  case SpikeReflection:
194  case BackScattering:
195  trackInformation->IncReflections();
196  fExpectedNextStatus=StepTooSmall;
197  break;
198  default:
199  break;
200  }
201  if(thePostPV->GetName()=="sphere")
202  trackInformation->AddTrackStatusFlag(hitSphere);
203  }
204  }
205 
206  if(fRecorder)fRecorder->RecordStep(theStep);
207 }