FDsimG4FilterSD.cc
Go to the documentation of this file.
1 #include "FDsimG4FilterSD.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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
14 FDsimG4FilterSD::FDsimG4FilterSD(G4String name):G4VSensitiveDetector(name)
15 {
16  collectionName.insert("FilterHits");
17 
18  G4cerr << "Created " << name << " Sensitive Detector." << G4endl ;
19 }
20 
21 
22 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
23 
25 {;}
26 
27 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
28 
29 void FDsimG4FilterSD::Initialize(G4HCofThisEvent* HCE)
30 {
31 
32  static int HCID1 = -1;
33 
34  // SensitiveDetectorName e collectionName sao data members de G4VSensitiveDetector
35 
37  new FDsimG4OpticalHitsCollection(SensitiveDetectorName,collectionName[0]);
38 
39  if(HCID1<0)
40  { HCID1 = GetCollectionID(0); }
41  HCE->AddHitsCollection(HCID1,FilterHitsCollection);
42 
43 }
44 
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46 
47 G4bool FDsimG4FilterSD::ProcessHits(G4Step* aStep,G4TouchableHistory*){
48 
49  G4TouchableHistory* theTouchable
50  = (G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable());
51 
52 
53  // Get Material
54 
55  G4String thisVolume = theTouchable->GetVolume()->GetName() ;
56  G4String particleName = aStep->GetTrack()->GetDefinition()->GetParticleName();
57 
58 
59 #ifdef VERBOSE_USE
60  G4cerr << "Volume name: " << thisVolume << G4endl ;;
61 #endif
62 
63 
64  G4int HitID ;
65 
66  G4VPhysicalVolume *fPhysVolume = theTouchable->GetVolume() ;
67  G4String VolumeName = fPhysVolume->GetName() ;
68  G4int copyNo = fPhysVolume->GetCopyNo();
69 
70  // G4cerr << "Volume name: " << VolumeName << " " << copyNo << G4endl ;
71 
72  if(particleName == "opticalphoton"){HitID=1;}else{HitID=10;}
73 
74  G4double fGlobalTOF = aStep->GetTrack()->GetGlobalTime() ;
75  G4double kineticEnergy = aStep->GetTrack()->GetKineticEnergy() ;
76  G4double Wavelength = h_Planck*c_light/kineticEnergy ;
77 
78  // G4StepPoint* preStepPoint ; // unused
79  G4StepPoint* posStepPoint ;
80  // preStepPoint = aStep->GetPreStepPoint(); // unused
81  posStepPoint = aStep->GetPostStepPoint();
82 
83  G4ThreeVector worldPos = posStepPoint->GetPosition();
84  G4ThreeVector localPos = theTouchable->GetHistory()->GetTopTransform().TransformPoint(worldPos);
85 
86  G4ThreeVector HitPosition = localPos;
87 
88 
89  G4ThreeVector HitDirection = posStepPoint->GetMomentumDirection() ;
90  G4ThreeVector SDposition = fPhysVolume->GetObjectTranslation() ;
91 
92 
93  if (HitDirection.z() < 0.){
94 
96 
97  Hit->SetEnergy(kineticEnergy);
98  Hit->SetWavelength(Wavelength);
99  Hit->SetTime(fGlobalTOF);
100  Hit->SetPosition(HitPosition);
101  Hit->SetWorldPosition(worldPos);
102  Hit->SetDirection(HitDirection);
103  Hit->SetPixelPosition(SDposition);
104  Hit->SetHitID(HitID);
105 
106  Hit->SetPMTid(copyNo) ;
107 
108 
109  FilterHitsCollection->insert(Hit);
110 
111  }
112 
113 #ifdef VERBOSE_USE
114 
115  G4cerr << fPhysVolume->GetObjectTranslation() << " "
116  << fPhysVolume->GetFrameTranslation() << " "
117  << fPhysVolume->GetTranslation() << G4endl ;
118 
119  G4cerr << fPhysVolume->GetMultiplicity() << " "
120  << fPhysVolume->GetCopyNo() << G4endl ;
121 
122 #endif
123  return true ;
124 
125 }
126 
127 
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
130 
131 void FDsimG4FilterSD::EndOfEvent(G4HCofThisEvent* HCE)
132 {
133  static G4int HCID = -1;
134  if(HCID<0)
135  {
136  HCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
137  }
138  HCE->AddHitsCollection(HCID,FilterHitsCollection);
139 
140 }
141 
142 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
143 
145 {;}
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
148 
150 {;}
151 
152 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
153 
155 {;}
156 
157 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4THitsCollection< FDsimG4OpticalHit > FDsimG4OpticalHitsCollection
FDsimG4OpticalHitsCollection * FilterHitsCollection
G4bool ProcessHits(G4Step *, G4TouchableHistory *)

, generated on Tue Sep 26 2023.