File indexing completed on 2025-02-23 09:22:47
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 "G4Electron.hh"
0039 #include "G4Gamma.hh"
0040 #include "G4ParticleDefinition.hh"
0041 #include "G4ParticleTable.hh"
0042 #include "G4Positron.hh"
0043 #include "G4SystemOfUnits.hh"
0044 #include "G4Track.hh"
0045 #include "G4UnitsTable.hh"
0046
0047 #include <iomanip>
0048
0049
0050
0051 Run::Run(DetectorConstruction* det)
0052 : G4Run(),
0053 fDetector(det),
0054 fParticle(nullptr),
0055 fEkin(0.),
0056 fChargedStep(0),
0057 fNeutralStep(0),
0058 fN_gamma(0),
0059 fN_elec(0),
0060 fN_pos(0)
0061 {
0062
0063
0064 for (G4int k = 0; k < kMaxAbsor; k++) {
0065 fSumEAbs[k] = fSum2EAbs[k] = fSumLAbs[k] = fSum2LAbs[k] = 0.;
0066 fEnergyDeposit[k].clear();
0067 }
0068 }
0069
0070
0071
0072 Run::~Run() {}
0073
0074
0075
0076 void Run::SetPrimary(G4ParticleDefinition* particle, G4double energy)
0077 {
0078 fParticle = particle;
0079 fEkin = energy;
0080 }
0081
0082
0083
0084 void Run::FillPerEvent(G4int kAbs, G4double EAbs, G4double LAbs)
0085 {
0086
0087
0088 fEnergyDeposit[kAbs].push_back(EAbs);
0089 fSumEAbs[kAbs] += EAbs;
0090 fSum2EAbs[kAbs] += EAbs * EAbs;
0091 fSumLAbs[kAbs] += LAbs;
0092 fSum2LAbs[kAbs] += LAbs * LAbs;
0093 }
0094
0095
0096
0097 void Run::AddChargedStep()
0098 {
0099 fChargedStep += 1.0;
0100 }
0101
0102
0103
0104 void Run::AddNeutralStep()
0105 {
0106 fNeutralStep += 1.0;
0107 }
0108
0109
0110
0111 void Run::AddSecondaryTrack(const G4Track* track)
0112 {
0113 const G4ParticleDefinition* d = track->GetDefinition();
0114 if (d == G4Gamma::Gamma()) {
0115 ++fN_gamma;
0116 }
0117 else if (d == G4Electron::Electron()) {
0118 ++fN_elec;
0119 }
0120 else if (d == G4Positron::Positron()) {
0121 ++fN_pos;
0122 }
0123 }
0124
0125
0126
0127 void Run::Merge(const G4Run* run)
0128 {
0129 const Run* localRun = static_cast<const Run*>(run);
0130
0131
0132 fParticle = localRun->fParticle;
0133 fEkin = localRun->fEkin;
0134
0135
0136
0137 for (G4int k = 0; k < kMaxAbsor; k++) {
0138 fSumEAbs[k] += localRun->fSumEAbs[k];
0139 fSum2EAbs[k] += localRun->fSum2EAbs[k];
0140 fSumLAbs[k] += localRun->fSumLAbs[k];
0141 fSum2LAbs[k] += localRun->fSum2LAbs[k];
0142 }
0143
0144 fChargedStep += localRun->fChargedStep;
0145 fNeutralStep += localRun->fNeutralStep;
0146
0147 fN_gamma += localRun->fN_gamma;
0148 fN_elec += localRun->fN_elec;
0149 fN_pos += localRun->fN_pos;
0150
0151 G4Run::Merge(run);
0152 }
0153
0154
0155
0156 void Run::EndOfRun()
0157 {
0158 G4int nEvt = numberOfEvent;
0159 G4double norm = G4double(nEvt);
0160 if (norm > 0) norm = 1. / norm;
0161 G4double qnorm = std::sqrt(norm);
0162
0163 fChargedStep *= norm;
0164 fNeutralStep *= norm;
0165
0166
0167
0168 G4double beamEnergy = fEkin;
0169 G4double sqbeam = std::sqrt(beamEnergy / GeV);
0170
0171 G4double MeanEAbs, MeanEAbs2, rmsEAbs, resolution, rmsres;
0172 G4double MeanLAbs, MeanLAbs2, rmsLAbs;
0173
0174 std::ios::fmtflags mode = G4cout.flags();
0175 G4int prec = G4cout.precision(2);
0176 G4cout << "\n------------------------------------------------------------\n";
0177 G4cout << std::setw(14) << "material" << std::setw(17) << "Edep RMS" << std::setw(33)
0178 << "sqrt(E0(GeV))*rmsE/Emean" << std::setw(23) << "total tracklen \n \n";
0179
0180 for (G4int k = 1; k <= fDetector->GetNbOfAbsor(); k++) {
0181 MeanEAbs = fSumEAbs[k] * norm;
0182 MeanEAbs2 = fSum2EAbs[k] * norm;
0183 rmsEAbs = std::sqrt(std::abs(MeanEAbs2 - MeanEAbs * MeanEAbs));
0184
0185 resolution = 100. * sqbeam * rmsEAbs / MeanEAbs;
0186 rmsres = resolution * qnorm;
0187
0188
0189 fSumEAbs[k] = MeanEAbs;
0190 fSum2EAbs[k] = rmsEAbs;
0191
0192 MeanLAbs = fSumLAbs[k] * norm;
0193 MeanLAbs2 = fSum2LAbs[k] * norm;
0194 rmsLAbs = std::sqrt(std::abs(MeanLAbs2 - MeanLAbs * MeanLAbs));
0195
0196
0197
0198 G4cout << std::setw(14) << fDetector->GetAbsorMaterial(k)->GetName() << ": "
0199 << std::setprecision(5) << std::setw(6) << G4BestUnit(MeanEAbs, "Energy") << " : "
0200 << std::setprecision(4) << std::setw(5) << G4BestUnit(rmsEAbs, "Energy") << std::setw(10)
0201 << resolution << " +- " << std::setw(5) << rmsres << " %" << std::setprecision(3)
0202 << std::setw(10) << G4BestUnit(MeanLAbs, "Length") << " +- " << std::setw(4)
0203 << G4BestUnit(rmsLAbs, "Length") << G4endl;
0204 }
0205 G4cout << "\n------------------------------------------------------------\n";
0206
0207 G4cout << " Beam particle " << fParticle->GetParticleName()
0208 << " E = " << G4BestUnit(beamEnergy, "Energy") << G4endl;
0209 G4cout << " Mean number of gamma " << (G4double)fN_gamma * norm << G4endl;
0210 G4cout << " Mean number of e- " << (G4double)fN_elec * norm << G4endl;
0211 G4cout << " Mean number of e+ " << (G4double)fN_pos * norm << G4endl;
0212 G4cout << std::setprecision(6) << " Mean number of charged steps " << fChargedStep << G4endl;
0213 G4cout << " Mean number of neutral steps " << fNeutralStep << G4endl;
0214 G4cout << "------------------------------------------------------------\n" << G4endl;
0215
0216 G4cout.setf(mode, std::ios::floatfield);
0217 G4cout.precision(prec);
0218 }
0219
0220