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