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