Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
examples/advanced/nanobeam/src/PhysicsList.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 // Please cite the following paper if you use this software
27 // Nucl.Instrum.Meth.B260:20-27, 2007
28 
29 #include "PhysicsList.hh"
30 #include "G4SystemOfUnits.hh"
31 
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
33 
35 {
37  fCutForGamma = defaultCutValue;
38  fCutForElectron = defaultCutValue;
39  fCutForPositron = defaultCutValue;
40  fCutForProton = defaultCutValue;
41 
42  SetVerboseLevel(1);
43 }
44 
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46 
48 {}
49 
50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
51 
53 {
54  // In this method, static member functions should be called
55  // for all particles which you want to use.
56  // This ensures that objects of these particle types will be
57  // created in the program.
58 
62 }
63 
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65 
67 {
68 
69  // gamma
71 
72  // optical photon
74 }
75  //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76 
78 {
79  // leptons
82 }
83 
84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85 
87 {
88  // barions
92 }
93 
94 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
95 
97 {
99  ConstructEM();
101 }
102 
103 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
104 
105 #include "G4PhotoElectricEffect.hh"
106 #include "G4ComptonScattering.hh"
107 #include "G4GammaConversion.hh"
108 
109 #include "G4eMultipleScattering.hh"
110 #include "G4eIonisation.hh"
111 #include "G4eBremsstrahlung.hh"
112 #include "G4eplusAnnihilation.hh"
113 
114 #include "G4MuMultipleScattering.hh"
115 #include "G4WentzelVIModel.hh"
116 
117 #include "G4MuIonisation.hh"
118 #include "G4MuBremsstrahlung.hh"
119 #include "G4MuPairProduction.hh"
120 #include "G4CoulombScattering.hh"
121 
122 #include "G4hMultipleScattering.hh"
123 #include "G4ionIonisation.hh"
124 #include "G4hIonisation.hh"
125 #include "G4hBremsstrahlung.hh"
126 #include "G4hPairProduction.hh"
127 
128 #include "G4StepLimiter.hh"
129 
130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
131 
132 void PhysicsList::ConstructEM()
133 {
134 
136 
137 // ****************************************************************
138 // Identical to G4EmStandardPhysics but added G4StepLimiter process
139 // ****************************************************************
140 
141  theParticleIterator->reset();
142 
143  while( (*theParticleIterator)() ){
144 
145  G4ParticleDefinition* particle = theParticleIterator->value();
146  G4String particleName = particle->GetParticleName();
147 
148  if (particleName == "gamma") {
149 
150  ph->RegisterProcess(new G4PhotoElectricEffect(), particle);
151  ph->RegisterProcess(new G4ComptonScattering(), particle);
152  ph->RegisterProcess(new G4GammaConversion(), particle);
153 
154  } else if (particleName == "e-") {
155 
156  ph->RegisterProcess(new G4eMultipleScattering(), particle);
157  ph->RegisterProcess(new G4eIonisation(), particle);
158  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
159 
160  } else if (particleName == "e+") {
161 
162  ph->RegisterProcess(new G4eMultipleScattering(), particle);
163  ph->RegisterProcess(new G4eIonisation(), particle);
164  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
165  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
166 
167  } else if( particleName == "mu+" ||
168  particleName == "mu-" ) {
169 
171  msc->AddEmModel(0, new G4WentzelVIModel());
172 
173  ph->RegisterProcess(msc, particle);
174  ph->RegisterProcess(new G4MuIonisation(), particle);
175  ph->RegisterProcess(new G4MuBremsstrahlung(), particle);
176  ph->RegisterProcess(new G4MuPairProduction(), particle);
177  ph->RegisterProcess(new G4CoulombScattering(), particle);
178 
179  } else if (particleName == "alpha" ||
180  particleName == "He3") {
181 
182  ph->RegisterProcess(new G4hMultipleScattering(), particle);
183  ph->RegisterProcess(new G4ionIonisation(), particle);
184 
185  } else if (particleName == "GenericIon") {
186 
187  ph->RegisterProcess(new G4hMultipleScattering(), particle);
188  ph->RegisterProcess(new G4ionIonisation(), particle);
189 
190  } else if (particleName == "proton") {
191  ph->RegisterProcess(new G4hMultipleScattering(), particle);
192  ph->RegisterProcess(new G4hIonisation(), particle);
193  ph->RegisterProcess(new G4hBremsstrahlung(), particle);
194  ph->RegisterProcess(new G4hPairProduction(), particle);
195 
196  ph->RegisterProcess(new G4StepLimiter(), particle);
197 
198  }
199  }
200 }
201 
202 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
203 
205 { }
206 
207 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
208 
210 {
211  if (verboseLevel >0){
212  G4cout << "PhysicsList::SetCuts:";
213  G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl;
214  }
215 
216  // set cut values for gamma at first and for e- second and next for e+,
217  // because some processes for e+/e- need cut values for gamma
218  SetCutValue(fCutForGamma, "gamma");
219  SetCutValue(fCutForElectron, "e-");
220  SetCutValue(fCutForPositron, "e+");
221 
222  // set cut values for proton and anti_proton before all other hadrons
223  // because some processes for hadrons need cut values for proton/anti_proton
224  SetCutValue(fCutForProton, "proton");
225  SetCutValue(fCutForProton, "anti_proton");
226 
227  //SetCutValueForOthers(defaultCutValue);
228 
230 }
231 
232 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
233 
235 {
236  if (verboseLevel >0){
237  G4cout << "PhysicsList::SetCuts:";
238  G4cout << "Gamma cut in energy: " << lowcut*MeV << " (MeV)" << G4endl;
239  }
240 
241  // G4Gamma::SetEnergyRange(lowcut,1e5);
242  SetGELowLimit(lowcut);
243 }
244 
245 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
246 
248 {
249  if (verboseLevel >0){
250 
251  G4cout << "PhysicsList::SetCuts:";
252  G4cout << "Electron cut in energy: " << lowcut*MeV << " (MeV)" << G4endl;
253  }
254 
255  // G4Electron::SetEnergyRange(lowcut,1e5);
256  SetGELowLimit(lowcut);
257 }
258 
259 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
260 
262 {
263  if (verboseLevel >0){
264  G4cout << "PhysicsList::SetGEPLowLimit:";
265  G4cout << "Gamma and Electron cut in energy: " << lowcut*MeV << " (MeV)" << G4endl;
266  }
267 
268  this->SetGELowLimit(lowcut);
269 
270  G4cerr << " SetGEPLowLimit : Uncertain whether setting Positron low limit " << G4endl;
271 }
272 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
273 
274 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
276 {
277  if (verboseLevel >0){
278  G4cout << "PhysicsList::SetGELowLimit:";
279  G4cout << "Gamma and Electron cut in energy: " << lowcut*MeV << " (MeV)" << G4endl;
280  }
281 
283 }
285 {
286  ResetCuts();
287  fCutForGamma = val;
288 }
289 
290 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
291 
293 {
294  // ResetCuts();
295  fCutForElectron = val;
296 }
297 
298 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
299 
301 {
302  // ResetCuts();
303  fCutForPositron = val;
304 }
305 
306 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
307 
309 {
310  //ResetCuts();
311  fCutForProton = val;
312 }
313 
314 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
void SetCutValue(G4double aCut, const G4String &pname)
void SetEnergyRange(G4double lowedge, G4double highedge)
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
static G4AntiProton * AntiProtonDefinition()
Definition: G4AntiProton.cc:88
const G4String & GetParticleName() const
void DumpCutValuesTable(G4int flag=1)
G4GLOB_DLL std::ostream G4cout
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
void SetVerboseLevel(G4int value)
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89
static G4ProductionCutsTable * GetProductionCutsTable()
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
static G4PhysicsListHelper * GetPhysicsListHelper()
void ResetCuts()
obsolete methods
#define G4endl
Definition: G4ios.hh:61
static G4OpticalPhoton * OpticalPhotonDefinition()
double G4double
Definition: G4Types.hh:76
int micrometer
Definition: hepunit.py:34
#define theParticleIterator
G4GLOB_DLL std::ostream G4cerr
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81