LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
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 //
28 //
29 //
30 // $Id: DirectAccess.cc 94302 2015-11-11 12:58:04Z 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"
51 #include "G4MuBremsstrahlungModel.hh"
52 #include "G4MuPairProductionModel.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 #include "G4DataVector.hh"
65 #include "G4NistManager.hh"
66 #include "G4ParticleTable.hh"
67 
68 int main() {
69 
70  G4UnitDefinition::BuildUnitsTable();
71 
72  G4ParticleDefinition* gamma = G4Gamma::Gamma();
73  G4ParticleDefinition* posit = G4Positron::Positron();
74  G4ParticleDefinition* elec = G4Electron::Electron();
75  G4ParticleDefinition* prot = G4Proton::Proton();
76  G4ParticleDefinition* muon = G4MuonPlus::MuonPlus();
77  G4ParticleTable* partTable = G4ParticleTable::GetParticleTable();
78  partTable->SetReadiness();
79 
80  G4DataVector cuts;
81  cuts.push_back(1*keV);
82 
83  // define materials
84  //
85  G4Material* material =
86  G4NistManager::Instance()->FindOrBuildMaterial("G4_Fe");
87 
88  G4cout << *(G4Material::GetMaterialTable()) << G4endl;
89 
90  G4MaterialCutsCouple* couple = new G4MaterialCutsCouple(material);
91  couple->SetIndex(0);
92 
93  // work only for simple materials
94  G4double Z = material->GetZ();
95  G4double A = material->GetA();
96 
97  // initialise gamma processes (models)
98  //
99  G4VEmModel* phot = new G4PEEffectFluoModel();
100  G4VEmModel* comp = new G4KleinNishinaCompton();
101  G4VEmModel* conv = new G4BetheHeitlerModel();
102  phot->Initialise(gamma, cuts);
103  comp->Initialise(gamma, cuts);
104  conv->Initialise(gamma, cuts);
105 
106  // valid pointer to a couple is needed for this model
107  phot->SetCurrentCouple(couple);
108 
109  // compute CrossSection per atom and MeanFreePath
110  //
111  G4double Emin = 1.01*MeV, Emax = 2.01*MeV, dE = 100*keV;
112 
113  G4cout << "\n #### Gamma : CrossSectionPerAtom and MeanFreePath for "
114  << material->GetName() << G4endl;
115  G4cout << "\n Energy \t PhotoElec \t Compton \t Conversion \t";
116  G4cout << "\t PhotoElec \t Compton \t Conversion" << G4endl;
117 
118  for (G4double Energy = Emin; Energy <= Emax; Energy += dE) {
119  G4cout << "\n " << G4BestUnit (Energy, "Energy")
120  << "\t"
121  << G4BestUnit (phot->ComputeCrossSectionPerAtom(gamma,Energy,Z),"Surface")
122  << "\t"
123  << G4BestUnit (comp->ComputeCrossSectionPerAtom(gamma,Energy,Z),"Surface")
124  << "\t"
125  << G4BestUnit (conv->ComputeCrossSectionPerAtom(gamma,Energy,Z),"Surface")
126  << "\t \t"
127  << G4BestUnit (phot->ComputeMeanFreePath(gamma,Energy,material),"Length")
128  << "\t"
129  << G4BestUnit (comp->ComputeMeanFreePath(gamma,Energy,material),"Length")
130  << "\t"
131  << G4BestUnit (conv->ComputeMeanFreePath(gamma,Energy,material),"Length");
132  }
133 
134  G4cout << G4endl;
135 
136  // initialise positron annihilation (model)
137  //
138  G4VEmModel* anni = new G4eeToTwoGammaModel();
139  anni->Initialise(posit, cuts);
140 
141  // compute CrossSection per atom and MeanFreePath
142  //
143  Emin = 1.01*MeV; Emax = 2.01*MeV; dE = 100*keV;
144 
145  G4cout << "\n #### e+ annihilation : CrossSectionPerAtom and MeanFreePath"
146  << " for " << material->GetName() << G4endl;
147  G4cout << "\n Energy \t e+ annihil \t";
148  G4cout << "\t e+ annihil" << G4endl;
149 
150  for (G4double Energy = Emin; Energy <= Emax; Energy += dE) {
151  G4cout << "\n " << G4BestUnit (Energy, "Energy")
152  << "\t"
153  << G4BestUnit (anni->ComputeCrossSectionPerAtom(posit,Energy,Z),"Surface")
154  << "\t \t"
155  << G4BestUnit (anni->ComputeMeanFreePath(posit,Energy,material),"Length");
156  }
157 
158  G4cout << G4endl;
159 
160  // initialise electron processes (models)
161  //
162  G4VEmModel* ioni = new G4MollerBhabhaModel();
163  G4VEmModel* brem = new G4SeltzerBergerModel();
164  ioni->Initialise(elec, cuts);
165  brem->Initialise(elec, cuts);
166 
167  // compute CrossSection per atom and MeanFreePath
168  //
169  Emin = 1.01*MeV; Emax = 101.01*MeV; dE = 10*MeV;
170  G4double Ecut = 100*keV;
171 
172  G4cout << "\n ####electron: CrossSection, MeanFreePath and StoppingPower"
173  << " for " << material->GetName()
174  << ";\tEnergy cut = " << G4BestUnit (Ecut, "Energy") << G4endl;
175 
176  G4cout << "\n Energy \t ionization \t bremsstra \t";
177  G4cout << "\t ionization \t bremsstra \t";
178  G4cout << "\t ionization \t bremsstra" << G4endl;
179 
180  for (G4double Energy = Emin; Energy <= Emax; Energy += dE) {
181  G4cout << "\n " << G4BestUnit (Energy, "Energy")
182  << "\t"
183  << G4BestUnit (ioni->ComputeCrossSectionPerAtom(elec,Energy,Z,A,Ecut),
184  "Surface")
185  << "\t"
186  << G4BestUnit (brem->ComputeCrossSectionPerAtom(elec,Energy,Z,A,Ecut),
187  "Surface")
188  << "\t \t"
189  << G4BestUnit (ioni->ComputeMeanFreePath(elec,Energy,material,Ecut),
190  "Length")
191  << "\t"
192  << G4BestUnit (brem->ComputeMeanFreePath(elec,Energy,material,Ecut),
193  "Length")
194  << "\t \t"
195  << G4BestUnit (ioni->ComputeDEDXPerVolume(material,elec,Energy,Ecut),
196  "Energy/Length")
197  << "\t"
198  << G4BestUnit (brem->ComputeDEDXPerVolume(material,elec,Energy,Ecut),
199  "Energy/Length");
200  }
201 
202  G4cout << G4endl;
203 
204  // initialise proton processes (models)
205  //
206  ioni = new G4BetheBlochModel();
207  ioni->Initialise(prot, cuts);
208 
209  // compute CrossSection per atom and MeanFreePath
210  //
211  Emin = 1.01*MeV; Emax = 102.01*MeV; dE = 10*MeV;
212  Ecut = 100*keV;
213 
214  G4cout << "\n #### proton : CrossSection, MeanFreePath and StoppingPower"
215  << " for " << material->GetName()
216  << ";\tEnergy cut = " << G4BestUnit (Ecut, "Energy") << G4endl;
217 
218  G4cout << "\n Energy \t ionization \t";
219  G4cout << "\t ionization \t";
220  G4cout << "\t ionization" << G4endl;
221 
222  for (G4double Energy = Emin; Energy <= Emax; Energy += dE) {
223  G4cout << "\n " << G4BestUnit (Energy, "Energy")
224  << "\t"
225  << G4BestUnit (ioni->ComputeCrossSectionPerAtom(prot,Energy,Z,A,Ecut),
226  "Surface")
227  << "\t \t"
228  << G4BestUnit (ioni->ComputeMeanFreePath(prot,Energy,material,Ecut),
229  "Length")
230  << "\t \t"
231  << G4BestUnit (ioni->ComputeDEDXPerVolume(material,prot,Energy,Ecut),
232  "Energy/Length");
233  }
234 
235  G4cout << G4endl;
236 
237  // low energy : Bragg Model
238  ioni = new G4BraggModel(prot);
239  ioni->Initialise(prot, cuts);
240 
241  // compute CrossSection per atom and MeanFreePath
242  //
243  Emin = 1.1*keV; Emax = 2.01*MeV; dE = 300*keV;
244  Ecut = 10*keV;
245 
246  G4cout << "\n #### proton : low energy model (Bragg) "
247  << ";\tEnergy cut = " << G4BestUnit (Ecut, "Energy") << G4endl;
248 
249  G4cout << "\n Energy \t ionization \t";
250  G4cout << "\t ionization \t";
251  G4cout << "\t ionization" << G4endl;
252 
253  for (G4double Energy = Emin; Energy <= Emax; Energy += dE) {
254  G4cout << "\n " << G4BestUnit (Energy, "Energy")
255  << "\t"
256  << G4BestUnit (ioni->ComputeCrossSectionPerAtom(prot,Energy,Z,A,Ecut),
257  "Surface")
258  << "\t \t"
259  << G4BestUnit (ioni->ComputeMeanFreePath(prot,Energy,material,Ecut),
260  "Length")
261  << "\t \t"
262  << G4BestUnit (ioni->ComputeDEDXPerVolume(material,prot,Energy,Ecut),
263  "Energy/Length");
264  }
265 
266  G4cout << G4endl;
267 
268  // initialise muon processes (models)
269  //
270  ioni = new G4MuBetheBlochModel();
271  brem = new G4MuBremsstrahlungModel();
272  G4VEmModel* pair = new G4MuPairProductionModel();
273  ioni->Initialise(muon, cuts);
274  brem->Initialise(muon, cuts);
275  pair->Initialise(muon, cuts);
276 
277  // compute CrossSection per atom and MeanFreePath
278  //
279  Emin = 1.01*GeV; Emax = 101.01*GeV; dE = 10*GeV;
280  Ecut = 10*MeV;
281 
282  G4cout << "\n ####muon: CrossSection and MeanFreePath for "
283  << material->GetName()
284  << ";\tEnergy cut = " << G4BestUnit (Ecut, "Energy") << G4endl;
285 
286  G4cout << "\n Energy \t ionization \t bremsstra \t pair_prod \t";
287  G4cout << "\t ionization \t bremsstra \t pair_prod" << G4endl;
288 
289  for (G4double Energy = Emin; Energy <= Emax; Energy += dE) {
290  G4cout << "\n " << G4BestUnit (Energy, "Energy")
291  << "\t"
292  << G4BestUnit (ioni->ComputeCrossSectionPerAtom(muon,Energy,Z,A,Ecut),
293  "Surface")
294  << "\t"
295  << G4BestUnit (brem->ComputeCrossSectionPerAtom(muon,Energy,Z,A,Ecut),
296  "Surface")
297  << "\t"
298  << G4BestUnit (pair->ComputeCrossSectionPerAtom(muon,Energy,Z,A,Ecut),
299  "Surface")
300  << "\t \t"
301  << G4BestUnit (ioni->ComputeMeanFreePath(muon,Energy,material,Ecut),
302  "Length")
303  << "\t"
304  << G4BestUnit (brem->ComputeMeanFreePath(muon,Energy,material,Ecut),
305  "Length")
306  << "\t"
307  << G4BestUnit (pair->ComputeMeanFreePath(muon,Energy,material,Ecut),
308  "Length");
309  }
310 
311  G4cout << G4endl;
312 
313  G4cout << "\n ####muon: StoppingPower for "
314  << material->GetName()
315  << ";\tEnergy cut = " << G4BestUnit (Ecut, "Energy") << G4endl;
316 
317  G4cout << "\n Energy \t ionization \t bremsstra \t pair_prod \t" << G4endl;
318 
319  for (G4double Energy = Emin; Energy <= Emax; Energy += dE) {
320  G4cout << "\n " << G4BestUnit (Energy, "Energy")
321  << "\t"
322  << G4BestUnit (ioni->ComputeDEDXPerVolume(material,muon,Energy,Ecut),
323  "Energy/Length")
324  << "\t"
325  << G4BestUnit (brem->ComputeDEDXPerVolume(material,muon,Energy,Ecut),
326  "Energy/Length")
327  << "\t"
328  << G4BestUnit (pair->ComputeDEDXPerVolume(material,muon,Energy,Ecut),
329  "Energy/Length");
330  }
331 
332  G4cout << G4endl;
333  return EXIT_SUCCESS;
334 }
335 
336 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
int main()
Definition: DirectAccess.cc:68
Float_t Z
Definition: plot.C:39
General LArSoft Utilities.