Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-12-16 09:30:17

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 /// \file electromagnetic/TestEm13/src/Run.cc
0027 /// \brief Implementation of the Run class
0028 //
0029 //
0030 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0031 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0046 
0047 Run::Run(DetectorConstruction* det) : fDetector(det) {}
0048 
0049 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0050 
0051 void Run::SetPrimary(G4ParticleDefinition* particle, G4double energy)
0052 {
0053   fParticle = particle;
0054   fEkin = energy;
0055 }
0056 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0070 
0071 void Run::Merge(const G4Run* run)
0072 {
0073   const Run* localRun = static_cast<const Run*>(run);
0074 
0075   // pass information about primary particle
0076   fParticle = localRun->fParticle;
0077   fEkin = localRun->fEkin;
0078 
0079   // map: processes count
0080   std::map<G4String, G4int>::const_iterator it;
0081   for (it = localRun->fProcCounter.begin(); it != localRun->fProcCounter.end(); ++it) {
0082     G4String procName = it->first;
0083     G4int localCount = it->second;
0084     if (fProcCounter.find(procName) == fProcCounter.end()) {
0085       fProcCounter[procName] = localCount;
0086     }
0087     else {
0088       fProcCounter[procName] += localCount;
0089     }
0090   }
0091 
0092   G4Run::Merge(run);
0093 }
0094 
0095 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0096 
0097 void Run::EndOfRun()
0098 {
0099   G4int prec = 5;
0100   G4int dfprec = G4cout.precision(prec);
0101 
0102   // run condition
0103   //
0104   G4String partName = fParticle->GetParticleName();
0105   G4Material* material = fDetector->GetMaterial();
0106   G4double density = material->GetDensity();
0107   G4double tickness = fDetector->GetSize();
0108 
0109   G4cout << "\n ======================== run summary ======================\n";
0110   G4cout << "\n The run is: " << numberOfEvent << " " << partName << " of "
0111          << G4BestUnit(fEkin, "Energy") << " through " << G4BestUnit(tickness, "Length") << " of "
0112          << material->GetName() << " (density: " << G4BestUnit(density, "Volumic Mass") << ")"
0113          << G4endl;
0114 
0115   // frequency of processes
0116   G4int totalCount = 0;
0117   G4int survive = 0;
0118   G4cout << "\n Process calls frequency --->";
0119   std::map<G4String, G4int>::iterator it;
0120   for (it = fProcCounter.begin(); it != fProcCounter.end(); it++) {
0121     G4String procName = it->first;
0122     G4int count = it->second;
0123     totalCount += count;
0124     G4cout << "\t" << procName << " = " << count;
0125     if (procName == "Transportation") survive = count;
0126   }
0127   G4cout << G4endl;
0128 
0129   if (totalCount == 0) {
0130     G4cout.precision(dfprec);
0131     return;
0132   };
0133   G4double ratio = double(survive) / totalCount;
0134 
0135   G4cout << "\n Nb of incident particles unaltered after " << G4BestUnit(tickness, "Length")
0136          << " of " << material->GetName() << " : " << survive << " over " << totalCount
0137          << " incident particles."
0138          << "  Ratio = " << 100 * ratio << " %" << G4endl;
0139 
0140   if (ratio == 0.) return;
0141 
0142   // compute cross section and related quantities
0143   //
0144   G4double CrossSection = -std::log(ratio) / tickness;
0145   G4double massicCS = CrossSection / density;
0146 
0147   G4cout << " ---> CrossSection per volume:\t" << CrossSection * cm << " cm^-1 "
0148          << "\tCrossSection per mass: " << G4BestUnit(massicCS, "Surface/Mass") << G4endl;
0149 
0150   // check cross section from G4EmCalculator
0151   //
0152   G4cout << "\n Verification from G4EmCalculator: \n";
0153   G4EmCalculator emCalculator;
0154   G4double sumc = 0.0;
0155   for (it = fProcCounter.begin(); it != fProcCounter.end(); it++) {
0156     G4String procName = it->first;
0157     G4double massSigma =
0158       emCalculator.GetCrossSectionPerVolume(fEkin, fParticle, procName, material) / density;
0159     if (fParticle == G4Gamma::Gamma())
0160       massSigma =
0161         emCalculator.ComputeCrossSectionPerVolume(fEkin, fParticle, procName, material) / density;
0162     sumc += massSigma;
0163     if (procName != "Transportation")
0164       G4cout << "\t" << procName << "= " << G4BestUnit(massSigma, "Surface/Mass");
0165   }
0166   G4cout << "\ttotal= " << G4BestUnit(sumc, "Surface/Mass") << G4endl;
0167 
0168   // expected ratio of transmitted particles
0169   G4double Ratio = std::exp(-sumc * density * tickness);
0170   G4cout << "\tExpected ratio of transmitted particles= " << 100 * Ratio << " %" << G4endl;
0171 
0172   // remove all contents in fProcCounter
0173   fProcCounter.clear();
0174 
0175   // restore default format
0176   G4cout.precision(dfprec);
0177 }
0178 
0179 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......