Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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#include "G4NuclearStopping.hh"
68
69// standard models
73#include "G4UrbanMscModel.hh"
78#include "G4Generator2BS.hh"
79#include "G4BraggModel.hh"
80#include "G4BraggIonModel.hh"
81#include "G4BetheBlochModel.hh"
82#include "G4DummyModel.hh"
83
84// DNA models
108
109static const G4double lowEnergyRPWBA = 100*CLHEP::MeV;
110static const G4double lowEnergyMSC = 1*CLHEP::MeV;
111static const G4double lowEnergyProtonIoni = 2*CLHEP::MeV;
112static const G4double highEnergyMillerGrean = 0.5*CLHEP::MeV;
113static const G4double highEnergyChargeExchange = 100*CLHEP::MeV;
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
116
118{
119 // standard particles
121
122 // DNA ions
123 G4DNAGenericIonsManager* genericIonsManager
125 genericIonsManager->GetIon("alpha+");
126 genericIonsManager->GetIon("helium");
127 genericIonsManager->GetIon("hydrogen");
128}
129
130//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
131
132void
134 const G4double emin_proton,
135 const G4double emin_alpha,
136 const G4double emin_ion,
137 const G4EmDNAMscModelType mscType,
138 const G4bool)
139{
142 const G4double emax = param->MaxKinEnergy();
144
145 // gamma
147
148 // photoelectric effect - Livermore model
149 auto thePEEffect = new G4PhotoElectricEffect();
150 thePEEffect->SetEmModel(new G4LivermorePhotoElectricModel());
151 ph->RegisterProcess(thePEEffect, gamma);
152
153 // Compton scattering - Klein-Nishina
154 auto theComptonScattering = new G4ComptonScattering();
155 theComptonScattering->SetEmModel(new G4KleinNishinaModel());
156 auto cModel = new G4LowEPComptonModel();
157 cModel->SetHighEnergyLimit(20*CLHEP::MeV);
158 theComptonScattering->AddEmModel(0, cModel);
159 ph->RegisterProcess(theComptonScattering, gamma);
160
161 // gamma conversion - 5D model
162 auto theGammaConversion = new G4GammaConversion();
163 ph->RegisterProcess(theGammaConversion, gamma);
164
165 // Rayleigh scattering - Livermore model
166 auto theRayleigh = new G4RayleighScattering();
167 ph->RegisterProcess(theRayleigh, gamma);
168
169 // electron
170 if(emin_elec < emax) {
172 auto msc_el = new G4eMultipleScattering();
173 G4VMscModel* msc_model_el;
174 if(mscType == dnaWVI) {
175 msc_model_el = new G4LowEWentzelVIModel();
176 } else if(mscType == dnaGS) {
177 msc_model_el = new G4GoudsmitSaundersonMscModel();
178 } else {
179 msc_model_el = new G4UrbanMscModel();
180 }
181 msc_model_el->SetActivationLowEnergyLimit(lowEnergyMSC);
182 msc_el->SetEmModel(msc_model_el);
183 ph->RegisterProcess(msc_el, elec);
184
185 auto ioni = new G4eIonisation();
186 auto mb_el = new G4MollerBhabhaModel();
187 mb_el->SetActivationLowEnergyLimit(emin_elec);
188 ioni->SetEmModel(mb_el);
189 ph->RegisterProcess(ioni, elec);
190
191 auto brem = new G4eBremsstrahlung();
192 auto sb_el = new G4SeltzerBergerModel();
193 sb_el->SetActivationLowEnergyLimit(emin_elec);
194 sb_el->SetHighEnergyLimit(emax);
195 sb_el->SetAngularDistribution(new G4Generator2BS());
196 brem->SetEmModel(sb_el);
197 ph->RegisterProcess(brem, elec);
198 }
199
200 // positron
202 auto msc_pos = new G4eMultipleScattering();
203 G4VMscModel* msc_model_pos;
204 if(mscType == dnaWVI) {
205 msc_model_pos = new G4LowEWentzelVIModel();
206 } else if(mscType == dnaGS) {
207 msc_model_pos = new G4GoudsmitSaundersonMscModel();
208 } else {
209 msc_model_pos = new G4UrbanMscModel();
210 }
211 msc_pos->SetEmModel(msc_model_pos);
212 ph->RegisterProcess(msc_pos, posi);
213 ph->RegisterProcess(new G4eIonisation(), posi);
214
215 auto brem = new G4eBremsstrahlung();
216 auto sb = new G4SeltzerBergerModel();
217 sb->SetHighEnergyLimit(emax);
218 sb->SetAngularDistribution(new G4Generator2BS());
219 brem->SetEmModel(sb);
220 ph->RegisterProcess(brem, posi);
221 ph->RegisterProcess(new G4eplusAnnihilation(), posi);
222
223 // proton
224 if(emin_proton < emax) {
226 StandardHadronPhysics(part, lowEnergyMSC, emin_proton, emax,
227 mscType, false);
228 }
229
230 // GenericIon
231 if(emin_ion < emax) {
233 StandardHadronPhysics(ion, lowEnergyMSC, emin_ion, emax,
234 dnaUrban, true);
235 }
236
237 // alpha
238 if(emin_alpha < emax) {
240 StandardHadronPhysics(part, lowEnergyMSC, emin_alpha, emax,
241 dnaUrban, true);
242
243 // alpha+
244 G4DNAGenericIonsManager* genericIonsManager
246 part = genericIonsManager->GetIon("alpha+");
247 StandardHadronPhysics(part, lowEnergyMSC, emin_alpha, emax,
248 dnaUrban, false);
249 }
250 // list of main standard particles
251 const std::vector<G4int> chargedParticles = {
252 13, -13, 211, -211, 321, -321, -2212,
253 1000010020, 1000010030, 1000020030
254 };
255 auto msc = new G4hMultipleScattering();
256 msc->SetEmModel(new G4WentzelVIModel());
257 G4EmBuilder::ConstructBasicEmPhysics(msc, chargedParticles);
258}
259
260//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
261
262void G4EmDNABuilder::StandardHadronPhysics(G4ParticleDefinition* part,
263 const G4double lowELimitForMSC,
264 const G4double lowELimitForIoni,
265 const G4double maxEnergy,
266 const G4EmDNAMscModelType mscType,
267 const G4bool isIon)
268{
271 G4VMscModel* msc_model = nullptr;
272 if(mscType == dnaWVI) {
273 msc_model = new G4LowEWentzelVIModel();
274 } else {
275 msc_model = new G4UrbanMscModel();
276 }
277 msc_model->SetActivationLowEnergyLimit(lowELimitForMSC);
278 msc_model->SetLowEnergyLimit(lowELimitForMSC);
279 msc_model->SetHighEnergyLimit(maxEnergy);
280 msc->SetEmModel(msc_model);
281 ph->RegisterProcess(msc, part);
282
283 G4VEnergyLossProcess* ioni = nullptr;
284 G4VEmModel* mod1 = nullptr;
285 if(isIon) {
286 ioni = new G4ionIonisation();
287 mod1 = new G4BraggIonModel();
288 } else {
289 ioni = new G4hIonisation();
290 mod1 = new G4BraggModel();
291 }
292 G4double eth = lowEnergyProtonIoni*part->GetPDGMass()/CLHEP::proton_mass_c2;
293 mod1->SetActivationLowEnergyLimit(lowELimitForIoni);
294 mod1->SetHighEnergyLimit(eth);
295 ioni->SetEmModel(mod1);
296
297 G4VEmModel* mod2 = new G4BetheBlochModel();
298 mod2->SetActivationLowEnergyLimit(lowELimitForIoni);
299 mod2->SetLowEnergyLimit(eth);
300 ioni->SetEmModel(mod2);
301
302 ph->RegisterProcess(ioni, part);
303}
304
305//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
306
307void
309 const G4int opt,
310 const G4bool fast,
311 const G4bool stationary,
312 const G4Region* reg)
313{
315
316 // limit of the Emfietzoglou models
317 G4double emaxE = 0.0;
318 // limit of the elastic and solvation models
319 G4double emaxT = 7.4*CLHEP::eV;
320 // limit for CPA100 models
321 G4double emaxCPA100 = 250*CLHEP::keV;
322 if (4 == opt) {
323 emaxE = 10.*CLHEP::keV;
324 emaxT = 10.*CLHEP::eV;
325 } else if(6 <= opt) {
326 emaxT = 11.*CLHEP::eV;
327 }
328
329 // *** Solvation ***
332 therm->SetHighEnergyLimit(emaxT);
333 pSolvation->AddEmModel(-1, therm, reg);
334
335 // *** Elastic scattering ***
336 auto pElasticProcess = FindOrBuildElastic(part, "e-_G4DNAElastic");
337 G4VEmModel* elast;
338 G4VEmModel* elast2 = nullptr;
339 if(4 == opt) {
341 } else if(6 <= opt) {
342 auto mod = new G4DNACPA100ElasticModel();
343 mod->SelectStationary(stationary);
344 elast = mod;
345 elast2 = new G4DNAChampionElasticModel();
346 } else {
347 elast = new G4DNAChampionElasticModel();
348 }
349 elast->SetLowEnergyLimit(emaxT);
350 elast->SetHighEnergyLimit(lowEnergyMSC);
351 pElasticProcess->AddEmModel(-2, elast, reg);
352
353 if(nullptr != elast2) {
354 elast->SetHighEnergyLimit(emaxCPA100);
355 elast2->SetLowEnergyLimit(emaxCPA100);
356 elast2->SetHighEnergyLimit(lowEnergyMSC);
357 pElasticProcess->AddEmModel(-3, elast2, reg);
358 }
359
360 // *** Excitation ***
361 auto theDNAExc = FindOrBuildExcitation(part, "e-_G4DNAExcitation");
362 if (emaxE > 0.0) {
363 auto modE = new G4DNAEmfietzoglouExcitationModel();
364 theDNAExc->AddEmModel(-1, modE, reg);
365 modE->SelectStationary(stationary);
366 modE->SetHighEnergyLimit(emaxE);
367 }
368 G4VEmModel* modB;
369 G4VEmModel* modB2 = nullptr;
370 if(6 == opt) {
371 auto mod = new G4DNACPA100ExcitationModel();
372 mod->SelectStationary(stationary);
373 modB = mod;
374 auto mod1 = new G4DNABornExcitationModel();
375 mod1->SelectStationary(stationary);
376 modB2 = mod1;
377 } else {
378 auto mod = new G4DNABornExcitationModel();
379 mod->SelectStationary(stationary);
380 modB = mod;
381 }
382 modB->SetLowEnergyLimit(emaxE);
383 modB->SetHighEnergyLimit(emaxDNA);
384 theDNAExc->AddEmModel(-2, modB, reg);
385 if(nullptr != modB2) {
386 modB->SetHighEnergyLimit(emaxCPA100);
387 modB2->SetLowEnergyLimit(emaxCPA100);
388 modB2->SetHighEnergyLimit(emaxDNA);
389 theDNAExc->AddEmModel(-3, modB2, reg);
390 }
391
392 // *** Ionisation ***
393 auto theDNAIoni = FindOrBuildIonisation(part, "e-_G4DNAIonisation");
394 if(emaxE > 0.0) {
395 auto modE = new G4DNAEmfietzoglouIonisationModel();
396 theDNAIoni->AddEmModel(-1, modE, reg);
397 modE->SelectFasterComputation(fast);
398 modE->SelectStationary(stationary);
399 modE->SetHighEnergyLimit(emaxE);
400 }
401 G4VEmModel* modI;
402 G4VEmModel* modI2 = nullptr;
403 if (6 == opt) {
404 auto mod = new G4DNACPA100IonisationModel();
405 mod->SelectStationary(stationary);
406 mod->SelectFasterComputation(fast);
407 modI = mod;
408 modI2 = new G4DNABornIonisationModel1();
409 } else {
410 modI = new G4DNABornIonisationModel1();
411 }
412 modI->SetLowEnergyLimit(emaxE);
413 modI->SetHighEnergyLimit(emaxDNA);
414 theDNAIoni->AddEmModel(-2, modI, reg);
415 if(nullptr != modI2) {
416 modI->SetHighEnergyLimit(emaxCPA100);
417 modI2->SetLowEnergyLimit(emaxCPA100);
418 modI2->SetHighEnergyLimit(emaxDNA);
419 theDNAIoni->AddEmModel(-3, modI2, reg);
420 }
421
422 if(4 != opt && 6 != opt) {
423 // *** Vibrational excitation ***
424 auto theDNAVibExc = FindOrBuildVibExcitation(part, "e-_G4DNAVibExcitation");
425 auto modS = new G4DNASancheExcitationModel();
426 theDNAVibExc->AddEmModel(-1, modS, reg);
427 modS->SelectStationary(stationary);
428
429 // *** Attachment ***
430 auto theDNAAttach = FindOrBuildAttachment(part, "e-_G4DNAAttachment");
431 auto modM = new G4DNAMeltonAttachmentModel();
432 theDNAAttach->AddEmModel(-1, modM, reg);
433 modM->SelectStationary(stationary);
434 }
435}
436
437//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
438
439void
441 const G4double emaxIonDNA,
442 const G4int opt,
443 const G4bool fast,
444 const G4bool stationary,
445 const G4Region* reg)
446{
448 const G4double emax = param->MaxKinEnergy();
450 G4double e2DNA = (8 == opt) ? std::min(lowEnergyRPWBA, emax) : e1DNA;
451
452 // *** Elastic scattering ***
453 auto pElasticProcess = FindOrBuildElastic(part, "proton_G4DNAElastic");
454 auto modE = new G4DNAIonElasticModel();
455 modE->SetHighEnergyLimit(lowEnergyMSC);
456 modE->SelectStationary(stationary);
457 pElasticProcess->AddEmModel(-1, modE, reg);
458
459 // *** Excitation ***
460 auto theDNAExc = FindOrBuildExcitation(part, "proton_G4DNAExcitation");
461 auto modMGE = new G4DNAMillerGreenExcitationModel();
462 modMGE->SetHighEnergyLimit(e2DNA);
463 modMGE->SelectStationary(stationary);
464 theDNAExc->AddEmModel(-1, modMGE, reg);
465
466 if(e2DNA < lowEnergyRPWBA) {
467 auto modB = new G4DNABornExcitationModel();
468 modB->SelectStationary(stationary);
469 modB->SetLowEnergyLimit(e2DNA);
470 modB->SetHighEnergyLimit(lowEnergyRPWBA);
471 theDNAExc->AddEmModel(-2, modB, reg);
472 }
473 if(lowEnergyRPWBA < emaxIonDNA) {
474 auto modC = new G4DNARPWBAExcitationModel();
475 modC->SelectStationary(stationary);
476 modC->SetLowEnergyLimit(lowEnergyRPWBA);
477 modC->SetHighEnergyLimit(emaxIonDNA);
478 theDNAExc->AddEmModel(-3, modC, reg);
479 }
480
481 // *** Ionisation ***
482 auto theDNAIoni = FindOrBuildIonisation(part, "proton_G4DNAIonisation");
483 G4VEmModel* modRI;
484 if (8 == opt) {
486 } else {
488 }
489 modRI->SetHighEnergyLimit(e2DNA);
490 theDNAIoni->AddEmModel(-1, modRI, reg);
491
492 if (e2DNA < lowEnergyRPWBA) {
494 modI->SetLowEnergyLimit(e2DNA);
495 modI->SetHighEnergyLimit(lowEnergyRPWBA);
496 theDNAIoni->AddEmModel(-2, modI, reg);
497 }
498 if (lowEnergyRPWBA < emaxIonDNA) {
499 auto modJ = new G4DNARPWBAIonisationModel();
500 modJ->SelectFasterComputation(fast);
501 modJ->SelectStationary(stationary);
502 modJ->SetLowEnergyLimit(lowEnergyRPWBA);
503 modJ->SetHighEnergyLimit(emaxIonDNA);
504 theDNAIoni->AddEmModel(-3, modJ, reg);
505 }
506
507 // *** Charge decrease ***
508 auto theDNAChargeDecreaseProcess =
509 FindOrBuildChargeDecrease(part, "proton_G4DNAChargeDecrease");
510 auto modDCD = new G4DNADingfelderChargeDecreaseModel();
511 modDCD->SelectStationary(stationary);
512 modDCD->SetLowEnergyLimit(0.0);
513 modDCD->SetHighEnergyLimit(highEnergyChargeExchange);
514 theDNAChargeDecreaseProcess->AddEmModel(-1, modDCD, reg);
515
516 // *** Tracking cut ***
517 G4double cut = (8 == opt) ? 0.05*CLHEP::keV : 1*CLHEP::keV;
518 FindOrBuildCapture(cut, part);
519}
520
521//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
522
523void
525 const G4int opt,
526 const G4Region* reg)
527{
529
530 // *** Ionisation ***
531 auto theDNAIoni = FindOrBuildIonisation(part, "GenericIon_G4DNAIonisation");
533 mod->SetHighEnergyLimit(emaxIonDNA);
534 theDNAIoni->AddEmModel(-1, mod, reg);
535
536 // *** NIEL ***
537 FindOrBuildNuclearStopping(part, CLHEP::MeV);
538
539 // *** Tracking cut ***
540 G4double cut = (8 == opt) ? 0.05*CLHEP::keV : 1*CLHEP::keV;
541 FindOrBuildCapture(cut, part);
542}
543
544//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
545
546void
548 const G4int charge,
549 const G4int opt,
550 const G4double emaxIonDNA,
551 const G4bool,
552 const G4bool stationary,
553 const G4Region* reg)
554{
555 const G4String& name = part->GetParticleName();
556 G4double elim1 = emaxIonDNA;
557 G4double elim2 = emaxIonDNA;
558 if (part->GetParticleName() == "hydrogen") {
559 elim1 = highEnergyMillerGrean;
560 elim2 = highEnergyChargeExchange;
561 }
562
563 // *** Elastic ***
564 auto theDNAElastic = FindOrBuildElastic(part, name + "_G4DNAElastic");
565 auto modEI = new G4DNAIonElasticModel();
566 modEI->SelectStationary(stationary);
567 modEI->SetHighEnergyLimit(lowEnergyMSC);
568 theDNAElastic->AddEmModel(-1, modEI, reg);
569
570 // *** Excitation ***
571 auto theDNAExc = FindOrBuildExcitation(part, name + "_G4DNAExcitation");
572 auto modMGE = new G4DNAMillerGreenExcitationModel();
573 modMGE->SelectStationary(stationary);
574 modMGE->SetLowEnergyLimit(0.0);
575 modMGE->SetHighEnergyLimit(elim1);
576 theDNAExc->AddEmModel(-1, modMGE, reg);
577
578 // *** Ionisation ***
579 auto theDNAIoni = FindOrBuildIonisation(part, name + "_G4DNAIonisation");
580 G4VEmModel* modRI;
581 if (8 == opt) {
583 } else {
585 }
586 modRI->SetHighEnergyLimit(elim2);
587 theDNAIoni->AddEmModel(-2, modRI, reg);
588
589 // *** Charge increase ***
590 if(2 > charge) {
591 auto theDNAChargeIncrease =
592 FindOrBuildChargeIncrease(part, name + "_G4DNAChargeIncrease");
593 auto modDCI = new G4DNADingfelderChargeIncreaseModel();
594 modDCI->SelectStationary(stationary);
595 modDCI->SetLowEnergyLimit(0.0);
596 modDCI->SetHighEnergyLimit(elim2);
597 theDNAChargeIncrease->AddEmModel(-1, modDCI, reg);
598 }
599
600 // *** Charge decrease ***
601 if(0 < charge) {
602 auto theDNAChargeDecrease =
603 FindOrBuildChargeDecrease(part, name + "_G4DNAChargeDecrease");
604 auto modDCD = new G4DNADingfelderChargeDecreaseModel();
605 modDCD->SelectStationary(stationary);
606 modDCD->SetLowEnergyLimit(0.0);
607 modDCD->SetHighEnergyLimit(elim2);
608 theDNAChargeDecrease->AddEmModel(-1, modDCD, reg);
609 }
610
611 // *** Tracking cut ***
612 G4double cut = (8 == opt) ? 0.05*CLHEP::keV : 1*CLHEP::keV;
613 FindOrBuildCapture(cut, part);
614}
615
616//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
617
619{
620 auto elec = G4Electron::Electron();
622 G4DNAElectronSolvation* ptr = dynamic_cast<G4DNAElectronSolvation*>(p);
623 if(nullptr == ptr) {
624 ptr = new G4DNAElectronSolvation("e-_G4DNAElectronSolvation");
626 ph->RegisterProcess(ptr, elec);
627 ptr->SetEmModel(new G4DummyModel());
628 }
629 return ptr;
630}
631
632//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
633
636 const G4String& name)
637{
639 G4DNAElastic* ptr = dynamic_cast<G4DNAElastic*>(p);
640 if(nullptr == ptr) {
641 ptr = new G4DNAElastic(name);
643 ph->RegisterProcess(ptr, part);
644 ptr->SetEmModel(new G4DummyModel());
645 }
646 return ptr;
647}
648
649//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
650
653 const G4String& name)
654{
656 G4DNAExcitation* ptr = dynamic_cast<G4DNAExcitation*>(p);
657 if(nullptr == ptr) {
658 ptr = new G4DNAExcitation(name);
660 ph->RegisterProcess(ptr, part);
661 ptr->SetEmModel(new G4DummyModel());
662 }
663 return ptr;
664}
665
666//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
667
670 const G4String& name)
671{
673 G4DNAVibExcitation* ptr = dynamic_cast<G4DNAVibExcitation*>(p);
674 if(nullptr == ptr) {
675 ptr = new G4DNAVibExcitation(name);
677 ph->RegisterProcess(ptr, part);
678 ptr->SetEmModel(new G4DummyModel());
679 }
680 return ptr;
681}
682
683//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
684
687 const G4String& name)
688{
690 G4DNAIonisation* ptr = dynamic_cast<G4DNAIonisation*>(p);
691 if(nullptr == ptr) {
692 ptr = new G4DNAIonisation(name);
694 ph->RegisterProcess(ptr, part);
695 ptr->SetEmModel(new G4DummyModel());
696 }
697 return ptr;
698}
699
700//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
701
704 const G4String& name)
705{
707 G4DNAAttachment* ptr = dynamic_cast<G4DNAAttachment*>(p);
708 if(nullptr == ptr) {
709 ptr = new G4DNAAttachment(name);
711 ph->RegisterProcess(ptr, part);
712 ptr->SetEmModel(new G4DummyModel());
713 }
714 return ptr;
715}
716
717//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
718
721 const G4String& name)
722{
724 G4DNAChargeDecrease* ptr = dynamic_cast<G4DNAChargeDecrease*>(p);
725 if(nullptr == ptr) {
726 ptr = new G4DNAChargeDecrease(name);
728 ph->RegisterProcess(ptr, part);
729 ptr->SetEmModel(new G4DummyModel());
730 }
731 return ptr;
732}
733
734//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
735
738 const G4String& name)
739{
741 G4DNAChargeIncrease* ptr = dynamic_cast<G4DNAChargeIncrease*>(p);
742 if(nullptr == ptr) {
743 ptr = new G4DNAChargeIncrease(name);
745 ph->RegisterProcess(ptr, part);
746 ptr->SetEmModel(new G4DummyModel());
747 }
748 return ptr;
749}
750
751//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
752
755{
756 auto p = G4PhysListUtil::FindProcess(part, -1);
757 G4LowECapture* ptr = dynamic_cast<G4LowECapture*>(p);
758 if (nullptr == ptr) {
759 ptr = new G4LowECapture(elim);
760 auto mng = part->GetProcessManager();
761 mng->AddDiscreteProcess(ptr);
762 }
763 return ptr;
764}
765
766//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
767
769 const G4double elim)
770{
772 auto ptr = dynamic_cast<G4NuclearStopping*>(p);
773 if (nullptr == ptr) {
774 ptr = new G4NuclearStopping();
775 }
776 ptr->SetMaxKinEnergy(elim);
778 ph->RegisterProcess(ptr, part);
779}
780
781//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4DNABornExcitationModel1 G4DNABornExcitationModel
G4EmDNAMscModelType
@ dnaWVI
@ dnaUrban
@ dnaGS
@ fNuclearStopping
@ fLowEnergyVibrationalExcitation
@ fLowEnergyElectronSolvation
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
static G4Alpha * Alpha()
Definition G4Alpha.cc:83
static G4DNAGenericIonsManager * Instance()
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:91
static void ConstructMinimalEmSet()
static void PrepareEMPhysics()
static void ConstructBasicEmPhysics(G4hMultipleScattering *hmsc, const std::vector< G4int > &listHadrons)
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 FindOrBuildNuclearStopping(G4ParticleDefinition *part, const G4double elim)
static void ConstructDNAIonPhysics(const G4double emax, const G4int opt, 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:81
static G4GenericIon * GenericIon()
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:90
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
static G4Proton * Proton()
Definition G4Proton.cc:90
G4Region defines a region or a group of regions in the detector geometry setup, sharing properties as...
Definition G4Region.hh:90
void SetHighEnergyLimit(G4double)
void SetActivationLowEnergyLimit(G4double)
void SetLowEnergyLimit(G4double)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetEmModel(G4VMscModel *, G4int idx=0)