File indexing completed on 2025-02-23 09:22:35
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 #include "Par04ParallelFullSensitiveDetector.hh"
0027
0028 #include "Par04EventInformation.hh" // for Par04EventInformation
0029 #include "Par04Hit.hh" // for Par04Hit, Par04HitsCollection
0030
0031 #include "G4Event.hh" // for G4Event
0032 #include "G4EventManager.hh" // for G4EventManager
0033 #include "G4HCofThisEvent.hh" // for G4HCofThisEvent
0034 #include "G4SDManager.hh" // for G4SDManager
0035 #include "G4Step.hh" // for G4Step
0036 #include "G4Track.hh" // for G4Track
0037
0038 #include <CLHEP/Vector/Rotation.h> // for HepRotation
0039 #include <CLHEP/Vector/ThreeVector.h> // for Hep3Vector
0040 #include <G4CollectionNameVector.hh> // for G4CollectionNameVector
0041 #include <G4FastHit.hh> // for G4FastHit
0042 #include <G4FastTrack.hh> // for G4FastTrack
0043 #include <G4RotationMatrix.hh> // for G4RotationMatrix
0044 #include <G4StepPoint.hh> // for G4StepPoint
0045 #include <G4THitsCollection.hh> // for G4THitsCollection
0046 #include <G4ThreeVector.hh> // for G4ThreeVector
0047 #include <G4VSensitiveDetector.hh> // for G4VSensitiveDetector
0048 #include <G4VUserEventInformation.hh> // for G4VUserEventInformation
0049 #include <cmath> // for floor
0050 #include <cstddef> // for size_t
0051 #include <vector> // for vector
0052
0053
0054
0055 Par04ParallelFullSensitiveDetector::Par04ParallelFullSensitiveDetector(G4String aName)
0056 : G4VSensitiveDetector(aName)
0057 {
0058 collectionName.insert("physicalCellsFullSim");
0059 }
0060
0061
0062 Par04ParallelFullSensitiveDetector::Par04ParallelFullSensitiveDetector(G4String aName,
0063 G4int aNbOfLayers,
0064 G4int aNbOfSlices,
0065 G4int aNbOfRows)
0066 : G4VSensitiveDetector(aName),
0067 fNbOfLayers(aNbOfLayers),
0068 fNbOfSlices(aNbOfSlices),
0069 fNbOfRows(aNbOfRows)
0070 {
0071 collectionName.insert("physicalCellsFullSim");
0072 }
0073
0074
0075
0076 Par04ParallelFullSensitiveDetector::~Par04ParallelFullSensitiveDetector() = default;
0077
0078
0079
0080 void Par04ParallelFullSensitiveDetector::Initialize(G4HCofThisEvent* aHCE)
0081 {
0082 fHitsCollection = new Par04HitsCollection(SensitiveDetectorName, collectionName[0]);
0083 if (fHitCollectionID < 0) {
0084 fHitCollectionID = G4SDManager::GetSDMpointer()->GetCollectionID(fHitsCollection);
0085 }
0086 aHCE->AddHitsCollection(fHitCollectionID, fHitsCollection);
0087 }
0088
0089
0090
0091 G4bool Par04ParallelFullSensitiveDetector::ProcessHits(G4Step* aStep, G4TouchableHistory*)
0092 {
0093 G4double edep = aStep->GetTotalEnergyDeposit();
0094 if (edep == 0.) return true;
0095
0096 auto touchable = (G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable());
0097
0098 G4int layerNo = touchable->GetCopyNumber(0);
0099 G4int sliceNo = touchable->GetCopyNumber(1);
0100 G4int rowNo = touchable->GetCopyNumber(2);
0101
0102 G4int hitID = fNbOfLayers * fNbOfSlices * rowNo + fNbOfLayers * sliceNo + layerNo;
0103 if (layerNo >= fNbOfLayers) {
0104 G4cout << "ERROR, problem with Layer IDs: " << layerNo << " > " << fNbOfLayers << G4endl;
0105 return false;
0106 }
0107 if (sliceNo >= fNbOfSlices) {
0108 G4cout << "ERROR, problem with Slice IDs: " << sliceNo << " >= " << fNbOfSlices << G4endl;
0109 return false;
0110 }
0111 if (rowNo >= fNbOfRows) {
0112 G4cout << "ERROR, problem with Row IDs: " << rowNo << " >= " << fNbOfRows << G4endl;
0113 return false;
0114 }
0115 auto hit = fHitsMap[hitID].get();
0116 if (hit == nullptr) {
0117 fHitsMap[hitID] = std::unique_ptr<Par04Hit>(new Par04Hit());
0118 hit = fHitsMap[hitID].get();
0119 hit->SetPhiId(sliceNo);
0120 hit->SetRhoId(layerNo);
0121 hit->SetZid(rowNo);
0122 }
0123
0124
0125 hit->AddEdep(edep);
0126
0127 hit->AddNdep(1);
0128
0129
0130
0131 if (hit->GetTime() == -1 || hit->GetTime() > aStep->GetTrack()->GetGlobalTime())
0132 hit->SetTime(aStep->GetTrack()->GetGlobalTime());
0133
0134
0135 hit->SetType(2);
0136
0137 return true;
0138 }
0139
0140
0141 void Par04ParallelFullSensitiveDetector::EndOfEvent(G4HCofThisEvent*)
0142 {
0143 for (const auto& hits : fHitsMap) {
0144 fHitsCollection->insert(new Par04Hit(*hits.second.get()));
0145 }
0146 fHitsMap.clear();
0147 }