Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
DicomPhysicsList.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 // $Id$
30 
31 #include "G4ParticleDefinition.hh"
32 #include "G4ProcessManager.hh"
33 #include "G4ParticleTypes.hh"
34 #include "G4StepLimiter.hh"
35 #include "G4BaryonConstructor.hh"
36 #include "G4IonConstructor.hh"
37 #include "G4MesonConstructor.hh"
38 #include "G4SystemOfUnits.hh"
39 
40 #include "DicomPhysicsList.hh"
41 
42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
44 {
46  fCutForGamma = defaultCutValue;
47  fCutForElectron = defaultCutValue;
48  fCutForPositron = defaultCutValue;
49  SetVerboseLevel(1);
50 }
51 
52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54 {}
55 
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
58 {
62 }
63 
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66 {
67  // gamma
69 
70  // optical photon
72 }
73 
74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76 {
77  // leptons
80 }
81 
82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
84 {
85  // baryons
86  G4BaryonConstructor bConstructor;
87  bConstructor.ConstructParticle();
88 
89  G4IonConstructor iConstructor;
90  iConstructor.ConstructParticle();
91 
92  G4MesonConstructor mConstructor;
93  mConstructor.ConstructParticle();
94 }
95 
96 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
98 {
100  ConstructEM();
101  ConstructHad();
103 }
104 
105 // *** Processes and models
106 
107 // gamma
108 
109 #include "G4PhotoElectricEffect.hh"
111 
112 #include "G4ComptonScattering.hh"
114 
115 #include "G4GammaConversion.hh"
117 
118 #include "G4RayleighScattering.hh"
120 
121 // e-
122 
123 #include "G4eMultipleScattering.hh"
124 #include "G4UniversalFluctuation.hh"
125 
126 #include "G4eIonisation.hh"
128 
129 #include "G4eBremsstrahlung.hh"
131 
132 // e+
133 
134 #include "G4eplusAnnihilation.hh"
135 
136 // mu
137 
138 #include "G4MuIonisation.hh"
139 #include "G4MuBremsstrahlung.hh"
140 #include "G4MuPairProduction.hh"
141 
142 // hadrons
143 
144 #include "G4hMultipleScattering.hh"
145 #include "G4MscStepLimitType.hh"
146 
147 #include "G4hBremsstrahlung.hh"
148 #include "G4hPairProduction.hh"
149 
150 #include "G4hIonisation.hh"
151 #include "G4ionIonisation.hh"
153 
154 //
155 
156 #include "G4LossTableManager.hh"
157 #include "G4EmProcessOptions.hh"
158 
159 
160 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
162 {
164 
165  while( (*theParticleIterator)() ){
166 
168  G4ProcessManager* pmanager = particle->GetProcessManager();
169  G4String particleName = particle->GetParticleName();
170 
171  if (particleName == "gamma") {
172 
173 
174  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
175  G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel = new G4LivermorePhotoElectricModel();
176  thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
177  pmanager->AddDiscreteProcess(thePhotoElectricEffect);
178 
179  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
180  G4LivermoreComptonModel* theLivermoreComptonModel = new G4LivermoreComptonModel();
181  theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
182  pmanager->AddDiscreteProcess(theComptonScattering);
183 
184  G4GammaConversion* theGammaConversion = new G4GammaConversion();
185  G4LivermoreGammaConversionModel* theLivermoreGammaConversionModel = new G4LivermoreGammaConversionModel();
186  theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
187  pmanager->AddDiscreteProcess(theGammaConversion);
188 
189  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
190  G4LivermoreRayleighModel* theRayleighModel = new G4LivermoreRayleighModel();
191  theRayleigh->AddEmModel(0, theRayleighModel);
192  pmanager->AddDiscreteProcess(theRayleigh);
193 
194  pmanager->AddProcess(new G4StepLimiter(), -1, -1, 5);
195 
196  } else if (particleName == "e-") {
197 
200  pmanager->AddProcess(msc, -1, 1, 1);
201 
202  // Ionisation
203  G4eIonisation* eIoni = new G4eIonisation();
205  eIoni->SetStepFunction(0.2, 100*um); //
206  pmanager->AddProcess(eIoni, -1, 2, 2);
207 
208  // Bremsstrahlung
209  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
211  pmanager->AddProcess(eBrem, -1,-3, 3);
212 
213  pmanager->AddProcess(new G4StepLimiter(), -1, -1, 4);
214 
215  } else if (particleName == "e+") {
216 
217  // Identical to G4EmStandardPhysics_option3
218 
221  pmanager->AddProcess(msc, -1, 1, 1);
222 
223  G4eIonisation* eIoni = new G4eIonisation();
224  eIoni->SetStepFunction(0.2, 100*um);
225  pmanager->AddProcess(eIoni, -1, 2, 2);
226 
227  pmanager->AddProcess(new G4eBremsstrahlung, -1,-3, 3);
228 
229  pmanager->AddProcess(new G4eplusAnnihilation,0,-1, 4);
230 
231  pmanager->AddProcess(new G4StepLimiter(), -1, -1, 5);
232 
233  } else if (particleName == "GenericIon") {
234 
235  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
236 
237  G4ionIonisation* ionIoni = new G4ionIonisation();
238  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
239  ionIoni->SetStepFunction(0.1, 20*um);
240  pmanager->AddProcess(ionIoni, -1, 2, 2);
241 
242  pmanager->AddProcess(new G4StepLimiter(), -1, -1, 3);
243 
244  } else if (particleName == "alpha" ||
245  particleName == "He3" ) {
246 
247  // Identical to G4EmStandardPhysics_option3
248 
249  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
250 
251  G4ionIonisation* ionIoni = new G4ionIonisation();
252  ionIoni->SetStepFunction(0.1, 20*um);
253  pmanager->AddProcess(ionIoni, -1, 2, 2);
254 
255  pmanager->AddProcess(new G4StepLimiter(), -1, -1, 3);
256  }
257 
258  //
259 
260  }
261 }
262 
263 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
264 #include "G4HadronElasticProcess.hh"
265 #include "G4LElastic.hh"
266 
269 #include "G4TripathiCrossSection.hh"
270 #include "G4IonsShenCrossSection.hh"
271 #include "G4LEAlphaInelastic.hh"
272 
274 {
275 
276  G4HadronElasticProcess * theElasticProcess = new G4HadronElasticProcess;
277  theElasticProcess->RegisterMe( new G4LElastic() );
278 
280  while( (*theParticleIterator)() )
281  {
283  G4ProcessManager* pManager = particle->GetProcessManager();
284 
285  if (particle->GetParticleName() == "alpha")
286  {
287 
288  // INELASTIC SCATTERING
289  // Binary Cascade
291  theBC -> SetMinEnergy(80.*MeV);
292  theBC -> SetMaxEnergy(40.*GeV);
293 
294  // TRIPATHI CROSS SECTION
295  // Implementation of formulas in analogy to NASA technical paper 3621 by
296  // Tripathi, et al. Cross-sections for ion ion scattering
297  G4TripathiCrossSection* TripathiCrossSection = new G4TripathiCrossSection;
298 
299  // IONS SHEN CROSS SECTION
300  // Implementation of formulas
301  // Shen et al. Nuc. Phys. A 491 130 (1989)
302  // Total Reaction Cross Section for Heavy-Ion Collisions
304 
305  // Final state production model for Alpha inelastic scattering below 20 GeV
306  G4LEAlphaInelastic* theAIModel = new G4LEAlphaInelastic;
307  theAIModel -> SetMaxEnergy(100.*MeV);
308 
310  theIPalpha->AddDataSet(TripathiCrossSection);
311  theIPalpha->AddDataSet(aShen);
312 
313  // Register the Alpha Inelastic and Binary Cascade Model
314  theIPalpha->RegisterMe(theAIModel);
315  theIPalpha->RegisterMe(theBC);
316 
317  // Activate the alpha inelastic scattering using the alpha inelastic and binary cascade model
318  pManager -> AddDiscreteProcess(theIPalpha);
319 
320  // Activate the Hadron Elastic Process
321  pManager -> AddDiscreteProcess(theElasticProcess);
322 
323  }
324  }
325 
326 }
327 
328 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
330 { }
331 
332 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
334 {
335  if (verboseLevel >0){
336  G4cout << "DicomPhysicsList::SetCuts:";
337  G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl;
338  }
339 
340  // set cut values for gamma at first and for e- second and next for e+,
341  // because some processes for e+/e- need cut values for gamma
342  SetCutValue(fCutForGamma, "gamma");
343  SetCutValue(fCutForElectron, "e-");
344  SetCutValue(fCutForPositron, "e+");
345 
347 
348 }
349