FDsimG4MercedesSD.cc
Go to the documentation of this file.
1 #include "FDsimG4MercedesSD.hh"
2 #include "G4RunManager.hh"
3 #include "G4Material.hh"
4 #include "G4Step.hh"
5 #include "G4VTouchable.hh"
6 #include "G4TouchableHistory.hh"
7 #include "G4SDManager.hh"
8 #include "G4ios.hh"
9 #include <sstream>
10 
11 using namespace TelescopeSimulatorLX;
12 
13 FDsimG4MercedesSD::FDsimG4MercedesSD(G4String name):G4VSensitiveDetector(name)
14 {
15  collectionName.insert("MercedesHits");
16 }
17 
18 
20 {
21 }
22 
23 
24 void FDsimG4MercedesSD::Initialize(G4HCofThisEvent* HCE)
25 {
26  static int HCID1 = -1;
28  new FDsimG4OpticalHitsCollection(SensitiveDetectorName,collectionName[0]);
29 
30  if (HCID1 < 0)
31  HCID1 = GetCollectionID(0);
32  HCE->AddHitsCollection(HCID1,MercedesHitsCollection);
33 }
34 
35 
36 G4bool FDsimG4MercedesSD::ProcessHits(G4Step* aStep,G4TouchableHistory*)
37 {
38  G4TouchableHistory* theTouchable
39  = (G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable());
40 
41  G4String thisVolume = theTouchable->GetVolume()->GetName();
42  G4String particleName = aStep->GetTrack()->GetDefinition()->GetParticleName();
43  G4int HitID;
44  G4VPhysicalVolume* fPhysVolume = theTouchable->GetVolume();
45  G4String VolumeName = fPhysVolume->GetName();
46  G4int copyNo = fPhysVolume->GetCopyNo();
47 
48  if (particleName == "opticalphoton")
49  HitID=1;
50  else
51  HitID=10;
52 
53  G4double globalTOF = aStep->GetTrack()->GetGlobalTime();
54  G4double kineticEnergy = aStep->GetTrack()->GetKineticEnergy();
55  G4double wavelength = h_Planck*c_light/kineticEnergy;
56  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
57  G4ThreeVector worldPos = preStepPoint->GetPosition();
58  G4ThreeVector localPos =
59  theTouchable->GetHistory()->GetTopTransform().TransformPoint(worldPos);
60  G4ThreeVector HitPosition = localPos;
61  G4ThreeVector HitDirection = aStep->GetPreStepPoint()->GetMomentumDirection();
62  G4ThreeVector SDposition = fPhysVolume->GetObjectTranslation();
63 
65 
66  Hit->SetEnergy(kineticEnergy);
67  Hit->SetWavelength(wavelength);
68  Hit->SetTime(globalTOF);
69  Hit->SetPosition(HitPosition);
70  Hit->SetWorldPosition(worldPos);
71  Hit->SetDirection(HitDirection);
72  Hit->SetPixelPosition(SDposition);
73  Hit->SetHitID(HitID);
74  Hit->SetPMTid(copyNo);
75 
76  MercedesHitsCollection->insert(Hit);
77 
78  return true;
79  }
80 
81 
82 void FDsimG4MercedesSD::EndOfEvent(G4HCofThisEvent* HCE)
83 {
84  static G4int HCID = -1;
85  if(HCID < 0)
86  HCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
87  HCE->AddHitsCollection(HCID,MercedesHitsCollection);
88 }
89 
90 
92 {
93 }
94 
95 
97 {
98 }
99 
100 
102 {
103 }
G4bool ProcessHits(G4Step *, G4TouchableHistory *)
G4THitsCollection< FDsimG4OpticalHit > FDsimG4OpticalHitsCollection
FDsimG4OpticalHitsCollection * MercedesHitsCollection

, generated on Tue Sep 26 2023.