Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4EmLivermorePhysics.cc
이 파일의 문서화 페이지로 가기
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 109526 2018-04-30 07:11:52Z gcosmo $
27 
28 #include "G4EmLivermorePhysics.hh"
29 #include "G4ParticleDefinition.hh"
30 #include "G4SystemOfUnits.hh"
31 #include "G4ParticleTable.hh"
32 
33 // *** Processes and models
34 
35 // gamma
36 #include "G4PhotoElectricEffect.hh"
38 
39 #include "G4ComptonScattering.hh"
41 
42 #include "G4GammaConversion.hh"
44 
45 #include "G4RayleighScattering.hh"
47 
48 // e+-
49 #include "G4eMultipleScattering.hh"
51 
52 #include "G4eIonisation.hh"
54 
55 #include "G4eBremsstrahlung.hh"
57 #include "G4Generator2BS.hh"
58 #include "G4SeltzerBergerModel.hh"
59 
60 // e+
61 #include "G4eplusAnnihilation.hh"
62 
63 // mu+-
65 #include "G4MuIonisation.hh"
66 #include "G4MuBremsstrahlung.hh"
67 #include "G4MuPairProduction.hh"
68 
73 
74 // hadrons
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"
85 #include "G4NuclearStopping.hh"
86 
87 // msc models
88 #include "G4UrbanMscModel.hh"
89 #include "G4WentzelVIModel.hh"
91 #include "G4CoulombScattering.hh"
93 
94 // interfaces
95 #include "G4LossTableManager.hh"
96 #include "G4UAtomicDeexcitation.hh"
97 #include "G4EmParameters.hh"
98 
99 // particles
100 
101 #include "G4Gamma.hh"
102 #include "G4Electron.hh"
103 #include "G4Positron.hh"
104 #include "G4MuonPlus.hh"
105 #include "G4MuonMinus.hh"
106 #include "G4PionPlus.hh"
107 #include "G4PionMinus.hh"
108 #include "G4KaonPlus.hh"
109 #include "G4KaonMinus.hh"
110 #include "G4Proton.hh"
111 #include "G4AntiProton.hh"
112 #include "G4Deuteron.hh"
113 #include "G4Triton.hh"
114 #include "G4He3.hh"
115 #include "G4Alpha.hh"
116 #include "G4GenericIon.hh"
117 
118 //
119 #include "G4PhysicsListHelper.hh"
120 #include "G4BuilderType.hh"
121 #include "G4EmModelActivator.hh"
122 
123 // factory
125 //
127 
128 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
129 
131  : G4VPhysicsConstructor("G4EmLivermore"), verbose(ver)
132 {
134  param->SetDefaults();
135  param->SetVerbose(verbose);
136  param->SetMinEnergy(100*eV);
137  param->SetLowestElectronEnergy(100*eV);
138  param->SetNumberOfBinsPerDecade(20);
140  param->SetMscRangeFactor(0.02);
141  param->SetMuHadLateralDisplacement(true);
143  param->SetFluo(true);
145 }
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148 
150 {}
151 
152 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
153 
155 {
156  // gamma
157  G4Gamma::Gamma();
158 
159  // leptons
164 
165  // mesons
170 
171  // baryons
174 
175  // ions
178  G4He3::He3();
179  G4Alpha::Alpha();
181 }
182 
183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
184 
186 {
187  if(verbose > 1) {
188  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
189  }
191 
192  // muon & hadron bremsstrahlung and pair production
201 
202  // muon & hadron multiple scattering
204  mumsc->SetEmModel(new G4WentzelVIModel());
205 
206  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
207 
208  // high energy limit for e+- scattering models
209  G4double highEnergyLimit = G4EmParameters::Instance()->MscEnergyLimit();
210 
211  // nuclear stopping
212  G4NuclearStopping* pnuc = new G4NuclearStopping();
213 
214  // Add Livermore EM Processes
216  for(const auto& particleName : partList.PartNames()) {
217  G4ParticleDefinition* particle = table->FindParticle(particleName);
218  if (!particle) { continue; }
219  if (particleName == "gamma") {
220 
221  // photoelectric effect - Livermore model only
222  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
223  thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
224  ph->RegisterProcess(thePhotoElectricEffect, particle);
225 
226  // Compton scattering - Livermore model only
227  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
228  theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
229  ph->RegisterProcess(theComptonScattering, particle);
230 
231  // gamma conversion - Livermore model below 80 GeV
232  G4GammaConversion* theGammaConversion = new G4GammaConversion();
233  theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
234  ph->RegisterProcess(theGammaConversion, particle);
235 
236  // default Rayleigh scattering is Livermore
237  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
238  ph->RegisterProcess(theRayleigh, particle);
239 
240  } else if (particleName == "e-") {
241 
242  // multiple and single scattering
244  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
245  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
246  msc1->SetHighEnergyLimit(highEnergyLimit);
247  msc2->SetLowEnergyLimit(highEnergyLimit);
248  msc->SetEmModel(msc1);
249  msc->SetEmModel(msc2);
250 
253  ss->SetEmModel(ssm);
254  ss->SetMinKinEnergy(highEnergyLimit);
255  ssm->SetLowEnergyLimit(highEnergyLimit);
256  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
257 
258  // Ionisation - Livermore should be used only for low energies
259  G4eIonisation* eIoni = new G4eIonisation();
260  G4LivermoreIonisationModel* theIoniLivermore = new
262  theIoniLivermore->SetHighEnergyLimit(0.1*MeV);
263  eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
264  eIoni->SetStepFunction(0.2, 100*um); //
265 
266  // Bremsstrahlung
267  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
268  G4VEmModel* theBremLivermore = new G4LivermoreBremsstrahlungModel();
269  theBremLivermore->SetHighEnergyLimit(1*GeV);
270  theBremLivermore->SetAngularDistribution(new G4Generator2BS());
271  eBrem->SetEmModel(theBremLivermore);
272 
273  // register processes
274  ph->RegisterProcess(msc, particle);
275  ph->RegisterProcess(eIoni, particle);
276  ph->RegisterProcess(eBrem, particle);
277  ph->RegisterProcess(ss, particle);
278 
279  } else if (particleName == "e+") {
280 
281  // multiple and single scattering
283  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
284  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
285  msc1->SetHighEnergyLimit(highEnergyLimit);
286  msc2->SetLowEnergyLimit(highEnergyLimit);
287  msc->SetEmModel(msc1);
288  msc->SetEmModel(msc2);
289 
292  ss->SetEmModel(ssm);
293  ss->SetMinKinEnergy(highEnergyLimit);
294  ssm->SetLowEnergyLimit(highEnergyLimit);
295  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
296 
297  // ionisation
298  G4eIonisation* eIoni = new G4eIonisation();
299  eIoni->SetStepFunction(0.2, 100*um);
300 
301  // Bremsstrahlung from standard
302  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
303  G4VEmModel* mod = new G4SeltzerBergerModel();
305  eBrem->SetEmModel(mod);
306 
307  // register processes
308  ph->RegisterProcess(msc, particle);
309  ph->RegisterProcess(eIoni, particle);
310  ph->RegisterProcess(eBrem, particle);
311  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
312  ph->RegisterProcess(ss, particle);
313 
314  } else if (particleName == "mu+" ||
315  particleName == "mu-" ) {
316 
317  G4MuIonisation* muIoni = new G4MuIonisation();
318  muIoni->SetStepFunction(0.2, 50*um);
319 
320  ph->RegisterProcess(mumsc, particle);
321  ph->RegisterProcess(muIoni, particle);
322  ph->RegisterProcess(mub, particle);
323  ph->RegisterProcess(mup, particle);
324  ph->RegisterProcess(new G4CoulombScattering(), particle);
325 
326  } else if (particleName == "alpha" ||
327  particleName == "He3" ) {
328 
330  G4ionIonisation* ionIoni = new G4ionIonisation();
331  ionIoni->SetStepFunction(0.1, 10*um);
332 
333  ph->RegisterProcess(msc, particle);
334  ph->RegisterProcess(ionIoni, particle);
335  ph->RegisterProcess(pnuc, particle);
336 
337  } else if (particleName == "GenericIon") {
338 
339  G4ionIonisation* ionIoni = new G4ionIonisation();
340  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
341  ionIoni->SetStepFunction(0.1, 1*um);
342 
343  ph->RegisterProcess(hmsc, particle);
344  ph->RegisterProcess(ionIoni, particle);
345  ph->RegisterProcess(pnuc, particle);
346 
347  } else if (particleName == "pi+" ||
348  particleName == "pi-" ) {
349 
351  G4hIonisation* hIoni = new G4hIonisation();
352  hIoni->SetStepFunction(0.2, 50*um);
353 
354  ph->RegisterProcess(pimsc, particle);
355  ph->RegisterProcess(hIoni, particle);
356  ph->RegisterProcess(pib, particle);
357  ph->RegisterProcess(pip, particle);
358 
359  } else if (particleName == "kaon+" ||
360  particleName == "kaon-" ) {
361 
363  G4hIonisation* hIoni = new G4hIonisation();
364  hIoni->SetStepFunction(0.2, 50*um);
365 
366  ph->RegisterProcess(kmsc, particle);
367  ph->RegisterProcess(hIoni, particle);
368  ph->RegisterProcess(kb, particle);
369  ph->RegisterProcess(kp, particle);
370 
371  } else if (particleName == "proton" ||
372  particleName == "anti_proton") {
373 
375  pmsc->SetEmModel(new G4WentzelVIModel());
376 
377  G4hIonisation* hIoni = new G4hIonisation();
378  hIoni->SetStepFunction(0.2, 50*um);
379 
380  ph->RegisterProcess(pmsc, particle);
381  ph->RegisterProcess(hIoni, particle);
382  ph->RegisterProcess(pb, particle);
383  ph->RegisterProcess(pp, particle);
384  ph->RegisterProcess(pnuc, particle);
385 
386  } else if (particleName == "B+" ||
387  particleName == "B-" ||
388  particleName == "D+" ||
389  particleName == "D-" ||
390  particleName == "Ds+" ||
391  particleName == "Ds-" ||
392  particleName == "anti_He3" ||
393  particleName == "anti_alpha" ||
394  particleName == "anti_deuteron" ||
395  particleName == "anti_lambda_c+" ||
396  particleName == "anti_omega-" ||
397  particleName == "anti_sigma_c+" ||
398  particleName == "anti_sigma_c++" ||
399  particleName == "anti_sigma+" ||
400  particleName == "anti_sigma-" ||
401  particleName == "anti_triton" ||
402  particleName == "anti_xi_c+" ||
403  particleName == "anti_xi-" ||
404  particleName == "deuteron" ||
405  particleName == "lambda_c+" ||
406  particleName == "omega-" ||
407  particleName == "sigma_c+" ||
408  particleName == "sigma_c++" ||
409  particleName == "sigma+" ||
410  particleName == "sigma-" ||
411  particleName == "tau+" ||
412  particleName == "tau-" ||
413  particleName == "triton" ||
414  particleName == "xi_c+" ||
415  particleName == "xi-" ) {
416 
417  ph->RegisterProcess(hmsc, particle);
418  ph->RegisterProcess(new G4hIonisation(), particle);
419  ph->RegisterProcess(pnuc, particle);
420  }
421  }
422 
423  // Nuclear stopping
424  pnuc->SetMaxKinEnergy(MeV);
425 
426  // Deexcitation
427  //
430 
432 }
433 
434 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)
static G4He3 * He3()
Definition: G4He3.cc:94
const std::vector< G4String > & PartNames() const
void SetStepFunction(G4double v1, G4double v2, G4bool lock=true)
void SetEmModel(G4VEmModel *, G4int index=0)
static G4ParticleTable * GetParticleTable()
static constexpr double MeV
Definition: G4SIunits.hh:214
void SetMaxKinEnergy(G4double e)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:707
void SetMinKinEnergy(G4double e)
#define G4endl
Definition: G4ios.hh:61
void SetMuHadLateralDisplacement(G4bool val)
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
void SetVerbose(G4int val)
G4EmLivermorePhysics(G4int ver=1, const G4String &name="")
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlusDefinition()
Definition: G4PionPlus.cc:93
void ActivateAngularGeneratorForIonisation(G4bool val)
static constexpr double um
Definition: G4SIunits.hh:113
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetMinEnergy(G4double val)
double G4double
Definition: G4Types.hh:76
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:728
void SetEmModel(G4VEmModel *, G4int index=0)
void SetLowestElectronEnergy(G4double val)
void SetMscRangeFactor(G4double val)
const G4String & GetPhysicsName() const
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
void SetFluo(G4bool val)
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4Positron * Positron()
Definition: G4Positron.cc:94
static constexpr double eV
Definition: G4SIunits.hh:215
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:582
void SetNumberOfBinsPerDecade(G4int val)
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:93
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
void SetMscStepLimitType(G4MscStepLimitType val)
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
static G4PhysicsListHelper * GetPhysicsListHelper()
int G4int
Definition: G4Types.hh:78
static G4Triton * Triton()
Definition: G4Triton.cc:95
void SetEmModel(G4VMscModel *, size_t index=0)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:714
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
G4GLOB_DLL std::ostream G4cout
static G4LossTableManager * Instance()
void SetAtomDeexcitation(G4VAtomDeexcitation *)
virtual void ConstructParticle()
static G4KaonPlus * KaonPlusDefinition()
Definition: G4KaonPlus.cc:108
G4double MscEnergyLimit() const
static constexpr double GeV
Definition: G4SIunits.hh:217
static G4EmParameters * Instance()