File indexing completed on 2025-10-31 08:24:23
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 #include "RE05CalorimeterSD.hh"
0032 
0033 #include "RE05CalorimeterHit.hh"
0034 
0035 #include "G4LogicalVolume.hh"
0036 #include "G4ParticleDefinition.hh"
0037 #include "G4Step.hh"
0038 #include "G4SystemOfUnits.hh"
0039 #include "G4TouchableHistory.hh"
0040 #include "G4Track.hh"
0041 #include "G4VPhysicalVolume.hh"
0042 #include "G4VTouchable.hh"
0043 #include "G4ios.hh"
0044 
0045 
0046 
0047 RE05CalorimeterSD::RE05CalorimeterSD(G4String name)
0048   : G4VSensitiveDetector(name), fNumberOfCellsInZ(20), fNumberOfCellsInPhi(48)
0049 {
0050   G4String HCname;
0051   collectionName.insert(HCname = "calCollection");
0052 }
0053 
0054 
0055 
0056 RE05CalorimeterSD::~RE05CalorimeterSD() {}
0057 
0058 
0059 
0060 void RE05CalorimeterSD::Initialize(G4HCofThisEvent*)
0061 {
0062   fCalCollection = new RE05CalorimeterHitsCollection(SensitiveDetectorName, collectionName[0]);
0063   for (G4int j = 0; j < fNumberOfCellsInZ; j++)
0064     for (G4int k = 0; k < fNumberOfCellsInPhi; k++) {
0065       fCellID[j][k] = -1;
0066     }
0067   verboseLevel = 0;
0068 }
0069 
0070 
0071 
0072 G4bool RE05CalorimeterSD::ProcessHits(G4Step* aStep, G4TouchableHistory*)
0073 {
0074   
0075   
0076   
0077   
0078 
0079   G4double edep = aStep->GetTotalEnergyDeposit();
0080   if (verboseLevel > 1) G4cout << "Next step edep(MeV) = " << edep / MeV << G4endl;
0081   if (edep == 0.) return false;
0082 
0083   const G4VTouchable* touchable = aStep->GetPreStepPoint()->GetTouchable();
0084   G4int copyIDinZ = touchable->GetReplicaNumber(0);
0085   G4int copyIDinPhi = touchable->GetReplicaNumber(1);
0086 
0087   if (fCellID[copyIDinZ][copyIDinPhi] == -1) {
0088     RE05CalorimeterHit* calHit =
0089       new RE05CalorimeterHit(touchable->GetVolume()->GetLogicalVolume(), copyIDinZ, copyIDinPhi);
0090     calHit->SetEdep(edep);
0091     G4AffineTransform aTrans = touchable->GetHistory()->GetTopTransform();
0092     aTrans.Invert();
0093     calHit->SetPos(aTrans.NetTranslation());
0094     calHit->SetRot(aTrans.NetRotation());
0095     G4int icell = fCalCollection->insert(calHit);
0096     fCellID[copyIDinZ][copyIDinPhi] = icell - 1;
0097     if (verboseLevel > 0) {
0098       G4cout << " New Calorimeter Hit on fCellID " << copyIDinZ << " " << copyIDinPhi << G4endl;
0099     }
0100   }
0101   else {
0102     (*fCalCollection)[fCellID[copyIDinZ][copyIDinPhi]]->AddEdep(edep);
0103     if (verboseLevel > 0) {
0104       G4cout << " Energy added to fCellID " << copyIDinZ << " " << copyIDinPhi << G4endl;
0105     }
0106   }
0107 
0108   return true;
0109 }
0110 
0111 
0112 
0113 void RE05CalorimeterSD::EndOfEvent(G4HCofThisEvent* HCE)
0114 {
0115   static G4int HCID = -1;
0116   if (HCID < 0) {
0117     HCID = GetCollectionID(0);
0118   }
0119   HCE->AddHitsCollection(HCID, fCalCollection);
0120 }
0121 
0122 
0123 
0124 void RE05CalorimeterSD::clear() {}
0125 
0126 
0127 
0128 void RE05CalorimeterSD::DrawAll() {}
0129 
0130 
0131 
0132 void RE05CalorimeterSD::PrintAll() {}
0133 
0134