File indexing completed on 2026-04-15 07:41:59
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 #include "F02CalorimeterSD.hh"
0030
0031 #include "F02CalorHit.hh"
0032 #include "F02DetectorConstruction.hh"
0033
0034 #include "G4SDManager.hh"
0035 #include "G4Step.hh"
0036 #include "G4TouchableHistory.hh"
0037 #include "G4VPhysicalVolume.hh"
0038 #include "G4VTouchable.hh"
0039
0040
0041
0042 F02CalorimeterSD::F02CalorimeterSD(G4String name, F02DetectorConstruction* det)
0043 : G4VSensitiveDetector(name), fDetector(det), fHitID(new G4int[500])
0044 {
0045 collectionName.insert("CalCollection");
0046 }
0047
0048
0049
0050 F02CalorimeterSD::~F02CalorimeterSD()
0051 {
0052 delete[] fHitID;
0053 }
0054
0055
0056
0057 void F02CalorimeterSD::Initialize(G4HCofThisEvent*)
0058 {
0059 fCalCollection = new F02CalorHitsCollection(SensitiveDetectorName, collectionName[0]);
0060 for (G4int j = 0; j < 1; j++) {
0061 fHitID[j] = -1;
0062 };
0063 }
0064
0065
0066
0067 G4bool F02CalorimeterSD::ProcessHits(G4Step* step, G4TouchableHistory*)
0068 {
0069 G4double edep = step->GetTotalEnergyDeposit();
0070 G4double stepl = 0.;
0071
0072 stepl = step->GetStepLength();
0073
0074 if ((edep == 0.) && (stepl == 0.)) return false;
0075
0076 auto theTouchable = (G4TouchableHistory*)(step->GetPreStepPoint()->GetTouchable());
0077
0078 G4VPhysicalVolume* physVol = theTouchable->GetVolume();
0079
0080 G4int number = 0;
0081 if (fHitID[number] == -1) {
0082 auto calHit = new F02CalorHit();
0083 if (physVol == fDetector->GetAbsorber()) calHit->AddAbs(edep, stepl);
0084 fHitID[number] = fCalCollection->insert(calHit) - 1;
0085 if (verboseLevel > 0) G4cout << " New Calorimeter Hit on F02: " << number << G4endl;
0086 }
0087 else {
0088 if (physVol == fDetector->GetAbsorber()) (*fCalCollection)[fHitID[number]]->AddAbs(edep, stepl);
0089 if (verboseLevel > 0) G4cout << " Energy added to F02: " << number << G4endl;
0090 }
0091
0092 return true;
0093 }
0094
0095
0096
0097 void F02CalorimeterSD::EndOfEvent(G4HCofThisEvent* hce)
0098 {
0099 static G4int hcID = -1;
0100 if (hcID < 0) {
0101 hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
0102 }
0103 hce->AddHitsCollection(hcID, fCalCollection);
0104 }
0105
0106