source: trunk/source/physics_lists/builders/src/G4EmPenelopePhysics.cc@ 1350

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

update ti head

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