File indexing completed on 2025-02-23 09:22:32
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
0032 #include "Par01CalorimeterSD.hh"
0033
0034 #include "Par01CalorimeterHit.hh"
0035
0036 #include "G4LogicalVolume.hh"
0037 #include "G4ParticleDefinition.hh"
0038 #include "G4SDManager.hh"
0039 #include "G4Track.hh"
0040 #include "G4VPhysicalVolume.hh"
0041 #include "G4ios.hh"
0042
0043
0044
0045 Par01CalorimeterSD::Par01CalorimeterSD(G4String name, G4int nCells, G4String colName)
0046 : G4VSensitiveDetector(name), fNumberOfCells(nCells), fHCID(-1)
0047 {
0048 G4String HCname;
0049 collectionName.insert(HCname = colName);
0050 fCellID = new G4int[fNumberOfCells];
0051 }
0052
0053
0054
0055 Par01CalorimeterSD::~Par01CalorimeterSD()
0056 {
0057 delete[] fCellID;
0058 }
0059
0060
0061
0062 void Par01CalorimeterSD::Initialize(G4HCofThisEvent*)
0063 {
0064 fCalCollection = new Par01CalorimeterHitsCollection(SensitiveDetectorName, collectionName[0]);
0065 for (G4int j = 0; j < fNumberOfCells; j++) {
0066 fCellID[j] = -1;
0067 }
0068 }
0069
0070
0071
0072 G4bool Par01CalorimeterSD::ProcessHits(G4Step* aStep, G4TouchableHistory*)
0073 {
0074 G4double edep = aStep->GetTotalEnergyDeposit();
0075 if (edep <= 0.) return false;
0076
0077 auto hist = (G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable());
0078 const G4VPhysicalVolume* physVol = hist->GetVolume();
0079 G4int copyID = hist->GetReplicaNumber();
0080
0081 if (fCellID[copyID] == -1) {
0082 Par01CalorimeterHit* calHit = new Par01CalorimeterHit(physVol->GetLogicalVolume());
0083 calHit->SetEdep(edep);
0084 G4AffineTransform aTrans = hist->GetHistory()->GetTopTransform();
0085 aTrans.Invert();
0086 calHit->SetPos(aTrans.NetTranslation());
0087 calHit->SetRot(aTrans.NetRotation());
0088 G4int icell = fCalCollection->insert(calHit);
0089 fCellID[copyID] = icell - 1;
0090 if (verboseLevel > 0) {
0091 G4cout << " New Calorimeter Hit on CellID " << copyID << G4endl;
0092 }
0093 }
0094 else {
0095 (*fCalCollection)[fCellID[copyID]]->AddEdep(edep);
0096 if (verboseLevel > 0) {
0097 G4cout << " Energy added to CellID " << copyID << G4endl;
0098 }
0099 }
0100
0101 return true;
0102 }
0103
0104
0105
0106 void Par01CalorimeterSD::EndOfEvent(G4HCofThisEvent* HCE)
0107 {
0108 if (fHCID < 0) {
0109 fHCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
0110 }
0111 HCE->AddHitsCollection(fHCID, fCalCollection);
0112 }
0113
0114
0115
0116 void Par01CalorimeterSD::clear() {}
0117
0118
0119
0120 void Par01CalorimeterSD::DrawAll() {}
0121
0122
0123
0124 void Par01CalorimeterSD::PrintAll() {}