Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmDNABuilder Class Reference

#include <G4EmDNABuilder.hh>

Static Public Member Functions

static void ConstructDNAParticles ()
 
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 void ConstructDNAElectronPhysics (const G4double emaxDNA, const G4int opt, const G4bool fast, const G4bool stationary, const G4Region *reg=nullptr)
 
static void ConstructDNAProtonPhysics (const G4double e1DNA, const G4double emaxDNA, const G4int opt, const G4bool fast, const G4bool stationary, const G4Region *reg=nullptr)
 
static void ConstructDNAIonPhysics (const G4double emax, const G4bool stationary, const G4Region *reg=nullptr)
 
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 G4DNAElectronSolvationFindOrBuildElectronSolvation ()
 
static G4DNAElasticFindOrBuildElastic (G4ParticleDefinition *part, const G4String &name)
 
static G4DNAExcitationFindOrBuildExcitation (G4ParticleDefinition *part, const G4String &name)
 
static G4DNAVibExcitationFindOrBuildVibExcitation (G4ParticleDefinition *part, const G4String &name)
 
static G4DNAIonisationFindOrBuildIonisation (G4ParticleDefinition *part, const G4String &name)
 
static G4DNAAttachmentFindOrBuildAttachment (G4ParticleDefinition *part, const G4String &name)
 
static G4DNAChargeDecreaseFindOrBuildChargeDecrease (G4ParticleDefinition *part, const G4String &name)
 
static G4DNAChargeIncreaseFindOrBuildChargeIncrease (G4ParticleDefinition *part, const G4String &name)
 
static G4LowECaptureFindOrBuildCapture (const G4double elim, G4ParticleDefinition *part)
 
static void FindOrBuildNuclearStopping (G4ParticleDefinition *part, const G4double elim)
 

Detailed Description

Definition at line 59 of file G4EmDNABuilder.hh.

Member Function Documentation

◆ ConstructDNAElectronPhysics()

void G4EmDNABuilder::ConstructDNAElectronPhysics ( const G4double emaxDNA,
const G4int opt,
const G4bool fast,
const G4bool stationary,
const G4Region * reg = nullptr )
static

Definition at line 303 of file G4EmDNABuilder.cc.

308{
309 G4ParticleDefinition* part = G4Electron::Electron();
310
311 // limit of the Emfietzoglou models
312 G4double emaxE = 0.0;
313 // limit of the elastic and solvation models
314 G4double emaxT = 7.4*CLHEP::eV;
315 // limit for CPA100 models
316 G4double emaxCPA100 = 250*CLHEP::keV;
317 if(4 == opt) {
318 emaxE = 10.*CLHEP::keV;
319 emaxT = 10.*CLHEP::eV;
320 } else if(5 < opt) {
321 emaxT = 11.*CLHEP::eV;
322 }
323
324 // *** Solvation ***
325 G4DNAElectronSolvation* pSolvation = FindOrBuildElectronSolvation();
327 therm->SetHighEnergyLimit(emaxT);
328 pSolvation->AddEmModel(-1, therm, reg);
329
330 // *** Elastic scattering ***
331 auto pElasticProcess = FindOrBuildElastic(part, "e-_G4DNAElastic");
332 G4VEmModel* elast = nullptr;
333 G4VEmModel* elast2 = nullptr;
334 if(4 == opt) {
335 elast = new G4DNAUeharaScreenedRutherfordElasticModel();
336 } else if(5 < opt) {
337 auto mod = new G4DNACPA100ElasticModel();
338 mod->SelectStationary(stationary);
339 elast = mod;
340 elast2 = new G4DNAChampionElasticModel();
341 } else {
342 elast = new G4DNAChampionElasticModel();
343 }
344 elast->SetHighEnergyLimit(lowEnergyMSC);
345 pElasticProcess->AddEmModel(-2, elast, reg);
346
347 if(nullptr != elast2) {
348 elast->SetHighEnergyLimit(emaxCPA100);
349 elast2->SetLowEnergyLimit(emaxCPA100);
350 elast2->SetHighEnergyLimit(lowEnergyMSC);
351 pElasticProcess->AddEmModel(-3, elast2, reg);
352 }
353
354 // *** Excitation ***
355 auto theDNAExc = FindOrBuildExcitation(part, "e-_G4DNAExcitation");
356 if(emaxE > 0.0) {
357 auto modE = new G4DNAEmfietzoglouExcitationModel();
358 theDNAExc->AddEmModel(-1, modE, reg);
359 modE->SelectStationary(stationary);
360 modE->SetHighEnergyLimit(emaxE);
361 }
362 G4VEmModel* modB = nullptr;
363 G4VEmModel* modB2 = nullptr;
364 if(6 == opt) {
365 auto mod = new G4DNACPA100ExcitationModel();
366 mod->SelectStationary(stationary);
367 modB = mod;
368 auto mod1 = new G4DNABornExcitationModel();
369 mod1->SelectStationary(stationary);
370 modB2 = mod1;
371 } else {
372 auto mod = new G4DNABornExcitationModel();
373 mod->SelectStationary(stationary);
374 modB = mod;
375 }
376 modB->SetLowEnergyLimit(emaxE);
377 modB->SetHighEnergyLimit(emaxDNA);
378 theDNAExc->AddEmModel(-2, modB, reg);
379 if(nullptr != modB2) {
380 modB->SetHighEnergyLimit(emaxCPA100);
381 modB2->SetLowEnergyLimit(emaxCPA100);
382 modB2->SetHighEnergyLimit(emaxDNA);
383 theDNAExc->AddEmModel(-3, modB2, reg);
384 }
385
386 // *** Ionisation ***
387 auto theDNAIoni = FindOrBuildIonisation(part, "e-_G4DNAIonisation");
388 if(emaxE > 0.0) {
389 auto modE = new G4DNAEmfietzoglouIonisationModel();
390 theDNAIoni->AddEmModel(-1, modE, reg);
391 modE->SelectFasterComputation(fast);
392 modE->SelectStationary(stationary);
393 modE->SetHighEnergyLimit(emaxE);
394 }
395 G4VEmModel* modI = nullptr;
396 G4VEmModel* modI2 = nullptr;
397 if(6 == opt) {
398 auto mod = new G4DNACPA100IonisationModel();
399 mod->SelectStationary(stationary);
400 mod->SelectFasterComputation(fast);
401 modI = mod;
402 auto mod1 = new G4DNABornIonisationModel();
403 mod1->SelectStationary(stationary);
404 modI2 = mod1;
405 } else {
406 auto mod = new G4DNABornIonisationModel1();
407 mod->SelectStationary(stationary);
408 mod->SelectFasterComputation(fast);
409 modI = mod;
410 }
411 modI->SetLowEnergyLimit(emaxE);
412 modI->SetHighEnergyLimit(emaxDNA);
413 theDNAIoni->AddEmModel(-2, modI, reg);
414 if(nullptr != modI2) {
415 modI->SetHighEnergyLimit(emaxCPA100);
416 modI2->SetLowEnergyLimit(emaxCPA100);
417 modI2->SetHighEnergyLimit(emaxDNA);
418 theDNAIoni->AddEmModel(-3, modI2, reg);
419 }
420
421 if(4 != opt && 6 != opt) {
422 // *** Vibrational excitation ***
423 auto theDNAVibExc = FindOrBuildVibExcitation(part, "e-_G4DNAVibExcitation");
424 auto modS = new G4DNASancheExcitationModel();
425 theDNAVibExc->AddEmModel(-1, modS, reg);
426 modS->SelectStationary(stationary);
427
428 // *** Attachment ***
429 auto theDNAAttach = FindOrBuildAttachment(part, "e-_G4DNAAttachment");
430 auto modM = new G4DNAMeltonAttachmentModel();
431 theDNAAttach->AddEmModel(-1, modM, reg);
432 modM->SelectStationary(stationary);
433 }
434}
G4DNABornExcitationModel1 G4DNABornExcitationModel
#define G4DNABornIonisationModel
double G4double
Definition G4Types.hh:83
static G4VEmModel * GetMacroDefinedModel()
One step thermalization model can be chosen via macro using /process/dna/e-SolvationSubType Ritchie19...
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4DNAExcitation * FindOrBuildExcitation(G4ParticleDefinition *part, const G4String &name)
static G4DNAElastic * FindOrBuildElastic(G4ParticleDefinition *part, const G4String &name)
static G4DNAElectronSolvation * FindOrBuildElectronSolvation()
static G4DNAVibExcitation * FindOrBuildVibExcitation(G4ParticleDefinition *part, const G4String &name)
static G4DNAAttachment * FindOrBuildAttachment(G4ParticleDefinition *part, const G4String &name)
static G4DNAIonisation * FindOrBuildIonisation(G4ParticleDefinition *part, const G4String &name)
void SetHighEnergyLimit(G4double)
void SetLowEnergyLimit(G4double)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)

Referenced by G4EmDNAPhysics::ConstructProcess(), G4EmDNAPhysics_option2::ConstructProcess(), G4EmDNAPhysics_option4::ConstructProcess(), G4EmDNAPhysics_option6::ConstructProcess(), G4EmDNAPhysics_option8::ConstructProcess(), and G4EmDNAPhysicsActivator::ConstructProcess().

◆ ConstructDNAIonPhysics()

void G4EmDNABuilder::ConstructDNAIonPhysics ( const G4double emax,
const G4bool stationary,
const G4Region * reg = nullptr )
static

Definition at line 527 of file G4EmDNABuilder.cc.

530{
531 G4ParticleDefinition* part = G4GenericIon::GenericIon();
532
533 // *** Ionisation ***
534 auto theDNAIoni = FindOrBuildIonisation(part, "GenericIon_G4DNAIonisation");
535 auto mod = new G4DNARuddIonisationExtendedModel();
536 mod->SelectStationary(stationary);
537 mod->SetHighEnergyLimit(emaxIonDNA);
538 theDNAIoni->AddEmModel(-1, mod, reg);
539
540 // *** NIEL ***
541 FindOrBuildNuclearStopping(part, CLHEP::MeV);
542
543 // *** Tracking cut ***
544 FindOrBuildCapture(1*CLHEP::keV, part);
545}
static G4LowECapture * FindOrBuildCapture(const G4double elim, G4ParticleDefinition *part)
static void FindOrBuildNuclearStopping(G4ParticleDefinition *part, const G4double elim)
static G4GenericIon * GenericIon()

Referenced by G4EmDNAPhysics::ConstructProcess(), G4EmDNAPhysics_option2::ConstructProcess(), G4EmDNAPhysics_option4::ConstructProcess(), G4EmDNAPhysics_option6::ConstructProcess(), G4EmDNAPhysics_option8::ConstructProcess(), and G4EmDNAPhysicsActivator::ConstructProcess().

◆ ConstructDNALightIonPhysics()

void G4EmDNABuilder::ConstructDNALightIonPhysics ( G4ParticleDefinition * part,
const G4int charge,
const G4int opt,
const G4double emax,
const G4bool fast,
const G4bool stationary,
const G4Region * reg = nullptr )
static

Definition at line 550 of file G4EmDNABuilder.cc.

557{
558 G4EmParameters* param = G4EmParameters::Instance();
559 const G4double emax = param->MaxKinEnergy();
560 const G4String& name = part->GetParticleName();
561
562 // *** Elastic ***
563 auto theDNAElastic = FindOrBuildElastic(part, name + "_G4DNAElastic");
564 auto modEI = new G4DNAIonElasticModel();
565 modEI->SelectStationary(stationary);
566 modEI->SetHighEnergyLimit(lowEnergyMSC);
567 theDNAElastic->AddEmModel(-1, modEI, reg);
568
569 // *** Excitation ***
570 auto theDNAExc = FindOrBuildExcitation(part, name + "_G4DNAExcitation");
571 auto modMGE = new G4DNAMillerGreenExcitationModel();
572 modMGE->SelectStationary(stationary);
573 modMGE->SetLowEnergyLimit(0.0);
574 modMGE->SetHighEnergyLimit(emaxIonDNA);
575 theDNAExc->AddEmModel(-1, modMGE, reg);
576
577 // *** Ionisation ***
578 auto theDNAIoni = FindOrBuildIonisation(part, name + "_G4DNAIonisation");
579 G4VEmModel* modRI = nullptr;
580 if(2 == opt) {
581 auto mod = new G4DNARuddIonisationExtendedModel();
582 mod->SelectStationary(stationary);
583 modRI = mod;
584 } else {
585 auto mod = new G4DNARuddIonisationModel();
586 mod->SelectStationary(stationary);
587 modRI = mod;
588 }
589 modRI->SetHighEnergyLimit(emaxIonDNA);
590 theDNAIoni->AddEmModel(-1, modRI, reg);
591
592 // *** Charge increase ***
593 if(2 > charge) {
594 auto theDNAChargeIncrease =
595 FindOrBuildChargeIncrease(part, name + "_G4DNAChargeIncrease");
596 auto modDCI = new G4DNADingfelderChargeIncreaseModel();
597 modDCI->SelectStationary(stationary);
598 modDCI->SetLowEnergyLimit(0.0);
599 modDCI->SetHighEnergyLimit(emax);
600 theDNAChargeIncrease->AddEmModel(-1, modDCI, reg);
601 }
602
603 // *** Charge decrease ***
604 if(0 < charge) {
605 auto theDNAChargeDecrease =
606 FindOrBuildChargeDecrease(part, name + "_G4DNAChargeDecrease");
607 auto modDCD = new G4DNADingfelderChargeDecreaseModel();
608 modDCD->SelectStationary(stationary);
609 modDCD->SetLowEnergyLimit(0.0);
610 modDCD->SetHighEnergyLimit(emax);
611 theDNAChargeDecrease->AddEmModel(-1, modDCD, reg);
612 }
613
614 // *** Tracking cut ***
615 FindOrBuildCapture(1*CLHEP::keV, part);
616}
static G4DNAChargeDecrease * FindOrBuildChargeDecrease(G4ParticleDefinition *part, const G4String &name)
static G4DNAChargeIncrease * FindOrBuildChargeIncrease(G4ParticleDefinition *part, const G4String &name)
static G4EmParameters * Instance()
G4double MaxKinEnergy() const
const G4String & GetParticleName() const
const char * name(G4int ptype)

Referenced by G4EmDNAPhysics::ConstructProcess(), G4EmDNAPhysics_option2::ConstructProcess(), G4EmDNAPhysics_option4::ConstructProcess(), G4EmDNAPhysics_option6::ConstructProcess(), G4EmDNAPhysics_option8::ConstructProcess(), and G4EmDNAPhysicsActivator::ConstructProcess().

◆ ConstructDNAParticles()

void G4EmDNABuilder::ConstructDNAParticles ( )
static

Definition at line 112 of file G4EmDNABuilder.cc.

113{
114 // standard particles
116
117 // DNA ions
118 G4DNAGenericIonsManager* genericIonsManager
120 genericIonsManager->GetIon("alpha+");
121 genericIonsManager->GetIon("helium");
122 genericIonsManager->GetIon("hydrogen");
123}
static G4DNAGenericIonsManager * Instance()
G4ParticleDefinition * GetIon(const G4String &name)
static void ConstructMinimalEmSet()

Referenced by G4EmDNAPhysics::ConstructParticle(), and G4EmDNAPhysicsActivator::ConstructParticle().

◆ ConstructDNAProtonPhysics()

void G4EmDNABuilder::ConstructDNAProtonPhysics ( const G4double e1DNA,
const G4double emaxDNA,
const G4int opt,
const G4bool fast,
const G4bool stationary,
const G4Region * reg = nullptr )
static

Definition at line 439 of file G4EmDNABuilder.cc.

445{
446 G4EmParameters* param = G4EmParameters::Instance();
447 const G4double emax = param->MaxKinEnergy();
448 G4ParticleDefinition* part = G4Proton::Proton();
449
450 // *** Elastic scattering ***
451 auto pElasticProcess = FindOrBuildElastic(part, "proton_G4DNAElastic");
452 auto modE = new G4DNAIonElasticModel();
453 modE->SetHighEnergyLimit(lowEnergyMSC);
454 modE->SelectStationary(stationary);
455 pElasticProcess->AddEmModel(-1, modE, reg);
456
457 // *** Excitation ***
458 G4double e2DNA = std::min(e1DNA, lowEnergyRPWBA);
459 auto theDNAExc = FindOrBuildExcitation(part, "proton_G4DNAExcitation");
460 auto modMGE = new G4DNAMillerGreenExcitationModel();
461 modMGE->SetHighEnergyLimit(e2DNA);
462 modMGE->SelectStationary(stationary);
463 theDNAExc->AddEmModel(-1, modMGE, reg);
464
465 if(e2DNA < lowEnergyRPWBA) {
466 auto modB = new G4DNABornExcitationModel();
467 modB->SelectStationary(stationary);
468 modB->SetLowEnergyLimit(e2DNA);
469 modB->SetHighEnergyLimit(lowEnergyRPWBA);
470 theDNAExc->AddEmModel(-2, modB, reg);
471 }
472 if(lowEnergyRPWBA < emaxIonDNA) {
473 auto modC = new G4DNARPWBAExcitationModel();
474 modC->SelectStationary(stationary);
475 modC->SetLowEnergyLimit(lowEnergyRPWBA);
476 modC->SetHighEnergyLimit(emaxIonDNA);
477 theDNAExc->AddEmModel(-3, modC, reg);
478 }
479
480 // *** Ionisation ***
481 auto theDNAIoni = FindOrBuildIonisation(part, "proton_G4DNAIonisation");
482 G4VEmModel* modRI = nullptr;
483 if(2 == opt) {
484 auto mod = new G4DNARuddIonisationExtendedModel();
485 mod->SelectStationary(stationary);
486 modRI = mod;
487 } else {
488 auto mod = new G4DNARuddIonisationModel();
489 mod->SelectStationary(stationary);
490 modRI = mod;
491 }
492 modRI->SetHighEnergyLimit(e1DNA);
493 theDNAIoni->AddEmModel(-1, modRI, reg);
494
495 if(e2DNA < lowEnergyRPWBA) {
496 auto modI = new G4DNABornIonisationModel1();
497 modI->SelectFasterComputation(fast);
498 modI->SelectStationary(stationary);
499 modI->SetLowEnergyLimit(e2DNA);
500 modI->SetHighEnergyLimit(lowEnergyRPWBA);
501 theDNAIoni->AddEmModel(-2, modI, reg);
502 }
503 if(lowEnergyRPWBA < emaxIonDNA) {
504 auto modJ = new G4DNARPWBAIonisationModel();
505 modJ->SelectFasterComputation(fast);
506 modJ->SelectStationary(stationary);
507 modJ->SetLowEnergyLimit(lowEnergyRPWBA);
508 modJ->SetHighEnergyLimit(emaxIonDNA);
509 theDNAIoni->AddEmModel(-3, modJ, reg);
510 }
511
512 // *** Charge decrease ***
513 auto theDNAChargeDecreaseProcess =
514 FindOrBuildChargeDecrease(part, "proton_G4DNAChargeDecrease");
515 auto modDCD = new G4DNADingfelderChargeDecreaseModel();
516 modDCD->SelectStationary(stationary);
517 modDCD->SetLowEnergyLimit(0.0);
518 modDCD->SetHighEnergyLimit(emax);
519 theDNAChargeDecreaseProcess->AddEmModel(-1, modDCD, reg);
520
521 FindOrBuildCapture(0.1*CLHEP::keV, part);
522}
static G4Proton * Proton()
Definition G4Proton.cc:90

Referenced by G4EmDNAPhysics::ConstructProcess(), G4EmDNAPhysics_option2::ConstructProcess(), G4EmDNAPhysics_option4::ConstructProcess(), G4EmDNAPhysics_option6::ConstructProcess(), G4EmDNAPhysics_option8::ConstructProcess(), and G4EmDNAPhysicsActivator::ConstructProcess().

◆ ConstructStandardEmPhysics()

void G4EmDNABuilder::ConstructStandardEmPhysics ( const G4double emin_electron,
const G4double emin_proton,
const G4double emin_alpha,
const G4double emin_ion,
const G4EmDNAMscModelType mscType,
const G4bool fast )
static

Definition at line 128 of file G4EmDNABuilder.cc.

134{
135 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
136 G4EmParameters* param = G4EmParameters::Instance();
137 const G4double emax = param->MaxKinEnergy();
139
140 // gamma
141 G4ParticleDefinition* gamma = G4Gamma::Gamma();
142
143 // photoelectric effect - Livermore model
144 auto thePEEffect = new G4PhotoElectricEffect();
145 thePEEffect->SetEmModel(new G4LivermorePhotoElectricModel());
146 ph->RegisterProcess(thePEEffect, gamma);
147
148 // Compton scattering - Klein-Nishina
149 auto theComptonScattering = new G4ComptonScattering();
150 theComptonScattering->SetEmModel(new G4KleinNishinaModel());
151 auto cModel = new G4LowEPComptonModel();
152 cModel->SetHighEnergyLimit(20*CLHEP::MeV);
153 theComptonScattering->AddEmModel(0, cModel);
154 ph->RegisterProcess(theComptonScattering, gamma);
155
156 // gamma conversion - 5D model
157 auto theGammaConversion = new G4GammaConversion();
158 ph->RegisterProcess(theGammaConversion, gamma);
159
160 // Rayleigh scattering - Livermore model
161 auto theRayleigh = new G4RayleighScattering();
162 ph->RegisterProcess(theRayleigh, gamma);
163
164 // electron
165 if(emin_elec < emax) {
166 G4ParticleDefinition* elec = G4Electron::Electron();
167 auto msc_el = new G4eMultipleScattering();
168 G4VMscModel* msc_model_el;
169 if(mscType == dnaWVI) {
170 msc_model_el = new G4LowEWentzelVIModel();
171 } else if(mscType == dnaGS) {
172 msc_model_el = new G4GoudsmitSaundersonMscModel();
173 } else {
174 msc_model_el = new G4UrbanMscModel();
175 }
176 msc_model_el->SetActivationLowEnergyLimit(lowEnergyMSC);
177 msc_el->SetEmModel(msc_model_el);
178 ph->RegisterProcess(msc_el, elec);
179
180 auto ioni = new G4eIonisation();
181 auto mb_el = new G4MollerBhabhaModel();
182 mb_el->SetActivationLowEnergyLimit(emin_elec);
183 ioni->SetEmModel(mb_el);
184 ph->RegisterProcess(ioni, elec);
185
186 auto brem = new G4eBremsstrahlung();
187 auto sb_el = new G4SeltzerBergerModel();
188 sb_el->SetActivationLowEnergyLimit(emin_elec);
189 sb_el->SetHighEnergyLimit(emax);
190 sb_el->SetAngularDistribution(new G4Generator2BS());
191 brem->SetEmModel(sb_el);
192 ph->RegisterProcess(brem, elec);
193 }
194
195 // positron
196 G4ParticleDefinition* posi = G4Positron::Positron();
197 auto msc_pos = new G4eMultipleScattering();
198 G4VMscModel* msc_model_pos;
199 if(mscType == dnaWVI) {
200 msc_model_pos = new G4LowEWentzelVIModel();
201 } else if(mscType == dnaGS) {
202 msc_model_pos = new G4GoudsmitSaundersonMscModel();
203 } else {
204 msc_model_pos = new G4UrbanMscModel();
205 }
206 msc_pos->SetEmModel(msc_model_pos);
207 ph->RegisterProcess(msc_pos, posi);
208 ph->RegisterProcess(new G4eIonisation(), posi);
209
210 auto brem = new G4eBremsstrahlung();
211 auto sb = new G4SeltzerBergerModel();
212 sb->SetHighEnergyLimit(emax);
213 sb->SetAngularDistribution(new G4Generator2BS());
214 brem->SetEmModel(sb);
215 ph->RegisterProcess(brem, posi);
216 ph->RegisterProcess(new G4eplusAnnihilation(), posi);
217
218 // proton
219 if(emin_proton < emax) {
220 G4ParticleDefinition* part = G4Proton::Proton();
221 StandardHadronPhysics(part, lowEnergyMSC, emin_proton, emax,
222 mscType, false);
223 }
224
225 // GenericIon
226 if(emin_ion < emax) {
227 G4ParticleDefinition* ion = G4GenericIon::GenericIon();
228 StandardHadronPhysics(ion, lowEnergyMSC, emin_ion, emax,
229 dnaUrban, true);
230 }
231
232 // alpha
233 if(emin_alpha < emax) {
234 G4ParticleDefinition* part = G4Alpha::Alpha();
235 StandardHadronPhysics(part, lowEnergyMSC, emin_alpha, emax,
236 dnaUrban, true);
237
238 // alpha+
239 G4DNAGenericIonsManager* genericIonsManager
241 part = genericIonsManager->GetIon("alpha+");
242 StandardHadronPhysics(part, lowEnergyMSC, emin_alpha, emax,
243 dnaUrban, false);
244 }
245 // list of main standard particles
246 const std::vector<G4int> chargedParticles = {
247 13, -13, 211, -211, 321, -321, -2212,
248 1000010020, 1000010030, 1000020030
249 };
250 auto msc = new G4hMultipleScattering();
251 msc->SetEmModel(new G4WentzelVIModel());
252 G4EmBuilder::ConstructBasicEmPhysics(msc, chargedParticles);
253}
@ dnaWVI
@ dnaUrban
@ dnaGS
static G4Alpha * Alpha()
Definition G4Alpha.cc:83
static void PrepareEMPhysics()
static void ConstructBasicEmPhysics(G4hMultipleScattering *hmsc, const std::vector< G4int > &listHadrons)
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Positron * Positron()
Definition G4Positron.cc:90
void SetActivationLowEnergyLimit(G4double)

Referenced by G4EmDNAPhysics::ConstructProcess(), G4EmDNAPhysics_option2::ConstructProcess(), G4EmDNAPhysics_option4::ConstructProcess(), G4EmDNAPhysics_option6::ConstructProcess(), and G4EmDNAPhysics_option8::ConstructProcess().

◆ FindOrBuildAttachment()

G4DNAAttachment * G4EmDNABuilder::FindOrBuildAttachment ( G4ParticleDefinition * part,
const G4String & name )
static

Definition at line 705 of file G4EmDNABuilder.cc.

707{
709 G4DNAAttachment* ptr = dynamic_cast<G4DNAAttachment*>(p);
710 if(nullptr == ptr) {
711 ptr = new G4DNAAttachment(name);
712 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
713 ph->RegisterProcess(ptr, part);
714 ptr->SetEmModel(new G4DummyModel());
715 }
716 return ptr;
717}
static G4VProcess * FindProcess(const G4ParticleDefinition *, G4int subtype)
void SetEmModel(G4VEmModel *, G4int index=0)

Referenced by ConstructDNAElectronPhysics().

◆ FindOrBuildCapture()

G4LowECapture * G4EmDNABuilder::FindOrBuildCapture ( const G4double elim,
G4ParticleDefinition * part )
static

Definition at line 756 of file G4EmDNABuilder.cc.

757{
758 auto p = G4PhysListUtil::FindProcess(part, -1);
759 G4LowECapture* ptr = dynamic_cast<G4LowECapture*>(p);
760 if (nullptr == ptr) {
761 ptr = new G4LowECapture(elim);
762 auto mng = part->GetProcessManager();
763 mng->AddDiscreteProcess(ptr);
764 }
765 return ptr;
766}
G4ProcessManager * GetProcessManager() const
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)

Referenced by ConstructDNAIonPhysics(), ConstructDNALightIonPhysics(), and ConstructDNAProtonPhysics().

◆ FindOrBuildChargeDecrease()

G4DNAChargeDecrease * G4EmDNABuilder::FindOrBuildChargeDecrease ( G4ParticleDefinition * part,
const G4String & name )
static

Definition at line 722 of file G4EmDNABuilder.cc.

724{
726 G4DNAChargeDecrease* ptr = dynamic_cast<G4DNAChargeDecrease*>(p);
727 if(nullptr == ptr) {
728 ptr = new G4DNAChargeDecrease(name);
729 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
730 ph->RegisterProcess(ptr, part);
731 ptr->SetEmModel(new G4DummyModel());
732 }
733 return ptr;
734}

Referenced by ConstructDNALightIonPhysics(), and ConstructDNAProtonPhysics().

◆ FindOrBuildChargeIncrease()

G4DNAChargeIncrease * G4EmDNABuilder::FindOrBuildChargeIncrease ( G4ParticleDefinition * part,
const G4String & name )
static

Definition at line 739 of file G4EmDNABuilder.cc.

741{
743 G4DNAChargeIncrease* ptr = dynamic_cast<G4DNAChargeIncrease*>(p);
744 if(nullptr == ptr) {
745 ptr = new G4DNAChargeIncrease(name);
746 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
747 ph->RegisterProcess(ptr, part);
748 ptr->SetEmModel(new G4DummyModel());
749 }
750 return ptr;
751}

Referenced by ConstructDNALightIonPhysics().

◆ FindOrBuildElastic()

G4DNAElastic * G4EmDNABuilder::FindOrBuildElastic ( G4ParticleDefinition * part,
const G4String & name )
static

Definition at line 637 of file G4EmDNABuilder.cc.

639{
641 G4DNAElastic* ptr = dynamic_cast<G4DNAElastic*>(p);
642 if(nullptr == ptr) {
643 ptr = new G4DNAElastic(name);
644 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
645 ph->RegisterProcess(ptr, part);
646 ptr->SetEmModel(new G4DummyModel());
647 }
648 return ptr;
649}

Referenced by ConstructDNAElectronPhysics(), ConstructDNALightIonPhysics(), and ConstructDNAProtonPhysics().

◆ FindOrBuildElectronSolvation()

G4DNAElectronSolvation * G4EmDNABuilder::FindOrBuildElectronSolvation ( )
static

Definition at line 620 of file G4EmDNABuilder.cc.

621{
622 auto elec = G4Electron::Electron();
624 G4DNAElectronSolvation* ptr = dynamic_cast<G4DNAElectronSolvation*>(p);
625 if(nullptr == ptr) {
626 ptr = new G4DNAElectronSolvation("e-_G4DNAElectronSolvation");
627 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
628 ph->RegisterProcess(ptr, elec);
629 ptr->SetEmModel(new G4DummyModel());
630 }
631 return ptr;
632}
@ fLowEnergyElectronSolvation

Referenced by ConstructDNAElectronPhysics().

◆ FindOrBuildExcitation()

G4DNAExcitation * G4EmDNABuilder::FindOrBuildExcitation ( G4ParticleDefinition * part,
const G4String & name )
static

Definition at line 654 of file G4EmDNABuilder.cc.

656{
658 G4DNAExcitation* ptr = dynamic_cast<G4DNAExcitation*>(p);
659 if(nullptr == ptr) {
660 ptr = new G4DNAExcitation(name);
661 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
662 ph->RegisterProcess(ptr, part);
663 ptr->SetEmModel(new G4DummyModel());
664 }
665 return ptr;
666}

Referenced by ConstructDNAElectronPhysics(), ConstructDNALightIonPhysics(), and ConstructDNAProtonPhysics().

◆ FindOrBuildIonisation()

G4DNAIonisation * G4EmDNABuilder::FindOrBuildIonisation ( G4ParticleDefinition * part,
const G4String & name )
static

Definition at line 688 of file G4EmDNABuilder.cc.

690{
692 G4DNAIonisation* ptr = dynamic_cast<G4DNAIonisation*>(p);
693 if(nullptr == ptr) {
694 ptr = new G4DNAIonisation(name);
695 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
696 ph->RegisterProcess(ptr, part);
697 ptr->SetEmModel(new G4DummyModel());
698 }
699 return ptr;
700}

Referenced by ConstructDNAElectronPhysics(), ConstructDNAIonPhysics(), ConstructDNALightIonPhysics(), and ConstructDNAProtonPhysics().

◆ FindOrBuildNuclearStopping()

void G4EmDNABuilder::FindOrBuildNuclearStopping ( G4ParticleDefinition * part,
const G4double elim )
static

Definition at line 770 of file G4EmDNABuilder.cc.

772{
774 auto ptr = dynamic_cast<G4NuclearStopping*>(p);
775 if (nullptr == ptr) {
776 ptr = new G4NuclearStopping();
777 }
778 ptr->SetMaxKinEnergy(elim);
780 ph->RegisterProcess(ptr, part);
781}
@ fNuclearStopping

Referenced by ConstructDNAIonPhysics().

◆ FindOrBuildVibExcitation()

G4DNAVibExcitation * G4EmDNABuilder::FindOrBuildVibExcitation ( G4ParticleDefinition * part,
const G4String & name )
static

Definition at line 671 of file G4EmDNABuilder.cc.

673{
675 G4DNAVibExcitation* ptr = dynamic_cast<G4DNAVibExcitation*>(p);
676 if(nullptr == ptr) {
677 ptr = new G4DNAVibExcitation(name);
678 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
679 ph->RegisterProcess(ptr, part);
680 ptr->SetEmModel(new G4DummyModel());
681 }
682 return ptr;
683}
@ fLowEnergyVibrationalExcitation

Referenced by ConstructDNAElectronPhysics().


The documentation for this class was generated from the following files: