Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-02-23 09:20:55

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/TestEm17/src/RunAction.cc
0027 /// \brief Implementation of the RunAction class
0028 //
0029 //
0030 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0031 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0032 
0033 #include "RunAction.hh"
0034 
0035 #include "DetectorConstruction.hh"
0036 #include "HistoManager.hh"
0037 #include "MuCrossSections.hh"
0038 #include "PrimaryGeneratorAction.hh"
0039 
0040 #include "G4EmCalculator.hh"
0041 #include "G4PhysicalConstants.hh"
0042 #include "G4ProductionCutsTable.hh"
0043 #include "G4Run.hh"
0044 #include "G4RunManager.hh"
0045 #include "G4SystemOfUnits.hh"
0046 #include "G4UnitsTable.hh"
0047 #include "Randomize.hh"
0048 
0049 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0050 
0051 RunAction::RunAction(DetectorConstruction* det, PrimaryGeneratorAction* prim, HistoManager* HistM)
0052   : G4UserRunAction(), fDetector(det), fPrimary(prim), fProcCounter(0), fHistoManager(HistM)
0053 {
0054   fMucs = new MuCrossSections();
0055 }
0056 
0057 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0058 
0059 RunAction::~RunAction()
0060 {
0061   delete fMucs;
0062 }
0063 
0064 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0065 
0066 void RunAction::BeginOfRunAction(const G4Run* aRun)
0067 {
0068   G4cout << "### Run " << aRun->GetRunID() << " start." << G4endl;
0069 
0070   // save Rndm status
0071   CLHEP::HepRandom::showEngineStatus();
0072 
0073   fProcCounter = new ProcessesCount();
0074   fHistoManager->Book();
0075 }
0076 
0077 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0078 
0079 void RunAction::CountProcesses(const G4String& procName)
0080 {
0081   // does the process  already encounted ?
0082   size_t n = fProcCounter->size();
0083   for (size_t i = 0; i < n; ++i) {
0084     if ((*fProcCounter)[i]->GetName() == procName) {
0085       (*fProcCounter)[i]->Count();
0086       return;
0087     }
0088   }
0089   OneProcessCount* count = new OneProcessCount(procName);
0090   count->Count();
0091   fProcCounter->push_back(count);
0092 }
0093 
0094 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0095 
0096 void RunAction::EndOfRunAction(const G4Run* aRun)
0097 {
0098   G4int NbOfEvents = aRun->GetNumberOfEvent();
0099   if (NbOfEvents == 0) return;
0100 
0101   //  std::ios::fmtflags mode = G4cout.flags();
0102   G4int prec = G4cout.precision(2);
0103 
0104   const G4Material* material = fDetector->GetMaterial();
0105   G4double length = fDetector->GetSize();
0106   G4double density = material->GetDensity();
0107 
0108   G4String particle = fPrimary->GetParticleGun()->GetParticleDefinition()->GetParticleName();
0109   G4double energy = fPrimary->GetParticleGun()->GetParticleEnergy();
0110 
0111   G4cout << "\n The run consists of " << NbOfEvents << " " << particle << " of "
0112          << G4BestUnit(energy, "Energy") << " through " << G4BestUnit(length, "Length") << " of "
0113          << material->GetName() << " (density: " << G4BestUnit(density, "Volumic Mass") << ")"
0114          << G4endl;
0115 
0116   // total number of process calls
0117   G4double countTot = 0.;
0118   G4cout << "\n Number of process calls --->";
0119   for (size_t i = 0; i < fProcCounter->size(); ++i) {
0120     G4String procName = (*fProcCounter)[i]->GetName();
0121     if (procName != "Transportation") {
0122       G4int count = (*fProcCounter)[i]->GetCounter();
0123       G4cout << "\t" << procName << " : " << count;
0124       countTot += count;
0125     }
0126   }
0127 
0128   // compute totalCrossSection, meanFreePath and massicCrossSection
0129   //
0130   G4double totalCrossSection = countTot / (NbOfEvents * length);
0131   G4double MeanFreePath = 1. / totalCrossSection;
0132   G4double massCrossSection = totalCrossSection / density;
0133 
0134   G4cout.precision(5);
0135   G4cout << "\n Simulation: "
0136          << "total CrossSection = " << totalCrossSection * cm << " /cm"
0137          << "\t MeanFreePath = " << G4BestUnit(MeanFreePath, "Length")
0138          << "\t massicCrossSection = " << massCrossSection * g / cm2 << " cm2/g" << G4endl;
0139 
0140   // compute theoretical predictions
0141   //
0142   if (particle == "mu+" || particle == "mu-") {
0143     totalCrossSection = 0.;
0144     for (size_t i = 0; i < fProcCounter->size(); ++i) {
0145       G4String procName = (*fProcCounter)[i]->GetName();
0146       if (procName != "Transportation") {
0147         totalCrossSection += ComputeTheory(procName, NbOfEvents);
0148         FillCrossSectionHisto(procName, NbOfEvents);
0149       }
0150     }
0151 
0152     MeanFreePath = 1. / totalCrossSection;
0153     massCrossSection = totalCrossSection / density;
0154 
0155     G4cout << " Theory:     "
0156            << "total CrossSection = " << totalCrossSection * cm << " /cm"
0157            << "\t MeanFreePath = " << G4BestUnit(MeanFreePath, "Length")
0158            << "\t massicCrossSection = " << massCrossSection * g / cm2 << " cm2/g" << G4endl;
0159   }
0160 
0161   //  G4cout.setf(mode,std::ios::floatfield);
0162   G4cout.precision(prec);
0163 
0164   // delete and remove all contents in fProcCounter
0165   size_t n = fProcCounter->size();
0166   for (size_t i = 0; i < n; ++i) {
0167     delete (*fProcCounter)[i];
0168   }
0169   delete fProcCounter;
0170 
0171   fHistoManager->Save();
0172 
0173   // show Rndm status
0174   // CLHEP::HepRandom::showEngineStatus();
0175 }
0176 
0177 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0178 
0179 G4double RunAction::ComputeTheory(const G4String& process, G4int NbOfMu)
0180 {
0181   const G4Material* material = fDetector->GetMaterial();
0182   G4double ekin = fPrimary->GetParticleGun()->GetParticleEnergy();
0183   G4double particleMass = fPrimary->GetParticleGun()->GetParticleDefinition()->GetPDGMass();
0184 
0185   G4int id = 0;
0186   G4double cut = 1.e-10 * ekin;
0187   if (process == "muIoni") {
0188     id = 11;
0189     cut = GetEnergyCut(material, 1);
0190   }
0191   else if (process == "muPairProd") {
0192     id = 12;
0193     cut = 2 * (GetEnergyCut(material, 1) + electron_mass_c2);
0194   }
0195   else if (process == "muBrems") {
0196     id = 13;
0197     cut = GetEnergyCut(material, 0);
0198   }
0199   else if (process == "muonNuclear") {
0200     id = 14;
0201     cut = 100 * MeV;
0202   }
0203   else if (process == "muToMuonPairProd") {
0204     id = 18;
0205     cut = 2 * particleMass;
0206   }
0207   if (id == 0) {
0208     return 0.;
0209   }
0210 
0211   G4int nbOfBins = 100;
0212   // G4double binMin = -10.;
0213   G4double binMin = std::log10(cut / ekin);
0214   G4double binMax = 0.;
0215   G4double binWidth = (binMax - binMin) / G4double(nbOfBins);
0216 
0217   // create histo for theoretical crossSections, with same bining as simulation
0218   //
0219   G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
0220 
0221   G4H1* histoTh = 0;
0222   if (fHistoManager->HistoExist(id)) {
0223     histoTh = analysisManager->GetH1(fHistoManager->GetHistoID(id));
0224     nbOfBins = fHistoManager->GetNbins(id);
0225     binMin = fHistoManager->GetVmin(id);
0226     binMax = fHistoManager->GetVmax(id);
0227     binWidth = fHistoManager->GetBinWidth(id);
0228   }
0229 
0230   // compute and plot differential crossSection, as function of energy transfert.
0231   // compute and return integrated crossSection for a given process.
0232   //(note: to compare with simulation, the integrated crossSection is function
0233   //        of the energy cut.)
0234   //
0235   G4double lgeps, etransf, sigmaE, dsigma;
0236   G4double sigmaTot = 0.;
0237   const G4double ln10 = std::log(10.);
0238   G4double length = fDetector->GetSize();
0239 
0240   // G4cout << "MU: " << process << " E= " << ekin
0241   //        <<"  binMin= " << binMin << " binW= " << binWidth << G4endl;
0242 
0243   for (G4int ibin = 0; ibin < nbOfBins; ibin++) {
0244     lgeps = binMin + (ibin + 0.5) * binWidth;
0245     etransf = ekin * std::pow(10., lgeps);
0246     sigmaE = fMucs->CR_Macroscopic(process, material, ekin, etransf);
0247     dsigma = sigmaE * etransf * binWidth * ln10;
0248     if (etransf > cut) sigmaTot += dsigma;
0249     if (histoTh) {
0250       G4double NbProcess = NbOfMu * length * dsigma;
0251       histoTh->fill(lgeps, NbProcess);
0252     }
0253   }
0254 
0255   // return integrated crossSection
0256   //
0257   return sigmaTot;
0258 }
0259 
0260 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0261 
0262 void RunAction::FillCrossSectionHisto(const G4String& process, G4int)
0263 {
0264   const G4Material* material = fDetector->GetMaterial();
0265   G4double ekin = fPrimary->GetParticleGun()->GetParticleEnergy();
0266   G4ParticleDefinition* particle = fPrimary->GetParticleGun()->GetParticleDefinition();
0267   G4double particleMass = particle->GetPDGMass();
0268 
0269   G4EmCalculator emCal;
0270 
0271   G4int id = 0;
0272   G4double cut = 1.e-10 * ekin;
0273   if (process == "muIoni") {
0274     id = 21;
0275     cut = GetEnergyCut(material, 1);
0276   }
0277   else if (process == "muPairProd") {
0278     id = 22;
0279     cut = 2 * (GetEnergyCut(material, 1) + electron_mass_c2);
0280   }
0281   else if (process == "muBrems") {
0282     id = 23;
0283     cut = GetEnergyCut(material, 0);
0284   }
0285   else if (process == "muonNuclear") {
0286     id = 24;
0287     cut = 100 * MeV;
0288   }
0289   else if (process == "muToMuonPairProd") {
0290     id = 28;
0291     cut = 2 * particleMass;
0292   }
0293   if (id == 0) {
0294     return;
0295   }
0296 
0297   G4int nbOfBins = 100;
0298   G4double binMin = cut;
0299   G4double binMax = ekin;
0300   G4double binWidth = (binMax - binMin) / G4double(nbOfBins);
0301 
0302   G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
0303 
0304   G4H1* histoTh = 0;
0305   if (fHistoManager->HistoExist(id)) {
0306     histoTh = analysisManager->GetH1(fHistoManager->GetHistoID(id));
0307     nbOfBins = fHistoManager->GetNbins(id);
0308     binMin = fHistoManager->GetVmin(id);
0309     binMax = fHistoManager->GetVmax(id);
0310     binWidth = fHistoManager->GetBinWidth(id);
0311   }
0312 
0313   G4double sigma, primaryEnergy;
0314 
0315   for (G4int ibin = 0; ibin < nbOfBins; ibin++) {
0316     primaryEnergy = binMin + (ibin + 0.5) * binWidth;
0317     sigma = emCal.GetCrossSectionPerVolume(primaryEnergy, particle, process, material);
0318     if (histoTh) {
0319       histoTh->fill(primaryEnergy, sigma);
0320     }
0321   }
0322 }
0323 
0324 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0325 
0326 G4double RunAction::GetEnergyCut(const G4Material* material, G4int idParticle)
0327 {
0328   G4ProductionCutsTable* table = G4ProductionCutsTable::GetProductionCutsTable();
0329 
0330   size_t index = 0;
0331   while ((table->GetMaterialCutsCouple(index)->GetMaterial() != material)
0332          && (index < table->GetTableSize()))
0333     index++;
0334 
0335   return (*(table->GetEnergyCutsVector(idParticle)))[index];
0336 }
0337 
0338 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......