Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4EmDNABuilder.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//
27// Geant4 class G4EmDNABuilder
28//
29// Author V.Ivanchenko 22.05.2020
30//
31
32#include "G4EmDNABuilder.hh"
33#include "G4SystemOfUnits.hh"
34
35// particles
36#include "G4Gamma.hh"
37#include "G4Electron.hh"
38#include "G4Positron.hh"
39#include "G4Proton.hh"
40#include "G4Alpha.hh"
41#include "G4GenericIon.hh"
44
45// utilities
46#include "G4SystemOfUnits.hh"
47#include "G4EmParameters.hh"
48#include "G4EmBuilder.hh"
51#include "G4PhysListUtil.hh"
52#include "G4ProcessManager.hh"
53#include "G4Region.hh"
54
55// standard processes
57#include "G4GammaConversion.hh"
62#include "G4eIonisation.hh"
63#include "G4hIonisation.hh"
64#include "G4ionIonisation.hh"
65#include "G4eBremsstrahlung.hh"
67
68// standard models
72#include "G4UrbanMscModel.hh"
77#include "G4Generator2BS.hh"
78#include "G4BraggModel.hh"
79#include "G4BraggIonModel.hh"
80#include "G4BetheBlochModel.hh"
81#include "G4DummyModel.hh"
82
83// DNA models
104
105static const G4double lowEnergyRPWBA = 100*CLHEP::MeV;
106static const G4double lowEnergyMSC = 1*CLHEP::MeV;
107static const G4double lowEnergyProtonIoni = 2*CLHEP::MeV;
108
109//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
110
112{
113 // standard particles
115
116 // DNA ions
117 G4DNAGenericIonsManager* genericIonsManager
119 genericIonsManager->GetIon("alpha+");
120 genericIonsManager->GetIon("helium");
121 genericIonsManager->GetIon("hydrogen");
122}
123
124//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125
126void
128 const G4double emin_proton,
129 const G4double emin_alpha,
130 const G4double emin_ion,
131 const G4EmDNAMscModelType mscType,
132 const G4bool)
133{
136 const G4double emax = param->MaxKinEnergy();
138
139 // gamma
141
142 // photoelectric effect - Livermore model
143 auto thePEEffect = new G4PhotoElectricEffect();
144 thePEEffect->SetEmModel(new G4LivermorePhotoElectricModel());
145 ph->RegisterProcess(thePEEffect, gamma);
146
147 // Compton scattering - Klein-Nishina
148 auto theComptonScattering = new G4ComptonScattering();
149 theComptonScattering->SetEmModel(new G4KleinNishinaModel());
150 auto cModel = new G4LowEPComptonModel();
151 cModel->SetHighEnergyLimit(20*CLHEP::MeV);
152 theComptonScattering->AddEmModel(0, cModel);
153 ph->RegisterProcess(theComptonScattering, gamma);
154
155 // gamma conversion - 5D model
156 auto theGammaConversion = new G4GammaConversion();
157 ph->RegisterProcess(theGammaConversion, gamma);
158
159 // Rayleigh scattering - Livermore model
160 auto theRayleigh = new G4RayleighScattering();
161 ph->RegisterProcess(theRayleigh, gamma);
162
163 // electron
164 if(emin_elec < emax) {
166 auto msc_el = new G4eMultipleScattering();
167 G4VMscModel* msc_model_el;
168 if(mscType == dnaWVI) {
169 msc_model_el = new G4LowEWentzelVIModel();
170 } else if(mscType == dnaGS) {
171 msc_model_el = new G4GoudsmitSaundersonMscModel();
172 } else {
173 msc_model_el = new G4UrbanMscModel();
174 }
175 msc_model_el->SetActivationLowEnergyLimit(lowEnergyMSC);
176 msc_el->SetEmModel(msc_model_el);
177 ph->RegisterProcess(msc_el, elec);
178
179 auto ioni = new G4eIonisation();
180 auto mb_el = new G4MollerBhabhaModel();
181 mb_el->SetActivationLowEnergyLimit(emin_elec);
182 ioni->SetEmModel(mb_el);
183 ph->RegisterProcess(ioni, elec);
184
185 auto brem = new G4eBremsstrahlung();
186 auto sb_el = new G4SeltzerBergerModel();
187 sb_el->SetActivationLowEnergyLimit(emin_elec);
188 sb_el->SetHighEnergyLimit(emax);
189 sb_el->SetAngularDistribution(new G4Generator2BS());
190 brem->SetEmModel(sb_el);
191 ph->RegisterProcess(brem, elec);
192 }
193
194 // positron
196 auto msc_pos = new G4eMultipleScattering();
197 G4VMscModel* msc_model_pos;
198 if(mscType == dnaWVI) {
199 msc_model_pos = new G4LowEWentzelVIModel();
200 } else if(mscType == dnaGS) {
201 msc_model_pos = new G4GoudsmitSaundersonMscModel();
202 } else {
203 msc_model_pos = new G4UrbanMscModel();
204 }
205 msc_pos->SetEmModel(msc_model_pos);
206 ph->RegisterProcess(msc_pos, posi);
207 ph->RegisterProcess(new G4eIonisation(), posi);
208
209 auto brem = new G4eBremsstrahlung();
210 auto sb = new G4SeltzerBergerModel();
211 sb->SetHighEnergyLimit(emax);
212 sb->SetAngularDistribution(new G4Generator2BS());
213 brem->SetEmModel(sb);
214 ph->RegisterProcess(brem, posi);
215 ph->RegisterProcess(new G4eplusAnnihilation(), posi);
216
217 // proton
218 if(emin_proton < emax) {
220 StandardHadronPhysics(part, lowEnergyMSC, emin_proton, emax,
221 mscType, false);
222 }
223
224 // GenericIon
225 if(emin_ion < emax) {
227 StandardHadronPhysics(ion, lowEnergyMSC, emin_ion, emax,
228 dnaUrban, true);
229 }
230
231 // alpha
232 if(emin_alpha < emax) {
234 StandardHadronPhysics(part, lowEnergyMSC, emin_alpha, emax,
235 dnaUrban, true);
236
237 // alpha+
238 G4DNAGenericIonsManager* genericIonsManager
240 part = genericIonsManager->GetIon("alpha+");
241 StandardHadronPhysics(part, lowEnergyMSC, emin_alpha, emax,
242 dnaUrban, false);
243 }
244 // list of main standard particles
245 const std::vector<G4int> chargedParticles = {
246 13, -13, 211, -211, 321, -321, -2212,
247 1000010020, 1000010030, 1000020030
248 };
249 auto msc = new G4hMultipleScattering();
250 msc->SetEmModel(new G4WentzelVIModel());
251 G4EmBuilder::ConstructBasicEmPhysics(msc, chargedParticles);
252}
253
254//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
255
256void G4EmDNABuilder::StandardHadronPhysics(G4ParticleDefinition* part,
257 const G4double lowELimitForMSC,
258 const G4double lowELimitForIoni,
259 const G4double maxEnergy,
260 const G4EmDNAMscModelType mscType,
261 const G4bool isIon)
262{
265 G4VMscModel* msc_model = nullptr;
266 if(mscType == dnaWVI) {
267 msc_model = new G4LowEWentzelVIModel();
268 } else {
269 msc_model = new G4UrbanMscModel();
270 }
271 msc_model->SetActivationLowEnergyLimit(lowELimitForMSC);
272 msc_model->SetLowEnergyLimit(lowELimitForMSC);
273 msc_model->SetHighEnergyLimit(maxEnergy);
274 msc->SetEmModel(msc_model);
275 ph->RegisterProcess(msc, part);
276
277 G4VEnergyLossProcess* ioni = nullptr;
278 G4VEmModel* mod1 = nullptr;
279 if(isIon) {
280 ioni = new G4ionIonisation();
281 mod1 = new G4BraggIonModel();
282 } else {
283 ioni = new G4hIonisation();
284 mod1 = new G4BraggModel();
285 }
286 G4double eth = lowEnergyProtonIoni*part->GetPDGMass()/CLHEP::proton_mass_c2;
287 mod1->SetActivationLowEnergyLimit(lowELimitForIoni);
288 mod1->SetHighEnergyLimit(eth);
289 ioni->SetEmModel(mod1);
290
291 G4VEmModel* mod2 = new G4BetheBlochModel();
292 mod2->SetActivationLowEnergyLimit(lowELimitForIoni);
293 mod2->SetLowEnergyLimit(eth);
294 ioni->SetEmModel(mod2);
295
296 ph->RegisterProcess(ioni, part);
297}
298
299//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
300
301void
303 const G4int opt,
304 const G4bool fast,
305 const G4bool stationary,
306 const G4Region* reg)
307{
309
310 // limit of the Emfietzoglou models
311 G4double emaxE = 0.0;
312 // limit of the elastic and solvation models
313 G4double emaxT = 7.4*CLHEP::eV;
314 // limit for CPA100 models
315 G4double emaxCPA100 = 250*CLHEP::keV;
316 if(4 == opt) {
317 emaxE = 10.*CLHEP::keV;
318 emaxT = 10.*CLHEP::eV;
319 } else if(5 < opt) {
320 emaxT = 11.*CLHEP::eV;
321 }
322
323 // *** Solvation ***
326 therm->SetHighEnergyLimit(emaxT);
327 pSolvation->AddEmModel(-1, therm, reg);
328
329 // *** Elastic scattering ***
330 auto pElasticProcess = FindOrBuildElastic(part, "e-_G4DNAElastic");
331 G4VEmModel* elast = nullptr;
332 G4VEmModel* elast2 = nullptr;
333 if(4 == opt) {
335 } else if(5 < opt) {
336 auto mod = new G4DNACPA100ElasticModel();
337 mod->SelectStationary(stationary);
338 elast = mod;
339 elast2 = new G4DNAChampionElasticModel();
340 } else {
341 elast = new G4DNAChampionElasticModel();
342 }
343 elast->SetHighEnergyLimit(lowEnergyMSC);
344 pElasticProcess->AddEmModel(-2, elast, reg);
345
346 if(nullptr != elast2) {
347 elast->SetHighEnergyLimit(emaxCPA100);
348 elast2->SetLowEnergyLimit(emaxCPA100);
349 elast2->SetHighEnergyLimit(lowEnergyMSC);
350 pElasticProcess->AddEmModel(-3, elast2, reg);
351 }
352
353 // *** Excitation ***
354 auto theDNAExc = FindOrBuildExcitation(part, "e-_G4DNAExcitation");
355 if(emaxE > 0.0) {
356 auto modE = new G4DNAEmfietzoglouExcitationModel();
357 theDNAExc->AddEmModel(-1, modE, reg);
358 modE->SelectStationary(stationary);
359 modE->SetHighEnergyLimit(emaxE);
360 }
361 G4VEmModel* modB = nullptr;
362 G4VEmModel* modB2 = nullptr;
363 if(6 == opt) {
364 auto mod = new G4DNACPA100ExcitationModel();
365 mod->SelectStationary(stationary);
366 modB = mod;
367 auto mod1 = new G4DNABornExcitationModel();
368 mod1->SelectStationary(stationary);
369 modB2 = mod1;
370 } else {
371 auto mod = new G4DNABornExcitationModel();
372 mod->SelectStationary(stationary);
373 modB = mod;
374 }
375 modB->SetLowEnergyLimit(emaxE);
376 modB->SetHighEnergyLimit(emaxDNA);
377 theDNAExc->AddEmModel(-2, modB, reg);
378 if(nullptr != modB2) {
379 modB->SetHighEnergyLimit(emaxCPA100);
380 modB2->SetLowEnergyLimit(emaxCPA100);
381 modB2->SetHighEnergyLimit(emaxDNA);
382 theDNAExc->AddEmModel(-3, modB2, reg);
383 }
384
385 // *** Ionisation ***
386 auto theDNAIoni = FindOrBuildIonisation(part, "e-_G4DNAIonisation");
387 if(emaxE > 0.0) {
388 auto modE = new G4DNAEmfietzoglouIonisationModel();
389 theDNAIoni->AddEmModel(-1, modE, reg);
390 modE->SelectFasterComputation(fast);
391 modE->SelectStationary(stationary);
392 modE->SetHighEnergyLimit(emaxE);
393 }
394 G4VEmModel* modI = nullptr;
395 G4VEmModel* modI2 = nullptr;
396 if(6 == opt) {
397 auto mod = new G4DNACPA100IonisationModel();
398 mod->SelectStationary(stationary);
399 mod->SelectFasterComputation(fast);
400 modI = mod;
401 auto mod1 = new G4DNABornIonisationModel();
402 mod1->SelectStationary(stationary);
403 modI2 = mod1;
404 } else {
405 auto mod = new G4DNABornIonisationModel1();
406 mod->SelectStationary(stationary);
407 mod->SelectFasterComputation(fast);
408 modI = mod;
409 }
410 modI->SetLowEnergyLimit(emaxE);
411 modI->SetHighEnergyLimit(emaxDNA);
412 theDNAIoni->AddEmModel(-2, modI, reg);
413 if(nullptr != modI2) {
414 modI->SetHighEnergyLimit(emaxCPA100);
415 modI2->SetLowEnergyLimit(emaxCPA100);
416 modI2->SetHighEnergyLimit(emaxDNA);
417 theDNAIoni->AddEmModel(-3, modI2, reg);
418 }
419
420 if(4 != opt && 6 != opt) {
421 // *** Vibrational excitation ***
422 auto theDNAVibExc = FindOrBuildVibExcitation(part, "e-_G4DNAVibExcitation");
423 auto modS = new G4DNASancheExcitationModel();
424 theDNAVibExc->AddEmModel(-1, modS, reg);
425 modS->SelectStationary(stationary);
426
427 // *** Attachment ***
428 auto theDNAAttach = FindOrBuildAttachment(part, "e-_G4DNAAttachment");
429 auto modM = new G4DNAMeltonAttachmentModel();
430 theDNAAttach->AddEmModel(-1, modM, reg);
431 modM->SelectStationary(stationary);
432 }
433}
434
435//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
436
437void
439 const G4double emaxIonDNA,
440 const G4int opt,
441 const G4bool fast,
442 const G4bool stationary,
443 const G4Region* reg)
444{
446 const G4double emax = param->MaxKinEnergy();
448
449 // *** Elastic scattering ***
450 auto pElasticProcess = FindOrBuildElastic(part, "proton_G4DNAElastic");
451 auto modE = new G4DNAIonElasticModel();
452 modE->SetHighEnergyLimit(lowEnergyMSC);
453 modE->SelectStationary(stationary);
454 pElasticProcess->AddEmModel(-1, modE, reg);
455
456 // *** Excitation ***
457 G4double e2DNA = std::min(e1DNA, lowEnergyRPWBA);
458 auto theDNAExc = FindOrBuildExcitation(part, "proton_G4DNAExcitation");
459 auto modMGE = new G4DNAMillerGreenExcitationModel();
460 modMGE->SetHighEnergyLimit(e2DNA);
461 modMGE->SelectStationary(stationary);
462 theDNAExc->AddEmModel(-1, modMGE, reg);
463
464 if(e2DNA < lowEnergyRPWBA) {
465 auto modB = new G4DNABornExcitationModel();
466 modB->SelectStationary(stationary);
467 modB->SetLowEnergyLimit(e2DNA);
468 modB->SetHighEnergyLimit(lowEnergyRPWBA);
469 theDNAExc->AddEmModel(-2, modB, reg);
470 }
471 if(lowEnergyRPWBA < emaxIonDNA) {
472 auto modC = new G4DNARPWBAExcitationModel();
473 modC->SelectStationary(stationary);
474 modC->SetLowEnergyLimit(lowEnergyRPWBA);
475 modC->SetHighEnergyLimit(emaxIonDNA);
476 theDNAExc->AddEmModel(-3, modC, reg);
477 }
478
479 // *** Ionisation ***
480 auto theDNAIoni = FindOrBuildIonisation(part, "proton_G4DNAIonisation");
481 G4VEmModel* modRI = nullptr;
482 if(2 == opt) {
483 auto mod = new G4DNARuddIonisationExtendedModel();
484 mod->SelectStationary(stationary);
485 modRI = mod;
486 } else {
487 auto mod = new G4DNARuddIonisationModel();
488 mod->SelectStationary(stationary);
489 modRI = mod;
490 }
491 modRI->SetHighEnergyLimit(e1DNA);
492 theDNAIoni->AddEmModel(-1, modRI, reg);
493
494 if(e2DNA < lowEnergyRPWBA) {
495 auto modI = new G4DNABornIonisationModel1();
496 modI->SelectFasterComputation(fast);
497 modI->SelectStationary(stationary);
498 modI->SetLowEnergyLimit(e2DNA);
499 modI->SetHighEnergyLimit(lowEnergyRPWBA);
500 theDNAIoni->AddEmModel(-2, modI, reg);
501 }
502 if(lowEnergyRPWBA < emaxIonDNA) {
503 auto modJ = new G4DNARPWBAIonisationModel();
504 modJ->SelectFasterComputation(fast);
505 modJ->SelectStationary(stationary);
506 modJ->SetLowEnergyLimit(lowEnergyRPWBA);
507 modJ->SetHighEnergyLimit(emaxIonDNA);
508 theDNAIoni->AddEmModel(-3, modJ, reg);
509 }
510
511 // *** Charge decrease ***
512 auto theDNAChargeDecreaseProcess =
513 FindOrBuildChargeDecrease(part, "proton_G4DNAChargeDecrease");
514 auto modDCD = new G4DNADingfelderChargeDecreaseModel();
515 modDCD->SelectStationary(stationary);
516 modDCD->SetLowEnergyLimit(0.0);
517 modDCD->SetHighEnergyLimit(emax);
518 theDNAChargeDecreaseProcess->AddEmModel(-1, modDCD, reg);
519
520 FindOrBuildCapture(0.1*CLHEP::keV, part);
521}
522
523//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
524
525void
527 const G4bool stationary,
528 const G4Region* reg)
529{
531
532 // *** Ionisation ***
533 auto theDNAIoni = FindOrBuildIonisation(part, "GenericIon_G4DNAIonisation");
534 auto mod = new G4DNARuddIonisationExtendedModel();
535 mod->SelectStationary(stationary);
536 mod->SetHighEnergyLimit(emaxIonDNA);
537 theDNAIoni->AddEmModel(-1, mod, reg);
538
539 FindOrBuildCapture(25.0*CLHEP::keV, part);
540}
541
542//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
543
544void
546 const G4int charge,
547 const G4int opt,
548 const G4double emaxIonDNA,
549 const G4bool,
550 const G4bool stationary,
551 const G4Region* reg)
552{
554 const G4double emax = param->MaxKinEnergy();
555 const G4String& name = part->GetParticleName();
556
557 // *** Elastic ***
558 auto theDNAElastic = FindOrBuildElastic(part, name + "_G4DNAElastic");
559 auto modEI = new G4DNAIonElasticModel();
560 modEI->SelectStationary(stationary);
561 modEI->SetHighEnergyLimit(lowEnergyMSC);
562 theDNAElastic->AddEmModel(-1, modEI, reg);
563
564 // *** Excitation ***
565 auto theDNAExc = FindOrBuildExcitation(part, name + "_G4DNAExcitation");
566 auto modMGE = new G4DNAMillerGreenExcitationModel();
567 modMGE->SelectStationary(stationary);
568 modMGE->SetLowEnergyLimit(0.0);
569 modMGE->SetHighEnergyLimit(emaxIonDNA);
570 theDNAExc->AddEmModel(-1, modMGE, reg);
571
572 // *** Ionisation ***
573 auto theDNAIoni = FindOrBuildIonisation(part, name + "_G4DNAIonisation");
574 G4VEmModel* modRI = nullptr;
575 if(2 == opt) {
576 auto mod = new G4DNARuddIonisationExtendedModel();
577 mod->SelectStationary(stationary);
578 modRI = mod;
579 } else {
580 auto mod = new G4DNARuddIonisationModel();
581 mod->SelectStationary(stationary);
582 modRI = mod;
583 }
584 modRI->SetHighEnergyLimit(emaxIonDNA);
585 theDNAIoni->AddEmModel(-1, modRI, reg);
586
587 // *** Charge increase ***
588 if(2 > charge) {
589 auto theDNAChargeIncrease =
590 FindOrBuildChargeIncrease(part, name + "_G4DNAChargeIncrease");
591 auto modDCI = new G4DNADingfelderChargeIncreaseModel();
592 modDCI->SelectStationary(stationary);
593 modDCI->SetLowEnergyLimit(0.0);
594 modDCI->SetHighEnergyLimit(emax);
595 theDNAChargeIncrease->AddEmModel(-1, modDCI, reg);
596 }
597
598 // *** Charge decrease ***
599 if(0 < charge) {
600 auto theDNAChargeDecrease =
601 FindOrBuildChargeDecrease(part, name + "_G4DNAChargeDecrease");
602 auto modDCD = new G4DNADingfelderChargeDecreaseModel();
603 modDCD->SelectStationary(stationary);
604 modDCD->SetLowEnergyLimit(0.0);
605 modDCD->SetHighEnergyLimit(emax);
606 theDNAChargeDecrease->AddEmModel(-1, modDCD, reg);
607 }
608 FindOrBuildCapture(1*CLHEP::keV, part);
609}
610
611//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
612
614{
615 auto elec = G4Electron::Electron();
617 G4DNAElectronSolvation* ptr = dynamic_cast<G4DNAElectronSolvation*>(p);
618 if(nullptr == ptr) {
619 ptr = new G4DNAElectronSolvation("e-_G4DNAElectronSolvation");
621 ph->RegisterProcess(ptr, elec);
622 ptr->SetEmModel(new G4DummyModel());
623 }
624 return ptr;
625}
626
627//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
628
631 const G4String& name)
632{
634 G4DNAElastic* ptr = dynamic_cast<G4DNAElastic*>(p);
635 if(nullptr == ptr) {
636 ptr = new G4DNAElastic(name);
638 ph->RegisterProcess(ptr, part);
639 ptr->SetEmModel(new G4DummyModel());
640 }
641 return ptr;
642}
643
644//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
645
648 const G4String& name)
649{
651 G4DNAExcitation* ptr = dynamic_cast<G4DNAExcitation*>(p);
652 if(nullptr == ptr) {
653 ptr = new G4DNAExcitation(name);
655 ph->RegisterProcess(ptr, part);
656 ptr->SetEmModel(new G4DummyModel());
657 }
658 return ptr;
659}
660
661//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
662
665 const G4String& name)
666{
668 G4DNAVibExcitation* ptr = dynamic_cast<G4DNAVibExcitation*>(p);
669 if(nullptr == ptr) {
670 ptr = new G4DNAVibExcitation(name);
672 ph->RegisterProcess(ptr, part);
673 ptr->SetEmModel(new G4DummyModel());
674 }
675 return ptr;
676}
677
678//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
679
682 const G4String& name)
683{
685 G4DNAIonisation* ptr = dynamic_cast<G4DNAIonisation*>(p);
686 if(nullptr == ptr) {
687 ptr = new G4DNAIonisation(name);
689 ph->RegisterProcess(ptr, part);
690 ptr->SetEmModel(new G4DummyModel());
691 }
692 return ptr;
693}
694
695//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
696
699 const G4String& name)
700{
702 G4DNAAttachment* ptr = dynamic_cast<G4DNAAttachment*>(p);
703 if(nullptr == ptr) {
704 ptr = new G4DNAAttachment(name);
706 ph->RegisterProcess(ptr, part);
707 ptr->SetEmModel(new G4DummyModel());
708 }
709 return ptr;
710}
711
712//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
713
716 const G4String& name)
717{
719 G4DNAChargeDecrease* ptr = dynamic_cast<G4DNAChargeDecrease*>(p);
720 if(nullptr == ptr) {
721 ptr = new G4DNAChargeDecrease(name);
723 ph->RegisterProcess(ptr, part);
724 ptr->SetEmModel(new G4DummyModel());
725 }
726 return ptr;
727}
728
729//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
730
733 const G4String& name)
734{
736 G4DNAChargeIncrease* ptr = dynamic_cast<G4DNAChargeIncrease*>(p);
737 if(nullptr == ptr) {
738 ptr = new G4DNAChargeIncrease(name);
740 ph->RegisterProcess(ptr, part);
741 ptr->SetEmModel(new G4DummyModel());
742 }
743 return ptr;
744}
745
746//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
747
750{
751 auto p = G4PhysListUtil::FindProcess(part, 66);
752 G4LowECapture* ptr = dynamic_cast<G4LowECapture*>(p);
753 if(nullptr == ptr) {
754 ptr = new G4LowECapture(elim);
755 auto mng = part->GetProcessManager();
756 mng->AddDiscreteProcess(ptr);
757 }
758 return ptr;
759}
760
761//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4DNABornExcitationModel1 G4DNABornExcitationModel
#define G4DNABornIonisationModel
G4EmDNAMscModelType
@ dnaWVI
@ dnaUrban
@ dnaGS
@ fLowEnergyChargeIncrease
@ fLowEnergyVibrationalExcitation
@ fLowEnergyElectronSolvation
@ fLowEnergyChargeDecrease
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
static G4DNAGenericIonsManager * Instance(void)
G4ParticleDefinition * GetIon(const G4String &name)
static G4VEmModel * GetMacroDefinedModel()
One step thermalization model can be chosen via macro using /process/dna/e-SolvationSubType Ritchie19...
static G4Electron * Electron()
Definition: G4Electron.cc:93
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:360
static void PrepareEMPhysics()
Definition: G4EmBuilder.cc:399
static void ConstructBasicEmPhysics(G4hMultipleScattering *hmsc, const std::vector< G4int > &listHadrons)
Definition: G4EmBuilder.cc:105
static G4LowECapture * FindOrBuildCapture(const G4double elim, G4ParticleDefinition *part)
static void ConstructDNALightIonPhysics(G4ParticleDefinition *part, const G4int charge, const G4int opt, const G4double emax, const G4bool fast, const G4bool stationary, const G4Region *reg=nullptr)
static G4DNAChargeDecrease * FindOrBuildChargeDecrease(G4ParticleDefinition *part, const G4String &name)
static G4DNAExcitation * FindOrBuildExcitation(G4ParticleDefinition *part, const G4String &name)
static G4DNAElastic * FindOrBuildElastic(G4ParticleDefinition *part, const G4String &name)
static G4DNAChargeIncrease * FindOrBuildChargeIncrease(G4ParticleDefinition *part, const G4String &name)
static void ConstructDNAParticles()
static void ConstructDNAIonPhysics(const G4double emax, const G4bool stationary, const G4Region *reg=nullptr)
static G4DNAElectronSolvation * FindOrBuildElectronSolvation()
static void ConstructDNAProtonPhysics(const G4double e1DNA, const G4double emaxDNA, const G4int opt, const G4bool fast, const G4bool stationary, const G4Region *reg=nullptr)
static G4DNAVibExcitation * FindOrBuildVibExcitation(G4ParticleDefinition *part, const G4String &name)
static void ConstructStandardEmPhysics(const G4double emin_electron, const G4double emin_proton, const G4double emin_alpha, const G4double emin_ion, const G4EmDNAMscModelType mscType, const G4bool fast)
static G4DNAAttachment * FindOrBuildAttachment(G4ParticleDefinition *part, const G4String &name)
static G4DNAIonisation * FindOrBuildIonisation(G4ParticleDefinition *part, const G4String &name)
static void ConstructDNAElectronPhysics(const G4double emaxDNA, const G4int opt, const G4bool fast, const G4bool stationary, const G4Region *reg=nullptr)
static G4EmParameters * Instance()
G4double MaxKinEnergy() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
static G4VProcess * FindProcess(const G4ParticleDefinition *, G4int subtype)
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Positron * Positron()
Definition: G4Positron.cc:93
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
static G4Proton * Proton()
Definition: G4Proton.cc:92
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetEmModel(G4VMscModel *, G4int idx=0)