QwGeant4
QweakSimCerenkov_PMTSD.cc
Go to the documentation of this file.
1 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
2 
3 // geant4 includes
4 #include "G4OpticalPhoton.hh"
5 
6 // user includes
10 #include "QweakSimTrajectory.hh"
11 
14 
15 
16 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
18 :G4VSensitiveDetector(name)
19 {
20  collectionName.insert("PMTHitCollection");
22  myUserInfo = userInfo;
23 }
24 
25 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
27 {
28 
29 }
30 
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 void QweakSimCerenkov_PMTSD::Initialize(G4HCofThisEvent* HCE)
33 {
34  CerenkovDetector_PMTHitsCollection = new QweakSimCerenkovDetector_PMTHitsCollection(SensitiveDetectorName,collectionName[0]);
35 
37  { CerenkovDetectorPMT_CollectionID = G4SDManager::GetSDMpointer()->GetCollectionID(CerenkovDetector_PMTHitsCollection); }
39 }
40 
41 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
42 G4bool QweakSimCerenkov_PMTSD::ProcessHits_constStep(const G4Step* aStep,G4TouchableHistory* /*ROhist*/)
43 {
44 
45  if (aStep->GetTrack()->GetDefinition() != G4OpticalPhoton::OpticalPhotonDefinition())
46  return false;
47 
48  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
49  G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
50 // G4VPhysicalVolume* thePrePV = preStepPoint->GetPhysicalVolume();
51 // G4VPhysicalVolume* thePostPV = postStepPoint->GetPhysicalVolume();
52 
53  G4TouchableHistory* theTouchable = (G4TouchableHistory*)(preStepPoint->GetTouchable());
54 // G4VPhysicalVolume* physVol = theTouchable->GetVolume();
55 
56  if (postStepPoint->GetStepStatus() != fGeomBoundary) return false;
57 
58  //if( physVol->GetName().compare("Cathode_Physical")!=0 ) return false;
59 
60  if (preStepPoint->GetStepStatus() != fGeomBoundary) return false; // Entering Geometry
61 
62  G4double currentPhotonEnergy = aStep->GetTrack()->GetTotalEnergy();
63  G4double currentHitTime = aStep->GetTrack()->GetGlobalTime();
64  G4int MotherCopyNo = theTouchable->GetVolume(1)->GetCopyNo();
66  G4int MotherReplicaNo2 = theTouchable->GetReplicaNumber(3); // Several MotherVolumes
67 
68 // std::cout<<"store hit info: detector "<<MotherReplicaNo2<<" PMTID "<<MotherCopyNo<<std::endl;
69  aHit->StoreDetectorID(MotherReplicaNo2); // which octant (number to be converted)
70  aHit->StorePMTID(MotherCopyNo); // left or right pmt
71 
72  aHit->StorePhotonEnergy(currentPhotonEnergy);
73  aHit->StoreHitTime(currentHitTime);
74  G4int hitCount = CerenkovDetector_PMTHitsCollection->insert(aHit);
75  aHit->StoreHitID(hitCount);
76  //aHit->SetHitValid(False);
77  myUserInfo->SetCurrentPMTHit(aHit,MotherCopyNo);
78 
79  return true;
80 }
81 
82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83 
84 void QweakSimCerenkov_PMTSD::EndOfEvent(G4HCofThisEvent* )
85 {
86 
87 }
88 
89 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
91 {
92  G4TrajectoryContainer* trajectoryContainer =
93  G4RunManager::GetRunManager()->GetCurrentEvent()->GetTrajectoryContainer();
94 
95  if(trajectoryContainer==0)
96  {
97  G4cout << " Could not find the TrajectoryContainer" << G4endl;
98  return 0;
99  }
100 // TrajectoryVector* vect = container->GetVector();
101 //
102 // (G4VTrajectory**) tr = vect->begin();
103 //
104 // while(tr != vect->end())
105 // {
106 // QweakSimTrajectory* tr1 = (QweakSimTrajectory*)(*tr);
107 // if(tr1->GetTrackID()==parentID) return tr1;
108 // tr++;
109 // }
110 
111 G4int i = 0;
112 
113 G4int n_trajectories = trajectoryContainer->entries();
114 
115 while(i != n_trajectories)
116 {
117 
118  QweakSimTrajectory* trj = (QweakSimTrajectory*)((*trajectoryContainer)[i]);
119  G4cout << " Current charge of possible parent track = " <<trj->GetCharge() << "--- Parent ID = "
120  << trj->GetTrackID() << G4endl;
121  if(trj->GetTrackID()==parentID) return trj;
122  i++;
123 }
124 
125 return 0;
126 }
127 
QweakSimTrajectory * GetParentTrajectory(G4int parentID)
G4THitsCollection< QweakSimCerenkov_PMTHit > QweakSimCerenkovDetector_PMTHitsCollection
Stores the information about the various tracks.
void EndOfEvent(G4HCofThisEvent *HCE)
void StorePhotonEnergy(G4double eng)
QweakSimUserInformation * myUserInfo
G4int GetTrackID() const
void StoreDetectorID(G4int detector_ID)
void SetCurrentPMTHit(QweakSimCerenkov_PMTHit *hit, G4int side)
G4double GetCharge() const
QweakSimCerenkovDetector_PMTHitsCollection * CerenkovDetector_PMTHitsCollection
QweakSimCerenkov_PMTSD(G4String name, QweakSimUserInformation *myUserInfo)
G4bool ProcessHits_constStep(const G4Step *aStep, G4TouchableHistory *ROhist)
void StoreHitTime(G4double time)
void Initialize(G4HCofThisEvent *HCE)