source: trunk/examples/extended/optical/LXe/src/LXeSteppingAction.cc @ 807

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

update

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