source: trunk/source/physics_lists/builders/src/G4EmLivermorePhysics.cc@ 1344

Last change on this file since 1344 was 1340, checked in by garnier, 15 years ago

update ti head

File size: 12.8 KB
Line 
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// $Id: G4EmLivermorePhysics.cc,v 1.12 2010/10/10 15:18:34 vnivanch Exp $
27// GEANT4 tag $Name: phys-lists-V09-03-34 $
28
29#include "G4EmLivermorePhysics.hh"
30
31#include "G4ParticleDefinition.hh"
32#include "G4ProcessManager.hh"
33
34// *** Processes and models
35
36// gamma
37
38#include "G4PhotoElectricEffect.hh"
39#include "G4LivermorePhotoElectricModel.hh"
40
41#include "G4ComptonScattering.hh"
42#include "G4LivermoreComptonModel.hh"
43
44#include "G4GammaConversion.hh"
45#include "G4LivermoreGammaConversionModel.hh"
46
47#include "G4RayleighScattering.hh"
48#include "G4LivermoreRayleighModel.hh"
49
50// e-
51
52#include "G4eMultipleScattering.hh"
53#include "G4UniversalFluctuation.hh"
54
55#include "G4eIonisation.hh"
56#include "G4LivermoreIonisationModel.hh"
57
58#include "G4eBremsstrahlung.hh"
59#include "G4LivermoreBremsstrahlungModel.hh"
60
61// e+
62
63#include "G4eplusAnnihilation.hh"
64
65// mu
66
67#include "G4MuMultipleScattering.hh"
68#include "G4MuIonisation.hh"
69#include "G4MuBremsstrahlung.hh"
70#include "G4MuPairProduction.hh"
71
72// hadrons
73
74#include "G4hMultipleScattering.hh"
75#include "G4MscStepLimitType.hh"
76
77#include "G4hBremsstrahlung.hh"
78#include "G4hPairProduction.hh"
79
80#include "G4hIonisation.hh"
81#include "G4ionIonisation.hh"
82#include "G4alphaIonisation.hh"
83#include "G4IonParametrisedLossModel.hh"
84#include "G4NuclearStopping.hh"
85
86// msc models
87#include "G4UrbanMscModel93.hh"
88#include "G4WentzelVIModel.hh"
89#include "G4GoudsmitSaundersonMscModel.hh"
90#include "G4CoulombScattering.hh"
91
92//
93
94#include "G4LossTableManager.hh"
95#include "G4EmProcessOptions.hh"
96
97// particles
98
99#include "G4Gamma.hh"
100#include "G4Electron.hh"
101#include "G4Positron.hh"
102#include "G4MuonPlus.hh"
103#include "G4MuonMinus.hh"
104#include "G4PionPlus.hh"
105#include "G4PionMinus.hh"
106#include "G4KaonPlus.hh"
107#include "G4KaonMinus.hh"
108#include "G4Proton.hh"
109#include "G4AntiProton.hh"
110#include "G4Deuteron.hh"
111#include "G4Triton.hh"
112#include "G4He3.hh"
113#include "G4Alpha.hh"
114#include "G4GenericIon.hh"
115
116//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
117
118G4EmLivermorePhysics::G4EmLivermorePhysics(G4int ver)
119 : G4VPhysicsConstructor("G4EmLivermorePhysics"), verbose(ver)
120{
121 G4LossTableManager::Instance();
122}
123
124//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125
126G4EmLivermorePhysics::G4EmLivermorePhysics(G4int ver, const G4String&)
127 : G4VPhysicsConstructor("G4EmLivermorePhysics"), verbose(ver)
128{
129 G4LossTableManager::Instance();
130}
131
132//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
133
134G4EmLivermorePhysics::~G4EmLivermorePhysics()
135{}
136
137//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
138
139void G4EmLivermorePhysics::ConstructParticle()
140{
141// gamma
142 G4Gamma::Gamma();
143
144// leptons
145 G4Electron::Electron();
146 G4Positron::Positron();
147 G4MuonPlus::MuonPlus();
148 G4MuonMinus::MuonMinus();
149
150// mesons
151 G4PionPlus::PionPlusDefinition();
152 G4PionMinus::PionMinusDefinition();
153 G4KaonPlus::KaonPlusDefinition();
154 G4KaonMinus::KaonMinusDefinition();
155
156// baryons
157 G4Proton::Proton();
158 G4AntiProton::AntiProton();
159
160// ions
161 G4Deuteron::Deuteron();
162 G4Triton::Triton();
163 G4He3::He3();
164 G4Alpha::Alpha();
165 G4GenericIon::GenericIonDefinition();
166}
167
168//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
169
170void G4EmLivermorePhysics::ConstructProcess()
171{
172 // Add Livermore EM Processes
173
174 theParticleIterator->reset();
175
176 while( (*theParticleIterator)() ){
177
178 G4ParticleDefinition* particle = theParticleIterator->value();
179 G4ProcessManager* pmanager = particle->GetProcessManager();
180 G4String particleName = particle->GetParticleName();
181
182 if(verbose > 1)
183 G4cout << "### " << GetPhysicsName() << " instantiates for "
184 << particleName << G4endl;
185
186 //Applicability range for Livermore models
187 //for higher energies, the Standard models are used
188 G4double LivermoreHighEnergyLimit = GeV;
189
190 if (particleName == "gamma") {
191
192 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
193 G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel =
194 new G4LivermorePhotoElectricModel();
195 theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
196 thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
197 pmanager->AddDiscreteProcess(thePhotoElectricEffect);
198
199 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
200 G4LivermoreComptonModel* theLivermoreComptonModel =
201 new G4LivermoreComptonModel();
202 theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
203 theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
204 pmanager->AddDiscreteProcess(theComptonScattering);
205
206 G4GammaConversion* theGammaConversion = new G4GammaConversion();
207 G4LivermoreGammaConversionModel* theLivermoreGammaConversionModel =
208 new G4LivermoreGammaConversionModel();
209 theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
210 theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
211 pmanager->AddDiscreteProcess(theGammaConversion);
212
213 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
214 G4LivermoreRayleighModel* theRayleighModel = new G4LivermoreRayleighModel();
215 theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
216 theRayleigh->AddEmModel(0, theRayleighModel);
217 pmanager->AddDiscreteProcess(theRayleigh);
218
219 } else if (particleName == "e-") {
220
221 G4eMultipleScattering* msc = new G4eMultipleScattering();
222 //msc->AddEmModel(0, new G4UrbanMscModel93());
223 msc->AddEmModel(0, new G4GoudsmitSaundersonMscModel());
224 msc->SetStepLimitType(fUseDistanceToBoundary);
225 pmanager->AddProcess(msc, -1, 1, 1);
226
227 // Ionisation
228 G4eIonisation* eIoni = new G4eIonisation();
229 G4LivermoreIonisationModel* theIoniLivermore = new
230 G4LivermoreIonisationModel();
231 theIoniLivermore->SetHighEnergyLimit(1*MeV);
232 eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
233 eIoni->SetStepFunction(0.2, 100*um); //
234 pmanager->AddProcess(eIoni, -1, 2, 2);
235
236 // Bremsstrahlung
237 G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
238 G4LivermoreBremsstrahlungModel* theBremLivermore = new
239 G4LivermoreBremsstrahlungModel();
240 theBremLivermore->SetHighEnergyLimit(LivermoreHighEnergyLimit);
241 eBrem->AddEmModel(0, theBremLivermore);
242 pmanager->AddProcess(eBrem, -1,-3, 3);
243
244 } else if (particleName == "e+") {
245
246 // Identical to G4EmStandardPhysics_option3
247
248 G4eMultipleScattering* msc = new G4eMultipleScattering();
249 //msc->AddEmModel(0, new G4UrbanMscModel93());
250 msc->AddEmModel(0, new G4GoudsmitSaundersonMscModel());
251 msc->SetStepLimitType(fUseDistanceToBoundary);
252 pmanager->AddProcess(msc, -1, 1, 1);
253
254 G4eIonisation* eIoni = new G4eIonisation();
255 eIoni->SetStepFunction(0.2, 100*um);
256
257 pmanager->AddProcess(eIoni, -1, 2, 2);
258 pmanager->AddProcess(new G4eBremsstrahlung, -1,-3, 3);
259 pmanager->AddProcess(new G4eplusAnnihilation,0,-1, 4);
260
261 } else if (particleName == "mu+" ||
262 particleName == "mu-" ) {
263
264 // Identical to G4EmStandardPhysics_option3
265
266 G4MuMultipleScattering* msc = new G4MuMultipleScattering();
267 msc->AddEmModel(0, new G4WentzelVIModel());
268 pmanager->AddProcess(msc, -1, 1, 1);
269
270 G4MuIonisation* muIoni = new G4MuIonisation();
271 muIoni->SetStepFunction(0.2, 50*um);
272
273 pmanager->AddProcess(muIoni, -1, 2, 2);
274 pmanager->AddProcess(new G4MuBremsstrahlung, -1,-3, 3);
275 pmanager->AddProcess(new G4MuPairProduction, -1,-4, 4);
276 pmanager->AddDiscreteProcess(new G4CoulombScattering());
277
278 } else if (particleName == "GenericIon") {
279
280 pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
281
282 G4ionIonisation* ionIoni = new G4ionIonisation();
283 ionIoni->SetEmModel(new G4IonParametrisedLossModel());
284 ionIoni->SetStepFunction(0.1, 10*um);
285 pmanager->AddProcess(ionIoni, -1, 2, 2);
286 pmanager->AddProcess(new G4NuclearStopping(), -1, 3,-1);
287
288 } else if (particleName == "alpha" ||
289 particleName == "He3" ) {
290
291 // Identical to G4EmStandardPhysics_option3
292
293 pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
294
295 G4ionIonisation* ionIoni = new G4ionIonisation();
296 ionIoni->SetStepFunction(0.1, 20*um);
297 pmanager->AddProcess(ionIoni, -1, 2, 2);
298 pmanager->AddProcess(new G4NuclearStopping(), -1, 3,-1);
299
300 } else if (particleName == "pi+" ||
301 particleName == "pi-" ||
302 particleName == "kaon+" ||
303 particleName == "kaon-" ||
304 particleName == "proton" ) {
305
306 // Identical to G4EmStandardPhysics_option3
307
308 pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
309
310 G4hIonisation* hIoni = new G4hIonisation();
311 hIoni->SetStepFunction(0.2, 50*um);
312
313 pmanager->AddProcess(hIoni, -1, 2, 2);
314 pmanager->AddProcess(new G4hBremsstrahlung, -1,-3, 3);
315 pmanager->AddProcess(new G4hPairProduction, -1,-4, 4);
316
317 } else if (particleName == "B+" ||
318 particleName == "B-" ||
319 particleName == "D+" ||
320 particleName == "D-" ||
321 particleName == "Ds+" ||
322 particleName == "Ds-" ||
323 particleName == "anti_He3" ||
324 particleName == "anti_alpha" ||
325 particleName == "anti_deuteron" ||
326 particleName == "anti_lambda_c+" ||
327 particleName == "anti_omega-" ||
328 particleName == "anti_proton" ||
329 particleName == "anti_sigma_c+" ||
330 particleName == "anti_sigma_c++" ||
331 particleName == "anti_sigma+" ||
332 particleName == "anti_sigma-" ||
333 particleName == "anti_triton" ||
334 particleName == "anti_xi_c+" ||
335 particleName == "anti_xi-" ||
336 particleName == "deuteron" ||
337 particleName == "lambda_c+" ||
338 particleName == "omega-" ||
339 particleName == "sigma_c+" ||
340 particleName == "sigma_c++" ||
341 particleName == "sigma+" ||
342 particleName == "sigma-" ||
343 particleName == "tau+" ||
344 particleName == "tau-" ||
345 particleName == "triton" ||
346 particleName == "xi_c+" ||
347 particleName == "xi-" ) {
348
349 // Identical to G4EmStandardPhysics_option3
350
351 pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
352 pmanager->AddProcess(new G4hIonisation, -1, 2, 2);
353
354 }
355 }
356
357 // Em options
358 //
359 G4EmProcessOptions opt;
360 opt.SetVerbose(verbose);
361
362 // Multiple Coulomb scattering
363 //
364 //opt.SetMscStepLimitation(fUseDistanceToBoundary);
365 //opt.SetMscRangeFactor(0.02);
366
367 // Physics tables
368 //
369
370 opt.SetMinEnergy(100*eV);
371 opt.SetMaxEnergy(10*TeV);
372 opt.SetDEDXBinning(220);
373 opt.SetLambdaBinning(220);
374
375 //opt.SetSplineFlag(true);
376 opt.SetPolarAngleLimit(0.2);
377
378 // Ionization
379 //
380 //opt.SetSubCutoff(true);
381}
382
383//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.