Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4EmDNAPhysics_stationary.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 
27 
28 #include "G4SystemOfUnits.hh"
29 
31 
32 // *** Processes and models for Geant4-DNA
33 
34 #include "G4DNAElastic.hh"
37 #include "G4DNAIonElasticModel.hh"
38 
39 #include "G4DNAExcitation.hh"
40 #include "G4DNAAttachment.hh"
41 #include "G4DNAVibExcitation.hh"
42 #include "G4DNAIonisation.hh"
43 #include "G4DNAChargeDecrease.hh"
44 #include "G4DNAChargeIncrease.hh"
45 
46 // particles
47 
48 #include "G4Electron.hh"
49 #include "G4Proton.hh"
50 #include "G4GenericIon.hh"
51 
52 // Warning : the following is needed in order to use EM Physics builders
53 // e+
54 #include "G4Positron.hh"
55 #include "G4eMultipleScattering.hh"
56 #include "G4eIonisation.hh"
57 #include "G4eBremsstrahlung.hh"
58 #include "G4eplusAnnihilation.hh"
59 // gamma
60 #include "G4Gamma.hh"
61 #include "G4PhotoElectricEffect.hh"
63 #include "G4ComptonScattering.hh"
65 #include "G4GammaConversion.hh"
67 #include "G4RayleighScattering.hh"
69 
70 #include "G4EmParameters.hh"
71 // end of warning
72 
73 #include "G4LossTableManager.hh"
74 #include "G4UAtomicDeexcitation.hh"
75 #include "G4PhysicsListHelper.hh"
76 #include "G4BuilderType.hh"
77 
78 // factory
80 //
82 
83 
84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
85 
87  : G4VPhysicsConstructor("G4EmDNAPhysics_stationary"), verbose(ver)
88 {
90  param->SetDefaults();
91  param->SetFluo(true);
92  param->SetAuger(true);
93  param->SetAugerCascade(true);
94  param->SetDeexcitationIgnoreCut(true);
95 
97 }
98 
99 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
100 
102  : G4VPhysicsConstructor("G4EmDNAPhysics_stationary"), verbose(ver)
103 {
105  param->SetDefaults();
106  param->SetFluo(true);
107  param->SetAuger(true);
108  param->SetAugerCascade(true);
109  param->SetDeexcitationIgnoreCut(true);
110 
112 }
113 
114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
115 
117 {}
118 
119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
120 
122 {
123 // bosons
124  G4Gamma::Gamma();
125 
126 // leptons
129 
130 // baryons
132 
134 
135  G4DNAGenericIonsManager * genericIonsManager;
136  genericIonsManager=G4DNAGenericIonsManager::Instance();
137  genericIonsManager->GetIon("alpha++");
138  genericIonsManager->GetIon("alpha+");
139  genericIonsManager->GetIon("helium");
140  genericIonsManager->GetIon("hydrogen");
141 
142 }
143 
144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145 
147 {
148  if(verbose > 1) {
149  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
150  }
152 
153  auto myParticleIterator=GetParticleIterator();
154  myParticleIterator->reset();
155  while( (*myParticleIterator)() )
156  {
157  G4ParticleDefinition* particle = myParticleIterator->value();
158  G4String particleName = particle->GetParticleName();
159 
160  if (particleName == "e-") {
161 
162  // *** Elastic scattering (two alternative models available) ***
163 
164  G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
165  theDNAElasticProcess->SetEmModel(new G4DNAChampionElasticModel());
166 
167  // or alternative model
168  //theDNAElasticProcess
169  //->SetEmModel(new G4DNAScreenedRutherfordElasticModel());
170 
171  ph->RegisterProcess(theDNAElasticProcess, particle);
172 
173  // *** Excitation ***
174 
175  G4DNAExcitation* theDNAExcitationProcess =
176  new G4DNAExcitation("e-_G4DNAExcitation");
178  theDNAExcitationProcess->SetEmModel(modB);
179  modB->SelectStationary(true);
180  ph->RegisterProcess(theDNAExcitationProcess, particle);
181 
182  // *** Ionisation ***
183 
184  G4DNAIonisation* theDNAIonisationProcess =
185  new G4DNAIonisation("e-_G4DNAIonisation");
187  theDNAIonisationProcess->SetEmModel(modI);
188  modI->SelectStationary(true);
189  ph->RegisterProcess(theDNAIonisationProcess, particle);
190 
191  // *** Vibrational excitation ***
192 
193  G4DNAVibExcitation* theDNAVibExcitationProcess =
194  new G4DNAVibExcitation("e-_G4DNAVibExcitation");
196  theDNAVibExcitationProcess->SetEmModel(modS);
197  modS->SelectStationary(true);
198  ph->RegisterProcess(theDNAVibExcitationProcess, particle);
199 
200  // *** Attachment ***
201 
202  G4DNAAttachment* theDNAAttachmentProcess =
203  new G4DNAAttachment("e-_G4DNAAttachment");
205  theDNAAttachmentProcess->SetEmModel(modM);
206  modM->SelectStationary(true);
207  ph->RegisterProcess(theDNAAttachmentProcess, particle);
208 
209  } else if ( particleName == "proton" ) {
210 
211  // *** Elastic ***
212 
213  G4DNAElastic* theDNAElasticProcess =
214  new G4DNAElastic("proton_G4DNAElastic");
215  theDNAElasticProcess->SetEmModel(new G4DNAIonElasticModel());
216  ((G4DNAIonElasticModel*)(theDNAElasticProcess->EmModel()))
217  ->SelectStationary(true);
218  ph->RegisterProcess(theDNAElasticProcess, particle);
219 
220  // *** Excitation ***
221 
222  G4DNAExcitation* theDNAExcitationProcess =
223  new G4DNAExcitation("proton_G4DNAExcitation");
224 
225  theDNAExcitationProcess
227  theDNAExcitationProcess
229 
231  (theDNAExcitationProcess->EmModel(0)))->SetLowEnergyLimit(10*eV);
233  (theDNAExcitationProcess->EmModel(0)))->SetHighEnergyLimit(500*keV);
235  (theDNAExcitationProcess->EmModel(0)))->SelectStationary(true);
236 
238  (theDNAExcitationProcess->EmModel(1)))->SetLowEnergyLimit(500*keV);
240  (theDNAExcitationProcess->EmModel(1)))->SetHighEnergyLimit(100*MeV);
242  (theDNAExcitationProcess->EmModel(1)))->SelectStationary(true);
243 
244  ph->RegisterProcess(theDNAExcitationProcess, particle);
245 
246  // *** Ionisation ***
247 
248  G4DNAIonisation* theDNAIonisationProcess =
249  new G4DNAIonisation("proton_G4DNAIonisation");
250 
251  theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationModel);
252  theDNAIonisationProcess->SetEmModel(new G4DNABornIonisationModel);
253 
254  ((G4DNARuddIonisationModel*)(theDNAIonisationProcess->EmModel(0)))
255  ->SetLowEnergyLimit(0*eV);
256  ((G4DNARuddIonisationModel*)(theDNAIonisationProcess->EmModel(0)))
257  ->SetHighEnergyLimit(500*keV);
258  ((G4DNARuddIonisationModel*)(theDNAIonisationProcess->EmModel(0)))
259  ->SelectStationary(true);
260 
261  ((G4DNABornIonisationModel*)(theDNAIonisationProcess->EmModel(1)))
262  ->SetLowEnergyLimit(500*keV);
263  ((G4DNABornIonisationModel*)(theDNAIonisationProcess->EmModel(1)))
264  ->SetHighEnergyLimit(100*MeV);
265  ((G4DNABornIonisationModel*)(theDNAIonisationProcess->EmModel(1)))
266  ->SelectStationary(true);
267 
268  ph->RegisterProcess(theDNAIonisationProcess, particle);
269 
270  // *** Charge decrease ***
271 
272  G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
273  new G4DNAChargeDecrease("proton_G4DNAChargeDecrease");
274  theDNAChargeDecreaseProcess
277  (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
278  ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
279 
280  } else if ( particleName == "hydrogen" ) {
281 
282  // *** Elastic ***
283 
284  G4DNAElastic* theDNAElasticProcess =
285  new G4DNAElastic("hydrogen_G4DNAElastic");
286  theDNAElasticProcess->SetEmModel(new G4DNAIonElasticModel());
287  ((G4DNAIonElasticModel*)(theDNAElasticProcess->EmModel()))
288  ->SelectStationary(true);
289  ph->RegisterProcess(theDNAElasticProcess, particle);
290 
291  // *** Excitation ***
292 
293  G4DNAExcitation* theDNAExcitationProcess =
294  new G4DNAExcitation("hydrogen_G4DNAExcitation");
295  theDNAExcitationProcess
297  ((G4DNAMillerGreenExcitationModel*)(theDNAExcitationProcess->EmModel()))
298  ->SelectStationary(true);
299  ph->RegisterProcess(theDNAExcitationProcess, particle);
300 
301  // *** Ionisation ***
302 
303  G4DNAIonisation* theDNAIonisationProcess =
304  new G4DNAIonisation("hydrogen_G4DNAIonisation");
305  theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationModel());
306  ((G4DNARuddIonisationModel*)(theDNAIonisationProcess->EmModel()))
307  ->SelectStationary(true);
308  ph->RegisterProcess(theDNAIonisationProcess, particle);
309 
310  // *** Charge increase ***
311 
312  G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
313  new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease");
314  theDNAChargeIncreaseProcess
317  (theDNAChargeIncreaseProcess->EmModel()))
318  ->SelectStationary(true);
319  ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
320 
321  } else if ( particleName == "alpha" ) {
322 
323  // *** Elastic ***
324 
325  G4DNAElastic* theDNAElasticProcess =
326  new G4DNAElastic("alpha_G4DNAElastic");
327  theDNAElasticProcess->SetEmModel(new G4DNAIonElasticModel());
328  ((G4DNAIonElasticModel*)(theDNAElasticProcess->EmModel()))
329  ->SelectStationary(true);
330  ph->RegisterProcess(theDNAElasticProcess, particle);
331 
332  // *** Excitation ***
333 
334  G4DNAExcitation* theDNAExcitationProcess =
335  new G4DNAExcitation("alpha_G4DNAExcitation");
336  theDNAExcitationProcess->SetEmModel
338  ((G4DNAMillerGreenExcitationModel*)(theDNAExcitationProcess->EmModel()))
339  ->SelectStationary(true);
340  ph->RegisterProcess(theDNAExcitationProcess, particle);
341 
342  // *** Ionisation ***
343 
344  G4DNAIonisation* theDNAIonisationProcess =
345  new G4DNAIonisation("alpha_G4DNAIonisation");
346  theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationModel());
347  ((G4DNARuddIonisationModel*)(theDNAIonisationProcess->EmModel()))
348  ->SelectStationary(true);
349  ph->RegisterProcess(theDNAIonisationProcess, particle);
350 
351  // *** Charge decrease ***
352 
353  G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
354  new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease");
355  theDNAChargeDecreaseProcess->SetEmModel
358  (theDNAChargeDecreaseProcess->EmModel()))
359  ->SelectStationary(true);
360  ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
361 
362  } else if ( particleName == "alpha+" ) {
363 
364  // *** Elastic ***
365 
366  G4DNAElastic* theDNAElasticProcess =
367  new G4DNAElastic("alpha+_G4DNAElastic");
368  theDNAElasticProcess->SetEmModel(new G4DNAIonElasticModel());
369  ((G4DNAIonElasticModel*)(theDNAElasticProcess->EmModel()))
370  ->SelectStationary(true);
371  ph->RegisterProcess(theDNAElasticProcess, particle);
372 
373  // *** Excitation ***
374 
375  G4DNAExcitation* theDNAExcitationProcess =
376  new G4DNAExcitation("alpha+_G4DNAExcitation");
377  theDNAExcitationProcess->SetEmModel
380  (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
381  ph->RegisterProcess(theDNAExcitationProcess, particle);
382 
383  // *** Ionisation ***
384 
385  G4DNAIonisation* theDNAIonisationProcess =
386  new G4DNAIonisation("alpha+_G4DNAIonisation");
387  theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationModel());
388  ((G4DNARuddIonisationModel*)(theDNAIonisationProcess->EmModel()))
389  ->SelectStationary(true);
390  ph->RegisterProcess(theDNAIonisationProcess, particle);
391 
392  // *** Charge decrease ***
393 
394  G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
395  new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease");
396  theDNAChargeDecreaseProcess->SetEmModel
399  (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
400  ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
401 
402  // *** Charge increase ***
403 
404  G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
405  new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease");
406  theDNAChargeIncreaseProcess->SetEmModel
409  (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
410  ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
411 
412  } else if ( particleName == "helium" ) {
413 
414  // *** Elastic ***
415 
416  G4DNAElastic* theDNAElasticProcess =
417  new G4DNAElastic("helium_G4DNAElastic");
418  theDNAElasticProcess->SetEmModel
419  (new G4DNAIonElasticModel());
420  ((G4DNAIonElasticModel*)(theDNAElasticProcess->EmModel()))
421  ->SelectStationary(true);
422  ph->RegisterProcess(theDNAElasticProcess, particle);
423 
424  // *** Excitation ***
425 
426  G4DNAExcitation* theDNAExcitationProcess =
427  new G4DNAExcitation("helium_G4DNAExcitation");
428  theDNAExcitationProcess->SetEmModel
430  ((G4DNAMillerGreenExcitationModel*)(theDNAExcitationProcess->EmModel()))
431  ->SelectStationary(true);
432  ph->RegisterProcess(theDNAExcitationProcess, particle);
433 
434  // *** Ionisation ***
435 
436  G4DNAIonisation* theDNAIonisationProcess =
437  new G4DNAIonisation("helium_G4DNAIonisation");
438  theDNAIonisationProcess->SetEmModel
439  (new G4DNARuddIonisationModel());
440  ((G4DNARuddIonisationModel*)(theDNAIonisationProcess->EmModel()))
441  ->SelectStationary(true);
442  ph->RegisterProcess(theDNAIonisationProcess, particle);
443 
444  // *** Charge increase ***
445 
446  G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
447  new G4DNAChargeIncrease("helium_G4DNAChargeIncrease");
448  theDNAChargeIncreaseProcess->SetEmModel
451  (theDNAChargeIncreaseProcess->EmModel()))
452  ->SelectStationary(true);
453  ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
454 
455  } else if ( particleName == "GenericIon" ) {
456 
457  // *** Ionisation ***
458 
459  G4DNAIonisation* theDNAIonisationProcess =
460  new G4DNAIonisation("GenericIon_G4DNAIonisation");
461  theDNAIonisationProcess->SetEmModel(
463  ((G4DNARuddIonisationExtendedModel*)(theDNAIonisationProcess->EmModel()))
464  ->SelectStationary(true);
465  ph->RegisterProcess(theDNAIonisationProcess, particle);
466 
467  }
468 
469  // Warning : the following particles and processes are needed by EM Physics
470  // builders
471  // They are taken from the default Livermore Physics list
472  // These particles are currently not handled by Geant4-DNA
473 
474  // e+
475 
476  else if (particleName == "e+") {
477 
478  // Identical to G4EmStandardPhysics_stationary
479 
482  G4eIonisation* eIoni = new G4eIonisation();
483  eIoni->SetStepFunction(0.2, 100*um);
484 
485  ph->RegisterProcess(msc, particle);
486  ph->RegisterProcess(eIoni, particle);
487  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
488  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
489 
490  } else if (particleName == "gamma") {
491  // photoelectric effect - Livermore model only
492  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
493  thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
494  ph->RegisterProcess(thePhotoElectricEffect, particle);
495 
496  // Compton scattering - Livermore model only
497  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
498  theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
499  ph->RegisterProcess(theComptonScattering, particle);
500 
501  // gamma conversion - Livermore model below 80 GeV
502  G4GammaConversion* theGammaConversion = new G4GammaConversion();
503  theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
504  ph->RegisterProcess(theGammaConversion, particle);
505 
506  // default Rayleigh scattering is Livermore
507  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
508  ph->RegisterProcess(theRayleigh, particle);
509  }
510 
511  // Warning : end of particles and processes are needed by EM Physics build.
512 
513  }
514 
515  // Deexcitation
516  //
519 }
520 
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)
void SetDeexcitationIgnoreCut(G4bool val)
void SetStepFunction(G4double v1, G4double v2, G4bool lock=true)
static constexpr double MeV
Definition: G4SIunits.hh:214
static constexpr double keV
Definition: G4SIunits.hh:216
#define G4endl
Definition: G4ios.hh:61
const G4String & GetParticleName() const
G4VEmModel * EmModel(size_t index=0) const
G4DNABornExcitationModel1 G4DNABornExcitationModel
static G4Proton * Proton()
Definition: G4Proton.cc:93
static constexpr double um
Definition: G4SIunits.hh:113
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetEmModel(G4VEmModel *, G4int index=0)
const G4String & GetPhysicsName() const
void SetFluo(G4bool val)
static G4Positron * Positron()
Definition: G4Positron.cc:94
static constexpr double eV
Definition: G4SIunits.hh:215
void SetAuger(G4bool val)
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
void SetAugerCascade(G4bool val)
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
static G4PhysicsListHelper * GetPhysicsListHelper()
int G4int
Definition: G4Types.hh:78
G4ParticleDefinition * GetIon(const G4String &name)
void SetStepLimitType(G4MscStepLimitType val)
#define G4DNABornIonisationModel
static G4DNAGenericIonsManager * Instance(void)
G4GLOB_DLL std::ostream G4cout
static G4LossTableManager * Instance()
void SetAtomDeexcitation(G4VAtomDeexcitation *)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
static G4EmParameters * Instance()