File indexing completed on 2025-01-31 09:21:48
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043 #include "G4VProcess.hh"
0044 #include "G4OpticalPhoton.hh"
0045 #include "G4HCofThisEvent.hh"
0046 #include "G4Step.hh"
0047 #include "G4SDManager.hh"
0048 #include "ConfigurationManager.hh"
0049
0050 #include "PhotonSD.hh"
0051 #ifdef WITH_G4OPTICKS
0052 # include "G4Opticks.hh"
0053 # include "TrackInfo.hh"
0054 # include "OpticksGenstep.h"
0055 # include "OpticksFlags.hh"
0056 # include "G4OpticksHit.hh"
0057 #endif
0058
0059 PhotonSD::PhotonSD(G4String name)
0060 : G4VSensitiveDetector(name)
0061 {
0062 G4String HCname = name + "_HC";
0063 collectionName.insert(HCname);
0064 G4cout << collectionName.size() << " PhotonSD name: " << name
0065 << " collection Name: " << HCname << G4endl;
0066 fHCID = -1;
0067 verbose = ConfigurationManager::getInstance()->isEnable_verbose();
0068 }
0069
0070 void PhotonSD::Initialize(G4HCofThisEvent* hce)
0071 {
0072 fPhotonHitsCollection =
0073 new PhotonHitsCollection(SensitiveDetectorName, collectionName[0]);
0074 if(fHCID < 0)
0075 {
0076 if(verbose)
0077 G4cout << "PhotonSD::Initialize: " << SensitiveDetectorName << " "
0078 << collectionName[0] << G4endl;
0079 fHCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
0080 }
0081 hce->AddHitsCollection(fHCID, fPhotonHitsCollection);
0082 }
0083
0084 G4bool PhotonSD::ProcessHits(G4Step* aStep, G4TouchableHistory*)
0085 {
0086 G4Track* theTrack = aStep->GetTrack();
0087
0088 if(theTrack->GetDefinition() != G4OpticalPhoton::OpticalPhotonDefinition())
0089 {
0090 return false;
0091 }
0092 G4double theEdep = theTrack->GetTotalEnergy() / CLHEP::eV;
0093 const G4VProcess* thisProcess = theTrack->GetCreatorProcess();
0094 G4String processname;
0095 if(thisProcess != NULL)
0096 processname = thisProcess->GetProcessName();
0097 else
0098 processname = "No Process";
0099 unsigned theCreationProcessid;
0100 if(processname == "Cerenkov")
0101 {
0102 theCreationProcessid = 0;
0103 }
0104 else if(processname == "Scintillation")
0105 {
0106 theCreationProcessid = 1;
0107 }
0108 else
0109 {
0110 theCreationProcessid = -1;
0111 }
0112 PhotonHit* newHit = new PhotonHit(
0113 0, theCreationProcessid, etolambda(theEdep), theTrack->GetGlobalTime(),
0114 aStep->GetPostStepPoint()->GetPosition(),
0115 aStep->GetPostStepPoint()->GetMomentumDirection(),
0116 aStep->GetPostStepPoint()->GetPolarization());
0117 fPhotonHitsCollection->insert(newHit);
0118 theTrack->SetTrackStatus(fStopAndKill);
0119 return true;
0120 }
0121
0122 void PhotonSD::EndOfEvent(G4HCofThisEvent*)
0123 {
0124 if(verbose)
0125 {
0126 G4int NbHits = fPhotonHitsCollection->entries();
0127 G4cout << " PhotonSD::EndOfEvent Number of PhotonHits: " << NbHits
0128 << G4endl;
0129 }
0130 }
0131 #ifdef WITH_G4OPTICKS
0132 void PhotonSD::AddOpticksHits()
0133 {
0134 G4Opticks* g4ok = G4Opticks::Get();
0135 bool way_enabled = g4ok->isWayEnabled();
0136 unsigned num_hits = g4ok->getNumHit();
0137 if(verbose)
0138 G4cout << "PhotonSD::AddOpticksHits PhotonHits: " << num_hits << G4endl;
0139 G4OpticksHit hit;
0140 G4OpticksHitExtra hit_extra;
0141 G4OpticksHitExtra* hit_extra_ptr = way_enabled ? &hit_extra : NULL;
0142 for(unsigned i = 0; i < num_hits; i++)
0143 {
0144 g4ok->getHit(i, &hit, hit_extra_ptr);
0145 PhotonHit* newHit =
0146 new PhotonHit(i, 0, hit.wavelength, hit.time, hit.global_position,
0147 hit.global_direction, hit.global_polarization);
0148 fPhotonHitsCollection->insert(newHit);
0149 }
0150 if(verbose)
0151 G4cout << "AddOpticksHits size: " << fPhotonHitsCollection->entries()
0152 << G4endl;
0153 }
0154 #endif