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