Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4EmLivermorePolarizedPhysics.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// $Id$
27
30#include "G4SystemOfUnits.hh"
31
32// *** Processes and models
33
34// gamma
37
40
41#include "G4GammaConversion.hh"
43
46
47// e+-
50
51#include "G4eIonisation.hh"
53
54#include "G4eBremsstrahlung.hh"
56
57// e+
59
60// mu+-
62#include "G4MuIonisation.hh"
63#include "G4MuBremsstrahlung.hh"
64#include "G4MuPairProduction.hh"
65
70
71// hadrons
73#include "G4MscStepLimitType.hh"
74
75#include "G4hBremsstrahlung.hh"
76#include "G4hPairProduction.hh"
77
78#include "G4hIonisation.hh"
79#include "G4ionIonisation.hh"
80#include "G4alphaIonisation.hh"
82#include "G4NuclearStopping.hh"
83
84// msc models
85#include "G4UrbanMscModel95.hh"
86#include "G4WentzelVIModel.hh"
90
91// interfaces
92#include "G4LossTableManager.hh"
93#include "G4EmProcessOptions.hh"
95
96// particles
97#include "G4Gamma.hh"
98#include "G4Electron.hh"
99#include "G4Positron.hh"
100#include "G4MuonPlus.hh"
101#include "G4MuonMinus.hh"
102#include "G4PionPlus.hh"
103#include "G4PionMinus.hh"
104#include "G4KaonPlus.hh"
105#include "G4KaonMinus.hh"
106#include "G4Proton.hh"
107#include "G4AntiProton.hh"
108#include "G4Deuteron.hh"
109#include "G4Triton.hh"
110#include "G4He3.hh"
111#include "G4Alpha.hh"
112#include "G4GenericIon.hh"
113
114//
115#include "G4PhysicsListHelper.hh"
116#include "G4BuilderType.hh"
117
118// factory
120//
122
123
124//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125
127 : G4VPhysicsConstructor("G4EmLivermorePolarizedPhysics"), verbose(ver)
128{
131}
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134
136 : G4VPhysicsConstructor("G4EmLivermorePolarizedPhysics"), verbose(ver)
137{
140}
141
142//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
143
145{}
146
147//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148
150{
151// gamma
153
154// leptons
159
160// mesons
165
166// baryons
169
170// ions
173 G4He3::He3();
176}
177
178//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
179
181{
183
184 // muon & hadron bremsstrahlung and pair production
193
194 // muon & hadron multiple scattering
196 mumsc->AddEmModel(0, new G4WentzelVIModel());
198 pimsc->AddEmModel(0, new G4WentzelVIModel());
200 kmsc->AddEmModel(0, new G4WentzelVIModel());
202 pmsc->AddEmModel(0, new G4WentzelVIModel());
203 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
204
205 // high energy limit for e+- scattering models
206 G4double highEnergyLimit = 100*MeV;
207
208 // nuclear stopping
209 G4NuclearStopping* ionnuc = new G4NuclearStopping();
211
212 // Add Livermore EM Processes
214
215 while( (*theParticleIterator)() ){
216
218 G4String particleName = particle->GetParticleName();
219
220 if(verbose > 1)
221 G4cout << "### " << GetPhysicsName() << " instantiates for "
222 << particleName << G4endl;
223
224 //Applicability range for Livermore models
225 //for higher energies, the Standard models are used
226 G4double LivermoreHighEnergyLimit = GeV;
227
228 if (particleName == "gamma") {
229
230 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
232 theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
233 thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
234 ph->RegisterProcess(thePhotoElectricEffect, particle);
235
236 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
238 theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
239 theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
240 ph->RegisterProcess(theComptonScattering, particle);
241
242 G4GammaConversion* theGammaConversion = new G4GammaConversion();
244 theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
245 theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
246 ph->RegisterProcess(theGammaConversion, particle);
247
248 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
250 theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
251 theRayleigh->AddEmModel(0, theRayleighModel);
252 ph->RegisterProcess(theRayleigh, particle);
253
254 } else if (particleName == "e-") {
255
256 // multiple scattering
261 msc1->SetHighEnergyLimit(highEnergyLimit);
262 msc2->SetLowEnergyLimit(highEnergyLimit);
263 msc->AddEmModel(0, msc1);
264 msc->AddEmModel(0, msc2);
265
268 ss->SetEmModel(ssm, 1);
269 ss->SetMinKinEnergy(highEnergyLimit);
270 ssm->SetLowEnergyLimit(highEnergyLimit);
271 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
272 ph->RegisterProcess(msc, particle);
273 ph->RegisterProcess(ss, particle);
274
275 // Ionisation
276 G4eIonisation* eIoni = new G4eIonisation();
277 G4LivermoreIonisationModel* theIoniLivermore = new
279 theIoniLivermore->SetHighEnergyLimit(0.1*MeV);
280 eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
281 eIoni->SetStepFunction(0.2, 100*um); //
282 ph->RegisterProcess(eIoni, particle);
283
284 // Bremsstrahlung from standard
286 ph->RegisterProcess(eBrem, particle);
287
288 } else if (particleName == "e+") {
289
290 // multiple scattering
295 msc1->SetHighEnergyLimit(highEnergyLimit);
296 msc2->SetLowEnergyLimit(highEnergyLimit);
297 msc->AddEmModel(0, msc1);
298 msc->AddEmModel(0, msc2);
299
302 ss->SetEmModel(ssm, 1);
303 ss->SetMinKinEnergy(highEnergyLimit);
304 ssm->SetLowEnergyLimit(highEnergyLimit);
305 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
306
307 // Ionisation
308 G4eIonisation* eIoni = new G4eIonisation();
309 eIoni->SetStepFunction(0.2, 100*um);
310
311 ph->RegisterProcess(msc, particle);
312 ph->RegisterProcess(eIoni, particle);
313 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
314 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
315 ph->RegisterProcess(ss, particle);
316
317 } else if (particleName == "mu+" ||
318 particleName == "mu-" ) {
319
320 G4MuIonisation* muIoni = new G4MuIonisation();
321 muIoni->SetStepFunction(0.2, 50*um);
322
323 ph->RegisterProcess(mumsc, particle);
324 ph->RegisterProcess(muIoni, particle);
325 ph->RegisterProcess(mub, particle);
326 ph->RegisterProcess(mup, particle);
327 ph->RegisterProcess(new G4CoulombScattering(), particle);
328
329 } else if (particleName == "alpha" ||
330 particleName == "He3" ) {
331
332 // Identical to G4EmStandardPhysics_option3
333
335 G4ionIonisation* ionIoni = new G4ionIonisation();
336 ionIoni->SetStepFunction(0.1, 10*um);
337
338 ph->RegisterProcess(msc, particle);
339 ph->RegisterProcess(ionIoni, particle);
340 ph->RegisterProcess(ionnuc, particle);
341
342 } else if (particleName == "GenericIon") {
343
344 // Identical to G4EmStandardPhysics_option3
345
346 G4ionIonisation* ionIoni = new G4ionIonisation();
348 ionIoni->SetStepFunction(0.1, 1*um);
349
350 ph->RegisterProcess(hmsc, particle);
351 ph->RegisterProcess(ionIoni, particle);
352 ph->RegisterProcess(ionnuc, particle);
353
354 } else if (particleName == "pi+" ||
355 particleName == "pi-" ) {
356
357 //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
358 G4hIonisation* hIoni = new G4hIonisation();
359 hIoni->SetStepFunction(0.2, 50*um);
360
361 ph->RegisterProcess(pimsc, particle);
362 ph->RegisterProcess(hIoni, particle);
363 ph->RegisterProcess(pib, particle);
364 ph->RegisterProcess(pip, particle);
365
366 } else if (particleName == "kaon+" ||
367 particleName == "kaon-" ) {
368
369 //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
370 G4hIonisation* hIoni = new G4hIonisation();
371 hIoni->SetStepFunction(0.2, 50*um);
372
373 ph->RegisterProcess(kmsc, particle);
374 ph->RegisterProcess(hIoni, particle);
375 ph->RegisterProcess(kb, particle);
376 ph->RegisterProcess(kp, particle);
377
378 } else if (particleName == "proton" ||
379 particleName == "anti_proton") {
380
381 //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
382 G4hIonisation* hIoni = new G4hIonisation();
383 hIoni->SetStepFunction(0.2, 50*um);
384
385 ph->RegisterProcess(pmsc, particle);
386 ph->RegisterProcess(hIoni, particle);
387 ph->RegisterProcess(pb, particle);
388 ph->RegisterProcess(pp, particle);
389 ph->RegisterProcess(pnuc, particle);
390
391 } else if (particleName == "B+" ||
392 particleName == "B-" ||
393 particleName == "D+" ||
394 particleName == "D-" ||
395 particleName == "Ds+" ||
396 particleName == "Ds-" ||
397 particleName == "anti_He3" ||
398 particleName == "anti_alpha" ||
399 particleName == "anti_deuteron" ||
400 particleName == "anti_lambda_c+" ||
401 particleName == "anti_omega-" ||
402 particleName == "anti_sigma_c+" ||
403 particleName == "anti_sigma_c++" ||
404 particleName == "anti_sigma+" ||
405 particleName == "anti_sigma-" ||
406 particleName == "anti_triton" ||
407 particleName == "anti_xi_c+" ||
408 particleName == "anti_xi-" ||
409 particleName == "deuteron" ||
410 particleName == "lambda_c+" ||
411 particleName == "omega-" ||
412 particleName == "sigma_c+" ||
413 particleName == "sigma_c++" ||
414 particleName == "sigma+" ||
415 particleName == "sigma-" ||
416 particleName == "tau+" ||
417 particleName == "tau-" ||
418 particleName == "triton" ||
419 particleName == "xi_c+" ||
420 particleName == "xi-" ) {
421
422 // Identical to G4EmStandardPhysics_option3
423
424 ph->RegisterProcess(hmsc, particle);
425 ph->RegisterProcess(new G4hIonisation(), particle);
426 ph->RegisterProcess(pnuc, particle);
427 }
428 }
429
430 // Em options
431 //
433 opt.SetVerbose(verbose);
434
435 // Multiple Coulomb scattering
436 //
437 opt.SetPolarAngleLimit(CLHEP::pi);
438
439 // Physics tables
440 //
441
442 opt.SetMinEnergy(100*eV);
443 opt.SetMaxEnergy(10*TeV);
444 opt.SetDEDXBinning(220);
445 opt.SetLambdaBinning(220);
446
447 // Nuclear stopping
448 pnuc->SetMaxKinEnergy(MeV);
449
450 // Ionization
451 //
452 //opt.SetSubCutoff(true);
453
454 // Deexcitation
455 //
458 de->SetFluo(true);
459}
460
461//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
@ fUseDistanceToBoundary
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetMaxEnergy(G4double val)
void SetDEDXBinning(G4int val)
void SetLambdaBinning(G4int val)
void SetPolarAngleLimit(G4double val)
void SetVerbose(G4int val, const G4String &name="all")
void SetMinEnergy(G4double val)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:87
static G4He3 * He3()
Definition: G4He3.cc:94
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
static G4KaonPlus * KaonPlusDefinition()
Definition: G4KaonPlus.cc:108
void SetAtomDeexcitation(G4VAtomDeexcitation *)
static G4LossTableManager * Instance()
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
const G4String & GetParticleName() const
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:93
static G4PionPlus * PionPlusDefinition()
Definition: G4PionPlus.cc:93
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Triton * Triton()
Definition: G4Triton.cc:95
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:606
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
void SetMinKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=1)
void SetMaxKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=1)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
void SetStepFunction(G4double v1, G4double v2)
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
void SetStepLimitType(G4MscStepLimitType val)
const G4String & GetPhysicsName() const
G4ParticleTable::G4PTblDicIterator * theParticleIterator