Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmModelActivator.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//---------------------------------------------------------------------------
28//
29// ClassName: G4EmModelActivator
30//
31// Author: V.Ivanchenko 01.06.2015
32//
33// Organisation: G4AI
34// Contract: ESA contract 4000107387/12/NL/AT
35// Customer: ESA/ESTEC
36//
37// Modified:
38//
39//----------------------------------------------------------------------------
40//
41
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
44
45#include "G4EmModelActivator.hh"
46#include "G4EmParameters.hh"
48#include "G4ParticleTable.hh"
49#include "G4RegionStore.hh"
50#include "G4Region.hh"
52#include "G4VMscModel.hh"
53#include "G4LossTableManager.hh"
54#include "G4EmConfigurator.hh"
55#include "G4PAIModel.hh"
56#include "G4PAIPhotModel.hh"
57#include "G4Gamma.hh"
58#include "G4Electron.hh"
59#include "G4Positron.hh"
60#include "G4MuonPlus.hh"
61#include "G4MuonMinus.hh"
62#include "G4Proton.hh"
63#include "G4GenericIon.hh"
64#include "G4Alpha.hh"
65#include "G4ProcessManager.hh"
66#include "G4DummyModel.hh"
67#include "G4EmProcessSubType.hh"
69#include "G4EmParticleList.hh"
70
71#include "G4MicroElecElastic.hh"
73
76
77#include "G4BraggModel.hh"
78#include "G4BraggIonModel.hh"
79#include "G4BetheBlochModel.hh"
80#include "G4UrbanMscModel.hh"
82#include "G4IonFluctuations.hh"
84#include "G4LowECapture.hh"
88#include "G4ionIonisation.hh"
90
93#include "G4WentzelVIModel.hh"
96#include "G4UrbanMscModel.hh"
99
106
114
115#include "G4SystemOfUnits.hh"
116#include "G4Threading.hh"
117#include <vector>
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
120
122 : baseName(emphys)
123{
124 theParameters = G4EmParameters::Instance();
125
126 const std::vector<G4String>& regnamesPAI = theParameters->RegionsPAI();
127 if(regnamesPAI.size() > 0)
128 {
129 ActivatePAI();
130 }
131 const std::vector<G4String>& regnamesME = theParameters->RegionsMicroElec();
132 if(regnamesME.size() > 0)
133 {
134 ActivateMicroElec();
135 }
136 const std::vector<G4String>& regnamesMSC = theParameters->RegionsPhysics();
137 if(regnamesMSC.size() > 0)
138 {
139 ActivateEmOptions();
140 }
141}
142
143//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
144
145void G4EmModelActivator::ActivateEmOptions()
146{
147 const std::vector<G4String>& regnamesPhys = theParameters->RegionsPhysics();
148 G4int nreg = regnamesPhys.size();
149 if(0 == nreg) { return; }
150 G4int verbose = theParameters->Verbose() - 1;
151 if(verbose > 0) {
152 G4cout << "### G4EmModelActivator::ActivateEmOptions for " << nreg << " regions"
153 << G4endl;
154 }
155 const std::vector<G4String>& typesPhys = theParameters->TypesPhysics();
156
157 // start configuration of models
160 const G4ParticleDefinition* phot = G4Gamma::Gamma();
163 G4EmConfigurator* em_config = man->EmConfigurator();
164 G4VAtomDeexcitation* adeexc = man->AtomDeexcitation();
166 G4VEmModel* mod;
167
168 // high energy limit for low-energy e+- model of msc
169 G4double mscEnergyLimit = theParameters->MscEnergyLimit();
170
171 // high energy limit for Livermore and Penelope models
172 G4double highEnergyLimit = 1*GeV;
173
174 // general high energy limit
175 G4double highEnergy = theParameters->MaxKinEnergy();
176
177 for(G4int i=0; i<nreg; ++i) {
178 G4String reg = regnamesPhys[i];
179 if(verbose > 0) {
180 G4cout << i << ". region <" << reg << ">; type <" << typesPhys[i] << "> "
181 << G4endl;
182 }
183
184 if(baseName == typesPhys[i]) { continue; }
185
186 if("G4EmStandard" == typesPhys[i]) {
187 G4UrbanMscModel* msc = new G4UrbanMscModel();
188 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
189
190 msc = new G4UrbanMscModel();
191 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
192
193 } else if("G4EmStandard_opt1" == typesPhys[i] || "G4EmStandard_opt2" == typesPhys[i]) {
194 G4UrbanMscModel* msc = new G4UrbanMscModel();
195 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
196
197 msc = new G4UrbanMscModel();
198 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
199
200 } else if("G4EmStandard_opt3" == typesPhys[i]) {
201
202 G4DummyModel* dummy = new G4DummyModel();
203 G4UrbanMscModel* msc = new G4UrbanMscModel();
204 SetMscParameters(elec, msc, typesPhys[i]);
205 em_config->SetExtraEmModel("e-", "msc", msc, reg);
206 FindOrAddProcess(elec, "CoulombScat");
207 em_config->SetExtraEmModel("e-", "CoulombScat", dummy, reg);
208
209 msc = new G4UrbanMscModel();
210 SetMscParameters(posi, msc, typesPhys[i]);
211 em_config->SetExtraEmModel("e+", "msc", msc, reg);
212 FindOrAddProcess(posi, "CoulombScat");
213 em_config->SetExtraEmModel("e+", "CoulombScat", dummy, reg);
214
215 msc = new G4UrbanMscModel();
216 SetMscParameters(prot, msc, typesPhys[i]);
217 em_config->SetExtraEmModel("proton", "msc", msc, reg);
218 FindOrAddProcess(prot, "CoulombScat");
219 em_config->SetExtraEmModel("proton", "CoulombScat", dummy, reg);
220
221 theParameters->SetNumberOfBinsPerDecade(20);
223 theParameters->SetDeexActiveRegion(reg, true, false, false);
224 }
225 theParameters->DefineRegParamForDeex(adeexc);
226 FindOrAddProcess(phot, "Rayl");
227 mod = new G4LivermoreRayleighModel();
228 em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
229 FindOrAddProcess(phot, "phot");
231 em_config->SetExtraEmModel("gamma", "phot", mod, reg);
232 FindOrAddProcess(phot, "compt");
233 mod = new G4KleinNishinaModel();
234 em_config->SetExtraEmModel("gamma", "compt", mod, reg);
235
236 } else if("G4EmStandard_opt4" == typesPhys[i]) {
238 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
239
241 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
242
243 theParameters->SetNumberOfBinsPerDecade(20);
244 theParameters->SetUseMottCorrection(true);
246 theParameters->SetDeexActiveRegion(reg, true, false, false);
247 }
248 theParameters->DefineRegParamForDeex(adeexc);
249
250 FindOrAddProcess(phot, "Rayl");
251 mod = new G4LivermoreRayleighModel();
252 em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
253 FindOrAddProcess(phot, "compt");
254 mod = new G4KleinNishinaModel();
255 em_config->SetExtraEmModel("gamma", "compt", mod, reg);
256 mod = new G4LowEPComptonModel();
257 mod->SetHighEnergyLimit(20*MeV);
258 em_config->SetExtraEmModel("gamma", "compt", mod, reg);
259
260 } else if("G4EmStandardGS" == typesPhys[i]) {
262 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
263
265 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
266
267 } else if("G4EmStandardWVI" == typesPhys[i]) {
269 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
270
271 msc = new G4WentzelVIModel();
272 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
273
274 theParameters->SetMscThetaLimit(0.15);
275
277 theParameters->SetDeexActiveRegion(regnamesPhys[i], true, false, false);
278 }
279 theParameters->DefineRegParamForDeex(adeexc);
280
281 } else if("G4EmStandardSS" == typesPhys[i] &&
282 baseName != "G4EmStandard_opt3") {
283 G4EmParticleList emList;
284 for(const auto& particleName : emList.PartNames()) {
285 G4ParticleDefinition* particle = table->FindParticle(particleName);
286 if(particle && 0.0 != particle->GetPDGCharge()) {
287 FindOrAddProcess(particle, "CoulombScat");
289 sc->SetPolarAngleLimit(0.0);
290 sc->SetLocked(true);
291 em_config->SetExtraEmModel(particleName, "CoulombScat", sc, reg);
292 if(particleName == "mu+" || particleName == "mu-") {
293 em_config->SetExtraEmModel(particleName, "muMsc",
294 new G4DummyModel(), reg);
295 } else {
296 em_config->SetExtraEmModel(particleName, "msc",
297 new G4DummyModel(), reg);
298 }
299 }
300 }
302 theParameters->SetDeexActiveRegion(regnamesPhys[i], true, true, true);
303 }
304 theParameters->DefineRegParamForDeex(adeexc);
305
306 } else if("G4EmLivermore" == typesPhys[i]) {
307
309 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
310
312 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
313
315 em_config->SetExtraEmModel("gamma", "phot", mod, reg);
316 mod = new G4LivermoreComptonModel();
317 em_config->SetExtraEmModel("gamma", "compt", mod, reg);
319 em_config->SetExtraEmModel("gamma", "conv", mod, reg);
320
321 FindOrAddProcess(phot, "Rayl");
322 mod = new G4LivermoreRayleighModel();
323 em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
324
325 mod = new G4LivermoreIonisationModel();
327 em_config->SetExtraEmModel("e-", "eIoni", mod, reg, 0.0, 0.1*MeV, uf);
329 em_config->SetExtraEmModel("e-", "eBrem", mod, reg, 0.0, highEnergyLimit);
330
331 theParameters->SetNumberOfBinsPerDecade(20);
332 theParameters->SetUseMottCorrection(true);
334 theParameters->SetDeexActiveRegion(regnamesPhys[i], true, false, false);
335 }
336 theParameters->DefineRegParamForDeex(adeexc);
337
338 } else if("G4EmPenelope" == typesPhys[i]) {
339
341 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
342
344 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
345
347 em_config->SetExtraEmModel("gamma", "phot", mod, reg);
348 mod = new G4PenelopeComptonModel();
349 em_config->SetExtraEmModel("gamma", "compt", mod, reg);
351 em_config->SetExtraEmModel("gamma", "conv", mod, reg);
352
353 FindOrAddProcess(phot, "Rayl");
354 mod = new G4PenelopeRayleighModel();
355 em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
356
357 mod = new G4PenelopeIonisationModel();
359 em_config->SetExtraEmModel("e-", "eIoni", mod, reg, 0.0, highEnergyLimit, uf);
361 em_config->SetExtraEmModel("e-", "eBrem", mod, reg, 0.0, highEnergyLimit);
362
363 mod = new G4PenelopeIonisationModel();
364 uf = new G4UniversalFluctuation();
365 em_config->SetExtraEmModel("e+", "eIoni", mod, reg, 0.0, highEnergyLimit, uf);
367 em_config->SetExtraEmModel("e+", "eBrem", mod, reg, 0.0, highEnergyLimit);
368 mod = new G4PenelopeAnnihilationModel();
369 em_config->SetExtraEmModel("e+", "annihil", mod, reg, 0.0, highEnergyLimit);
370
371 theParameters->SetNumberOfBinsPerDecade(20);
372 theParameters->SetUseMottCorrection(true);
374 theParameters->SetDeexActiveRegion(regnamesPhys[i], true, false, false);
375 }
376 theParameters->DefineRegParamForDeex(adeexc);
377
378 } else {
379 if(verbose > 0 && G4Threading::IsMasterThread()) {
380 G4cout << "### G4EmModelActivator::ActivateEmOptions WARNING: \n"
381 << " EM Physics configuration name <" << typesPhys[i]
382 << "> is not known - ignored" << G4endl;
383 }
384 }
385 }
386}
387
388//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
389
390void G4EmModelActivator::ActivatePAI()
391{
392 const std::vector<G4String> regnamesPAI = theParameters->RegionsPAI();
393 G4int nreg = regnamesPAI.size();
394 if(0 == nreg) { return; }
395 G4int verbose = theParameters->Verbose() - 1;
396 if(verbose > 0) {
397 G4cout << "### G4EmModelActivator::ActivatePAI for " << nreg << " regions"
398 << G4endl;
399 }
400 const std::vector<G4String> particlesPAI = theParameters->ParticlesPAI();
401 const std::vector<G4String> typesPAI = theParameters->TypesPAI();
402
403 const std::vector<G4VEnergyLossProcess*>& v = G4LossTableManager::Instance()
405
407
413
414 for(G4int i = 0; i < nreg; ++i) {
415 const G4ParticleDefinition* p = nullptr;
416 if(particlesPAI[i] != "all") {
417 p = G4ParticleTable::GetParticleTable()->FindParticle(particlesPAI[i]);
418 if(!p) {
419 G4cout << "### WARNING: ActivatePAI::FindParticle fails to find "
420 << particlesPAI[i] << G4endl;
421 continue;
422 }
423 }
424 const G4Region* r = regionStore->GetRegion(regnamesPAI[i], false);
425 if(!r) {
426 G4cout << "### WARNING: ActivatePAI::GetRegion fails to find "
427 << regnamesPAI[i] << G4endl;
428 continue;
429 }
430
431 G4String name = "hIoni";
432 if(p == elec || p == posi)
433 { name = "eIoni"; }
434 else if (p == mupl || p == mumi)
435 { name = "muIoni"; }
436 else if (p == gion)
437 { name = "ionIoni"; }
438
439 for(auto proc : v) {
440
441 if(!proc->IsIonisationProcess()) { continue; }
442
443 G4String namep = proc->GetProcessName();
444 if(p) {
445 if(name != namep) { continue; }
446 } else {
447 if(namep != "hIoni" && namep != "muIoni" &&
448 namep != "eIoni" && namep != "ionIoni")
449 { continue; }
450 }
451 G4double emin = 50*CLHEP::keV;
452 if(namep == "eIoni") emin = 110*CLHEP::eV;
453 else if(namep == "muIoni") emin = 5*CLHEP::keV;
454
455 G4VEmModel* em = nullptr;
456 G4VEmFluctuationModel* fm = nullptr;
457 if(typesPAI[i] == "PAIphoton" || typesPAI[i] == "pai_photon") {
458 G4PAIPhotModel* mod = new G4PAIPhotModel(p,"PAIPhotModel");
459 em = mod;
460 fm = mod;
461 } else {
462 G4PAIModel* mod = new G4PAIModel(p,"PAIModel");
463 em = mod;
464 fm = mod;
465 }
466 // added PAI model above low energy limit
467 em->SetLowEnergyLimit(emin);
468 proc->AddEmModel(0, em, fm, r);
469
470 // added low energy model
471 if(namep == "eIoni") {
472 em = new G4MollerBhabhaModel();
473 fm = new G4UniversalFluctuation();
474 } else if(namep == "ionIoni") {
475 em = new G4BraggIonModel();
476 fm = new G4IonFluctuations();
477 } else {
478 em = new G4BraggModel();
479 fm = new G4UniversalFluctuation();
480 }
481 em->SetHighEnergyLimit(emin);
482 proc->AddEmModel(-1, em, fm, r);
483
484 if(verbose > 0) {
485 G4cout << "### G4EmModelActivator: add <" << typesPAI[i]
486 << "> model for " << particlesPAI[i]
487 << " in the " << regnamesPAI[i]
488 << " Emin(keV)= " << emin/CLHEP::keV << G4endl;
489 }
490 }
491 }
492}
493
494//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
495
496void G4EmModelActivator::ActivateMicroElec()
497{
498 const std::vector<G4String> regnamesME = theParameters->RegionsMicroElec();
499 G4int nreg = regnamesME.size();
500 if(0 == nreg)
501 {
502 return;
503 }
504 G4int verbose = theParameters->Verbose() - 1;
505 if(verbose > 0)
506 {
507 G4cout << "### G4EmModelActivator::ActivateMicroElec for " << nreg
508 << " regions" << G4endl;
509 }
511
515 G4ProcessManager* eman = elec->GetProcessManager();
516 G4ProcessManager* pman = prot->GetProcessManager();
517 G4ProcessManager* iman = gion->GetProcessManager();
518
519 G4bool emsc = HasMsc(eman);
520
521 // MicroElec elastic is not active in the world
522 G4MicroElecElastic* eElasProc =
523 new G4MicroElecElastic("e-G4MicroElecElastic");
524 eman->AddDiscreteProcess(eElasProc);
525
526 G4MicroElecInelastic* eInelProc =
527 new G4MicroElecInelastic("e-G4MicroElecInelastic");
528 eman->AddDiscreteProcess(eInelProc);
529
530 G4MicroElecInelastic* pInelProc =
531 new G4MicroElecInelastic("p_G4MicroElecInelastic");
532 pman->AddDiscreteProcess(pInelProc);
533
534 G4MicroElecInelastic* iInelProc =
535 new G4MicroElecInelastic("ion_G4MicroElecInelastic");
536 iman->AddDiscreteProcess(iInelProc);
537
538 // start configuration of models
539 G4EmConfigurator* em_config = man->EmConfigurator();
540 G4VEmModel* mod;
541
542 // limits for MicroElec applicability
543 G4double elowest = 16.7 * eV;
544 G4double elimel = 100 * MeV;
545 G4double elimin = 9 * MeV;
546 G4double pmin = 50 * keV;
547 G4double pmax = 99.9 * MeV;
548
549 G4LowECapture* ecap = new G4LowECapture(elowest);
550 eman->AddDiscreteProcess(ecap);
551
552 for(G4int i = 0; i < nreg; ++i)
553 {
554
555 G4String reg = regnamesME[i];
556 G4cout << "### MicroElec models are activated for G4Region " << reg
557 << G4endl
558 << " Energy limits for e- elastic: " << elowest/eV << " eV - "
559 << elimel/MeV << " MeV"
560 << G4endl
561 << " Energy limits for e- inelastic: " << elowest/eV << " eV - "
562 << elimin/MeV << " MeV"
563 << G4endl
564 << " Energy limits for hadrons/ions: " << pmin/MeV << " MeV - "
565 << pmax/MeV << " MeV"
566 << G4endl;
567
568 // e-
569 if(emsc)
570 {
571 G4UrbanMscModel* msc = new G4UrbanMscModel();
572 msc->SetActivationLowEnergyLimit(elimel);
573 em_config->SetExtraEmModel("e-", "msc", msc, reg);
574 }
575 else
576 {
577 mod = new G4DummyModel();
578 em_config->SetExtraEmModel("e-", "CoulombScat", mod, reg, 0.0, elimel);
579 }
580
581 mod = new G4MicroElecElasticModel();
582 em_config->SetExtraEmModel("e-",
583 "e-G4MicroElecElastic",
584 mod,
585 reg,
586 elowest,
587 elimin);
588
589 mod = new G4MollerBhabhaModel();
590 mod->SetActivationLowEnergyLimit(elimin);
591 em_config->SetExtraEmModel("e-",
592 "eIoni",
593 mod,
594 reg,
595 0.0,
596 10 * TeV,
598
599 mod = new G4MicroElecInelasticModel();
600 em_config->SetExtraEmModel("e-",
601 "e-G4MicroElecInelastic",
602 mod,
603 reg,
604 elowest,
605 elimin);
606
607 // proton
608 mod = new G4BraggModel();
610 em_config->SetExtraEmModel("proton",
611 "hIoni",
612 mod,
613 reg,
614 0.0,
615 2 * MeV,
617
618 mod = new G4BetheBlochModel();
620 em_config->SetExtraEmModel("proton",
621 "hIoni",
622 mod,
623 reg,
624 2 * MeV,
625 10 * TeV,
627
628 mod = new G4MicroElecInelasticModel();
629 em_config->SetExtraEmModel("proton",
630 "p_G4MicroElecInelastic",
631 mod,
632 reg,
633 pmin,
634 pmax);
635
636 // ions
637 mod = new G4BraggIonModel();
639 em_config->SetExtraEmModel("GenericIon",
640 "ionIoni",
641 mod,
642 reg,
643 0.0,
644 2 * MeV,
645 new G4IonFluctuations());
646
647 mod = new G4BetheBlochModel();
649 em_config->SetExtraEmModel("GenericIon",
650 "ionIoni",
651 mod,
652 reg,
653 2 * MeV,
654 10 * TeV,
655 new G4IonFluctuations());
656
657 mod = new G4MicroElecInelasticModel();
658 em_config->SetExtraEmModel("GenericIon",
659 "ion_G4MicroElecInelastic",
660 mod,
661 reg,
662 pmin,
663 pmax);
664 }
665}
666
667//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
668
669G4bool G4EmModelActivator::HasMsc(G4ProcessManager* pm) const
670{
671 G4bool res = false;
673 G4int nproc = pm->GetProcessListLength();
674 for(G4int i = 0; i < nproc; ++i)
675 {
676 if(((*pv)[i])->GetProcessSubType() == fMultipleScattering)
677 {
678 res = true;
679 break;
680 }
681 }
682 return res;
683}
684
685//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
686
687void G4EmModelActivator::AddStandardScattering(const G4ParticleDefinition* part,
688 G4EmConfigurator* em_config,
689 G4VMscModel* mscmod,
690 const G4String& reg,
691 G4double e1, G4double e2,
692 const G4String& type)
693{
694 G4String pname = part->GetParticleName();
695
696 // low-energy msc model
697 SetMscParameters(part, mscmod, type);
698 em_config->SetExtraEmModel(pname, "msc", mscmod, reg, 0.0, e1);
699
700 // high energy msc model
702 SetMscParameters(part, msc, type);
703 em_config->SetExtraEmModel(pname, "msc", msc, reg, e1, e2);
704
705 // high energy single scattering model
706 FindOrAddProcess(part, "CoulombScat");
709 mod->SetLocked(true);
710 em_config->SetExtraEmModel(pname, "CoulombScat", mod, reg, 0.0, e2);
711}
712
713//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
714
715void G4EmModelActivator::SetMscParameters(const G4ParticleDefinition* part,
716 G4VMscModel* msc, const G4String& phys)
717{
718 if(part == G4Electron::Electron() || part == G4Positron::Positron()) {
719 if(phys == "G4EmStandard_opt1" || phys == "G4EmStandard_opt2") {
720 msc->SetRangeFactor(0.2);
722 } else if(phys == "G4EmStandard_opt3") {
724 } else if(phys == "G4EmStandard_opt4" || phys == "G4EmLivermore" || phys == "G4EmPenelope") {
725 msc->SetRangeFactor(0.08);
727 msc->SetSkin(3);
728 } else if(phys == "G4EmStandardGS") {
729 msc->SetRangeFactor(0.06);
730 }
731 } else {
732 if(phys != "G4EmStandard" && phys != "G4EmStandard_opt1" && phys != "G4EmStandard_opt2") {
733 msc->SetLateralDisplasmentFlag(true);
734 }
735 }
736 msc->SetLocked(true);
737}
738
739//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
740
741void G4EmModelActivator::FindOrAddProcess(const G4ParticleDefinition* part,
742 const G4String& name)
743{
746 G4int nproc = pm->GetProcessListLength();
747 for(G4int i = 0; i<nproc; ++i) {
748 if(((*pv)[i])->GetProcessName() == name) { return; }
749 }
750 if(name == "CoulombScat") {
752 cs->SetEmModel(new G4DummyModel());
753 pm->AddDiscreteProcess(cs);
754 } else if(name == "Rayl") {
756 rs->SetEmModel(new G4DummyModel());
757 pm->AddDiscreteProcess(rs);
758 }
759}
760
761//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ fMultipleScattering
@ fMinimal
@ fUseSafetyPlus
@ fUseDistanceToBoundary
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4Electron * Electron()
Definition: G4Electron.cc:93
void SetExtraEmModel(const G4String &particleName, const G4String &processName, G4VEmModel *, const G4String &regionName="", G4double emin=0.0, G4double emax=DBL_MAX, G4VEmFluctuationModel *fm=nullptr)
G4EmModelActivator(const G4String &emphys="")
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
const std::vector< G4String > & TypesPhysics() const
void SetMscThetaLimit(G4double val)
G4double MscEnergyLimit() const
const std::vector< G4String > & RegionsPAI() const
void SetDeexActiveRegion(const G4String &region, G4bool fdeex, G4bool fauger, G4bool fpixe)
G4int Verbose() const
const std::vector< G4String > & ParticlesPAI() const
const std::vector< G4String > & RegionsPhysics() const
void DefineRegParamForDeex(G4VAtomDeexcitation *) const
G4double MaxKinEnergy() const
void SetUseMottCorrection(G4bool val)
const std::vector< G4String > & RegionsMicroElec() const
const std::vector< G4String > & TypesPAI() const
const std::vector< G4String > & PartNames() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
static G4LossTableManager * Instance()
G4EmConfigurator * EmConfigurator()
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
G4VAtomDeexcitation * AtomDeexcitation()
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:99
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:98
G4ProcessManager * GetProcessManager() const
G4double GetPDGCharge() const
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
static G4Positron * Positron()
Definition: G4Positron.cc:93
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4int GetProcessListLength() const
G4ProcessVector * GetProcessList() const
static G4Proton * Proton()
Definition: G4Proton.cc:92
static G4RegionStore * GetInstance()
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:792
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:778
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
void SetActivationHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:771
void SetLocked(G4bool)
Definition: G4VEmModel.hh:874
void SetEmModel(G4VEmModel *, G4int index=0)
void SetSkin(G4double)
Definition: G4VMscModel.hh:220
void SetRangeFactor(G4double)
Definition: G4VMscModel.hh:227
void SetLateralDisplasmentFlag(G4bool val)
Definition: G4VMscModel.hh:213
void SetStepLimitType(G4MscStepLimitType)
Definition: G4VMscModel.hh:255
const char * name(G4int ptype)
G4bool IsMasterThread()
Definition: G4Threading.cc:124