File indexing completed on 2025-02-23 09:20:53
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
0033 #include "Run.hh"
0034
0035 #include "DetectorConstruction.hh"
0036 #include "PrimaryGeneratorAction.hh"
0037
0038 #include "G4EmCalculator.hh"
0039 #include "G4Gamma.hh"
0040 #include "G4SystemOfUnits.hh"
0041 #include "G4UnitsTable.hh"
0042
0043 #include <iomanip>
0044
0045
0046
0047 Run::Run(DetectorConstruction* det) : fDetector(det) {}
0048
0049
0050
0051 void Run::SetPrimary(G4ParticleDefinition* particle, G4double energy)
0052 {
0053 fParticle = particle;
0054 fEkin = energy;
0055 }
0056
0057
0058 void Run::CountProcesses(G4String procName)
0059 {
0060 std::map<G4String, G4int>::iterator it = fProcCounter.find(procName);
0061 if (it == fProcCounter.end()) {
0062 fProcCounter[procName] = 1;
0063 }
0064 else {
0065 fProcCounter[procName]++;
0066 }
0067 }
0068
0069
0070
0071 void Run::SumTrack(G4double track)
0072 {
0073 fTotalCount++;
0074 fSumTrack += track;
0075 fSumTrack2 += track * track;
0076 }
0077
0078
0079
0080 void Run::SumeTransf(G4double energy)
0081 {
0082 fEnTransfer += energy;
0083 }
0084
0085
0086
0087 void Run::Merge(const G4Run* run)
0088 {
0089 const Run* localRun = static_cast<const Run*>(run);
0090
0091
0092 fParticle = localRun->fParticle;
0093 fEkin = localRun->fEkin;
0094
0095
0096 std::map<G4String, G4int>::const_iterator it;
0097 for (it = localRun->fProcCounter.begin(); it != localRun->fProcCounter.end(); ++it) {
0098 G4String procName = it->first;
0099 G4int localCount = it->second;
0100 if (fProcCounter.find(procName) == fProcCounter.end()) {
0101 fProcCounter[procName] = localCount;
0102 }
0103 else {
0104 fProcCounter[procName] += localCount;
0105 }
0106 }
0107
0108 fTotalCount += localRun->fTotalCount;
0109 fSumTrack += localRun->fSumTrack;
0110 fSumTrack2 += localRun->fSumTrack2;
0111 fEnTransfer += localRun->fEnTransfer;
0112
0113 G4Run::Merge(run);
0114 }
0115
0116
0117
0118 void Run::EndOfRun()
0119 {
0120 G4int prec = 5;
0121 G4int dfprec = G4cout.precision(prec);
0122
0123
0124
0125 G4String partName = fParticle->GetParticleName();
0126 G4Material* material = fDetector->GetMaterial();
0127 G4double density = material->GetDensity();
0128 G4double tickness = fDetector->GetSize();
0129
0130 G4cout << "\n ======================== run summary ======================\n";
0131 G4cout << "\n The run is: " << numberOfEvent << " " << partName << " of "
0132 << G4BestUnit(fEkin, "Energy") << " through " << G4BestUnit(tickness, "Length") << " of "
0133 << material->GetName() << " (density: " << G4BestUnit(density, "Volumic Mass") << ")"
0134 << G4endl;
0135
0136
0137 G4int survive = 0;
0138 G4cout << "\n Process calls frequency --->";
0139 std::map<G4String, G4int>::iterator it;
0140 for (it = fProcCounter.begin(); it != fProcCounter.end(); it++) {
0141 G4String procName = it->first;
0142 G4int count = it->second;
0143 G4cout << "\t" << procName << " = " << count;
0144 if (procName == "Transportation") survive = count;
0145 }
0146
0147 if (survive > 0) {
0148 G4cout << "\n\n Nb of incident particles surviving after "
0149 << G4BestUnit(fDetector->GetSize(), "Length") << " of " << material->GetName() << " : "
0150 << survive << G4endl;
0151 }
0152
0153 if (fTotalCount == 0) fTotalCount = 1;
0154
0155
0156
0157 G4double MeanFreePath = fSumTrack / fTotalCount;
0158 G4double MeanTrack2 = fSumTrack2 / fTotalCount;
0159 G4double rms = std::sqrt(std::fabs(MeanTrack2 - MeanFreePath * MeanFreePath));
0160 G4double CrossSection = 1. / MeanFreePath;
0161 G4double massicMFP = MeanFreePath * density;
0162 G4double massicCS = 1. / massicMFP;
0163
0164 G4cout << "\n\n MeanFreePath:\t" << G4BestUnit(MeanFreePath, "Length") << " +- "
0165 << G4BestUnit(rms, "Length") << "\tmassic: " << G4BestUnit(massicMFP, "Mass/Surface")
0166 << "\n CrossSection:\t" << CrossSection * cm << " cm^-1 "
0167 << "\t\t\tmassic: " << G4BestUnit(massicCS, "Surface/Mass") << G4endl;
0168
0169
0170
0171 G4double MeanTransfer = fEnTransfer / fTotalCount;
0172 G4double massTransfCoef = massicCS * MeanTransfer / fEkin;
0173
0174 G4cout << "\n mean energy of charged secondaries: " << G4BestUnit(MeanTransfer, "Energy")
0175 << "\n ---> mass_energy_transfer coef: " << G4BestUnit(massTransfCoef, "Surface/Mass")
0176 << G4endl;
0177
0178
0179
0180 G4cout << "\n Verification : "
0181 << "crossSections from G4EmCalculator \n";
0182
0183 G4EmCalculator emCalculator;
0184 G4double sumc = 0.0;
0185 for (it = fProcCounter.begin(); it != fProcCounter.end(); it++) {
0186 G4String procName = it->first;
0187 G4double massSigma =
0188 emCalculator.GetCrossSectionPerVolume(fEkin, fParticle, procName, material) / density;
0189 if (fParticle == G4Gamma::Gamma())
0190 massSigma =
0191 emCalculator.ComputeCrossSectionPerVolume(fEkin, fParticle, procName, material) / density;
0192 sumc += massSigma;
0193 G4cout << " " << procName << "= " << G4BestUnit(massSigma, "Surface/Mass");
0194 }
0195 G4cout << " total= " << G4BestUnit(sumc, "Surface/Mass") << G4endl;
0196
0197
0198 fProcCounter.clear();
0199
0200
0201 G4cout.precision(dfprec);
0202 }
0203
0204