FDsimG4PMT_SD.cc
Go to the documentation of this file.
1 #include "FDsimG4PMT_SD.hh"
3 #include "G4RunManager.hh"
4 #include "G4Material.hh"
5 #include "G4Step.hh"
6 #include "G4VTouchable.hh"
7 #include "G4TouchableHistory.hh"
8 #include "G4SDManager.hh"
9 #include "G4ios.hh"
10 
11 #include <sstream>
12 
13 using namespace TelescopeSimulatorLX;
14 
15 FDsimG4PMT_SD::FDsimG4PMT_SD(G4String name):G4VSensitiveDetector(name)
16 {
17  collectionName.insert("OpticalHits");
18 }
19 
20 
22 {
23 }
24 
25 
26 void FDsimG4PMT_SD::Initialize(G4HCofThisEvent* HCE)
27 {
28  static int HCID1 = -1;
30  new FDsimG4OpticalHitsCollection(SensitiveDetectorName,collectionName[0]);
31 
32  if (HCID1 < 0)
33  HCID1 = GetCollectionID(0);
34  HCE->AddHitsCollection(HCID1,fOpticalHitsCollection);
35 }
36 
37 
38 G4bool FDsimG4PMT_SD::ProcessHits(G4Step* aStep,G4TouchableHistory*)
39 {
40  G4bool status = false;
41  G4String particleName = aStep->GetTrack()->GetDefinition()->GetParticleName();
42 
43  G4TouchableHistory* theTouchable
44  = (G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable());
45 
47  static FDsimG4OpBoundaryProcess* boundary = NULL;
48 
49  if (particleName == "opticalphoton"){
50  if (!boundary){
51  G4ProcessManager* pm
52  = aStep->GetTrack()->GetDefinition()->GetProcessManager();
53  G4int nprocesses = pm->GetProcessListLength();
54  G4ProcessVector* pv = pm->GetProcessList();
55  G4int i;
56  for( i=0;i<nprocesses;i++){
57  if ((*pv)[i]->GetProcessName() == "OpBoundary"){
58  boundary = (FDsimG4OpBoundaryProcess*)(*pv)[i];
59  break;
60  }
61  }
62  }
63 
64  boundaryStatus=boundary->GetStatus();
65 
66  //Check to see if the partcile was actually at a boundary
67  G4StepPoint* posStepPoint = aStep->GetPostStepPoint();
68 
69  if (posStepPoint->GetStepStatus() == fGeomBoundary) {
70  switch(boundaryStatus){
71  case Absorption:
72  break;
73  case Detection:
74  status = ProcessOpticalHit(aStep,theTouchable);
75  break;
76  default:
77  break;
78  }
79  }
80  }
81  else
82  status = true;
83 
84  return status;
85 }
86 
87 
88 G4bool FDsimG4PMT_SD::ProcessOpticalHit(G4Step* aStep, G4TouchableHistory* theTouchable)
89 {
90  G4Track* track = aStep->GetTrack();
91  G4String particleName = track->GetDefinition()->GetParticleName();
92 
93  G4int HitID;
94  G4VPhysicalVolume* fPhysVolume = theTouchable->GetVolume();
95  G4String VolumeName = fPhysVolume->GetName();
96  G4int copyNo = fPhysVolume->GetCopyNo();
97 
98  if (particleName == "opticalphoton")
99  HitID = 1;
100  else
101  HitID = 10;
102 
103  G4double fGlobalTOF = track->GetGlobalTime();
104  G4double kineticEnergy = track->GetKineticEnergy();
105  G4double Wavelength = h_Planck*c_light/kineticEnergy;
106  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
107  G4ThreeVector worldPos = preStepPoint->GetPosition();
108  G4ThreeVector HitPosition =
109  theTouchable->GetHistory()->GetTopTransform().TransformPoint(worldPos);
110  G4ThreeVector HitDirection = aStep->GetPreStepPoint()->GetMomentumDirection();
111  G4ThreeVector SDposition = fPhysVolume->GetObjectTranslation();
112  G4double weight = track->GetWeight();
113 
114  FDsimG4OpticalHit* OpticalHit = new FDsimG4OpticalHit;
115 
116  OpticalHit->SetEnergy(kineticEnergy);
117  OpticalHit->SetWavelength(Wavelength);
118  OpticalHit->SetTime(fGlobalTOF);
119  OpticalHit->SetPosition(HitPosition);
120  OpticalHit->SetWorldPosition(worldPos);
121  OpticalHit->SetDirection(HitDirection);
122  OpticalHit->SetPixelPosition(SDposition);
123  OpticalHit->SetHitID(HitID);
124  OpticalHit->SetPMTid(copyNo);
125  OpticalHit->SetWeight(weight);
126 
127  fOpticalHitsCollection->insert(OpticalHit);
128 
129  return true;
130 }
131 
132 
133 void FDsimG4PMT_SD::EndOfEvent(G4HCofThisEvent* HCE)
134 {
135  static G4int HCID = -1;
136  if (HCID < 0)
137  HCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
138 
139  HCE->AddHitsCollection(HCID,fOpticalHitsCollection);
140 }
141 
142 
144 {
145 }
146 
147 
149 {
150 }
151 
152 
154 {
155 }
G4THitsCollection< FDsimG4OpticalHit > FDsimG4OpticalHitsCollection
void EndOfEvent(G4HCofThisEvent *)
G4bool ProcessHits(G4Step *, G4TouchableHistory *)
FDsimG4OpticalHitsCollection * fOpticalHitsCollection
FDsimG4OpBoundaryProcessStatus GetStatus() const
G4bool ProcessOpticalHit(G4Step *, G4TouchableHistory *)
void Initialize(G4HCofThisEvent *)

, generated on Tue Sep 26 2023.