EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PHG4FPbScRegionSteppingAction.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file PHG4FPbScRegionSteppingAction.cc
2 
3 #include "PHG4FPbScDetector.h"
4 
5 #include <g4main/PHG4Hit.h>
7 #include <g4main/PHG4Hitv1.h>
8 #include <g4main/PHG4Shower.h>
10 
11 #include <phool/getClass.h>
12 
13 #include <Geant4/G4Step.hh>
14 #include <Geant4/G4StepPoint.hh> // for G4StepPoint
15 #include <Geant4/G4StepStatus.hh> // for fGeomBoundary, fUndefined
16 #include <Geant4/G4SystemOfUnits.hh> // for cm, nanosecond, GeV
17 #include <Geant4/G4ThreeVector.hh> // for G4ThreeVector
18 #include <Geant4/G4TouchableHandle.hh> // for G4TouchableHandle
19 #include <Geant4/G4Track.hh> // for G4Track
20 #include <Geant4/G4Types.hh> // for G4double
21 #include <Geant4/G4VTouchable.hh> // for G4VTouchable
22 #include <Geant4/G4VUserTrackInformation.hh> // for G4VUserTrackInformation
23 
24 #include <iostream>
25 #include <string> // for operator+, operator<<
26 
27 class G4VPhysicalVolume;
28 
29 using namespace std;
30 //____________________________________________________________________________..
32  : detector_(detector)
33  , hits_(nullptr)
34  , hit(nullptr)
35 {
36 }
37 
38 //____________________________________________________________________________..
40 {
41  // get volume of the current step
42  G4VPhysicalVolume* volume = aStep->GetPreStepPoint()->GetTouchableHandle()->GetVolume();
43 
44  // collect energy and track length step by step
45  G4double edep = aStep->GetTotalEnergyDeposit() / GeV;
46  if (edep == 0.) return;
47 
48  const G4Track* aTrack = aStep->GetTrack();
49 
50  // make sure we are in a volume
51  if (detector_->isInScintillator(volume))
52  {
53  int layer_id = 0;
54  G4StepPoint* prePoint = aStep->GetPreStepPoint();
55  G4StepPoint* postPoint = aStep->GetPostStepPoint();
56  cout << "track id " << aTrack->GetTrackID() << endl;
57  cout << "time prepoint: " << prePoint->GetGlobalTime() << endl;
58  cout << "time postpoint: " << postPoint->GetGlobalTime() << endl;
59  switch (prePoint->GetStepStatus())
60  {
61  case fGeomBoundary:
62  case fUndefined:
63  hit = new PHG4Hitv1();
64  //here we set the entrance values in cm
65  hit->set_x(0, prePoint->GetPosition().x() / cm);
66  hit->set_y(0, prePoint->GetPosition().y() / cm);
67  hit->set_z(0, prePoint->GetPosition().z() / cm);
68  // time in ns
69  hit->set_t(0, prePoint->GetGlobalTime() / nanosecond);
70  //set the track ID
71  {
72  hit->set_trkid(aTrack->GetTrackID());
73  if (G4VUserTrackInformation* p = aTrack->GetUserInformation())
74  {
75  if (PHG4TrackUserInfoV1* pp = dynamic_cast<PHG4TrackUserInfoV1*>(p))
76  {
77  hit->set_trkid(pp->GetUserTrackId());
78  hit->set_shower_id(pp->GetShower()->get_id());
79  }
80  }
81  }
82 
83  //set the initial energy deposit
84  hit->set_edep(0);
85 
86  // Now add the hit
87  hits_->AddHit(layer_id, hit);
88 
89  {
90  if (G4VUserTrackInformation* p = aTrack->GetUserInformation())
91  {
92  if (PHG4TrackUserInfoV1* pp = dynamic_cast<PHG4TrackUserInfoV1*>(p))
93  {
94  pp->GetShower()->add_g4hit_id(hits_->GetID(), hit->get_hit_id());
95  }
96  }
97  }
98 
99  break;
100  default:
101  break;
102  }
103  // here we just update the exit values, it will be overwritten
104  // for every step until we leave the volume or the particle
105  // ceases to exist
106  hit->set_x(1, postPoint->GetPosition().x() / cm);
107  hit->set_y(1, postPoint->GetPosition().y() / cm);
108  hit->set_z(1, postPoint->GetPosition().z() / cm);
109 
110  hit->set_t(1, postPoint->GetGlobalTime() / nanosecond);
111  //sum up the energy to get total deposited
112  hit->set_edep(hit->get_edep() + edep);
113 
114  // hit->identify();
115  // return true to indicate the hit was used
116  return;
117  }
118  else
119  {
120  return;
121  }
122 }
123 
124 //____________________________________________________________________________..
126 {
127  string hitnodename = "G4HIT_" + detector_->GetName();
128  //now look for the map and grab a pointer to it.
129  hits_ = findNode::getClass<PHG4HitContainer>(topNode, hitnodename.c_str());
130 
131  // if we do not find the node we need to make it.
132  if (!hits_)
133  {
134  std::cout << "PHG4FPbScRegionSteppingAction::SetTopNode - unable to find " << hitnodename << std::endl;
135  }
136 }