FDsimG4FocalSurfaceSD.cc
Go to the documentation of this file.
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 FDsimG4FocalSurfaceSD::FDsimG4FocalSurfaceSD(G4String name):G4VSensitiveDetector(name)
15 {
16  collectionName.insert("FocalSurfaceHits");
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 FDsimG4FocalSurfaceSD::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,OpticalHitsCollection);
42 
43 }
44 
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46 
47 G4bool FDsimG4FocalSurfaceSD::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  if(particleName == "opticalphoton"){HitID=1;}else{HitID=10;}
71 
72 
73  G4double fGlobalTOF = aStep->GetTrack()->GetGlobalTime() ;
74  G4double kineticEnergy = aStep->GetTrack()->GetKineticEnergy() ;
75  G4double Wavelength = h_Planck*c_light/kineticEnergy ;
76 
77 
78  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
79  G4ThreeVector worldPos = preStepPoint->GetPosition();
80  G4ThreeVector localPos = theTouchable->GetHistory()->GetTopTransform().TransformPoint(worldPos);
81 
82  G4ThreeVector HitPosition = localPos;
83 
84 
85  G4ThreeVector HitDirection = aStep->GetPreStepPoint()->GetMomentumDirection() ;
86  G4ThreeVector SDposition = fPhysVolume->GetObjectTranslation() ;
87 
88  if (HitDirection.z() > 0.){
89 
91 
92  Hit->SetEnergy(kineticEnergy);
93  Hit->SetWavelength(Wavelength);
94  Hit->SetTime(fGlobalTOF);
95  Hit->SetPosition(HitPosition);
96  Hit->SetWorldPosition(worldPos);
97  Hit->SetDirection(HitDirection);
98  Hit->SetPixelPosition(SDposition);
99  Hit->SetHitID(HitID);
100 
101  Hit->SetPMTid(copyNo) ;
102 
103  OpticalHitsCollection->insert(Hit);
104 
105  }
106 
107 #ifdef VERBOSE_USE
108 
109  G4cerr << fPhysVolume->GetObjectTranslation() << " "
110  << fPhysVolume->GetFrameTranslation() << " "
111  << fPhysVolume->GetTranslation() << G4endl ;
112 
113  G4cerr << fPhysVolume->GetMultiplicity() << " "
114  << fPhysVolume->GetCopyNo() << G4endl ;
115 
116 #endif
117  return true ;
118 
119 }
120 
121 
122 
123 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
124 
125 void FDsimG4FocalSurfaceSD::EndOfEvent(G4HCofThisEvent* HCE)
126 {
127  static G4int HCID = -1;
128  if(HCID<0)
129  {
130  HCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
131  }
132  HCE->AddHitsCollection(HCID,OpticalHitsCollection);
133 
134 }
135 
136 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
137 
139 {;}
140 
141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
142 
144 {;}
145 
146 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
147 
149 {;}
150 
151 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4THitsCollection< FDsimG4OpticalHit > FDsimG4OpticalHitsCollection
G4bool ProcessHits(G4Step *, G4TouchableHistory *)
FDsimG4OpticalHitsCollection * OpticalHitsCollection

, generated on Tue Sep 26 2023.