Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 /// \file hadronic/Hadr00/src/HistoManager.cc 27 /// \brief Implementation of the HistoManager class 28 // 29 // 30 //--------------------------------------------------------------------------- 31 // 32 // ClassName: HistoManager 33 // 34 // 35 // Author: V.Ivanchenko 30/01/01 36 // 37 // Modified: 38 // 04.06.2006 Adoptation of hadr01 (V.Ivanchenko) 39 // 40 //---------------------------------------------------------------------------- 41 // 42 43 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 45 46 #include "HistoManager.hh" 47 48 #include "HistoManagerMessenger.hh" 49 50 #include "G4HadronicProcessStore.hh" 51 #include "G4Neutron.hh" 52 #include "G4NistManager.hh" 53 #include "G4NucleiProperties.hh" 54 #include "G4ParticleDefinition.hh" 55 #include "G4ParticleTable.hh" 56 #include "G4StableIsotopes.hh" 57 #include "G4SystemOfUnits.hh" 58 #include "G4UnitsTable.hh" 59 #include "G4ios.hh" 60 #include "globals.hh" 61 62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 63 64 HistoManager::HistoManager() 65 { 66 fAnalysisManager = 0; 67 fHistoName = "hadr00"; 68 69 fNeutron = G4Neutron::Neutron(); 70 fMessenger = new HistoManagerMessenger(this); 71 fVerbose = 1; 72 73 fParticleName = "proton"; 74 fElementName = "Al"; 75 76 fTargetMaterial = 0; 77 78 fMinKinEnergy = 0.1 * MeV; 79 fMaxKinEnergy = 10 * TeV; 80 fMinMomentum = 1 * MeV; 81 fMaxMomentum = 10 * TeV; 82 83 fBinsE = 800; 84 fBinsP = 700; 85 } 86 87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 88 89 HistoManager::~HistoManager() 90 { 91 delete fMessenger; 92 } 93 94 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 95 96 void HistoManager::BeginOfRun() 97 { 98 G4double p1 = std::log10(fMinMomentum / GeV); 99 G4double p2 = std::log10(fMaxMomentum / GeV); 100 G4double e1 = std::log10(fMinKinEnergy / MeV); 101 G4double e2 = std::log10(fMaxKinEnergy / MeV); 102 103 // G4cout<<"e1= "<<e1<<" e2= "<<e2<<" p1= "<<p1<<" p2= "<<p2<<G4endl; 104 105 fAnalysisManager = G4AnalysisManager::Instance(); 106 fAnalysisManager->OpenFile(fHistoName + ".root"); 107 fAnalysisManager->SetFirstHistoId(1); 108 109 fAnalysisManager->CreateH1("h1", "Elastic cross section (barn) as a functions of log10(p/GeV)", 110 fBinsP, p1, p2); 111 fAnalysisManager->CreateH1("h2", "Elastic cross section (barn) as a functions of log10(E/MeV)", 112 fBinsE, e1, e2); 113 fAnalysisManager->CreateH1("h3", "Inelastic cross section (barn) as a functions of log10(p/GeV)", 114 fBinsP, p1, p2); 115 fAnalysisManager->CreateH1("h4", "Inelastic cross section (barn) as a functions of log10(E/MeV)", 116 fBinsE, e1, e2); 117 fAnalysisManager->CreateH1("h5", "Capture cross section (barn) as a functions of log10(E/MeV)", 118 fBinsE, e1, e2); 119 fAnalysisManager->CreateH1("h6", "Fission cross section (barn) as a functions of log10(E/MeV)", 120 fBinsE, e1, e2); 121 fAnalysisManager->CreateH1( 122 "h7", "Charge exchange cross section (barn) as a functions of log10(E/MeV)", fBinsE, e1, e2); 123 fAnalysisManager->CreateH1("h8", "Total cross section (barn) as a functions of log10(E/MeV)", 124 fBinsE, e1, e2); 125 fAnalysisManager->CreateH1( 126 "h9", "Inelastic cross section per volume as a functions of log10(E/MeV)", fBinsE, e1, e2); 127 fAnalysisManager->CreateH1( 128 "h10", "Elastic cross section per volume as a functions of log10(E/MeV)", fBinsE, e1, e2); 129 } 130 131 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 132 133 void HistoManager::EndOfRun() 134 { 135 if (fVerbose > 0) { 136 G4cout << "HistoManager: End of run actions are started" << G4endl; 137 } 138 139 const G4Element* elm = G4NistManager::Instance()->FindOrBuildElement(fElementName); 140 const G4Material* mat = G4NistManager::Instance()->FindOrBuildMaterial("G4_" + fElementName); 141 const G4ParticleDefinition* particle = 142 G4ParticleTable::GetParticleTable()->FindParticle(fParticleName); 143 144 G4cout << "### Fill Cross Sections for " << fParticleName << " off " << fElementName << G4endl; 145 if (fVerbose > 0) { 146 G4cout << "-------------------------------------------------------------" << G4endl; 147 G4cout << " N E(MeV) Elastic(b) Inelastic(b)"; 148 if (particle == fNeutron) { 149 G4cout << " Capture(b) Fission(b)"; 150 } 151 G4cout << " Total(b)" << G4endl; 152 G4cout << "-------------------------------------------------------------" << G4endl; 153 } 154 if (!particle || !elm) { 155 G4cout << "HistoManager WARNING Particle or element undefined" << G4endl; 156 return; 157 } 158 159 G4int prec = G4cout.precision(); 160 G4cout.precision(4); 161 162 G4HadronicProcessStore* store = G4HadronicProcessStore::Instance(); 163 G4double mass = particle->GetPDGMass(); 164 165 // Build histograms 166 167 G4double p1 = std::log10(fMinMomentum / GeV); 168 G4double p2 = std::log10(fMaxMomentum / GeV); 169 G4double e1 = std::log10(fMinKinEnergy / MeV); 170 G4double e2 = std::log10(fMaxKinEnergy / MeV); 171 G4double de = (e2 - e1) / G4double(fBinsE); 172 G4double dp = (p2 - p1) / G4double(fBinsP); 173 174 G4double x = e1 - de * 0.5; 175 G4double e, p, xs, xtot; 176 G4int i; 177 178 G4double coeff = 1.0; 179 if (fTargetMaterial) { 180 coeff = fTargetMaterial->GetDensity() * cm2 / g; 181 } 182 183 for (i = 0; i < fBinsE; i++) { 184 x += de; 185 e = std::pow(10., x) * MeV; 186 if (fVerbose > 0) G4cout << std::setw(5) << i << std::setw(12) << e; 187 xs = store->GetElasticCrossSectionPerAtom(particle, e, elm, mat); 188 xtot = xs; 189 if (fVerbose > 0) G4cout << std::setw(12) << xs / barn; 190 fAnalysisManager->FillH1(2, x, xs / barn); 191 xs = store->GetInelasticCrossSectionPerAtom(particle, e, elm, mat); 192 xtot += xs; 193 if (fVerbose > 0) G4cout << " " << std::setw(12) << xs / barn; 194 fAnalysisManager->FillH1(4, x, xs / barn); 195 if (fTargetMaterial) { 196 xs = store->GetInelasticCrossSectionPerVolume(particle, e, fTargetMaterial); 197 fAnalysisManager->FillH1(9, x, xs / coeff); 198 xs = store->GetElasticCrossSectionPerVolume(particle, e, fTargetMaterial); 199 fAnalysisManager->FillH1(10, x, xs / coeff); 200 } 201 if (particle == fNeutron) { 202 xs = store->GetCaptureCrossSectionPerAtom(particle, e, elm, mat); 203 xtot += xs; 204 if (fVerbose > 0) G4cout << " " << std::setw(12) << xs / barn; 205 fAnalysisManager->FillH1(5, x, xs / barn); 206 xs = store->GetFissionCrossSectionPerAtom(particle, e, elm, mat); 207 xtot += xs; 208 if (fVerbose > 0) G4cout << " " << std::setw(12) << xs / barn; 209 fAnalysisManager->FillH1(6, x, xs / barn); 210 } 211 xs = store->GetChargeExchangeCrossSectionPerAtom(particle, e, elm, mat); 212 if (fVerbose > 0) G4cout << " " << std::setw(12) << xtot / barn << G4endl; 213 fAnalysisManager->FillH1(7, x, xs / barn); 214 fAnalysisManager->FillH1(8, x, xtot / barn); 215 } 216 217 x = p1 - dp * 0.5; 218 for (i = 0; i < fBinsP; i++) { 219 x += dp; 220 p = std::pow(10., x) * GeV; 221 e = std::sqrt(p * p + mass * mass) - mass; 222 xs = store->GetElasticCrossSectionPerAtom(particle, e, elm, mat); 223 fAnalysisManager->FillH1(1, x, xs / barn); 224 xs = store->GetInelasticCrossSectionPerAtom(particle, e, elm, mat); 225 fAnalysisManager->FillH1(3, x, xs / barn); 226 } 227 if (fVerbose > 0) { 228 G4cout << "-------------------------------------------------------------" << G4endl; 229 } 230 G4cout.precision(prec); 231 fAnalysisManager->Write(); 232 fAnalysisManager->CloseFile(); 233 fAnalysisManager->Clear(); 234 235 G4bool extra = true; 236 if (fTargetMaterial && extra) { 237 G4double E = 5 * GeV; 238 G4double cross = store->GetInelasticCrossSectionPerVolume(particle, E, fTargetMaterial); 239 if (cross <= 0.0) { 240 cross = 1.e-100; 241 } 242 G4cout << "### " << particle->GetParticleName() << " " << E / GeV << " GeV on " 243 << fTargetMaterial->GetName() 244 << " xs/X0= " << 1.0 / (cross * fTargetMaterial->GetRadlen()) << G4endl; 245 } 246 } 247 248 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 249 250 void HistoManager::SetVerbose(G4int val) 251 { 252 fVerbose = val; 253 } 254 255 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 256