File indexing completed on 2026-04-07 07:52:20
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 "Run.hh"
0030
0031 #include "G4Run.hh"
0032 #include "G4RunManager.hh"
0033 #include "G4SystemOfUnits.hh"
0034
0035
0036
0037 Run::Run()
0038 : G4Run(),
0039 fNumEvents(0),
0040 fPrimaryParticleId(0),
0041 fPrimaryParticleEnergy(0.0),
0042 fPrimaryParticleDirection(G4ThreeVector(0.0, 0.0, 0.0)),
0043 fTrackerMaterialName(""),
0044 fEmCaloMaterialName(""),
0045 fHadCaloMaterialName(""),
0046 fCubicVolumeScoringTrackerShell(1.0),
0047 fCubicVolumeScoringEmCaloShell(1.0),
0048 fCubicVolumeScoringHadCaloShell(1.0)
0049 {
0050 fSteppingArray.fill(0.0);
0051 fTrackingArray1.fill(0);
0052 fTrackingArray2.fill(0.0);
0053 }
0054
0055
0056
0057 void Run::RecordEvent(const G4Event* anEvent)
0058 {
0059
0060
0061 G4int nEvt = anEvent->GetEventID();
0062 if (nEvt % 10 == 0) G4cout << " Event#=" << nEvt << G4endl;
0063 G4Run::RecordEvent(anEvent);
0064 }
0065
0066
0067
0068 void Run::Merge(const G4Run* aRun)
0069 {
0070
0071
0072 const Run* localRun = static_cast<const Run*>(aRun);
0073 fPrimaryParticleId = localRun->GetPrimaryParticleId();
0074 fPrimaryParticleEnergy = localRun->GetPrimaryParticleEnergy();
0075 fPrimaryParticleDirection = localRun->GetPrimaryParticleDirection();
0076 fTrackerMaterialName = localRun->GetTrackerMaterialName();
0077 fEmCaloMaterialName = localRun->GetEmCaloMaterialName();
0078 fHadCaloMaterialName = localRun->GetHadCaloMaterialName();
0079 fCubicVolumeScoringTrackerShell = localRun->GetCubicVolumeScoringTrackerShell();
0080 fCubicVolumeScoringEmCaloShell = localRun->GetCubicVolumeScoringEmCaloShell();
0081 fCubicVolumeScoringHadCaloShell = localRun->GetCubicVolumeScoringHadCaloShell();
0082 fNumEvents += localRun->GetNumberOfEvent();
0083 for (G4int i = 0; i < SteppingAction::fkNumberCombinations; ++i) {
0084 fSteppingArray[i] += localRun->GetSteppingArray()[i];
0085 }
0086 for (G4int i = 0; i < TrackingAction::fkNumberCombinations; ++i) {
0087 fTrackingArray1[i] += localRun->GetTrackingArray1()[i];
0088 fTrackingArray2[i] += localRun->GetTrackingArray2()[i];
0089 }
0090 G4Run::Merge(aRun);
0091 }
0092
0093
0094
0095 void Run::PrintInfo() const
0096 {
0097
0098
0099 const G4double floatingNumberOfEvents =
0100 std::max(1.0, fNumEvents > 0 ? fNumEvents * 1.0 : GetNumberOfEvent() * 1.0);
0101
0102
0103 const G4double conversionFactor = CLHEP::cm * CLHEP::cm;
0104 const G4double factorTracker =
0105 conversionFactor / (0.5 * fCubicVolumeScoringTrackerShell * floatingNumberOfEvents);
0106 const G4double factorEmCalo =
0107 conversionFactor / (0.5 * fCubicVolumeScoringEmCaloShell * floatingNumberOfEvents);
0108 const G4double factorHadCalo =
0109 conversionFactor / (0.5 * fCubicVolumeScoringHadCaloShell * floatingNumberOfEvents);
0110 G4cout << std::setprecision(6) << G4endl << G4endl
0111 << " =============== Run::PrintInfo() =============== \t RunID = " << GetRunID()
0112 << G4endl << " Primary particle PDG code = " << fPrimaryParticleId << G4endl
0113 << " Primary particle kinetic energy = " << fPrimaryParticleEnergy / CLHEP::GeV << " GeV"
0114 << G4endl << " Primary particle direction = " << fPrimaryParticleDirection << G4endl
0115 << " Tracker material = " << fTrackerMaterialName << G4endl
0116 << " EmCalo material = " << fEmCaloMaterialName << G4endl
0117 << " HadCalo material = " << fHadCaloMaterialName << G4endl
0118 << " Cubic-volume scoring tracker shell = " << fCubicVolumeScoringTrackerShell << " mm^3"
0119 << G4endl << " Cubic-volume scoring emCalo shell = " << fCubicVolumeScoringEmCaloShell
0120 << " mm^3" << G4endl
0121 << " Cubic-volume scoring hadCalo shell = " << fCubicVolumeScoringHadCaloShell << " mm^3"
0122 << G4endl << " Number of events = " << floatingNumberOfEvents << G4endl
0123 << " Conversion factor: fluence from mm^-2 to cm^-2 = " << conversionFactor << G4endl
0124 << " Particle fluence in unit of cm^-2 :" << G4endl;
0125 for (G4int i = 0; i < SteppingAction::fkNumberScoringShells; ++i) {
0126 G4double factor = factorTracker;
0127 if (i == 1)
0128 factor = factorEmCalo;
0129 else if (i == 2)
0130 factor = factorHadCalo;
0131 for (G4int j = 0; j < SteppingAction::fkNumberKinematicRegions; ++j) {
0132 for (G4int k = 0; k < SteppingAction::fkNumberScoringPositions; ++k) {
0133 for (G4int ll = 0; ll < SteppingAction::fkNumberParticleTypes; ++ll) {
0134 G4int index = SteppingAction::GetIndex(i, j, k, ll);
0135
0136
0137 G4cout << " case=" << std::setw(3) << index << " " << std::setw(12)
0138 << SteppingAction::fkArrayScoringShellNames[i] << " " << std::setw(12)
0139 << SteppingAction::fkArrayKinematicRegionNames[j] << " " << std::setw(12)
0140 << SteppingAction::fkArrayScoringPositionNames[k] << " " << std::setw(12)
0141 << SteppingAction::fkArrayParticleTypeNames[ll] << " " << std::setw(8)
0142 << factor * fSteppingArray[index] << G4endl;
0143 }
0144 }
0145 }
0146 }
0147 G4cout << " ------------------------------------------------------------- " << G4endl
0148 << " Extra information: particle production \t \t <N> <E_kin> <Sum_Ekin> [MeV]"
0149 << G4endl;
0150 const G4double normalization = 1.0 / floatingNumberOfEvents;
0151 for (G4int i = 0; i < TrackingAction::fkNumberScoringVolumes; ++i) {
0152 for (G4int j = 0; j < TrackingAction::fkNumberKinematicRegions; ++j) {
0153 for (G4int k = 0; k < TrackingAction::fkNumberParticleTypes; ++k) {
0154 G4int index = TrackingAction::GetIndex(i, j, k);
0155
0156 G4cout << " case=" << std::setw(3) << index << " " << std::setw(12)
0157 << TrackingAction::fkArrayScoringVolumeNames[i] << " " << std::setw(12)
0158 << TrackingAction::fkArrayKinematicRegionNames[j] << " " << std::setw(12)
0159 << TrackingAction::fkArrayParticleTypeNames[k] << " " << std::setw(8)
0160 << normalization * fTrackingArray1[index] << " " << std::setw(8)
0161 << (fTrackingArray1[index] > 0 ? fTrackingArray2[index] / fTrackingArray1[index]
0162 : 0.0)
0163 << " " << std::setw(8) << normalization * fTrackingArray2[index] << G4endl;
0164 }
0165 }
0166 }
0167 G4cout << " ============================================================= " << G4endl << G4endl;
0168 }
0169
0170
0171
0172 void Run::SetSteppingArray(
0173 const std::array<G4double, SteppingAction::fkNumberCombinations>& inputArray)
0174 {
0175 for (G4int i = 0; i < SteppingAction::fkNumberCombinations; ++i) {
0176 fSteppingArray[i] = inputArray[i];
0177 }
0178 }
0179
0180
0181
0182 void Run::SetTrackingArray1(
0183 const std::array<G4long, TrackingAction::fkNumberCombinations>& inputArray)
0184 {
0185 for (G4int i = 0; i < TrackingAction::fkNumberCombinations; ++i) {
0186 fTrackingArray1[i] = inputArray[i];
0187 }
0188 }
0189
0190
0191
0192 void Run::SetTrackingArray2(
0193 const std::array<G4double, TrackingAction::fkNumberCombinations>& inputArray)
0194 {
0195 for (G4int i = 0; i < TrackingAction::fkNumberCombinations; ++i) {
0196 fTrackingArray2[i] = inputArray[i];
0197 }
0198 }
0199
0200