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 "Par04EventAction.hh"
0027
0028 #include "Par04DetectorConstruction.hh" // for Par04DetectorConstruction
0029 #include "Par04Hit.hh" // for Par04Hit, Par04HitsCollection
0030 #include "Par04ParallelFullWorld.hh"
0031
0032 #include "G4AnalysisManager.hh" // for G4AnalysisManager
0033 #include "G4Event.hh" // for G4Event
0034 #include "G4EventManager.hh" // for G4EventManager
0035 #include "G4Exception.hh" // for G4Exception, G4ExceptionDesc...
0036 #include "G4ExceptionSeverity.hh" // for FatalException
0037 #include "G4GenericAnalysisManager.hh" // for G4GenericAnalysisManager
0038 #include "G4HCofThisEvent.hh" // for G4HCofThisEvent
0039 #include "G4PrimaryParticle.hh" // for G4PrimaryParticle
0040 #include "G4PrimaryVertex.hh" // for G4PrimaryVertex
0041 #include "G4SDManager.hh" // for G4SDManager
0042 #include "G4SystemOfUnits.hh" // for GeV
0043 #include "G4THitsCollection.hh" // for G4THitsCollection
0044 #include "G4ThreeVector.hh" // for G4ThreeVector
0045 #include "G4Timer.hh" // for G4Timer
0046 #include "G4UserEventAction.hh" // for G4UserEventAction
0047
0048 #include <CLHEP/Units/SystemOfUnits.h> // for GeV
0049 #include <CLHEP/Vector/ThreeVector.h> // for Hep3Vector
0050 #include <algorithm> // for max
0051 #include <cmath> // for log10
0052 #include <cstddef> // for size_t
0053 #include <ostream> // for basic_ostream::operator<<
0054
0055
0056
0057 Par04EventAction::Par04EventAction(Par04DetectorConstruction* aDetector,
0058 Par04ParallelFullWorld* aParallel)
0059 : G4UserEventAction(),
0060 fHitCollectionID(-1),
0061 fPhysicalFullHitCollectionID(-1),
0062 fPhysicalFastHitCollectionID(-1),
0063 fTimer(),
0064 fDetector(aDetector),
0065 fParallel(aParallel)
0066 {
0067 fCellNbRho = aDetector->GetMeshNbOfCells().x();
0068 fCellNbPhi = aDetector->GetMeshNbOfCells().y();
0069 fCellNbZ = aDetector->GetMeshNbOfCells().z();
0070 fCalEdep.reserve(fCellNbRho * fCellNbPhi * fCellNbZ);
0071 fCalRho.reserve(fCellNbRho * fCellNbPhi * fCellNbZ);
0072 fCalPhi.reserve(fCellNbRho * fCellNbPhi * fCellNbZ);
0073 fCalZ.reserve(fCellNbRho * fCellNbPhi * fCellNbZ);
0074 fCalPhysicalEdep.reserve(fPhysicalNbLayers * fPhysicalNbRows * fPhysicalNbSlices);
0075 fCalPhysicalLayer.reserve(fPhysicalNbLayers * fPhysicalNbRows * fPhysicalNbSlices);
0076 fCalPhysicalSlice.reserve(fPhysicalNbLayers * fPhysicalNbRows * fPhysicalNbSlices);
0077 fCalPhysicalRow.reserve(fPhysicalNbLayers * fPhysicalNbRows * fPhysicalNbSlices);
0078 }
0079
0080
0081
0082 Par04EventAction::~Par04EventAction() = default;
0083
0084
0085
0086 void Par04EventAction::BeginOfEventAction(const G4Event*)
0087 {
0088 StartTimer();
0089 }
0090
0091
0092
0093 void Par04EventAction::StartTimer()
0094 {
0095 fTimer.Start();
0096 }
0097
0098
0099
0100 void Par04EventAction::StopTimer()
0101 {
0102 fTimer.Stop();
0103 }
0104
0105
0106
0107 void Par04EventAction::EndOfEventAction(const G4Event* aEvent)
0108 {
0109 G4SDManager::GetSDMpointer()->GetHCtable();
0110 StopTimer();
0111
0112
0113 if (fHitCollectionID == -1) {
0114 fHitCollectionID = G4SDManager::GetSDMpointer()->GetCollectionID("hits");
0115 }
0116 if (fPhysicalFullHitCollectionID == -1) {
0117 fPhysicalFullHitCollectionID =
0118 G4SDManager::GetSDMpointer()->GetCollectionID("physicalCellsFullSim");
0119 }
0120 if (fPhysicalFastHitCollectionID == -1) {
0121 fPhysicalFastHitCollectionID =
0122 G4SDManager::GetSDMpointer()->GetCollectionID("physicalCellsFastSim");
0123 }
0124
0125 auto hitsCollection =
0126 static_cast<Par04HitsCollection*>(aEvent->GetHCofThisEvent()->GetHC(fHitCollectionID));
0127 auto physicalFullHitsCollection = static_cast<Par04HitsCollection*>(
0128 aEvent->GetHCofThisEvent()->GetHC(fPhysicalFullHitCollectionID));
0129 auto physicalFastHitsCollection = static_cast<Par04HitsCollection*>(
0130 aEvent->GetHCofThisEvent()->GetHC(fPhysicalFastHitCollectionID));
0131
0132 if (hitsCollection == nullptr) {
0133 G4ExceptionDescription msg;
0134 msg << "Cannot access hitsCollection ID " << fHitCollectionID;
0135 G4Exception("Par04EventAction::GetHitsCollection()", "MyCode0001", FatalException, msg);
0136 }
0137 if (physicalFullHitsCollection == nullptr) {
0138 G4ExceptionDescription msg;
0139 msg << "Cannot access physical full sim hitsCollection ID " << fPhysicalFullHitCollectionID;
0140 G4Exception("Par04EventAction::GetHitsCollection()", "MyCode0001", FatalException, msg);
0141 }
0142 if (physicalFastHitsCollection == nullptr) {
0143 G4ExceptionDescription msg;
0144 msg << "Cannot access physical fast sim hitsCollection ID " << fPhysicalFastHitCollectionID;
0145 G4Exception("Par04EventAction::GetHitsCollection()", "MyCode0001", FatalException, msg);
0146 }
0147
0148 auto analysisManager = G4AnalysisManager::Instance();
0149
0150 if (fCellSizeZ == 0) {
0151 fCellSizeZ = fDetector->GetMeshSizeOfCells().z();
0152 fCellSizePhi = fDetector->GetMeshSizeOfCells().y();
0153 fCellSizeRho = fDetector->GetMeshSizeOfCells().x();
0154 fCellNbRho = fDetector->GetMeshNbOfCells().x();
0155 fCellNbPhi = fDetector->GetMeshNbOfCells().y();
0156 fCellNbZ = fDetector->GetMeshNbOfCells().z();
0157 }
0158 if (fPhysicalNbLayers == 0) {
0159 fPhysicalNbLayers = fParallel->GetNbOfLayers();
0160 fPhysicalNbSlices = fParallel->GetNbOfSlices();
0161 fPhysicalNbRows = fParallel->GetNbOfRows();
0162 }
0163
0164
0165
0166 auto primaryVertex =
0167 G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetPrimaryVertex();
0168 auto primaryParticle = primaryVertex->GetPrimary(0);
0169 G4double primaryEnergy = primaryParticle->GetTotalEnergy();
0170
0171
0172 auto primaryDirection = primaryParticle->GetMomentumDirection();
0173 auto primaryEntrance =
0174 primaryVertex->GetPosition() - primaryVertex->GetPosition().z() * primaryDirection;
0175
0176
0177 fCalEdep.resize(fCellNbRho * fCellNbPhi * fCellNbZ, 0);
0178 fCalRho.resize(fCellNbRho * fCellNbPhi * fCellNbZ, 0);
0179 fCalPhi.resize(fCellNbRho * fCellNbPhi * fCellNbZ, 0);
0180 fCalZ.resize(fCellNbRho * fCellNbPhi * fCellNbZ, 0);
0181 fCalPhysicalEdep.resize(fPhysicalNbLayers * fPhysicalNbRows * fPhysicalNbSlices, 0);
0182 fCalPhysicalLayer.resize(fPhysicalNbLayers * fPhysicalNbRows * fPhysicalNbSlices, 0);
0183 fCalPhysicalSlice.resize(fPhysicalNbLayers * fPhysicalNbRows * fPhysicalNbSlices, 0);
0184 fCalPhysicalRow.resize(fPhysicalNbLayers * fPhysicalNbRows * fPhysicalNbSlices, 0);
0185
0186
0187 Par04Hit* hit = nullptr;
0188 G4double hitEn = 0;
0189 G4double totalEnergy = 0;
0190 G4int hitNum = 0;
0191 G4int totalNum = 0;
0192 G4int hitZ = -1;
0193 G4int hitRho = -1;
0194 G4int hitPhi = -1;
0195 G4int hitType = -1;
0196 G4int numNonZeroThresholdCells = 0;
0197 G4double tDistance = 0., rDistance = 0., phiDistance = 0.;
0198 G4double tFirstMoment = 0., tSecondMoment = 0.;
0199 G4double rFirstMoment = 0., rSecondMoment = 0.;
0200 G4double phiMean = 0.;
0201 for (size_t iHit = 0; iHit < hitsCollection->entries(); iHit++) {
0202 hit = static_cast<Par04Hit*>(hitsCollection->GetHit(iHit));
0203 hitZ = hit->GetZid();
0204 hitRho = hit->GetRhoId();
0205 hitPhi = hit->GetPhiId();
0206 hitEn = hit->GetEdep();
0207 hitNum = hit->GetNdep();
0208 hitType = hit->GetType();
0209 if (hitEn > 0) {
0210 totalEnergy += hitEn;
0211 totalNum += hitNum;
0212 tDistance = hitZ * fCellSizeZ;
0213 rDistance = hitRho * fCellSizeRho;
0214 phiDistance = hitPhi * fCellSizePhi;
0215 tFirstMoment += hitEn * tDistance;
0216 rFirstMoment += hitEn * rDistance;
0217 phiMean += hitEn * phiDistance;
0218 analysisManager->FillH1(4, tDistance, hitEn);
0219 analysisManager->FillH1(5, rDistance, hitEn);
0220 analysisManager->FillH1(10, hitType);
0221 if (hitEn > 0.0005) {
0222 fCalEdep[numNonZeroThresholdCells] = hitEn;
0223 fCalRho[numNonZeroThresholdCells] = hitRho;
0224 fCalPhi[numNonZeroThresholdCells] = hitPhi;
0225 fCalZ[numNonZeroThresholdCells] = hitZ;
0226 numNonZeroThresholdCells++;
0227 analysisManager->FillH1(13, std::log10(hitEn));
0228 analysisManager->FillH1(15, hitNum);
0229 }
0230 }
0231 }
0232 tFirstMoment /= totalEnergy;
0233 rFirstMoment /= totalEnergy;
0234 phiMean /= totalEnergy;
0235 analysisManager->FillH1(0, primaryEnergy / GeV);
0236 analysisManager->FillH1(1, totalEnergy / GeV);
0237 analysisManager->FillH1(2, totalEnergy / primaryEnergy);
0238 analysisManager->FillH1(3, fTimer.GetRealElapsed());
0239 analysisManager->FillH1(6, tFirstMoment);
0240 analysisManager->FillH1(7, rFirstMoment);
0241 analysisManager->FillH1(12, numNonZeroThresholdCells);
0242 analysisManager->FillH1(14, totalNum);
0243
0244 fCalEdep.resize(numNonZeroThresholdCells);
0245 fCalRho.resize(numNonZeroThresholdCells);
0246 fCalPhi.resize(numNonZeroThresholdCells);
0247 fCalZ.resize(numNonZeroThresholdCells);
0248 analysisManager->FillNtupleDColumn(0, 0, primaryEnergy);
0249 analysisManager->FillNtupleDColumn(0, 1, fTimer.GetRealElapsed());
0250
0251 for (size_t iHit = 0; iHit < hitsCollection->entries(); iHit++) {
0252 hit = static_cast<Par04Hit*>(hitsCollection->GetHit(iHit));
0253 hitEn = hit->GetEdep();
0254 hitZ = hit->GetZid();
0255 hitRho = hit->GetRhoId();
0256 hitPhi = hit->GetPhiId();
0257 if (hitEn > 0) {
0258 tDistance = hitZ * fCellSizeZ;
0259 rDistance = hitRho * fCellSizeRho;
0260 phiDistance = hitPhi * fCellSizePhi;
0261 tSecondMoment += hitEn * std::pow(tDistance - tFirstMoment, 2);
0262 rSecondMoment += hitEn * std::pow(rDistance - rFirstMoment, 2);
0263 analysisManager->FillH1(11, phiDistance - phiMean, hitEn);
0264 }
0265 }
0266 tSecondMoment /= totalEnergy;
0267 rSecondMoment /= totalEnergy;
0268 analysisManager->FillH1(8, tSecondMoment);
0269 analysisManager->FillH1(9, rSecondMoment);
0270
0271
0272 G4double totalPhysicalEnergy = 0;
0273 totalNum = 0;
0274 hitEn = 0;
0275 hitNum = 0;
0276 G4int hitLayer = -1;
0277 G4int hitRow = -1;
0278 G4int hitSlice = -1;
0279 numNonZeroThresholdCells = 0;
0280 for (size_t iHit = 0; iHit < physicalFullHitsCollection->entries(); iHit++) {
0281 hit = static_cast<Par04Hit*>(physicalFullHitsCollection->GetHit(iHit));
0282 hitLayer = hit->GetRhoId();
0283 hitRow = hit->GetZid();
0284 hitSlice = hit->GetPhiId();
0285 hitEn = hit->GetEdep();
0286 hitNum = hit->GetNdep();
0287 if (hitEn > 0) {
0288 totalPhysicalEnergy += hitEn;
0289 totalNum += hitNum;
0290 if (hitEn > 0.0005) {
0291 fCalPhysicalEdep[numNonZeroThresholdCells] = hitEn;
0292 fCalPhysicalLayer[numNonZeroThresholdCells] = hitLayer;
0293 fCalPhysicalRow[numNonZeroThresholdCells] = hitRow;
0294 fCalPhysicalSlice[numNonZeroThresholdCells] = hitSlice;
0295 numNonZeroThresholdCells++;
0296 analysisManager->FillH1(19, std::log10(hitEn));
0297 analysisManager->FillH1(21, hitNum);
0298 }
0299 }
0300 }
0301 for (size_t iHit = 0; iHit < physicalFastHitsCollection->entries(); iHit++) {
0302 hit = static_cast<Par04Hit*>(physicalFastHitsCollection->GetHit(iHit));
0303 hitLayer = hit->GetRhoId();
0304 hitRow = hit->GetZid();
0305 hitSlice = hit->GetPhiId();
0306 hitEn = hit->GetEdep();
0307 hitNum = hit->GetNdep();
0308 if (hitEn > 0) {
0309 totalPhysicalEnergy += hitEn;
0310 totalNum += hitNum;
0311 if (hitEn > 0.0005) {
0312 fCalPhysicalEdep[numNonZeroThresholdCells] = hitEn;
0313 fCalPhysicalLayer[numNonZeroThresholdCells] = hitLayer;
0314 fCalPhysicalRow[numNonZeroThresholdCells] = hitRow;
0315 fCalPhysicalSlice[numNonZeroThresholdCells] = hitSlice;
0316 numNonZeroThresholdCells++;
0317 analysisManager->FillH1(19, std::log10(hitEn));
0318 analysisManager->FillH1(21, hitNum);
0319 }
0320 }
0321 }
0322 analysisManager->FillH1(16, totalPhysicalEnergy / GeV);
0323 analysisManager->FillH1(17, totalPhysicalEnergy / primaryEnergy);
0324 analysisManager->FillH1(18, numNonZeroThresholdCells);
0325 analysisManager->FillH1(20, totalNum);
0326 fCalPhysicalEdep.resize(numNonZeroThresholdCells);
0327 fCalPhysicalLayer.resize(numNonZeroThresholdCells);
0328 fCalPhysicalSlice.resize(numNonZeroThresholdCells);
0329 fCalPhysicalRow.resize(numNonZeroThresholdCells);
0330 analysisManager->AddNtupleRow(0);
0331 analysisManager->AddNtupleRow(1);
0332 analysisManager->AddNtupleRow(2);
0333 }