Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4EmDNAPhysics_stationary_option2.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_option2"), 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 const G4String&)
103 : G4VPhysicsConstructor("G4EmDNAPhysics_stationary_option2"), verbose(ver)
104 {
106  param->SetDefaults();
107  param->SetFluo(true);
108  param->SetAuger(true);
109  param->SetAugerCascade(true);
110  param->SetDeexcitationIgnoreCut(true);
111 
113 }
114 
115 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
116 
118 {}
119 
120 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
121 
123 {
124 // bosons
125  G4Gamma::Gamma();
126 
127 // leptons
130 
131 // baryons
133 
135 
136  G4DNAGenericIonsManager * genericIonsManager;
137  genericIonsManager=G4DNAGenericIonsManager::Instance();
138  genericIonsManager->GetIon("alpha++");
139  genericIonsManager->GetIon("alpha+");
140  genericIonsManager->GetIon("helium");
141  genericIonsManager->GetIon("hydrogen");
142 
143 }
144 
145 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
146 
148 {
149  if(verbose > 1) {
150  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
151  }
153 
154  auto myParticleIterator=GetParticleIterator();
155  myParticleIterator->reset();
156  while( (*myParticleIterator)() )
157  {
158  G4ParticleDefinition* particle = myParticleIterator->value();
159  G4String particleName = particle->GetParticleName();
160 
161  if (particleName == "e-") {
162 
163  // *** Elastic scattering (two alternative models available) ***
164 
165  G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
166  theDNAElasticProcess->SetEmModel(new G4DNAChampionElasticModel());
167 
168  // or alternative model
169  //theDNAElasticProcess
170  //->SetEmModel(new G4DNAScreenedRutherfordElasticModel());
171 
172  ph->RegisterProcess(theDNAElasticProcess, particle);
173 
174  // *** Excitation ***
175 
176  G4DNAExcitation* theDNAExcitationProcess =
177  new G4DNAExcitation("e-_G4DNAExcitation");
178  theDNAExcitationProcess->SetEmModel(new G4DNABornExcitationModel());
179  ((G4DNABornExcitationModel*)(theDNAExcitationProcess->EmModel()))
180  ->SelectStationary(true);
181  ph->RegisterProcess(theDNAExcitationProcess, particle);
182 
183  // *** Ionisation ***
184 
185  G4DNAIonisation* theDNAIonisationProcess =
186  new G4DNAIonisation("e-_G4DNAIonisation");
187  theDNAIonisationProcess->SetEmModel(new G4DNABornIonisationModel());
188  ((G4DNABornIonisationModel*)(theDNAIonisationProcess->EmModel()))
189  ->SelectStationary(true);
190  //
191  ((G4DNABornIonisationModel*)(theDNAIonisationProcess->EmModel()))
192  ->SelectFasterComputation(true);
193  //
194  ph->RegisterProcess(theDNAIonisationProcess, particle);
195 
196  // *** Vibrational excitation ***
197 
198  G4DNAVibExcitation* theDNAVibExcitationProcess =
199  new G4DNAVibExcitation("e-_G4DNAVibExcitation");
200  theDNAVibExcitationProcess->SetEmModel(new G4DNASancheExcitationModel());
201  ((G4DNASancheExcitationModel*)(theDNAVibExcitationProcess->EmModel()))
202  ->SelectStationary(true);
203  ph->RegisterProcess(theDNAVibExcitationProcess, particle);
204 
205  // *** Attachment ***
206 
207  G4DNAAttachment* theDNAAttachmentProcess =
208  new G4DNAAttachment("e-_G4DNAAttachment");
209  theDNAAttachmentProcess->SetEmModel(new G4DNAMeltonAttachmentModel());
210  ((G4DNAMeltonAttachmentModel*)(theDNAAttachmentProcess->EmModel()))
211  ->SelectStationary(true);
212  ph->RegisterProcess(theDNAAttachmentProcess, particle);
213 
214  } else if ( particleName == "proton" ) {
215 
216  // *** Elastic ***
217 
218  G4DNAElastic* theDNAElasticProcess =
219  new G4DNAElastic("proton_G4DNAElastic");
220  theDNAElasticProcess->SetEmModel(new G4DNAIonElasticModel());
221  ((G4DNAIonElasticModel*)(theDNAElasticProcess->EmModel()))
222  ->SelectStationary(true);
223  ph->RegisterProcess(theDNAElasticProcess, particle);
224 
225  // *** Excitation ***
226 
227  G4DNAExcitation* theDNAExcitationProcess =
228  new G4DNAExcitation("proton_G4DNAExcitation");
229 
230  theDNAExcitationProcess->SetEmModel
232  theDNAExcitationProcess->SetEmModel
233  (new G4DNABornExcitationModel());
234 
235  ((G4DNAMillerGreenExcitationModel*)(theDNAExcitationProcess->EmModel()))
236  ->SetLowEnergyLimit(10*eV);
237  ((G4DNAMillerGreenExcitationModel*)(theDNAExcitationProcess->EmModel()))
238  ->SetHighEnergyLimit(500*keV);
239  ((G4DNAMillerGreenExcitationModel*)(theDNAExcitationProcess->EmModel()))
240  ->SelectStationary(true);
241 
242  ((G4DNABornExcitationModel*)(theDNAExcitationProcess->EmModel(1)))
243  ->SetLowEnergyLimit(500*keV);
244  ((G4DNABornExcitationModel*)(theDNAExcitationProcess->EmModel(1)))
245  ->SetHighEnergyLimit(100*MeV);
246  ((G4DNABornExcitationModel*)(theDNAExcitationProcess->EmModel(1)))
247  ->SelectStationary(true);
248 
249  ph->RegisterProcess(theDNAExcitationProcess, particle);
250 
251  // *** Ionisation ***
252 
253  G4DNAIonisation* theDNAIonisationProcess =
254  new G4DNAIonisation("proton_G4DNAIonisation");
255 
256  theDNAIonisationProcess->SetEmModel(
258  theDNAIonisationProcess->SetEmModel(
260 
262  (theDNAIonisationProcess->EmModel()))->SetLowEnergyLimit(0*eV);
264  (theDNAIonisationProcess->EmModel()))->SetHighEnergyLimit(500*keV);
266  (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
267 
269  (theDNAIonisationProcess->EmModel(1)))->SetLowEnergyLimit(500*keV);
271  (theDNAIonisationProcess->EmModel(1)))->SetHighEnergyLimit(100*MeV);
273  (theDNAIonisationProcess->EmModel(1)))->SelectStationary(true);
274  //
276  (theDNAIonisationProcess->EmModel(1)))->SelectFasterComputation(true);
277  //
278 
279  ph->RegisterProcess(theDNAIonisationProcess, particle);
280 
281  // *** Charge decrease ***
282 
283  G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
284  new G4DNAChargeDecrease("proton_G4DNAChargeDecrease");
285  theDNAChargeDecreaseProcess->SetEmModel(
288  (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
289  ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
290 
291  } else if ( particleName == "hydrogen" ) {
292 
293  // *** Elastic ***
294 
295  G4DNAElastic* theDNAElasticProcess =
296  new G4DNAElastic("hydrogen_G4DNAElastic");
297  theDNAElasticProcess->SetEmModel(new G4DNAIonElasticModel());
298  ((G4DNAIonElasticModel*)(theDNAElasticProcess->EmModel()))
299  ->SelectStationary(true);
300  ph->RegisterProcess(theDNAElasticProcess, particle);
301 
302  // *** Excitation ***
303 
304  G4DNAExcitation* theDNAExcitationProcess =
305  new G4DNAExcitation("hydrogen_G4DNAExcitation");
306  theDNAExcitationProcess->SetEmModel(
308  ((G4DNAMillerGreenExcitationModel*)(theDNAExcitationProcess->EmModel()))
309  ->SelectStationary(true);
310  ph->RegisterProcess(theDNAExcitationProcess, particle);
311 
312  // *** Ionisation ***
313 
314  G4DNAIonisation* theDNAIonisationProcess =
315  new G4DNAIonisation("hydrogen_G4DNAIonisation");
316  theDNAIonisationProcess->SetEmModel(
319  (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
320  ph->RegisterProcess(theDNAIonisationProcess, particle);
321 
322  // *** Charge increase ***
323 
324  G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
325  new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease");
326  theDNAChargeIncreaseProcess->SetEmModel(
329  (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
330  ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
331 
332  } else if ( particleName == "alpha" ) {
333 
334  // *** Elastic ***
335 
336  G4DNAElastic* theDNAElasticProcess =
337  new G4DNAElastic("alpha_G4DNAElastic");
338  theDNAElasticProcess->SetEmModel(new G4DNAIonElasticModel());
340  (theDNAElasticProcess->EmModel()))->SelectStationary(true);
341  ph->RegisterProcess(theDNAElasticProcess, particle);
342 
343  // *** Excitation ***
344 
345  G4DNAExcitation* theDNAExcitationProcess =
346  new G4DNAExcitation("alpha_G4DNAExcitation");
347  theDNAExcitationProcess->SetEmModel(
350  (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
351  ph->RegisterProcess(theDNAExcitationProcess, particle);
352 
353  // *** Ionisation ***
354 
355  G4DNAIonisation* theDNAIonisationProcess =
356  new G4DNAIonisation("alpha_G4DNAIonisation");
357  theDNAIonisationProcess->SetEmModel(
360  (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
361  ph->RegisterProcess(theDNAIonisationProcess, particle);
362 
363  // *** Charge decrease ***
364 
365  G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
366  new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease");
367  theDNAChargeDecreaseProcess->SetEmModel(
370  (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
371  ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
372 
373  } else if ( particleName == "alpha+" ) {
374 
375  // *** Elastic ***
376 
377  G4DNAElastic* theDNAElasticProcess =
378  new G4DNAElastic("alpha+_G4DNAElastic");
379  theDNAElasticProcess->SetEmModel(new G4DNAIonElasticModel());
381  (theDNAElasticProcess->EmModel()))->SelectStationary(true);
382  ph->RegisterProcess(theDNAElasticProcess, particle);
383 
384  // *** Excitation ***
385 
386  G4DNAExcitation* theDNAExcitationProcess =
387  new G4DNAExcitation("alpha+_G4DNAExcitation");
388  theDNAExcitationProcess->SetEmModel(
391  (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
392  ph->RegisterProcess(theDNAExcitationProcess, particle);
393 
394  // *** Ionisation ***
395 
396  G4DNAIonisation* theDNAIonisationProcess =
397  new G4DNAIonisation("alpha+_G4DNAIonisation");
398  theDNAIonisationProcess->SetEmModel(
401  (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
402  ph->RegisterProcess(theDNAIonisationProcess, particle);
403 
404  // *** Charge decrease ***
405 
406  G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
407  new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease");
408  theDNAChargeDecreaseProcess->SetEmModel(
411  (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
412  ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
413 
414  // *** Charge increase ***
415 
416  G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
417  new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease");
418  theDNAChargeIncreaseProcess->SetEmModel(
421  (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
422  ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
423 
424  } else if ( particleName == "helium" ) {
425 
426  // *** Elastic ***
427 
428  G4DNAElastic* theDNAElasticProcess =
429  new G4DNAElastic("helium_G4DNAElastic");
430  theDNAElasticProcess->SetEmModel(
431  new G4DNAIonElasticModel());
433  (theDNAElasticProcess->EmModel()))->SelectStationary(true);
434  ph->RegisterProcess(theDNAElasticProcess, particle);
435 
436  // *** Excitation ***
437 
438  G4DNAExcitation* theDNAExcitationProcess =
439  new G4DNAExcitation("helium_G4DNAExcitation");
440  theDNAExcitationProcess->SetEmModel(
443  (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
444  ph->RegisterProcess(theDNAExcitationProcess, particle);
445 
446  // *** Ionisation ***
447 
448  G4DNAIonisation* theDNAIonisationProcess =
449  new G4DNAIonisation("helium_G4DNAIonisation");
450  theDNAIonisationProcess->SetEmModel(
453  (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
454  ph->RegisterProcess(theDNAIonisationProcess, particle);
455 
456  // *** Charge increase ***
457 
458  G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
459  new G4DNAChargeIncrease("helium_G4DNAChargeIncrease");
460  theDNAChargeIncreaseProcess->SetEmModel(
463  (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
464  ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
465 
466  } else if ( particleName == "GenericIon" ) {
467 
468  // *** Ionisation ***
469 
470  G4DNAIonisation* theDNAIonisationProcess =
471  new G4DNAIonisation("GenericIon_G4DNAIonisation");
472  theDNAIonisationProcess->SetEmModel(
475  (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
476  ph->RegisterProcess(theDNAIonisationProcess, particle);
477 
478  }
479 
480  // Warning : the following particles and processes are needed by EM Physics
481  // builders
482  // They are taken from the default Livermore Physics list
483  // These particles are currently not handled by Geant4-DNA
484 
485  // e+
486 
487  else if (particleName == "e+") {
488 
489  // Identical to G4EmStandardPhysics_stationary
490 
493  G4eIonisation* eIoni = new G4eIonisation();
494  eIoni->SetStepFunction(0.2, 100*um);
495 
496  ph->RegisterProcess(msc, particle);
497  ph->RegisterProcess(eIoni, particle);
498  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
499  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
500 
501  } else if (particleName == "gamma") {
502 
503  // photoelectric effect - Livermore model only
504  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
505  thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
506  ph->RegisterProcess(thePhotoElectricEffect, particle);
507 
508  // Compton scattering - Livermore model only
509  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
510  theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
511  ph->RegisterProcess(theComptonScattering, particle);
512 
513  // gamma conversion - Livermore model below 80 GeV
514  G4GammaConversion* theGammaConversion = new G4GammaConversion();
515  theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
516  ph->RegisterProcess(theGammaConversion, particle);
517 
518  // default Rayleigh scattering is Livermore
519  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
520  ph->RegisterProcess(theRayleigh, particle);
521  }
522 
523  // Warning : end of particles and processes are needed by EM Physics build.
524 
525  }
526 
527  // Deexcitation
528  //
531  de->SetFluo(true);
532 }
533 
#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()