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 #include "Par03SensitiveDetector.hh" 27 28 #include "Par03Hit.hh" 29 30 #include "G4HCofThisEvent.hh" 31 #include "G4SDManager.hh" 32 #include "G4Step.hh" 33 #include "G4TouchableHistory.hh" 34 #include "G4Track.hh" 35 36 Par03SensitiveDetector::Par03SensitiveDetector(G4String aName) : G4VSensitiveDetector(aName) 37 { 38 collectionName.insert("hits"); 39 } 40 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 41 42 Par03SensitiveDetector::Par03SensitiveDetector(G4String aName, G4int aNumLayers, G4int aNumRho, 43 G4int aNumPhi) 44 : G4VSensitiveDetector(aName), fCellNoZ(aNumLayers), fCellNoRho(aNumRho), fCellNoPhi(aNumPhi) 45 { 46 collectionName.insert("hits"); 47 } 48 49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 50 51 Par03SensitiveDetector::~Par03SensitiveDetector() = default; 52 53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 54 55 void Par03SensitiveDetector::Initialize(G4HCofThisEvent* aHCE) 56 { 57 fHitsCollection = new Par03HitsCollection(SensitiveDetectorName, collectionName[0]); 58 if (fHitCollectionID < 0) { 59 fHitCollectionID = G4SDManager::GetSDMpointer()->GetCollectionID(fHitsCollection); 60 } 61 aHCE->AddHitsCollection(fHitCollectionID, fHitsCollection); 62 63 // fill calorimeter hits with zero energy deposition 64 for (G4int iphi = 0; iphi < fCellNoPhi; iphi++) 65 for (G4int irho = 0; irho < fCellNoRho; irho++) 66 for (G4int iz = 0; iz < fCellNoZ; iz++) { 67 auto hit = new Par03Hit(); 68 fHitsCollection->insert(hit); 69 } 70 } 71 72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 73 74 G4bool Par03SensitiveDetector::ProcessHits(G4Step* aStep, G4TouchableHistory*) 75 { 76 G4double edep = aStep->GetTotalEnergyDeposit(); 77 if (edep == 0.) return true; 78 79 auto aTouchable = (G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable()); 80 81 auto hit = RetrieveAndSetupHit(aTouchable); 82 83 // Add energy deposit from G4Step 84 hit->AddEdep(edep); 85 86 // Fill time information from G4Step 87 // If it's already filled, choose hit with earliest global time 88 if (hit->GetTime() == -1 || hit->GetTime() > aStep->GetTrack()->GetGlobalTime()) 89 hit->SetTime(aStep->GetTrack()->GetGlobalTime()); 90 91 // Set hit type to full simulation (only if hit is not already marked as fast 92 // sim) 93 if (hit->GetType() != 1) hit->SetType(0); 94 95 return true; 96 } 97 98 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 99 100 G4bool Par03SensitiveDetector::ProcessHits(const G4FastHit* aHit, const G4FastTrack* aTrack, 101 G4TouchableHistory* aTouchable) 102 { 103 G4double edep = aHit->GetEnergy(); 104 if (edep == 0.) return true; 105 106 auto hit = RetrieveAndSetupHit(aTouchable); 107 108 // Add energy deposit from G4FastHit 109 hit->AddEdep(edep); 110 111 // Fill time information from G4FastTrack 112 // If it's already filled, choose hit with earliest global time 113 if (hit->GetTime() == -1 || hit->GetTime() > aTrack->GetPrimaryTrack()->GetGlobalTime()) { 114 hit->SetTime(aTrack->GetPrimaryTrack()->GetGlobalTime()); 115 } 116 117 // Set hit type to fast simulation (even if hit was already marked as full 118 // sim, overwrite it) 119 hit->SetType(1); 120 121 return true; 122 } 123 124 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 125 126 Par03Hit* Par03SensitiveDetector::RetrieveAndSetupHit(G4TouchableHistory* aTouchable) 127 { 128 G4int rhoNo = aTouchable->GetCopyNumber(0); // cell 129 G4int phiNo = aTouchable->GetCopyNumber(1); // segment 130 G4int zNo = aTouchable->GetCopyNumber(2); // layer 131 132 std::size_t hitID = fCellNoRho * fCellNoZ * phiNo + fCellNoZ * rhoNo + zNo; 133 134 if (hitID >= fHitsCollection->entries()) { 135 G4Exception("Par03SensitiveDetector::RetrieveAndSetupHit()", "InvalidSetup", FatalException, 136 "Size of hit collection in Par03SensitiveDetector is smaller than the " 137 "number of cells created in Par03DetectorConstruction!"); 138 } 139 Par03Hit* hit = (*fHitsCollection)[hitID]; 140 141 if (hit->GetRhoId() < 0) { 142 hit->SetRhoId(rhoNo); 143 hit->SetPhiId(phiNo); 144 hit->SetZid(zNo); 145 hit->SetLogV(aTouchable->GetVolume(0)->GetLogicalVolume()); 146 G4AffineTransform transform = aTouchable->GetHistory()->GetTopTransform(); 147 hit->SetRot(transform.NetRotation()); 148 transform.Invert(); 149 hit->SetPos(transform.NetTranslation()); 150 } 151 return hit; 152 }