Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 // 27 /// \file optical/LXe/src/LXePMTSD.cc 28 /// \brief Implementation of the LXePMTSD class 29 // 30 // 31 #include "LXePMTSD.hh" 32 33 #include "LXeDetectorConstruction.hh" 34 #include "LXePMTHit.hh" 35 #include "LXeUserTrackInformation.hh" 36 37 #include "G4LogicalVolume.hh" 38 #include "G4ParticleDefinition.hh" 39 #include "G4ParticleTypes.hh" 40 #include "G4SDManager.hh" 41 #include "G4Step.hh" 42 #include "G4TouchableHistory.hh" 43 #include "G4Track.hh" 44 #include "G4VPhysicalVolume.hh" 45 #include "G4VTouchable.hh" 46 #include "G4ios.hh" 47 48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 49 50 LXePMTSD::LXePMTSD(G4String name) : G4VSensitiveDetector(name) 51 { 52 collectionName.insert("pmtHitCollection"); 53 } 54 55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 56 57 LXePMTSD::~LXePMTSD() 58 { 59 delete fPMTPositionsX; 60 delete fPMTPositionsY; 61 delete fPMTPositionsZ; 62 } 63 64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 65 66 void LXePMTSD::SetPmtPositions(const std::vector<G4ThreeVector>& positions) 67 { 68 for (size_t i = 0; i < positions.size(); ++i) { 69 if (fPMTPositionsX) fPMTPositionsX->push_back(positions[i].x()); 70 if (fPMTPositionsY) fPMTPositionsY->push_back(positions[i].y()); 71 if (fPMTPositionsZ) fPMTPositionsZ->push_back(positions[i].z()); 72 } 73 } 74 75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 76 77 void LXePMTSD::Initialize(G4HCofThisEvent* hitsCE) 78 { 79 fPMTHitCollection = new LXePMTHitsCollection(SensitiveDetectorName, collectionName[0]); 80 81 if (fHitCID < 0) { 82 fHitCID = G4SDManager::GetSDMpointer()->GetCollectionID(fPMTHitCollection); 83 } 84 hitsCE->AddHitsCollection(fHitCID, fPMTHitCollection); 85 } 86 87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 88 89 G4bool LXePMTSD::ProcessHits(G4Step*, G4TouchableHistory*) 90 { 91 return false; 92 } 93 94 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 95 96 // Generates a hit and uses the postStepPoint's mother volume replica number 97 // PostStepPoint because the hit is generated manually when the photon is 98 // absorbed by the photocathode 99 100 G4bool LXePMTSD::ProcessHits_boundary(const G4Step* aStep, G4TouchableHistory*) 101 { 102 // need to know if this is an optical photon 103 if (aStep->GetTrack()->GetDefinition() != G4OpticalPhoton::OpticalPhotonDefinition()) 104 return false; 105 106 // User replica number 1 since photocathode is a daughter volume 107 // to the pmt which was replicated 108 G4int pmtNumber = aStep->GetPostStepPoint()->GetTouchable()->GetReplicaNumber(1); 109 G4VPhysicalVolume* physVol = aStep->GetPostStepPoint()->GetTouchable()->GetVolume(1); 110 111 // Find the correct hit collection 112 size_t n = fPMTHitCollection->entries(); 113 LXePMTHit* hit = nullptr; 114 for (size_t i = 0; i < n; ++i) { 115 if ((*fPMTHitCollection)[i]->GetPMTNumber() == pmtNumber) { 116 hit = (*fPMTHitCollection)[i]; 117 break; 118 } 119 } 120 121 if (hit == nullptr) { // this pmt wasn't previously hit in this event 122 hit = new LXePMTHit(); // so create new hit 123 hit->SetPMTNumber(pmtNumber); 124 hit->SetPMTPhysVol(physVol); 125 fPMTHitCollection->insert(hit); 126 hit->SetPMTPos((*fPMTPositionsX)[pmtNumber], (*fPMTPositionsY)[pmtNumber], 127 (*fPMTPositionsZ)[pmtNumber]); 128 } 129 130 hit->IncPhotonCount(); // increment hit for the selected pmt 131 132 if (!LXeDetectorConstruction::GetSphereOn()) { 133 hit->SetDrawit(true); 134 // If the sphere is disabled then this hit is automaticaly drawn 135 } 136 else { // sphere enabled 137 auto trackInfo = (LXeUserTrackInformation*)aStep->GetTrack()->GetUserInformation(); 138 if (trackInfo->GetTrackStatus() & hitSphere) 139 // only draw this hit if the photon has hit the sphere first 140 hit->SetDrawit(true); 141 } 142 143 return true; 144 } 145