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 electromagnetic/TestEm0/DirectAccess.cc 27 /// \brief Main program of the electromagnetic/TestEm0 example 28 // 29 // 30 // 31 // ------------------------------------------------------------ 32 // 33 // To print cross sections per atom and mean free path for simple material 34 // 35 #include "G4BetheBlochModel.hh" 36 #include "G4BetheHeitlerModel.hh" 37 #include "G4BraggModel.hh" 38 #include "G4DataVector.hh" 39 #include "G4Electron.hh" 40 #include "G4Gamma.hh" 41 #include "G4KleinNishinaCompton.hh" 42 #include "G4Material.hh" 43 #include "G4MollerBhabhaModel.hh" 44 #include "G4MuBetheBlochModel.hh" 45 #include "G4MuBremsstrahlungModel.hh" 46 #include "G4MuPairProductionModel.hh" 47 #include "G4MuonPlus.hh" 48 #include "G4NistManager.hh" 49 #include "G4PEEffectFluoModel.hh" 50 #include "G4ParticleTable.hh" 51 #include "G4Positron.hh" 52 #include "G4Proton.hh" 53 #include "G4SeltzerBergerModel.hh" 54 #include "G4SystemOfUnits.hh" 55 #include "G4UnitsTable.hh" 56 #include "G4eeToTwoGammaModel.hh" 57 #include "globals.hh" 58 59 int main() 60 { 61 G4UnitDefinition::BuildUnitsTable(); 62 63 G4ParticleDefinition* gamma = G4Gamma::Gamma(); 64 G4ParticleDefinition* posit = G4Positron::Positron(); 65 G4ParticleDefinition* elec = G4Electron::Electron(); 66 G4ParticleDefinition* prot = G4Proton::Proton(); 67 G4ParticleDefinition* muon = G4MuonPlus::MuonPlus(); 68 G4ParticleTable* partTable = G4ParticleTable::GetParticleTable(); 69 partTable->SetReadiness(); 70 71 G4DataVector cuts; 72 cuts.push_back(1 * keV); 73 74 // define materials 75 // 76 G4Material* material = G4NistManager::Instance()->FindOrBuildMaterial("G4_Fe"); 77 78 G4cout << *(G4Material::GetMaterialTable()) << G4endl; 79 80 G4MaterialCutsCouple* couple = new G4MaterialCutsCouple(material); 81 couple->SetIndex(0); 82 83 // work only for simple materials 84 G4double Z = material->GetZ(); 85 G4double A = material->GetA(); 86 87 // initialise gamma processes (models) 88 // 89 G4VEmModel* phot = new G4PEEffectFluoModel(); 90 G4VEmModel* comp = new G4KleinNishinaCompton(); 91 G4VEmModel* conv = new G4BetheHeitlerModel(); 92 phot->Initialise(gamma, cuts); 93 comp->Initialise(gamma, cuts); 94 conv->Initialise(gamma, cuts); 95 96 // valid pointer to a couple is needed for this model 97 phot->SetCurrentCouple(couple); 98 99 // compute CrossSection per atom and MeanFreePath 100 // 101 G4double Emin = 1.01 * MeV, Emax = 2.01 * MeV, dE = 100 * keV; 102 103 G4cout << "\n #### Gamma : CrossSectionPerAtom and MeanFreePath for " << material->GetName() 104 << G4endl; 105 G4cout << "\n Energy \t PhotoElec \t Compton \t Conversion \t"; 106 G4cout << "\t PhotoElec \t Compton \t Conversion" << G4endl; 107 108 for (G4double Energy = Emin; Energy <= Emax; Energy += dE) { 109 G4cout << "\n " << G4BestUnit(Energy, "Energy") << "\t" 110 << G4BestUnit(phot->ComputeCrossSectionPerAtom(gamma, Energy, Z), "Surface") << "\t" 111 << G4BestUnit(comp->ComputeCrossSectionPerAtom(gamma, Energy, Z), "Surface") << "\t" 112 << G4BestUnit(conv->ComputeCrossSectionPerAtom(gamma, Energy, Z), "Surface") << "\t \t" 113 << G4BestUnit(phot->ComputeMeanFreePath(gamma, Energy, material), "Length") << "\t" 114 << G4BestUnit(comp->ComputeMeanFreePath(gamma, Energy, material), "Length") << "\t" 115 << G4BestUnit(conv->ComputeMeanFreePath(gamma, Energy, material), "Length"); 116 } 117 118 G4cout << G4endl; 119 120 // initialise positron annihilation (model) 121 // 122 G4VEmModel* anni = new G4eeToTwoGammaModel(); 123 anni->Initialise(posit, cuts); 124 125 // compute CrossSection per atom and MeanFreePath 126 // 127 Emin = 1.01 * MeV; 128 Emax = 2.01 * MeV; 129 dE = 100 * keV; 130 131 G4cout << "\n #### e+ annihilation : CrossSectionPerAtom and MeanFreePath" 132 << " for " << material->GetName() << G4endl; 133 G4cout << "\n Energy \t e+ annihil \t"; 134 G4cout << "\t e+ annihil" << G4endl; 135 136 for (G4double Energy = Emin; Energy <= Emax; Energy += dE) { 137 G4cout << "\n " << G4BestUnit(Energy, "Energy") << "\t" 138 << G4BestUnit(anni->ComputeCrossSectionPerAtom(posit, Energy, Z), "Surface") << "\t \t" 139 << G4BestUnit(anni->ComputeMeanFreePath(posit, Energy, material), "Length"); 140 } 141 142 G4cout << G4endl; 143 144 // initialise electron processes (models) 145 // 146 G4VEmModel* ioni = new G4MollerBhabhaModel(); 147 G4VEmModel* brem = new G4SeltzerBergerModel(); 148 ioni->Initialise(elec, cuts); 149 brem->Initialise(elec, cuts); 150 151 // compute CrossSection per atom and MeanFreePath 152 // 153 Emin = 1.01 * MeV; 154 Emax = 101.01 * MeV; 155 dE = 10 * MeV; 156 G4double Ecut = 100 * keV; 157 158 G4cout << "\n ####electron: CrossSection, MeanFreePath and StoppingPower" 159 << " for " << material->GetName() << ";\tEnergy cut = " << G4BestUnit(Ecut, "Energy") 160 << G4endl; 161 162 G4cout << "\n Energy \t ionization \t bremsstra \t"; 163 G4cout << "\t ionization \t bremsstra \t"; 164 G4cout << "\t ionization \t bremsstra" << G4endl; 165 166 for (G4double Energy = Emin; Energy <= Emax; Energy += dE) { 167 G4cout << "\n " << G4BestUnit(Energy, "Energy") << "\t" 168 << G4BestUnit(ioni->ComputeCrossSectionPerAtom(elec, Energy, Z, A, Ecut), "Surface") 169 << "\t" 170 << G4BestUnit(brem->ComputeCrossSectionPerAtom(elec, Energy, Z, A, Ecut), "Surface") 171 << "\t \t" 172 << G4BestUnit(ioni->ComputeMeanFreePath(elec, Energy, material, Ecut), "Length") << "\t" 173 << G4BestUnit(brem->ComputeMeanFreePath(elec, Energy, material, Ecut), "Length") 174 << "\t \t" 175 << G4BestUnit(ioni->ComputeDEDXPerVolume(material, elec, Energy, Ecut), "Energy/Length") 176 << "\t" 177 << G4BestUnit(brem->ComputeDEDXPerVolume(material, elec, Energy, Ecut), "Energy/Length"); 178 } 179 180 G4cout << G4endl; 181 182 // initialise proton processes (models) 183 // 184 ioni = new G4BetheBlochModel(); 185 ioni->Initialise(prot, cuts); 186 187 // compute CrossSection per atom and MeanFreePath 188 // 189 Emin = 1.01 * MeV; 190 Emax = 102.01 * MeV; 191 dE = 10 * MeV; 192 Ecut = 100 * keV; 193 194 G4cout << "\n #### proton : CrossSection, MeanFreePath and StoppingPower" 195 << " for " << material->GetName() << ";\tEnergy cut = " << G4BestUnit(Ecut, "Energy") 196 << G4endl; 197 198 G4cout << "\n Energy \t ionization \t"; 199 G4cout << "\t ionization \t"; 200 G4cout << "\t ionization" << G4endl; 201 202 for (G4double Energy = Emin; Energy <= Emax; Energy += dE) { 203 G4cout << "\n " << G4BestUnit(Energy, "Energy") << "\t" 204 << G4BestUnit(ioni->ComputeCrossSectionPerAtom(prot, Energy, Z, A, Ecut), "Surface") 205 << "\t \t" 206 << G4BestUnit(ioni->ComputeMeanFreePath(prot, Energy, material, Ecut), "Length") 207 << "\t \t" 208 << G4BestUnit(ioni->ComputeDEDXPerVolume(material, prot, Energy, Ecut), "Energy/Length"); 209 } 210 211 G4cout << G4endl; 212 213 // low energy : Bragg Model 214 ioni = new G4BraggModel(prot); 215 ioni->Initialise(prot, cuts); 216 217 // compute CrossSection per atom and MeanFreePath 218 // 219 Emin = 1.1 * keV; 220 Emax = 2.01 * MeV; 221 dE = 300 * keV; 222 Ecut = 10 * keV; 223 224 G4cout << "\n #### proton : low energy model (Bragg) " 225 << ";\tEnergy cut = " << G4BestUnit(Ecut, "Energy") << G4endl; 226 227 G4cout << "\n Energy \t ionization \t"; 228 G4cout << "\t ionization \t"; 229 G4cout << "\t ionization" << G4endl; 230 231 for (G4double Energy = Emin; Energy <= Emax; Energy += dE) { 232 G4cout << "\n " << G4BestUnit(Energy, "Energy") << "\t" 233 << G4BestUnit(ioni->ComputeCrossSectionPerAtom(prot, Energy, Z, A, Ecut), "Surface") 234 << "\t \t" 235 << G4BestUnit(ioni->ComputeMeanFreePath(prot, Energy, material, Ecut), "Length") 236 << "\t \t" 237 << G4BestUnit(ioni->ComputeDEDXPerVolume(material, prot, Energy, Ecut), "Energy/Length"); 238 } 239 240 G4cout << G4endl; 241 242 // initialise muon processes (models) 243 // 244 ioni = new G4MuBetheBlochModel(); 245 brem = new G4MuBremsstrahlungModel(); 246 G4VEmModel* pair = new G4MuPairProductionModel(); 247 ioni->Initialise(muon, cuts); 248 brem->Initialise(muon, cuts); 249 pair->Initialise(muon, cuts); 250 251 // compute CrossSection per atom and MeanFreePath 252 // 253 Emin = 1.01 * GeV; 254 Emax = 101.01 * GeV; 255 dE = 10 * GeV; 256 Ecut = 10 * MeV; 257 258 G4cout << "\n ####muon: CrossSection and MeanFreePath for " << material->GetName() 259 << ";\tEnergy cut = " << G4BestUnit(Ecut, "Energy") << G4endl; 260 261 G4cout << "\n Energy \t ionization \t bremsstra \t pair_prod \t"; 262 G4cout << "\t ionization \t bremsstra \t pair_prod" << G4endl; 263 264 for (G4double Energy = Emin; Energy <= Emax; Energy += dE) { 265 G4cout << "\n " << G4BestUnit(Energy, "Energy") << "\t" 266 << G4BestUnit(ioni->ComputeCrossSectionPerAtom(muon, Energy, Z, A, Ecut), "Surface") 267 << "\t" 268 << G4BestUnit(brem->ComputeCrossSectionPerAtom(muon, Energy, Z, A, Ecut), "Surface") 269 << "\t" 270 << G4BestUnit(pair->ComputeCrossSectionPerAtom(muon, Energy, Z, A, Ecut), "Surface") 271 << "\t \t" 272 << G4BestUnit(ioni->ComputeMeanFreePath(muon, Energy, material, Ecut), "Length") << "\t" 273 << G4BestUnit(brem->ComputeMeanFreePath(muon, Energy, material, Ecut), "Length") << "\t" 274 << G4BestUnit(pair->ComputeMeanFreePath(muon, Energy, material, Ecut), "Length"); 275 } 276 277 G4cout << G4endl; 278 279 G4cout << "\n ####muon: StoppingPower for " << material->GetName() 280 << ";\tEnergy cut = " << G4BestUnit(Ecut, "Energy") << G4endl; 281 282 G4cout << "\n Energy \t ionization \t bremsstra \t pair_prod \t" << G4endl; 283 284 for (G4double Energy = Emin; Energy <= Emax; Energy += dE) { 285 G4cout << "\n " << G4BestUnit(Energy, "Energy") << "\t" 286 << G4BestUnit(ioni->ComputeDEDXPerVolume(material, muon, Energy, Ecut), "Energy/Length") 287 << "\t" 288 << G4BestUnit(brem->ComputeDEDXPerVolume(material, muon, Energy, Ecut), "Energy/Length") 289 << "\t" 290 << G4BestUnit(pair->ComputeDEDXPerVolume(material, muon, Energy, Ecut), "Energy/Length"); 291 } 292 293 G4cout << G4endl; 294 return EXIT_SUCCESS; 295 } 296 297 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 298