Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4EmPenelopePhysics.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: G4EmPenelopePhysics.cc 109526 2018-04-30 07:11:52Z gcosmo $
27 
28 #include "G4EmPenelopePhysics.hh"
29 #include "G4ParticleDefinition.hh"
30 #include "G4SystemOfUnits.hh"
31 #include "G4ParticleTable.hh"
32 
33 // *** Processes and models
34 
35 // gamma
36 
37 #include "G4PhotoElectricEffect.hh"
39 
40 #include "G4ComptonScattering.hh"
42 
43 #include "G4GammaConversion.hh"
45 
46 #include "G4RayleighScattering.hh"
48 
49 // e- and e+
50 
51 #include "G4eMultipleScattering.hh"
53 
54 #include "G4eIonisation.hh"
56 
57 #include "G4eBremsstrahlung.hh"
59 
60 // e+ only
61 
62 #include "G4eplusAnnihilation.hh"
64 
65 // mu
66 
68 #include "G4MuIonisation.hh"
69 #include "G4MuBremsstrahlung.hh"
70 #include "G4MuPairProduction.hh"
71 
76 
77 // hadrons
78 
79 #include "G4hMultipleScattering.hh"
80 #include "G4MscStepLimitType.hh"
81 
82 #include "G4hBremsstrahlung.hh"
83 #include "G4hPairProduction.hh"
84 
85 #include "G4hIonisation.hh"
86 #include "G4ionIonisation.hh"
87 #include "G4alphaIonisation.hh"
89 #include "G4NuclearStopping.hh"
90 
91 // msc models
92 #include "G4UrbanMscModel.hh"
94 #include "G4WentzelVIModel.hh"
95 #include "G4CoulombScattering.hh"
97 //
98 
99 #include "G4LossTableManager.hh"
100 #include "G4VAtomDeexcitation.hh"
101 #include "G4UAtomicDeexcitation.hh"
102 
103 // particles
104 
105 #include "G4Gamma.hh"
106 #include "G4Electron.hh"
107 #include "G4Positron.hh"
108 #include "G4MuonPlus.hh"
109 #include "G4MuonMinus.hh"
110 #include "G4PionPlus.hh"
111 #include "G4PionMinus.hh"
112 #include "G4KaonPlus.hh"
113 #include "G4KaonMinus.hh"
114 #include "G4Proton.hh"
115 #include "G4AntiProton.hh"
116 #include "G4Deuteron.hh"
117 #include "G4Triton.hh"
118 #include "G4He3.hh"
119 #include "G4Alpha.hh"
120 #include "G4GenericIon.hh"
121 
122 //
123 #include "G4PhysicsListHelper.hh"
124 #include "G4BuilderType.hh"
125 #include "G4EmModelActivator.hh"
126 
127 // factory
129 //
131 
132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
133 
135  : G4VPhysicsConstructor("G4EmPenelope"), verbose(ver)
136 {
138  param->SetDefaults();
139  param->SetVerbose(verbose);
140  param->SetMinEnergy(100*eV);
141  param->SetLowestElectronEnergy(100*eV);
142  param->SetNumberOfBinsPerDecade(20);
143  param->SetMscRangeFactor(0.02);
145  param->SetMuHadLateralDisplacement(true);
146  param->SetFluo(true);
147  param->SetPIXEElectronCrossSectionModel("Penelope");
149 }
150 
151 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
152 
154 {}
155 
156 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
157 
159 {
160  // gamma
161  G4Gamma::Gamma();
162 
163  // leptons
168 
169  // mesons
174 
175  // baryons
178 
179  // ions
182  G4He3::He3();
183  G4Alpha::Alpha();
185 }
186 
187 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
188 
190 {
191  if(verbose > 1) {
192  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
193  }
195 
196  // muon & hadron bremsstrahlung and pair production
205 
206  // muon & hadron multiple scattering
208  mumsc->SetEmModel(new G4WentzelVIModel());
209  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
210 
211  // high energy limit for e+- scattering models
212  G4double highEnergyLimit = G4EmParameters::Instance()->MscEnergyLimit();
213 
214  // nuclear stopping
215  G4NuclearStopping* pnuc = new G4NuclearStopping();
216 
217  //Applicability range for Penelope models
218  //for higher energies, the Standard models are used
219  G4double PenelopeHighEnergyLimit = 1.0*GeV;
220 
221  // Add Penelope EM Processes
223  for(const auto& particleName : partList.PartNames()) {
224  G4ParticleDefinition* particle = table->FindParticle(particleName);
225  if (!particle) { continue; }
226  if (particleName == "gamma") {
227 
228  //Photo-electric effect
229  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
230  G4PenelopePhotoElectricModel* thePEPenelopeModel = new
232  thePEPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
233  thePhotoElectricEffect->SetEmModel(thePEPenelopeModel);
234  ph->RegisterProcess(thePhotoElectricEffect, particle);
235 
236  //Compton scattering
237  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
238  G4PenelopeComptonModel* theComptonPenelopeModel =
240  theComptonPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
241  theComptonScattering->SetEmModel(theComptonPenelopeModel);
242  ph->RegisterProcess(theComptonScattering, particle);
243 
244  //Gamma conversion
245  G4GammaConversion* theGammaConversion = new G4GammaConversion();
246  G4PenelopeGammaConversionModel* theGCPenelopeModel =
248  theGammaConversion->SetEmModel(theGCPenelopeModel);
249  ph->RegisterProcess(theGammaConversion, particle);
250 
251  //Rayleigh scattering
252  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
253  G4PenelopeRayleighModel* theRayleighPenelopeModel =
255  //theRayleighPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
256  theRayleigh->SetEmModel(theRayleighPenelopeModel);
257  ph->RegisterProcess(theRayleigh, particle);
258 
259  } else if (particleName == "e-") {
260 
261  // multiple scattering
263  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
264  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
265  msc1->SetHighEnergyLimit(highEnergyLimit);
266  msc2->SetLowEnergyLimit(highEnergyLimit);
267  msc->SetEmModel(msc1);
268  msc->SetEmModel(msc2);
269 
272  ss->SetEmModel(ssm);
273  ss->SetMinKinEnergy(highEnergyLimit);
274  ssm->SetLowEnergyLimit(highEnergyLimit);
275  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
276 
277  //Ionisation
278  G4eIonisation* eIoni = new G4eIonisation();
279  G4PenelopeIonisationModel* theIoniPenelope =
281  theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
282  eIoni->AddEmModel(0,theIoniPenelope,new G4UniversalFluctuation());
283  eIoni->SetStepFunction(0.2, 100*um); //
284 
285  //Bremsstrahlung
286  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
287  G4PenelopeBremsstrahlungModel* theBremPenelope = new
289  theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
290  eBrem->SetEmModel(theBremPenelope);
291 
292  // register processes
293  ph->RegisterProcess(msc, particle);
294  ph->RegisterProcess(eIoni, particle);
295  ph->RegisterProcess(eBrem, particle);
296  ph->RegisterProcess(ss, particle);
297 
298  } else if (particleName == "e+") {
299 
300  // multiple scattering
302  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
303  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
304  msc1->SetHighEnergyLimit(highEnergyLimit);
305  msc2->SetLowEnergyLimit(highEnergyLimit);
306  msc->SetEmModel(msc1);
307  msc->SetEmModel(msc2);
308 
311  ss->SetEmModel(ssm);
312  ss->SetMinKinEnergy(highEnergyLimit);
313  ssm->SetLowEnergyLimit(highEnergyLimit);
314  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
315 
316  //Ionisation
317  G4eIonisation* eIoni = new G4eIonisation();
318  G4PenelopeIonisationModel* theIoniPenelope =
320  theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
321  eIoni->AddEmModel(0,theIoniPenelope,new G4UniversalFluctuation());
322  eIoni->SetStepFunction(0.2, 100*um); //
323 
324  //Bremsstrahlung
325  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
326  G4PenelopeBremsstrahlungModel* theBremPenelope = new
328  theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
329  eBrem->SetEmModel(theBremPenelope);
330 
331  //Annihilation
333  G4PenelopeAnnihilationModel* theAnnPenelope = new
335  theAnnPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
336  eAnni->AddEmModel(0, theAnnPenelope);
337 
338  // register processes
339  ph->RegisterProcess(msc, particle);
340  ph->RegisterProcess(eIoni, particle);
341  ph->RegisterProcess(eBrem, particle);
342  ph->RegisterProcess(eAnni, particle);
343  ph->RegisterProcess(ss, particle);
344 
345  } else if (particleName == "mu+" ||
346  particleName == "mu-" ) {
347 
348  G4MuIonisation* muIoni = new G4MuIonisation();
349  muIoni->SetStepFunction(0.2, 50*um);
350 
351  ph->RegisterProcess(mumsc, particle);
352  ph->RegisterProcess(muIoni, particle);
353  ph->RegisterProcess(mub, particle);
354  ph->RegisterProcess(mup, particle);
355  ph->RegisterProcess(new G4CoulombScattering(), particle);
356 
357  } else if (particleName == "alpha" ||
358  particleName == "He3" ) {
359 
361  G4ionIonisation* ionIoni = new G4ionIonisation();
362  ionIoni->SetStepFunction(0.1, 10*um);
363 
364  ph->RegisterProcess(msc, particle);
365  ph->RegisterProcess(ionIoni, particle);
366  ph->RegisterProcess(pnuc, particle);
367 
368  } else if (particleName == "GenericIon") {
369 
370  G4ionIonisation* ionIoni = new G4ionIonisation();
371  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
372  ionIoni->SetStepFunction(0.1, 1*um);
373 
374  ph->RegisterProcess(hmsc, particle);
375  ph->RegisterProcess(ionIoni, particle);
376  ph->RegisterProcess(pnuc, particle);
377 
378  } else if (particleName == "pi+" ||
379  particleName == "pi-" ) {
380 
382  G4hIonisation* hIoni = new G4hIonisation();
383  hIoni->SetStepFunction(0.2, 50*um);
384 
385  ph->RegisterProcess(pimsc, particle);
386  ph->RegisterProcess(hIoni, particle);
387  ph->RegisterProcess(pib, particle);
388  ph->RegisterProcess(pip, particle);
389 
390  } else if (particleName == "kaon+" ||
391  particleName == "kaon-" ) {
392 
394  G4hIonisation* hIoni = new G4hIonisation();
395  hIoni->SetStepFunction(0.2, 50*um);
396 
397  ph->RegisterProcess(kmsc, particle);
398  ph->RegisterProcess(hIoni, particle);
399  ph->RegisterProcess(kb, particle);
400  ph->RegisterProcess(kp, particle);
401 
402  } else if (particleName == "proton" ||
403  particleName == "anti_proton") {
404 
406  G4hIonisation* hIoni = new G4hIonisation();
407  hIoni->SetStepFunction(0.2, 50*um);
408 
409  ph->RegisterProcess(pmsc, particle);
410  ph->RegisterProcess(hIoni, particle);
411  ph->RegisterProcess(pb, particle);
412  ph->RegisterProcess(pp, particle);
413  ph->RegisterProcess(pnuc, particle);
414 
415  } else if (particleName == "B+" ||
416  particleName == "B-" ||
417  particleName == "D+" ||
418  particleName == "D-" ||
419  particleName == "Ds+" ||
420  particleName == "Ds-" ||
421  particleName == "anti_He3" ||
422  particleName == "anti_alpha" ||
423  particleName == "anti_deuteron" ||
424  particleName == "anti_lambda_c+" ||
425  particleName == "anti_omega-" ||
426  particleName == "anti_sigma_c+" ||
427  particleName == "anti_sigma_c++" ||
428  particleName == "anti_sigma+" ||
429  particleName == "anti_sigma-" ||
430  particleName == "anti_triton" ||
431  particleName == "anti_xi_c+" ||
432  particleName == "anti_xi-" ||
433  particleName == "deuteron" ||
434  particleName == "lambda_c+" ||
435  particleName == "omega-" ||
436  particleName == "sigma_c+" ||
437  particleName == "sigma_c++" ||
438  particleName == "sigma+" ||
439  particleName == "sigma-" ||
440  particleName == "tau+" ||
441  particleName == "tau-" ||
442  particleName == "triton" ||
443  particleName == "xi_c+" ||
444  particleName == "xi-" ) {
445 
446  ph->RegisterProcess(hmsc, particle);
447  ph->RegisterProcess(new G4hIonisation(), particle);
448  }
449  }
450 
451  // Nuclear stopping
452  pnuc->SetMaxKinEnergy(MeV);
453 
454  // Deexcitation
455  //
456  G4VAtomDeexcitation* deexcitation = new G4UAtomicDeexcitation();
458 
460 }
461 
462 //....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
G4EmParticleList partList
void SetStepFunction(G4double v1, G4double v2, G4bool lock=true)
void SetEmModel(G4VEmModel *, G4int index=0)
virtual void ConstructParticle()
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)
void SetPIXEElectronCrossSectionModel(const G4String &)
#define G4endl
Definition: G4ios.hh:61
void SetMuHadLateralDisplacement(G4bool val)
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
void SetVerbose(G4int val)
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlusDefinition()
Definition: G4PionPlus.cc:93
virtual void ConstructProcess()
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)
G4EmPenelopePhysics(G4int ver=1, const G4String &name="")
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 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 *)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
static G4KaonPlus * KaonPlusDefinition()
Definition: G4KaonPlus.cc:108
G4double MscEnergyLimit() const
static constexpr double GeV
Definition: G4SIunits.hh:217
static G4EmParameters * Instance()