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

Last change on this file since 1331 was 1315, checked in by garnier, 15 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

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