EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PHG4RICHSteppingAction.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file PHG4RICHSteppingAction.cc
1 
10 
11 #include "PHG4RICHDetector.h"
12 #include "ePHENIXRICHConstruction.h" // for ePHENIXRICHConstruction
13 
14 #include <g4main/PHG4Hit.h>
16 #include <g4main/PHG4Hitv1.h>
17 #include <g4main/PHG4Shower.h>
19 
20 #include <phool/getClass.h>
21 
22 #include <Geant4/G4ExceptionSeverity.hh> // for FatalException
23 #include <Geant4/G4OpBoundaryProcess.hh> // for StepTooSmall, Undefined
24 #include <Geant4/G4OpticalPhoton.hh> // for G4OpticalPhoton
25 #include <Geant4/G4ParticleDefinition.hh> // for G4ParticleDefinition
26 #include <Geant4/G4ProcessManager.hh>
27 #include <Geant4/G4ProcessVector.hh> // for G4ProcessVector
28 #include <Geant4/G4Step.hh>
29 #include <Geant4/G4StepPoint.hh> // for G4StepPoint
30 #include <Geant4/G4StepStatus.hh> // for fGeomBoundary
31 #include <Geant4/G4String.hh> // for G4String
32 #include <Geant4/G4SystemOfUnits.hh> // for cm, nanosecond, GeV
33 #include <Geant4/G4ThreeVector.hh> // for G4ThreeVector
34 #include <Geant4/G4Track.hh> // for G4Track
35 #include <Geant4/G4Types.hh> // for G4int, G4double
36 #include <Geant4/G4VPhysicalVolume.hh> // for G4VPhysicalVolume
37 #include <Geant4/G4VProcess.hh> // for G4VProcess
38 #include <Geant4/G4VTouchable.hh> // for G4VTouchable
39 #include <Geant4/G4VUserTrackInformation.hh> // for G4VUserTrackInformation
40 #include <Geant4/globals.hh> // for G4Exception, G4Exceptio...
41 
42 #include <cassert> // for assert
43 #include <iostream> // for operator<<, basic_ostream
44 
45 using namespace std;
46 
48  : detector_(detector)
49  , hits_(nullptr)
50  , hit(nullptr)
51  , fExpectedNextStatus(Undefined)
52 {
53 }
54 
56 {
57  G4Track* theTrack = aStep->GetTrack();
58 
59  G4StepPoint* thePostPoint = aStep->GetPostStepPoint();
60  G4VPhysicalVolume* thePostPV = thePostPoint->GetPhysicalVolume();
61 
62  G4OpBoundaryProcessStatus boundaryStatus = Undefined;
63  static G4OpBoundaryProcess* boundary = nullptr;
64 
65  /* find the boundary process only once */
66  if (!boundary)
67  {
68  G4ProcessManager* pm = aStep->GetTrack()->GetDefinition()->GetProcessManager();
69  G4int nprocesses = pm->GetProcessListLength();
70  G4ProcessVector* pv = pm->GetProcessList();
71  G4int i;
72  for (i = 0; i < nprocesses; i++)
73  {
74  if ((*pv)[i]->GetProcessName() == "OpBoundary")
75  {
76  boundary = (G4OpBoundaryProcess*) (*pv)[i];
77  break;
78  }
79  }
80  }
81 
82  if (!thePostPV)
83  { //out of world
84  return;
85  }
86 
87  /* Optical photon only */
88  G4ParticleDefinition* particleType = theTrack->GetDefinition();
89  if (particleType == G4OpticalPhoton::OpticalPhotonDefinition())
90  {
91  /* Was the photon absorbed by the absorption process? */
92  if (thePostPoint->GetProcessDefinedStep()->GetProcessName() == "OpAbsorption")
93  {
94  }
95  assert(boundary != nullptr);
96  boundaryStatus = boundary->GetStatus();
97 
98  /*Check to see if the partcile was actually at a boundary
99  Otherwise the boundary status may not be valid
100  Prior to Geant4.6.0-p1 this would not have been enough to check */
101  if (thePostPoint->GetStepStatus() == fGeomBoundary)
102  {
103  if (fExpectedNextStatus == StepTooSmall)
104  {
105  if (boundaryStatus != StepTooSmall)
106  {
107  G4ExceptionDescription ed;
108  ed << "EicRichGemTbSteppingAction::UserSteppingAction(): "
109  << "No reallocation step after reflection!"
110  << std::endl;
111  G4Exception("EicRichGemTbSteppingAction::UserSteppingAction()", "EicRichGemTbExpl01",
112  FatalException, ed,
113  "Something is wrong with the surface normal or geometry");
114  }
115  }
117  switch (boundaryStatus)
118  {
119  case Absorption:
120  break;
121  case Detection: /*Note, this assumes that the volume causing detection
122  is the photocathode because it is the only one with
123  non-zero efficiency */
124  {
125  /* make sure the photon actually did hit the GEM stack */
126  if (thePostPV->GetName() == "RICHHBDGEMFrontCu1Physical")
127  {
128  MakeHit(aStep);
129  }
130  break;
131  }
132  case FresnelReflection:
133  case TotalInternalReflection:
134  case LambertianReflection:
135  case LobeReflection:
136  case SpikeReflection:
137  case BackScattering:
138  fExpectedNextStatus = StepTooSmall;
139  break;
140  default:
141  break;
142  }
143  }
144  }
145 
146  return;
147 }
148 
149 bool PHG4RICHSteppingAction::MakeHit(const G4Step* aStep)
150 {
151  // collect energy and track
152  G4double edep = aStep->GetTotalEnergyDeposit() / GeV;
153  const G4Track* aTrack = aStep->GetTrack();
154  const G4VTouchable* aTouch = aTrack->GetTouchable();
155  G4StepPoint* postPoint = aStep->GetPostStepPoint();
156 
157  // set sector number
158  int sector_id = -999;
159  bool sector_found = false;
160 
161  // Check if volume(1) is in sector volume, if not check volume(0)
162  if (detector_->ePHENIXRICHConstruction::is_in_sector(aTouch->GetVolume(1)) > -1)
163  {
164  sector_id = aTouch->GetCopyNumber(1);
165  sector_found = true;
166  }
167  else if (detector_->ePHENIXRICHConstruction::is_in_sector(aTouch->GetVolume()) > -1)
168  {
169  sector_id = aTouch->GetCopyNumber();
170  sector_found = true;
171  }
172 
173  if (!sector_found)
174  {
175  if (!aTouch->GetVolume(1) || !aTouch->GetVolume())
176  cout << "WARNING: Missing volumes for hit!" << endl;
177  else
178  cout << "WARNING: Photon hit volume is not the RICH readout plane volume!" << endl;
179  }
180 
181  hit = new PHG4Hitv1();
182 
183  //here we set the entrance values in cm
184  hit->set_x(0, postPoint->GetPosition().x() / cm);
185  hit->set_y(0, postPoint->GetPosition().y() / cm);
186  hit->set_z(0, postPoint->GetPosition().z() / cm);
187 
188  // time in ns
189  hit->set_t(0, postPoint->GetGlobalTime() / nanosecond);
190 
191  //same for exit values (photons absorbed/detected at boundary to post step volume)
192  hit->set_x(1, postPoint->GetPosition().x() / cm);
193  hit->set_y(1, postPoint->GetPosition().y() / cm);
194  hit->set_z(1, postPoint->GetPosition().z() / cm);
195 
196  hit->set_t(1, postPoint->GetGlobalTime() / nanosecond);
197 
198  //set the track ID
199  {
200  hit->set_trkid(aTrack->GetTrackID());
201  if (G4VUserTrackInformation* p = aTrack->GetUserInformation())
202  {
203  if (PHG4TrackUserInfoV1* pp = dynamic_cast<PHG4TrackUserInfoV1*>(p))
204  {
205  hit->set_trkid(pp->GetUserTrackId());
206  hit->set_shower_id(pp->GetShower()->get_id());
207 
208  pp->SetKeep(true); // we want to keep the track
209  }
210  }
211  }
212 
213  // set optical photon energy deposition
214  hit->set_edep(edep);
215 
216  // Now add the hit
217  hits_->AddHit(sector_id, hit);
218 
219  {
220  if (G4VUserTrackInformation* p = aTrack->GetUserInformation())
221  {
222  if (PHG4TrackUserInfoV1* pp = dynamic_cast<PHG4TrackUserInfoV1*>(p))
223  {
224  pp->GetShower()->add_g4hit_id(hits_->GetID(), hit->get_hit_id());
225  }
226  }
227  }
228 
229  // return true to indicate the hit was used
230  return true;
231 }
232 
234 {
235  //now look for the map and grab a pointer to it.
236  hits_ = findNode::getClass<PHG4HitContainer>(topNode, "G4HIT_RICH");
237 
238  // if we do not find the node we need to make it.
239  if (!hits_)
240  {
241  std::cout
242  << "PHG4RICHSteppingAction::SetTopNode - unable to find G4HIT_RICH"
243  << std::endl;
244  }
245 }