Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4EmDNAPhysicsActivator.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: G4EmDNAPhysicsActivator.cc 97499 2016-06-03 12:02:26Z matkara $
27 // add elastic scattering processes of proton, hydrogen, helium, alpha+, alpha++
28 
30 
31 #include "G4EmParameters.hh"
32 #include "G4ParticleDefinition.hh"
33 #include "G4ParticleTable.hh"
34 #include "G4RegionStore.hh"
35 #include "G4Region.hh"
36 #include "G4VEnergyLossProcess.hh"
37 #include "G4LossTableManager.hh"
38 #include "G4EmConfigurator.hh"
39 
40 #include "G4Gamma.hh"
41 #include "G4Electron.hh"
42 #include "G4Positron.hh"
43 #include "G4MuonPlus.hh"
44 #include "G4MuonMinus.hh"
45 #include "G4Proton.hh"
46 #include "G4GenericIon.hh"
47 #include "G4Alpha.hh"
48 
49 #include "G4ProcessManager.hh"
50 #include "G4DummyModel.hh"
51 #include "G4EmProcessSubType.hh"
52 #include "G4PhysicsListHelper.hh"
53 
54 #include "G4BraggModel.hh"
55 #include "G4BraggIonModel.hh"
56 #include "G4BetheBlochModel.hh"
57 #include "G4UrbanMscModel.hh"
58 #include "G4WentzelVIModel.hh"
59 #include "G4MollerBhabhaModel.hh"
60 #include "G4IonFluctuations.hh"
62 #include "G4LowECapture.hh"
63 #include "G4hMultipleScattering.hh"
66 #include "G4hIonisation.hh"
68 
69 // Processes and models for Geant4-DNA
71 
72 #include "G4DNAElastic.hh"
74 //#include "G4DNAScreenedRutherfordElasticModel.hh"
76 #include "G4DNAIonElasticModel.hh"
77 
78 #include "G4DNAExcitation.hh"
79 #include "G4DNAAttachment.hh"
80 #include "G4DNAVibExcitation.hh"
81 #include "G4DNAIonisation.hh"
82 #include "G4DNAChargeDecrease.hh"
83 #include "G4DNAChargeIncrease.hh"
84 
85 #include "G4SystemOfUnits.hh"
86 #include <vector>
87 
88 #include "G4DNAChemistryManager.hh"
92 
93 #include "G4Threading.hh"
94 
95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96 
98  : G4VPhysicsConstructor("G4EmDNAPhysicsActivator"), verbose(ver)
99 {
101 }
102 
103 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
104 
106 {}
107 
109 {
110  return (0 < verbose && G4Threading::IsMasterThread());
111 }
112 
113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
114 
116 {
117 // bosons
118  G4Gamma::Gamma();
119 
120 // leptons
123 
124 // baryons
126 
128  G4Alpha::Alpha();
129 
130  G4DNAGenericIonsManager * genericIonsManager;
131  genericIonsManager=G4DNAGenericIonsManager::Instance();
132  genericIonsManager->GetIon("alpha+");
133  genericIonsManager->GetIon("helium");
134  genericIonsManager->GetIon("hydrogen");
135  //genericIonsManager->GetIon("carbon");
136  //genericIonsManager->GetIon("nitrogen");
137  //genericIonsManager->GetIon("oxygen");
138  //genericIonsManager->GetIon("iron");
139 
140 }
141 
142 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
143 
145 {
146  const std::vector<G4String>& regnamesDNA = theParameters->RegionsDNA();
147  G4int nreg = regnamesDNA.size();
148  if(0 == nreg)
149  {
150  return;
151  }
152  const std::vector<G4String>& typesDNA = theParameters->TypesDNA();
153 
154  if(IsVerbose()) {
155  G4cout << "### G4EmDNAPhysicsActivator::ConstructProcess for " << nreg
156  << " regions; DNA physics type " << typesDNA[0] << G4endl;
157  }
158 
159  // list of particles
161  const G4ParticleDefinition* prot = G4Proton::Proton();
163 
164  G4DNAGenericIonsManager * genericIonsManager =
166  const G4ParticleDefinition* alpha2 = G4Alpha::Alpha();
167  const G4ParticleDefinition* alpha1 = genericIonsManager->GetIon("alpha+");
168  const G4ParticleDefinition* alpha0 = genericIonsManager->GetIon("helium");
169  const G4ParticleDefinition* h0 = genericIonsManager->GetIon("hydrogen");
170 
171  G4ProcessManager* eman = elec->GetProcessManager();
172  G4ProcessManager* pman = prot->GetProcessManager();
173  G4ProcessManager* iman = gion->GetProcessManager();
174  G4ProcessManager* a2man = alpha2->GetProcessManager();
175  G4ProcessManager* a1man = alpha1->GetProcessManager();
176  G4ProcessManager* a0man = alpha0->GetProcessManager();
177  G4ProcessManager* h0man = h0->GetProcessManager();
178 
179  // alpha+ standard processes
181  G4ParticleDefinition* alpha11 = const_cast<G4ParticleDefinition*>(alpha1);
182  ph->RegisterProcess(new G4hMultipleScattering(), alpha11);
183  ph->RegisterProcess(new G4hIonisation(), alpha11);
184 
185  G4bool emsc = HasMsc(eman);
186  G4bool pmsc = HasMsc(pman);
187  G4bool a2msc = HasMsc(a2man);
188  G4bool a1msc = HasMsc(a1man);
189  // G4bool imsc = HasMsc(iman);
190 
191  // processes are defined with dummy models for the world
192  // elastic scatetring
193  G4DNAElastic* theDNAeElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
194  theDNAeElasticProcess->SetEmModel(new G4DummyModel());
195  eman->AddDiscreteProcess(theDNAeElasticProcess);
196 
197  G4DNAElastic* theDNApElasticProcess =
198  new G4DNAElastic("proton_G4DNAElastic");
199  theDNApElasticProcess->SetEmModel(new G4DummyModel());
200  pman->AddDiscreteProcess(theDNApElasticProcess);
201 
202  G4DNAElastic* theDNAa2ElasticProcess =
203  new G4DNAElastic("alpha_G4DNAElastic");
204  theDNAa2ElasticProcess->SetEmModel(new G4DummyModel());
205  a2man->AddDiscreteProcess(theDNAa2ElasticProcess);
206 
207  G4DNAElastic* theDNAa1ElasticProcess =
208  new G4DNAElastic("alpha+_G4DNAElastic");
209  theDNAa1ElasticProcess->SetEmModel(new G4DummyModel());
210  a1man->AddDiscreteProcess(theDNAa1ElasticProcess);
211 
212  G4DNAElastic* theDNAa0ElasticProcess =
213  new G4DNAElastic("helium_G4DNAElastic");
214  theDNAa0ElasticProcess->SetEmModel(new G4DummyModel());
215  a0man->AddDiscreteProcess(theDNAa0ElasticProcess);
216 
217  G4DNAElastic* theDNAh0ElasticProcess =
218  new G4DNAElastic("hydrogen_G4DNAElastic");
219  theDNAh0ElasticProcess->SetEmModel(new G4DummyModel());
220  h0man->AddDiscreteProcess(theDNAh0ElasticProcess);
221 
222  // excitation
223  G4DNAExcitation* theDNAeExcProcess =
224  new G4DNAExcitation("e-_G4DNAExcitation");
225  theDNAeExcProcess->SetEmModel(new G4DummyModel());
226  eman->AddDiscreteProcess(theDNAeExcProcess);
227 
228  G4DNAExcitation* theDNApExcProcess =
229  new G4DNAExcitation("proton_G4DNAExcitation");
230  theDNApExcProcess->SetEmModel(new G4DummyModel());
231  pman->AddDiscreteProcess(theDNApExcProcess);
232 
233  G4DNAExcitation* theDNAa2ExcProcess =
234  new G4DNAExcitation("alpha_G4DNAExcitation");
235  theDNAa2ExcProcess->SetEmModel(new G4DummyModel());
236  a2man->AddDiscreteProcess(theDNAa2ExcProcess);
237 
238  G4DNAExcitation* theDNAa1ExcProcess =
239  new G4DNAExcitation("alpha+_G4DNAExcitation");
240  theDNAa1ExcProcess->SetEmModel(new G4DummyModel());
241  a1man->AddDiscreteProcess(theDNAa1ExcProcess);
242 
243  G4DNAExcitation* theDNAa0ExcProcess =
244  new G4DNAExcitation("helium_G4DNAExcitation");
245  theDNAa0ExcProcess->SetEmModel(new G4DummyModel());
246  a0man->AddDiscreteProcess(theDNAa0ExcProcess);
247 
248  G4DNAExcitation* theDNAh0ExcProcess =
249  new G4DNAExcitation("hydrogen_G4DNAExcitation");
250  theDNAh0ExcProcess->SetEmModel(new G4DummyModel());
251  h0man->AddDiscreteProcess(theDNAh0ExcProcess);
252 
253  // vibration excitation
254  G4DNAVibExcitation* theDNAeVibExcProcess =
255  new G4DNAVibExcitation("e-_G4DNAVibExcitation");
256  theDNAeVibExcProcess->SetEmModel(new G4DummyModel());
257  eman->AddDiscreteProcess(theDNAeVibExcProcess);
258 
259  // ionisation
260  G4DNAIonisation* theDNAeIoniProcess =
261  new G4DNAIonisation("e-_G4DNAIonisation");
262  theDNAeIoniProcess->SetEmModel(new G4DummyModel());
263  eman->AddDiscreteProcess(theDNAeIoniProcess);
264 
265  G4DNAIonisation* theDNApIoniProcess =
266  new G4DNAIonisation("proton_G4DNAIonisation");
267  theDNApIoniProcess->SetEmModel(new G4DummyModel());
268  pman->AddDiscreteProcess(theDNApIoniProcess);
269 
270  G4DNAIonisation* theDNAa2IoniProcess =
271  new G4DNAIonisation("alpha_G4DNAIonisation");
272  theDNAa2IoniProcess->SetEmModel(new G4DummyModel());
273  a2man->AddDiscreteProcess(theDNAa2IoniProcess);
274 
275  G4DNAIonisation* theDNAa1IoniProcess =
276  new G4DNAIonisation("alpha+_G4DNAIonisation");
277  theDNAa1IoniProcess->SetEmModel(new G4DummyModel());
278  a1man->AddDiscreteProcess(theDNAa1IoniProcess);
279 
280  G4DNAIonisation* theDNAa0IoniProcess =
281  new G4DNAIonisation("helium_G4DNAIonisation");
282  theDNAa0IoniProcess->SetEmModel(new G4DummyModel());
283  a0man->AddDiscreteProcess(theDNAa0IoniProcess);
284 
285  G4DNAIonisation* theDNAh0IoniProcess =
286  new G4DNAIonisation("hydrogen_G4DNAIonisation");
287  theDNAh0IoniProcess->SetEmModel(new G4DummyModel());
288  h0man->AddDiscreteProcess(theDNAh0IoniProcess);
289 
290  G4DNAIonisation* theDNAiIoniProcess =
291  new G4DNAIonisation("GenericIon_G4DNAIonisation");
292  theDNAiIoniProcess->SetEmModel(new G4DummyModel());
293  iman->AddDiscreteProcess(theDNAiIoniProcess);
294 
295  // attachment
296  G4DNAAttachment* theDNAAttachProcess =
297  new G4DNAAttachment("e-_G4DNAAttachment");
298  theDNAAttachProcess->SetEmModel(new G4DummyModel());
299  eman->AddDiscreteProcess(theDNAAttachProcess);
300 
301  // charge exchange
302  G4DNAChargeDecrease* theDNApChargeDecreaseProcess =
303  new G4DNAChargeDecrease("proton_G4DNAChargeDecrease");
304  theDNApChargeDecreaseProcess->SetEmModel(new G4DummyModel());
305  pman->AddDiscreteProcess(theDNApChargeDecreaseProcess);
306 
307  G4DNAChargeDecrease* theDNAa2ChargeDecreaseProcess =
308  new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease");
309  theDNAa2ChargeDecreaseProcess->SetEmModel(new G4DummyModel());
310  a2man->AddDiscreteProcess(theDNAa2ChargeDecreaseProcess);
311 
312  G4DNAChargeDecrease* theDNAa1ChargeDecreaseProcess =
313  new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease");
314  theDNAa1ChargeDecreaseProcess->SetEmModel(new G4DummyModel());
315  a1man->AddDiscreteProcess(theDNAa1ChargeDecreaseProcess);
316 
317  G4DNAChargeIncrease* theDNAa1ChargeIncreaseProcess =
318  new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease");
319  theDNAa1ChargeIncreaseProcess->SetEmModel(new G4DummyModel());
320  a1man->AddDiscreteProcess(theDNAa1ChargeIncreaseProcess);
321 
322  G4DNAChargeIncrease* theDNAa0ChargeIncreaseProcess =
323  new G4DNAChargeIncrease("helium_G4DNAChargeIncrease");
324  theDNAa0ChargeIncreaseProcess->SetEmModel(new G4DummyModel());
325  a0man->AddDiscreteProcess(theDNAa0ChargeIncreaseProcess);
326 
327  G4DNAChargeIncrease* theDNAh0ChargeIncreaseProcess =
328  new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease");
329  theDNAh0ChargeIncreaseProcess->SetEmModel(new G4DummyModel());
330  h0man->AddDiscreteProcess(theDNAh0ChargeIncreaseProcess);
331 
332  // limits for DNA model applicability
333  static const G4double elowest= 7.4 * eV;
334  static const G4double elimel = 1 * MeV;
335  static const G4double pminbb = 2 * MeV;
336  static const G4double pmin = 0.1 * keV;
337  static const G4double pmax = 100 * MeV;
338  static const G4double hemin = 1 * keV;
339  static const G4double ionmin = 0.5 * MeV;
340 
341  // low-energy capture
342  G4LowECapture* ecap = nullptr;
344  {
345  // Note: G4DNAElectronSolvation could also be used instead of G4LowECapture
346  ecap = new G4LowECapture(elowest);
347  // ecap->SetVerboseLevel(1);
348  eman->AddDiscreteProcess(ecap);
349  }
350  else
351  {
352  // When chemistry is activated: G4DNAElectronSolvation turns the electron
353  // to a solvated electron, otherwise it kills the electron at the
354  // corresponding high energy limit of the model
355  G4DNAElectronSolvation* solvatation =
356  new G4DNAElectronSolvation("e-_G4DNAElectronSolvation");
357  solvatation->AddEmModel(0, new G4DummyModel());
358  eman->AddDiscreteProcess(solvatation);
359  }
360 
361  G4LowECapture* pcap = new G4LowECapture(pmin);
362  pman->AddDiscreteProcess(pcap);
363  G4LowECapture* icap = new G4LowECapture(ionmin);
364  iman->AddDiscreteProcess(icap);
365  G4LowECapture* a2cap = new G4LowECapture(hemin);
366  a2man->AddDiscreteProcess(a2cap);
367  G4LowECapture* a1cap = new G4LowECapture(hemin);
368  a1man->AddDiscreteProcess(a1cap);
369  G4LowECapture* a0cap = new G4LowECapture(hemin);
370  a0man->AddDiscreteProcess(a0cap);
371  G4LowECapture* h0cap = new G4LowECapture(ionmin);
372  h0man->AddDiscreteProcess(h0cap);
373 
374  // loop over regions
375  for(G4int i = 0; i < nreg; ++i)
376  {
377  G4String reg = regnamesDNA[i];
378  if(IsVerbose())
379  {
380  G4cout << "### DNA models type " << typesDNA[i]
381  << " are activated for G4Region " << reg << G4endl;
382  }
383 
384  // type of DNA physics
385  G4int itype = 0;
386 
387  if(0 == itype) {
388  AddElectronModels0(reg, ecap, emsc, elowest, elimel);
389  AddProtonModels0(reg, pmsc, elimel, pminbb, pmax);
390  AddHeliumModels0(reg, a1msc, a2msc, elimel, pminbb, pmax);
391  AddGenericIonModels0(reg, pminbb);
392  DeactivateNuclearStopping(pman, elimel);
393  DeactivateNuclearStopping(a1man, elimel);
394  DeactivateNuclearStopping(a2man, elimel);
395  }
396  }
398 }
399 
400 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
401 
403  G4LowECapture* ecap, G4bool emsc,
404  G4double elowest, G4double elimel)
405 {
406  G4EmConfigurator* em_config =
408  G4VEmModel* mod;
409 
410  static const G4double elimin = 1 * MeV;
411  static const G4double elimvb = 100 * eV;
412  static const G4double elimat = 13 * eV;
413  static const G4double elim1 = 10 * keV;
414 
415  // for e- 100 MeV is a limit between different msc models
417 
418  if(emsc) {
419  G4UrbanMscModel* msc = new G4UrbanMscModel();
420  msc->SetActivationLowEnergyLimit(elimel);
421  G4double emaxmsc = std::min(100*MeV, emax);
422  em_config->SetExtraEmModel("e-", "msc", msc, reg, 0.0, emaxmsc);
423  } else {
424  mod = new G4eCoulombScatteringModel();
425  mod->SetActivationLowEnergyLimit(elimel);
426  em_config->SetExtraEmModel("e-", "CoulombScat", mod, reg, 0.0, emax);
427  }
428 
429  // cuts and solvation
430  if(ecap) {
431  ecap->AddRegion(reg);
432  } else {
433  // For this condition to work, the chemistry list has to be called
434  // before any standard EM physics list
436  em_config->SetExtraEmModel("e-", "e-_G4DNAElectronSolvation",
437  mod, reg, 0., elowest);
438  }
439 
440  // elastic
441  mod = new G4DNAChampionElasticModel();
442  em_config->SetExtraEmModel("e-", "e-_G4DNAElastic",
443  mod, reg, 0.0, elimel);
444  // ionisation
445  mod = new G4MollerBhabhaModel();
446  mod->SetActivationLowEnergyLimit(elimin);
447  em_config->SetExtraEmModel("e-", "eIoni",
448  mod, reg, 0.0, emax,
449  new G4UniversalFluctuation());
450 
451  mod = new G4DNABornIonisationModel();
452  em_config->SetExtraEmModel("e-", "e-_G4DNAIonisation",
453  mod, reg, elim1, elimin);
454 
456  em_config->SetExtraEmModel("e-", "e-_G4DNAIonisation",
457  mod, reg, 0.0, elim1);
458 
459  // exc
461  em_config->SetExtraEmModel("e-", "e-_G4DNAExcitation",
462  mod, reg, 0.0, elim1);
463 
464  mod = new G4DNABornExcitationModel();
465  em_config->SetExtraEmModel("e-", "e-_G4DNAExcitation",
466  mod, reg, elim1, elimin);
467 
468  mod = new G4DNASancheExcitationModel();
469  em_config->SetExtraEmModel("e-", "e-_G4DNAVibExcitation",
470  mod, reg, 0.0, elimvb);
471 
472  mod = new G4DNAMeltonAttachmentModel();
473  em_config->SetExtraEmModel("e-", "e-_G4DNAAttachment",
474  mod, reg, 0.0, elimat);
475 
476 }
477 
478 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
479 
481  G4bool pmsc, G4double elimel,
482  G4double pminbb, G4double pmax)
483 {
484  G4EmConfigurator* em_config =
486  G4VEmModel* mod;
487 
488  static const G4double gmmax = 500 * keV;
489 
491 
492  // proton
493 
494  // if SS physics list msc process does not exist
495  if(pmsc) {
496  G4WentzelVIModel* msc = new G4WentzelVIModel();
497  msc->SetActivationLowEnergyLimit(elimel);
498  em_config->SetExtraEmModel("proton", "msc", msc, reg, 0.0, emax);
499  }
500  // single scattering always applied
501  mod = new G4eCoulombScatteringModel();
502  mod->SetActivationLowEnergyLimit(elimel);
503  em_config->SetExtraEmModel("proton", "CoulombScat", mod, reg, 0.0, emax);
504 
505  mod = new G4BraggModel();
506  mod->SetActivationLowEnergyLimit(std::min(pminbb, pmax));
507  em_config->SetExtraEmModel("proton", "hIoni",
508  mod, reg, 0.0, pminbb,
509  new G4UniversalFluctuation());
510 
511  mod = new G4BetheBlochModel();
512  mod->SetActivationLowEnergyLimit(pmax);
513  em_config->SetExtraEmModel("proton", "hIoni",
514  mod, reg, pminbb, emax,
515  new G4UniversalFluctuation());
516 
517  mod = new G4DNARuddIonisationModel();
518  em_config->SetExtraEmModel("proton", "proton_G4DNAIonisation",
519  mod, reg, 0.0, gmmax);
520 
521  mod = new G4DNABornIonisationModel();
522  em_config->SetExtraEmModel("proton", "proton_G4DNAIonisation",
523  mod, reg, gmmax, pmax);
524 
526  em_config->SetExtraEmModel("proton", "proton_G4DNAExcitation",
527  mod, reg, 0.0, gmmax);
528 
529  mod = new G4DNABornExcitationModel();
530  em_config->SetExtraEmModel("proton", "proton_G4DNAExcitation",
531  mod, reg, gmmax, pmax);
532 
534  em_config->SetExtraEmModel("proton", "proton_G4DNAChargeDecrease",
535  mod, reg, 0.0, pmax);
536 
537  mod = new G4DNAIonElasticModel();
538  em_config->SetExtraEmModel("proton", "proton_G4DNAElastic",
539  mod, reg, 0.0, elimel);
540 
541  // hydrogen
542  mod = new G4DNARuddIonisationModel();
543  em_config->SetExtraEmModel("hydrogen", "hydrogen_G4DNAIonisation",
544  mod, reg, 0.0, pmax);
545 
547  em_config->SetExtraEmModel("hydrogen", "hydrogen_G4DNAExcitation",
548  mod, reg, 0.0, gmmax);
549 
551  em_config->SetExtraEmModel("hydrogen", "hydrogen_G4DNAChargeIncrease",
552  mod, reg, 0.0, pmax);
553 
554  mod = new G4DNAIonElasticModel();
555  em_config->SetExtraEmModel("hydrogen", "hydrogen_G4DNAElastic",
556  mod, reg, 0.0, elimel);
557 
558 }
559 
560 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
561 
563  G4double pminbb)
564 {
565  G4EmConfigurator* em_config =
567  G4VEmModel* mod;
568 
570  G4double iemax = std::min(10*MeV, emax);
571  //G4double iemin = 100*eV;
572 
573  mod = new G4BraggIonModel();
574  mod->SetActivationLowEnergyLimit(iemax);
575  em_config->SetExtraEmModel("GenericIon", "ionIoni",
576  mod, reg, 0.0, pminbb,
577  new G4IonFluctuations());
578 
579  mod = new G4BetheBlochModel();
580  mod->SetActivationLowEnergyLimit(iemax);
581  em_config->SetExtraEmModel("GenericIon", "ionIoni",
582  mod, reg, pminbb, emax,
583  new G4IonFluctuations());
584 
586  em_config->SetExtraEmModel("GenericIon", "GenericIon_G4DNAIonisation",
587  mod, reg, 0.0, iemax);
588 
589 }
590 
591 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
592 
594  G4bool a1msc,
595  G4bool a2msc, G4double elimel,
596  G4double pminbb, G4double)
597 {
598  G4EmConfigurator* em_config =
600  G4VEmModel* mod;
601 
602  static const G4double hemax = 400 * MeV;
603  static const G4double massRatio = G4Alpha::Alpha()->GetPDGMass()/CLHEP::proton_mass_c2;
604 
606  G4double pminbba = massRatio*pminbb;
607  if(IsVerbose()) {
608  G4cout << "AddHeliumModels0 for <" << reg << "> a1msc: " << a1msc <<" a2msc: " << a2msc
609  << " elimel= " << elimel << " pminbba= " << pminbba << G4endl;
610  }
611  // alpha++
612  if(elimel < emax) {
613  if(a2msc) {
614  G4UrbanMscModel* msc = new G4UrbanMscModel();
615  msc->SetActivationLowEnergyLimit(elimel);
616  em_config->SetExtraEmModel("alpha", "msc", msc, reg, 0.0, emax);
617  } else {
618  mod = new G4IonCoulombScatteringModel();
619  mod->SetActivationLowEnergyLimit(elimel);
620  em_config->SetExtraEmModel("alpha", "CoulombScat", mod, reg, 0.0, emax);
621  }
622  }
623 
624  mod = new G4BraggIonModel();
625  mod->SetActivationLowEnergyLimit(hemax/massRatio);
626  em_config->SetExtraEmModel("alpha", "ionIoni",
627  mod, reg, 0.0, pminbba,
628  new G4IonFluctuations());
629 
630  mod = new G4BetheBlochModel();
631  mod->SetActivationLowEnergyLimit(hemax/massRatio);
632  em_config->SetExtraEmModel("alpha", "ionIoni",
633  mod, reg, pminbba, emax,
634  new G4IonFluctuations());
635 
636  mod = new G4DNARuddIonisationModel();
637  em_config->SetExtraEmModel("alpha", "alpha_G4DNAIonisation",
638  mod, reg, 0.0, hemax);
639 
641  em_config->SetExtraEmModel("alpha", "alpha_G4DNAExcitation",
642  mod, reg, 0.0, hemax);
643 
645  em_config->SetExtraEmModel("alpha", "alpha_G4DNAChargeDecrease",
646  mod, reg, 0.0, hemax);
647 
648  mod = new G4DNAIonElasticModel();
649  em_config->SetExtraEmModel("alpha", "alpha_G4DNAElastic",
650  mod, reg, 0.0, elimel);
651 
652  // ---
653  // alpha+
654  if(elimel < emax) {
655  if(a1msc) {
656  G4UrbanMscModel* msc = new G4UrbanMscModel();
657  msc->SetActivationLowEnergyLimit(elimel);
658  em_config->SetExtraEmModel("alpha+", "msc", msc, reg, 0.0, emax);
659  } else {
660  mod = new G4IonCoulombScatteringModel();
661  mod->SetActivationLowEnergyLimit(elimel);
662  em_config->SetExtraEmModel("alpha+", "CoulombScat", mod, reg, 0.0, emax);
663  }
664  }
665 
666  mod = new G4BraggIonModel();
667  mod->SetActivationLowEnergyLimit(hemax/massRatio);
668  em_config->SetExtraEmModel("alpha+", "hIoni",
669  mod, reg, 0.0, pminbba,
670  new G4IonFluctuations());
671 
672  mod = new G4BetheBlochModel();
673  mod->SetActivationLowEnergyLimit(hemax/massRatio);
674  em_config->SetExtraEmModel("alpha+", "hIoni",
675  mod, reg, pminbba, emax,
676  new G4IonFluctuations());
677 
678  mod = new G4DNARuddIonisationModel();
679  em_config->SetExtraEmModel("alpha+", "alpha+_G4DNAIonisation",
680  mod, reg, 0.0, hemax);
681 
683  em_config->SetExtraEmModel("alpha+", "alpha+_G4DNAExcitation",
684  mod, reg, 0.0, hemax);
685 
687  em_config->SetExtraEmModel("alpha+", "alpha+_G4DNAChargeDecrease",
688  mod, reg, 0.0, hemax);
689 
691  em_config->SetExtraEmModel("alpha+", "alpha+_G4DNAChargeIncrease",
692  mod, reg, 0.0, hemax);
693 
694  mod = new G4DNAIonElasticModel();
695  em_config->SetExtraEmModel("alpha+", "alpha+_G4DNAElastic",
696  mod, reg, 0.0, elimel);
697 
698  // ---
699  // helium
700  mod = new G4DNARuddIonisationModel();
701  em_config->SetExtraEmModel("helium", "helium_G4DNAIonisation",
702  mod, reg, 0.0, hemax);
703 
705  em_config->SetExtraEmModel("helium", "helium_G4DNAExcitation",
706  mod, reg, 0.0, hemax);
707 
709  em_config->SetExtraEmModel("helium", "helium_G4DNAChargeIncrease",
710  mod, reg, 0.0, hemax);
711 
712  mod = new G4DNAIonElasticModel();
713  em_config->SetExtraEmModel("helium", "helium_G4DNAElastic",
714  mod, reg, 0.0, elimel);
715 }
716 
717 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
718 
720  G4double elimel)
721 {
722  G4ProcessVector* pv = pman->GetProcessList();
723  G4int nproc = pman->GetProcessListLength();
724  for(G4int i = 0; i < nproc; ++i) {
725  if(((*pv)[i])->GetProcessSubType() == fNuclearStopping) {
726  G4VEmProcess* proc = static_cast<G4VEmProcess*>((*pv)[i]);
727  if(proc) {
729  mod->SetActivationLowEnergyLimit(elimel);
730  proc->SetEmModel(mod);
731  }
732  break;
733  }
734  }
735 }
736 
737 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
738 
740 {
741  G4bool res = false;
742  G4ProcessVector* pv = pman->GetProcessList();
743  G4int nproc = pman->GetProcessListLength();
744  for(G4int i = 0; i < nproc; ++i)
745  {
746  if(((*pv)[i])->GetProcessSubType() == fMultipleScattering)
747  {
748  res = true;
749  break;
750  }
751  }
752  return res;
753 }
754 
755 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:93
static constexpr double MeV
Definition: G4SIunits.hh:214
G4int GetProcessListLength() const
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
void AddRegion(const G4String &)
static constexpr double keV
Definition: G4SIunits.hh:216
G4TDNAOneStepThermalizationModel< DNA::Penetration::Meesungnoen2002 > G4DNAOneStepThermalizationModel
#define G4endl
Definition: G4ios.hh:61
G4double MaxKinEnergy() const
G4bool IsMasterThread()
Definition: G4Threading.cc:130
G4DNABornExcitationModel1 G4DNABornExcitationModel
const std::vector< G4String > & TypesDNA() const
static const G4double emax
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetExtraEmModel(const G4String &particleName, const G4String &processName, G4VEmModel *, const G4String &regionName="", G4double emin=0.0, G4double emax=DBL_MAX, G4VEmFluctuationModel *fm=nullptr)
G4ProcessVector * GetProcessList() const
G4double GetPDGMass() const
void DeactivateNuclearStopping(G4ProcessManager *, G4double elimel)
void AddElectronModels0(const G4String &region, G4LowECapture *ecap, G4bool emsc, G4double elowest, G4double elimel)
void AddGenericIonModels0(const G4String &region, G4double pminbb)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static constexpr double proton_mass_c2
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
void AddHeliumModels0(const G4String &region, G4bool a1msc, G4bool a2msc, G4double elimel, G4double pminbb, G4double pmax)
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:728
void SetEmModel(G4VEmModel *, G4int index=0)
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4Positron * Positron()
Definition: G4Positron.cc:94
static constexpr double eV
Definition: G4SIunits.hh:215
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
G4bool HasMsc(G4ProcessManager *) const
static G4PhysicsListHelper * GetPhysicsListHelper()
int G4int
Definition: G4Types.hh:78
G4ParticleDefinition * GetIon(const G4String &name)
G4ProcessManager * GetProcessManager() const
#define G4DNABornIonisationModel
static G4DNAGenericIonsManager * Instance(void)
void AddProtonModels0(const G4String &region, G4bool pmsc, G4double elimel, G4double pminbb, G4double pmax)
G4GLOB_DLL std::ostream G4cout
static G4LossTableManager * Instance()
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
static const G4double reg
G4EmConfigurator * EmConfigurator()
const std::vector< G4String > & RegionsDNA() const
static G4EmParameters * Instance()
T min(const T t1, const T t2)
brief Return the smallest of the two arguments