Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-04-04 08:05:04

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 hadronic/Hadr00/src/HistoManager.cc
0027 /// \brief Implementation of the HistoManager class
0028 //
0029 //
0030 //---------------------------------------------------------------------------
0031 //
0032 // ClassName:   HistoManager
0033 //
0034 //
0035 // Author:      V.Ivanchenko 30/01/01
0036 //
0037 // Modified:
0038 // 04.06.2006 Adoptation of hadr01 (V.Ivanchenko)
0039 //
0040 //----------------------------------------------------------------------------
0041 //
0042 
0043 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0044 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0045 
0046 #include "HistoManager.hh"
0047 
0048 #include "HistoManagerMessenger.hh"
0049 
0050 #include "G4HadronicProcessStore.hh"
0051 #include "G4Neutron.hh"
0052 #include "G4NistManager.hh"
0053 #include "G4NucleiProperties.hh"
0054 #include "G4ParticleDefinition.hh"
0055 #include "G4ParticleTable.hh"
0056 #include "G4StableIsotopes.hh"
0057 #include "G4SystemOfUnits.hh"
0058 #include "G4UnitsTable.hh"
0059 #include "G4ios.hh"
0060 #include "globals.hh"
0061 
0062 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0063 
0064 HistoManager::HistoManager()
0065 {
0066   fAnalysisManager = 0;
0067   fHistoName = "hadr00";
0068 
0069   fNeutron = G4Neutron::Neutron();
0070   fMessenger = new HistoManagerMessenger(this);
0071   fVerbose = 1;
0072 
0073   fParticleName = "proton";
0074   fElementName = "Al";
0075 
0076   fTargetMaterial = 0;
0077 
0078   fMinKinEnergy = 0.1 * MeV;
0079   fMaxKinEnergy = 10 * TeV;
0080   fMinMomentum = 1 * MeV;
0081   fMaxMomentum = 10 * TeV;
0082 
0083   fBinsE = 800;
0084   fBinsP = 700;
0085 }
0086 
0087 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0088 
0089 HistoManager::~HistoManager()
0090 {
0091   delete fMessenger;
0092 }
0093 
0094 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0095 
0096 void HistoManager::BeginOfRun()
0097 {
0098   G4double p1 = std::log10(fMinMomentum / GeV);
0099   G4double p2 = std::log10(fMaxMomentum / GeV);
0100   G4double e1 = std::log10(fMinKinEnergy / MeV);
0101   G4double e2 = std::log10(fMaxKinEnergy / MeV);
0102 
0103   // G4cout<<"e1= "<<e1<<" e2= "<<e2<<" p1= "<<p1<<" p2= "<<p2<<G4endl;
0104 
0105   fAnalysisManager = G4AnalysisManager::Instance();
0106   fAnalysisManager->OpenFile(fHistoName + ".root");
0107   fAnalysisManager->SetFirstHistoId(1);
0108 
0109   fAnalysisManager->CreateH1("h1", "Elastic cross section (barn) as a functions of log10(p/GeV)",
0110                              fBinsP, p1, p2);
0111   fAnalysisManager->CreateH1("h2", "Elastic cross section (barn) as a functions of log10(E/MeV)",
0112                              fBinsE, e1, e2);
0113   fAnalysisManager->CreateH1("h3", "Inelastic cross section (barn) as a functions of log10(p/GeV)",
0114                              fBinsP, p1, p2);
0115   fAnalysisManager->CreateH1("h4", "Inelastic cross section (barn) as a functions of log10(E/MeV)",
0116                              fBinsE, e1, e2);
0117   fAnalysisManager->CreateH1("h5", "Capture cross section (barn) as a functions of log10(E/MeV)",
0118                              fBinsE, e1, e2);
0119   fAnalysisManager->CreateH1("h6", "Fission cross section (barn) as a functions of log10(E/MeV)",
0120                              fBinsE, e1, e2);
0121   fAnalysisManager->CreateH1(
0122     "h7", "Charge exchange cross section (barn) as a functions of log10(E/MeV)", fBinsE, e1, e2);
0123   fAnalysisManager->CreateH1("h8", "Total cross section (barn) as a functions of log10(E/MeV)",
0124                              fBinsE, e1, e2);
0125   fAnalysisManager->CreateH1(
0126     "h9", "Inelastic cross section per volume as a functions of log10(E/MeV)", fBinsE, e1, e2);
0127   fAnalysisManager->CreateH1(
0128     "h10", "Elastic cross section per volume as a functions of log10(E/MeV)", fBinsE, e1, e2);
0129 }
0130 
0131 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0132 
0133 void HistoManager::EndOfRun()
0134 {
0135   if (fVerbose > 0) {
0136     G4cout << "HistoManager: End of run actions are started" << G4endl;
0137   }
0138 
0139   const G4Element* elm = G4NistManager::Instance()->FindOrBuildElement(fElementName);
0140   const G4Material* mat = G4NistManager::Instance()->FindOrBuildMaterial("G4_" + fElementName);
0141   const G4ParticleDefinition* particle =
0142     G4ParticleTable::GetParticleTable()->FindParticle(fParticleName);
0143 
0144   G4cout << "### Fill Cross Sections for " << fParticleName << " off " << fElementName << G4endl;
0145   if (fVerbose > 0) {
0146     G4cout << "-------------------------------------------------------------" << G4endl;
0147     G4cout << "    N     E(MeV)   Elastic(b)   Inelastic(b)";
0148     if (particle == fNeutron) {
0149       G4cout << " Capture(b)   Fission(b)";
0150     }
0151     G4cout << "   Total(b)" << G4endl;
0152     G4cout << "-------------------------------------------------------------" << G4endl;
0153   }
0154   if (!particle || !elm) {
0155     G4cout << "HistoManager WARNING Particle or element undefined" << G4endl;
0156     return;
0157   }
0158 
0159   G4int prec = G4cout.precision();
0160   G4cout.precision(4);
0161 
0162   G4HadronicProcessStore* store = G4HadronicProcessStore::Instance();
0163   G4double mass = particle->GetPDGMass();
0164 
0165   // Build histograms
0166 
0167   G4double p1 = std::log10(fMinMomentum / GeV);
0168   G4double p2 = std::log10(fMaxMomentum / GeV);
0169   G4double e1 = std::log10(fMinKinEnergy / MeV);
0170   G4double e2 = std::log10(fMaxKinEnergy / MeV);
0171   G4double de = (e2 - e1) / G4double(fBinsE);
0172   G4double dp = (p2 - p1) / G4double(fBinsP);
0173 
0174   G4double x = e1 - de * 0.5;
0175   G4double e, p, xs, xtot;
0176   G4int i;
0177 
0178   G4double coeff = 1.0;
0179   if (fTargetMaterial) {
0180     coeff = fTargetMaterial->GetDensity() * cm2 / g;
0181   }
0182 
0183   for (i = 0; i < fBinsE; i++) {
0184     x += de;
0185     e = std::pow(10., x) * MeV;
0186     if (fVerbose > 0) G4cout << std::setw(5) << i << std::setw(12) << e;
0187     xs = store->GetElasticCrossSectionPerAtom(particle, e, elm, mat);
0188     xtot = xs;
0189     if (fVerbose > 0) G4cout << std::setw(12) << xs / barn;
0190     fAnalysisManager->FillH1(2, x, xs / barn);
0191     xs = store->GetInelasticCrossSectionPerAtom(particle, e, elm, mat);
0192     xtot += xs;
0193     if (fVerbose > 0) G4cout << " " << std::setw(12) << xs / barn;
0194     fAnalysisManager->FillH1(4, x, xs / barn);
0195     if (fTargetMaterial) {
0196       xs = store->GetInelasticCrossSectionPerVolume(particle, e, fTargetMaterial);
0197       fAnalysisManager->FillH1(9, x, xs / coeff);
0198       xs = store->GetElasticCrossSectionPerVolume(particle, e, fTargetMaterial);
0199       fAnalysisManager->FillH1(10, x, xs / coeff);
0200     }
0201     if (particle == fNeutron) {
0202       xs = store->GetCaptureCrossSectionPerAtom(particle, e, elm, mat);
0203       xtot += xs;
0204       if (fVerbose > 0) G4cout << " " << std::setw(12) << xs / barn;
0205       fAnalysisManager->FillH1(5, x, xs / barn);
0206       xs = store->GetFissionCrossSectionPerAtom(particle, e, elm, mat);
0207       xtot += xs;
0208       if (fVerbose > 0) G4cout << " " << std::setw(12) << xs / barn;
0209       fAnalysisManager->FillH1(6, x, xs / barn);
0210     }
0211     xs = store->GetChargeExchangeCrossSectionPerAtom(particle, e, elm, mat);
0212     if (fVerbose > 0) G4cout << " " << std::setw(12) << xtot / barn << G4endl;
0213     fAnalysisManager->FillH1(7, x, xs / barn);
0214     fAnalysisManager->FillH1(8, x, xtot / barn);
0215   }
0216 
0217   x = p1 - dp * 0.5;
0218   for (i = 0; i < fBinsP; i++) {
0219     x += dp;
0220     p = std::pow(10., x) * GeV;
0221     e = std::sqrt(p * p + mass * mass) - mass;
0222     xs = store->GetElasticCrossSectionPerAtom(particle, e, elm, mat);
0223     fAnalysisManager->FillH1(1, x, xs / barn);
0224     xs = store->GetInelasticCrossSectionPerAtom(particle, e, elm, mat);
0225     fAnalysisManager->FillH1(3, x, xs / barn);
0226   }
0227   if (fVerbose > 0) {
0228     G4cout << "-------------------------------------------------------------" << G4endl;
0229   }
0230   G4cout.precision(prec);
0231   fAnalysisManager->Write();
0232   fAnalysisManager->CloseFile();
0233   fAnalysisManager->Clear();
0234 
0235   G4bool extra = true;
0236   if (fTargetMaterial && extra) {
0237     G4double E = 5 * GeV;
0238     G4double cross = store->GetInelasticCrossSectionPerVolume(particle, E, fTargetMaterial);
0239     if (cross <= 0.0) {
0240       cross = 1.e-100;
0241     }
0242     G4cout << "### " << particle->GetParticleName() << " " << E / GeV << " GeV on "
0243            << fTargetMaterial->GetName()
0244            << " xs/X0= " << 1.0 / (cross * fTargetMaterial->GetRadlen()) << G4endl;
0245   }
0246 }
0247 
0248 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0249 
0250 void HistoManager::SetVerbose(G4int val)
0251 {
0252   fVerbose = val;
0253 }
0254 
0255 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......