Geant4-11
Public Member Functions | Private Member Functions | Private Attributes
G4EmModelActivator Class Reference

#include <G4EmModelActivator.hh>

Public Member Functions

 G4EmModelActivator (const G4EmModelActivator &)=delete
 
 G4EmModelActivator (const G4String &emphys="")
 
G4EmModelActivatoroperator= (const G4EmModelActivator &right)=delete
 

Private Member Functions

void ActivateEmOptions ()
 
void ActivateMicroElec ()
 
void ActivatePAI ()
 
void AddStandardScattering (const G4ParticleDefinition *, G4EmConfigurator *, G4VMscModel *, const G4String &, G4double, G4double, const G4String &)
 
void FindOrAddProcess (const G4ParticleDefinition *, const G4String &)
 
G4bool HasMsc (G4ProcessManager *) const
 
void SetMscParameters (const G4ParticleDefinition *, G4VMscModel *, const G4String &phys)
 

Private Attributes

G4String baseName
 
G4EmParameterstheParameters
 

Detailed Description

Definition at line 58 of file G4EmModelActivator.hh.

Constructor & Destructor Documentation

◆ G4EmModelActivator() [1/2]

G4EmModelActivator::G4EmModelActivator ( const G4String emphys = "")
explicit

Definition at line 122 of file G4EmModelActivator.cc.

123 : baseName(emphys)
124{
126
127 const std::vector<G4String>& regnamesPAI = theParameters->RegionsPAI();
128 if(regnamesPAI.size() > 0)
129 {
130 ActivatePAI();
131 }
132 const std::vector<G4String>& regnamesME = theParameters->RegionsMicroElec();
133 if(regnamesME.size() > 0)
134 {
136 }
137 const std::vector<G4String>& regnamesMSC = theParameters->RegionsPhysics();
138 if(regnamesMSC.size() > 0)
139 {
141 }
142}
G4EmParameters * theParameters
static G4EmParameters * Instance()
const std::vector< G4String > & RegionsPAI() const
const std::vector< G4String > & RegionsPhysics() const
const std::vector< G4String > & RegionsMicroElec() const

References ActivateEmOptions(), ActivateMicroElec(), ActivatePAI(), G4EmParameters::Instance(), G4EmParameters::RegionsMicroElec(), G4EmParameters::RegionsPAI(), G4EmParameters::RegionsPhysics(), and theParameters.

◆ G4EmModelActivator() [2/2]

G4EmModelActivator::G4EmModelActivator ( const G4EmModelActivator )
delete

Member Function Documentation

◆ ActivateEmOptions()

void G4EmModelActivator::ActivateEmOptions ( )
private

Definition at line 146 of file G4EmModelActivator.cc.

147{
148 const std::vector<G4String>& regnamesPhys = theParameters->RegionsPhysics();
149 G4int nreg = regnamesPhys.size();
150 if(0 == nreg) { return; }
151 G4int verbose = theParameters->Verbose() - 1;
152 if(verbose > 0) {
153 G4cout << "### G4EmModelActivator::ActivateEmOptions for " << nreg << " regions"
154 << G4endl;
155 }
156 const std::vector<G4String>& typesPhys = theParameters->TypesPhysics();
157
158 // start configuration of models
161 const G4ParticleDefinition* phot = G4Gamma::Gamma();
164 G4EmConfigurator* em_config = man->EmConfigurator();
165 G4VAtomDeexcitation* adeexc = man->AtomDeexcitation();
167 G4VEmModel* mod;
168
169 // high energy limit for low-energy e+- model of msc
170 G4double mscEnergyLimit = theParameters->MscEnergyLimit();
171
172 // high energy limit for Livermore and Penelope models
173 G4double highEnergyLimit = 1*GeV;
174
175 // general high energy limit
176 G4double highEnergy = theParameters->MaxKinEnergy();
177
178 for(G4int i=0; i<nreg; ++i) {
179 G4String reg = regnamesPhys[i];
180 if(verbose > 0) {
181 G4cout << i << ". region <" << reg << ">; type <" << typesPhys[i] << "> "
182 << G4endl;
183 }
184
185 if(baseName == typesPhys[i]) { continue; }
186
187 if("G4EmStandard" == typesPhys[i]) {
188 G4UrbanMscModel* msc = new G4UrbanMscModel();
189 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
190
191 msc = new G4UrbanMscModel();
192 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
193
194 } else if("G4EmStandard_opt1" == typesPhys[i] || "G4EmStandard_opt2" == typesPhys[i]) {
195 G4UrbanMscModel* msc = new G4UrbanMscModel();
196 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
197
198 msc = new G4UrbanMscModel();
199 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
200
201 } else if("G4EmStandard_opt3" == typesPhys[i]) {
202
203 G4DummyModel* dummy = new G4DummyModel();
204 G4UrbanMscModel* msc = new G4UrbanMscModel();
205 SetMscParameters(elec, msc, typesPhys[i]);
206 em_config->SetExtraEmModel("e-", "msc", msc, reg);
207 FindOrAddProcess(elec, "CoulombScat");
208 em_config->SetExtraEmModel("e-", "CoulombScat", dummy, reg);
209
210 msc = new G4UrbanMscModel();
211 SetMscParameters(posi, msc, typesPhys[i]);
212 em_config->SetExtraEmModel("e+", "msc", msc, reg);
213 FindOrAddProcess(posi, "CoulombScat");
214 em_config->SetExtraEmModel("e+", "CoulombScat", dummy, reg);
215
216 msc = new G4UrbanMscModel();
217 SetMscParameters(prot, msc, typesPhys[i]);
218 em_config->SetExtraEmModel("proton", "msc", msc, reg);
219 FindOrAddProcess(prot, "CoulombScat");
220 em_config->SetExtraEmModel("proton", "CoulombScat", dummy, reg);
221
224 theParameters->SetDeexActiveRegion(reg, true, false, false);
225 }
227 FindOrAddProcess(phot, "Rayl");
228 mod = new G4LivermoreRayleighModel();
229 em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
230 FindOrAddProcess(phot, "phot");
232 em_config->SetExtraEmModel("gamma", "phot", mod, reg);
233 FindOrAddProcess(phot, "compt");
234 mod = new G4KleinNishinaModel();
235 em_config->SetExtraEmModel("gamma", "compt", mod, reg);
236
237 } else if("G4EmStandard_opt4" == typesPhys[i]) {
239 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
240
242 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
243
247 theParameters->SetDeexActiveRegion(reg, true, false, false);
248 }
250
251 FindOrAddProcess(phot, "Rayl");
252 mod = new G4LivermoreRayleighModel();
253 em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
254 FindOrAddProcess(phot, "phot");
256 FindOrAddProcess(phot, "compt");
257 mod = new G4KleinNishinaModel();
258 em_config->SetExtraEmModel("gamma", "compt", mod, reg);
259 mod = new G4LowEPComptonModel();
260 mod->SetHighEnergyLimit(20*MeV);
261 em_config->SetExtraEmModel("gamma", "compt", mod, reg);
262 FindOrAddProcess(phot, "conv");
263 mod = new G4BetheHeitler5DModel();
264 em_config->SetExtraEmModel("gamma", "conv", mod, reg);
265
266 } else if("G4EmStandardGS" == typesPhys[i]) {
268 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
269
271 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
272
273 } else if("G4EmStandardWVI" == typesPhys[i]) {
275 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
276
277 msc = new G4WentzelVIModel();
278 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
279
281
283 theParameters->SetDeexActiveRegion(regnamesPhys[i], true, false, false);
284 }
286
287 } else if("G4EmStandardSS" == typesPhys[i] &&
288 baseName != "G4EmStandard_opt3") {
289 G4EmParticleList emList;
290 for(const auto& particleName : emList.PartNames()) {
291 G4ParticleDefinition* particle = table->FindParticle(particleName);
292 if(particle && 0.0 != particle->GetPDGCharge()) {
293 FindOrAddProcess(particle, "CoulombScat");
295 sc->SetPolarAngleLimit(0.0);
296 sc->SetLocked(true);
297 em_config->SetExtraEmModel(particleName, "CoulombScat", sc, reg);
298 if(particleName == "mu+" || particleName == "mu-") {
299 em_config->SetExtraEmModel(particleName, "muMsc",
300 new G4DummyModel(), reg);
301 } else {
302 em_config->SetExtraEmModel(particleName, "msc",
303 new G4DummyModel(), reg);
304 }
305 }
306 }
308 theParameters->SetDeexActiveRegion(regnamesPhys[i], true, true, true);
309 }
311
312 } else if("G4EmLivermore" == typesPhys[i]) {
313
315 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
316
318 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
319
321 em_config->SetExtraEmModel("gamma", "phot", mod, reg);
322 mod = new G4LivermoreComptonModel();
323 em_config->SetExtraEmModel("gamma", "compt", mod, reg);
325 em_config->SetExtraEmModel("gamma", "conv", mod, reg);
326
327 FindOrAddProcess(phot, "Rayl");
328 mod = new G4LivermoreRayleighModel();
329 em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
330
331 mod = new G4LivermoreIonisationModel();
333 em_config->SetExtraEmModel("e-", "eIoni", mod, reg, 0.0, 0.1*MeV, uf);
335 em_config->SetExtraEmModel("e-", "eBrem", mod, reg, 0.0, highEnergyLimit);
336
340 theParameters->SetDeexActiveRegion(regnamesPhys[i], true, false, false);
341 }
343
344 } else if("G4EmPenelope" == typesPhys[i]) {
345
347 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
348
350 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
351
353 em_config->SetExtraEmModel("gamma", "phot", mod, reg);
354 mod = new G4PenelopeComptonModel();
355 em_config->SetExtraEmModel("gamma", "compt", mod, reg);
357 em_config->SetExtraEmModel("gamma", "conv", mod, reg);
358
359 FindOrAddProcess(phot, "Rayl");
360 mod = new G4PenelopeRayleighModel();
361 em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
362
363 mod = new G4PenelopeIonisationModel();
365 em_config->SetExtraEmModel("e-", "eIoni", mod, reg, 0.0, highEnergyLimit, uf);
367 em_config->SetExtraEmModel("e-", "eBrem", mod, reg, 0.0, highEnergyLimit);
368
369 mod = new G4PenelopeIonisationModel();
370 uf = new G4UniversalFluctuation();
371 em_config->SetExtraEmModel("e+", "eIoni", mod, reg, 0.0, highEnergyLimit, uf);
373 em_config->SetExtraEmModel("e+", "eBrem", mod, reg, 0.0, highEnergyLimit);
374 mod = new G4PenelopeAnnihilationModel();
375 em_config->SetExtraEmModel("e+", "annihil", mod, reg, 0.0, highEnergyLimit);
376
380 theParameters->SetDeexActiveRegion(regnamesPhys[i], true, false, false);
381 }
383
384 } else {
385 if(verbose > 0 && G4Threading::IsMasterThread()) {
386 G4cout << "### G4EmModelActivator::ActivateEmOptions WARNING: \n"
387 << " EM Physics configuration name <" << typesPhys[i]
388 << "> is not known - ignored" << G4endl;
389 }
390 }
391 }
392}
static const G4double reg
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
double G4double
Definition: G4Types.hh:83
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)
void FindOrAddProcess(const G4ParticleDefinition *, const G4String &)
void AddStandardScattering(const G4ParticleDefinition *, G4EmConfigurator *, G4VMscModel *, const G4String &, G4double, G4double, const G4String &)
void SetMscParameters(const G4ParticleDefinition *, G4VMscModel *, const G4String &phys)
void SetNumberOfBinsPerDecade(G4int val)
const std::vector< G4String > & TypesPhysics() const
void SetMscThetaLimit(G4double val)
G4double MscEnergyLimit() const
void SetDeexActiveRegion(const G4String &region, G4bool fdeex, G4bool fauger, G4bool fpixe)
G4int Verbose() const
void DefineRegParamForDeex(G4VAtomDeexcitation *) const
G4double MaxKinEnergy() const
void SetUseMottCorrection(G4bool val)
const std::vector< G4String > & PartNames() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4LossTableManager * Instance()
G4EmConfigurator * EmConfigurator()
G4VAtomDeexcitation * AtomDeexcitation()
G4double GetPDGCharge() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
static G4Positron * Positron()
Definition: G4Positron.cc:93
static G4Proton * Proton()
Definition: G4Proton.cc:92
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:802
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
void SetLocked(G4bool)
Definition: G4VEmModel.hh:884
G4bool IsMasterThread()
Definition: G4Threading.cc:124

References AddStandardScattering(), G4LossTableManager::AtomDeexcitation(), baseName, G4EmParameters::DefineRegParamForDeex(), G4Electron::Electron(), G4LossTableManager::EmConfigurator(), FindOrAddProcess(), G4ParticleTable::FindParticle(), G4cout, G4endl, G4Gamma::Gamma(), G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetPDGCharge(), GeV, G4LossTableManager::Instance(), G4Threading::IsMasterThread(), G4EmParameters::MaxKinEnergy(), MeV, G4EmParameters::MscEnergyLimit(), G4EmParticleList::PartNames(), G4Positron::Positron(), G4Proton::Proton(), reg, G4EmParameters::RegionsPhysics(), G4EmParameters::SetDeexActiveRegion(), G4EmConfigurator::SetExtraEmModel(), G4VEmModel::SetHighEnergyLimit(), G4VEmModel::SetLocked(), SetMscParameters(), G4EmParameters::SetMscThetaLimit(), G4EmParameters::SetNumberOfBinsPerDecade(), G4VEmModel::SetPolarAngleLimit(), G4EmParameters::SetUseMottCorrection(), theParameters, G4EmParameters::TypesPhysics(), and G4EmParameters::Verbose().

Referenced by G4EmModelActivator().

◆ ActivateMicroElec()

void G4EmModelActivator::ActivateMicroElec ( )
private

Definition at line 502 of file G4EmModelActivator.cc.

503{
504 const std::vector<G4String> regnamesME = theParameters->RegionsMicroElec();
505 G4int nreg = regnamesME.size();
506 if(0 == nreg)
507 {
508 return;
509 }
510 G4int verbose = theParameters->Verbose() - 1;
511 if(verbose > 0)
512 {
513 G4cout << "### G4EmModelActivator::ActivateMicroElec for " << nreg
514 << " regions" << G4endl;
515 }
517
521 G4ProcessManager* eman = elec->GetProcessManager();
522 G4ProcessManager* pman = prot->GetProcessManager();
523 G4ProcessManager* iman = gion->GetProcessManager();
524
525 G4bool emsc = HasMsc(eman);
526
527 // MicroElec elastic is not active in the world
528 G4MicroElecElastic* eElasProc =
529 new G4MicroElecElastic("e-G4MicroElecElastic");
530 eman->AddDiscreteProcess(eElasProc);
531
532 G4MicroElecInelastic* eInelProc =
533 new G4MicroElecInelastic("e-G4MicroElecInelastic");
534 eman->AddDiscreteProcess(eInelProc);
535
536 G4MicroElecInelastic* pInelProc =
537 new G4MicroElecInelastic("p_G4MicroElecInelastic");
538 pman->AddDiscreteProcess(pInelProc);
539
540 G4MicroElecInelastic* iInelProc =
541 new G4MicroElecInelastic("ion_G4MicroElecInelastic");
542 iman->AddDiscreteProcess(iInelProc);
543
544 // start configuration of models
545 G4EmConfigurator* em_config = man->EmConfigurator();
546 G4VEmModel* mod;
547
548 // limits for MicroElec applicability
549 G4double elowest = 16.7 * eV;
550 G4double elimel = 100 * MeV;
551 G4double elimin = 9 * MeV;
552 G4double pmin = 50 * keV;
553 G4double pmax = 99.9 * MeV;
554
555 G4LowECapture* ecap = new G4LowECapture(elowest);
556 eman->AddDiscreteProcess(ecap);
557
558 for(G4int i = 0; i < nreg; ++i)
559 {
560
561 G4String reg = regnamesME[i];
562 G4cout << "### MicroElec models are activated for G4Region " << reg
563 << G4endl
564 << " Energy limits for e- elastic: " << elowest/eV << " eV - "
565 << elimel/MeV << " MeV"
566 << G4endl
567 << " Energy limits for e- inelastic: " << elowest/eV << " eV - "
568 << elimin/MeV << " MeV"
569 << G4endl
570 << " Energy limits for hadrons/ions: " << pmin/MeV << " MeV - "
571 << pmax/MeV << " MeV"
572 << G4endl;
573
574 // e-
575 if(emsc)
576 {
577 G4UrbanMscModel* msc = new G4UrbanMscModel();
578 msc->SetActivationLowEnergyLimit(elimel);
579 em_config->SetExtraEmModel("e-", "msc", msc, reg);
580 }
581 else
582 {
583 mod = new G4DummyModel();
584 em_config->SetExtraEmModel("e-", "CoulombScat", mod, reg, 0.0, elimel);
585 }
586
587 mod = new G4MicroElecElasticModel();
588 em_config->SetExtraEmModel("e-",
589 "e-G4MicroElecElastic",
590 mod,
591 reg,
592 elowest,
593 elimin);
594
595 mod = new G4MollerBhabhaModel();
596 mod->SetActivationLowEnergyLimit(elimin);
597 em_config->SetExtraEmModel("e-",
598 "eIoni",
599 mod,
600 reg,
601 0.0,
602 10 * TeV,
604
605 mod = new G4MicroElecInelasticModel();
606 em_config->SetExtraEmModel("e-",
607 "e-G4MicroElecInelastic",
608 mod,
609 reg,
610 elowest,
611 elimin);
612
613 // proton
614 mod = new G4BraggModel();
616 em_config->SetExtraEmModel("proton",
617 "hIoni",
618 mod,
619 reg,
620 0.0,
621 2 * MeV,
623
624 mod = new G4BetheBlochModel();
626 em_config->SetExtraEmModel("proton",
627 "hIoni",
628 mod,
629 reg,
630 2 * MeV,
631 10 * TeV,
633
634 mod = new G4MicroElecInelasticModel();
635 em_config->SetExtraEmModel("proton",
636 "p_G4MicroElecInelastic",
637 mod,
638 reg,
639 pmin,
640 pmax);
641
642 // ions
643 mod = new G4BraggIonModel();
645 em_config->SetExtraEmModel("GenericIon",
646 "ionIoni",
647 mod,
648 reg,
649 0.0,
650 2 * MeV,
651 new G4IonFluctuations());
652
653 mod = new G4BetheBlochModel();
655 em_config->SetExtraEmModel("GenericIon",
656 "ionIoni",
657 mod,
658 reg,
659 2 * MeV,
660 10 * TeV,
661 new G4IonFluctuations());
662
663 mod = new G4MicroElecInelasticModel();
664 em_config->SetExtraEmModel("GenericIon",
665 "ion_G4MicroElecInelastic",
666 mod,
667 reg,
668 pmin,
669 pmax);
670 }
671}
static constexpr double keV
Definition: G4SIunits.hh:202
static constexpr double eV
Definition: G4SIunits.hh:201
static constexpr double TeV
Definition: G4SIunits.hh:204
bool G4bool
Definition: G4Types.hh:86
G4bool HasMsc(G4ProcessManager *) const
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
G4ProcessManager * GetProcessManager() const
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:788
void SetActivationHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:781

References G4ProcessManager::AddDiscreteProcess(), G4Electron::Electron(), G4LossTableManager::EmConfigurator(), eV, G4cout, G4endl, G4GenericIon::GenericIon(), G4ParticleDefinition::GetProcessManager(), HasMsc(), G4LossTableManager::Instance(), keV, MeV, G4Proton::Proton(), reg, G4EmParameters::RegionsMicroElec(), G4VEmModel::SetActivationHighEnergyLimit(), G4VEmModel::SetActivationLowEnergyLimit(), G4EmConfigurator::SetExtraEmModel(), TeV, theParameters, and G4EmParameters::Verbose().

Referenced by G4EmModelActivator().

◆ ActivatePAI()

void G4EmModelActivator::ActivatePAI ( )
private

Definition at line 396 of file G4EmModelActivator.cc.

397{
398 const std::vector<G4String> regnamesPAI = theParameters->RegionsPAI();
399 G4int nreg = regnamesPAI.size();
400 if(0 == nreg) { return; }
401 G4int verbose = theParameters->Verbose() - 1;
402 if(verbose > 0) {
403 G4cout << "### G4EmModelActivator::ActivatePAI for " << nreg << " regions"
404 << G4endl;
405 }
406 const std::vector<G4String> particlesPAI = theParameters->ParticlesPAI();
407 const std::vector<G4String> typesPAI = theParameters->TypesPAI();
408
409 const std::vector<G4VEnergyLossProcess*>& v = G4LossTableManager::Instance()
411
413
419
420 for(G4int i = 0; i < nreg; ++i) {
421 const G4ParticleDefinition* p = nullptr;
422 if(particlesPAI[i] != "all") {
423 p = G4ParticleTable::GetParticleTable()->FindParticle(particlesPAI[i]);
424 if(!p) {
425 G4cout << "### WARNING: ActivatePAI::FindParticle fails to find "
426 << particlesPAI[i] << G4endl;
427 continue;
428 }
429 }
430 const G4Region* r = regionStore->GetRegion(regnamesPAI[i], false);
431 if(!r) {
432 G4cout << "### WARNING: ActivatePAI::GetRegion fails to find "
433 << regnamesPAI[i] << G4endl;
434 continue;
435 }
436
437 G4String name = "hIoni";
438 if(p == elec || p == posi)
439 { name = "eIoni"; }
440 else if (p == mupl || p == mumi)
441 { name = "muIoni"; }
442 else if (p == gion)
443 { name = "ionIoni"; }
444
445 for(auto proc : v) {
446
447 if(!proc->IsIonisationProcess()) { continue; }
448
449 G4String namep = proc->GetProcessName();
450 if(p) {
451 if(name != namep) { continue; }
452 } else {
453 if(namep != "hIoni" && namep != "muIoni" &&
454 namep != "eIoni" && namep != "ionIoni")
455 { continue; }
456 }
457 G4double emin = 50*CLHEP::keV;
458 if(namep == "eIoni") emin = 110*CLHEP::eV;
459 else if(namep == "muIoni") emin = 5*CLHEP::keV;
460
461 G4VEmModel* em = nullptr;
462 G4VEmFluctuationModel* fm = nullptr;
463 if(typesPAI[i] == "PAIphoton" || typesPAI[i] == "pai_photon") {
464 G4PAIPhotModel* mod = new G4PAIPhotModel(p,"PAIPhotModel");
465 em = mod;
466 fm = mod;
467 } else {
468 G4PAIModel* mod = new G4PAIModel(p,"PAIModel");
469 em = mod;
470 fm = mod;
471 }
472 // added PAI model above low energy limit
473 em->SetLowEnergyLimit(emin);
474 proc->AddEmModel(0, em, fm, r);
475
476 // added low energy model
477 if(namep == "eIoni") {
478 em = new G4MollerBhabhaModel();
479 fm = new G4UniversalFluctuation();
480 } else if(namep == "ionIoni") {
481 em = new G4BraggIonModel();
482 fm = new G4IonFluctuations();
483 } else {
484 em = new G4BraggModel();
485 fm = new G4UniversalFluctuation();
486 }
487 em->SetHighEnergyLimit(emin);
488 proc->AddEmModel(-1, em, fm, r);
489
490 if(verbose > 0) {
491 G4cout << "### G4EmModelActivator: add <" << typesPAI[i]
492 << "> model for " << particlesPAI[i]
493 << " in the " << regnamesPAI[i]
494 << " Emin(keV)= " << emin/CLHEP::keV << G4endl;
495 }
496 }
497 }
498}
const std::vector< G4String > & ParticlesPAI() const
const std::vector< G4String > & TypesPAI() const
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:99
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:98
static G4RegionStore * GetInstance()
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:774
static constexpr double keV
static constexpr double eV
const char * name(G4int ptype)

References G4Electron::Electron(), CLHEP::eV, G4ParticleTable::FindParticle(), G4cout, G4endl, G4GenericIon::GenericIon(), G4LossTableManager::GetEnergyLossProcessVector(), G4RegionStore::GetInstance(), G4ParticleTable::GetParticleTable(), G4RegionStore::GetRegion(), G4LossTableManager::Instance(), CLHEP::keV, G4MuonMinus::MuonMinus(), G4MuonPlus::MuonPlus(), G4InuclParticleNames::name(), G4EmParameters::ParticlesPAI(), G4Positron::Positron(), G4EmParameters::RegionsPAI(), G4VEmModel::SetHighEnergyLimit(), G4VEmModel::SetLowEnergyLimit(), theParameters, G4EmParameters::TypesPAI(), and G4EmParameters::Verbose().

Referenced by G4EmModelActivator().

◆ AddStandardScattering()

void G4EmModelActivator::AddStandardScattering ( const G4ParticleDefinition part,
G4EmConfigurator em_config,
G4VMscModel mscmod,
const G4String reg,
G4double  e1,
G4double  e2,
const G4String type 
)
private

Definition at line 693 of file G4EmModelActivator.cc.

699{
701
702 // low-energy msc model
703 SetMscParameters(part, mscmod, type);
704 em_config->SetExtraEmModel(pname, "msc", mscmod, reg, 0.0, e1);
705
706 // high energy msc model
708 SetMscParameters(part, msc, type);
709 em_config->SetExtraEmModel(pname, "msc", msc, reg, e1, e2);
710
711 // high energy single scattering model
712 FindOrAddProcess(part, "CoulombScat");
715 mod->SetLocked(true);
716 em_config->SetExtraEmModel(pname, "CoulombScat", mod, reg, 0.0, e2);
717}
static const G4double e1[44]
static const G4double e2[44]
const G4String & GetParticleName() const
string pname
Definition: eplot.py:33

References e1, e2, FindOrAddProcess(), G4ParticleDefinition::GetParticleName(), eplot::pname, reg, G4VEmModel::SetActivationLowEnergyLimit(), G4EmConfigurator::SetExtraEmModel(), G4VEmModel::SetLocked(), and SetMscParameters().

Referenced by ActivateEmOptions().

◆ FindOrAddProcess()

void G4EmModelActivator::FindOrAddProcess ( const G4ParticleDefinition part,
const G4String name 
)
private

Definition at line 747 of file G4EmModelActivator.cc.

749{
752 G4int nproc = pm->GetProcessListLength();
753 for(G4int i = 0; i<nproc; ++i) {
754 if(((*pv)[i])->GetProcessName() == name) { return; }
755 }
756 if(name == "CoulombScat") {
758 cs->SetEmModel(new G4DummyModel());
759 pm->AddDiscreteProcess(cs);
760 } else if(name == "Rayl") {
762 rs->SetEmModel(new G4DummyModel());
763 pm->AddDiscreteProcess(rs);
764 }
765}
G4int GetProcessListLength() const
G4ProcessVector * GetProcessList() const
void SetEmModel(G4VEmModel *, G4int index=0)

References G4ProcessManager::AddDiscreteProcess(), G4ProcessManager::GetProcessList(), G4ProcessManager::GetProcessListLength(), G4ParticleDefinition::GetProcessManager(), G4InuclParticleNames::name(), and G4VEmProcess::SetEmModel().

Referenced by ActivateEmOptions(), and AddStandardScattering().

◆ HasMsc()

G4bool G4EmModelActivator::HasMsc ( G4ProcessManager pm) const
private

Definition at line 675 of file G4EmModelActivator.cc.

676{
677 G4bool res = false;
679 G4int nproc = pm->GetProcessListLength();
680 for(G4int i = 0; i < nproc; ++i)
681 {
682 if(((*pv)[i])->GetProcessSubType() == fMultipleScattering)
683 {
684 res = true;
685 break;
686 }
687 }
688 return res;
689}
@ fMultipleScattering

References fMultipleScattering, G4ProcessManager::GetProcessList(), and G4ProcessManager::GetProcessListLength().

Referenced by ActivateMicroElec().

◆ operator=()

G4EmModelActivator & G4EmModelActivator::operator= ( const G4EmModelActivator right)
delete

◆ SetMscParameters()

void G4EmModelActivator::SetMscParameters ( const G4ParticleDefinition part,
G4VMscModel msc,
const G4String phys 
)
private

Definition at line 721 of file G4EmModelActivator.cc.

723{
724 if(part == G4Electron::Electron() || part == G4Positron::Positron()) {
725 if(phys == "G4EmStandard_opt1" || phys == "G4EmStandard_opt2") {
726 msc->SetRangeFactor(0.2);
728 } else if(phys == "G4EmStandard_opt3") {
730 } else if(phys == "G4EmStandard_opt4" || phys == "G4EmLivermore" || phys == "G4EmPenelope") {
731 msc->SetRangeFactor(0.08);
733 msc->SetSkin(3);
734 } else if(phys == "G4EmStandardGS") {
735 msc->SetRangeFactor(0.06);
736 }
737 } else {
738 if(phys != "G4EmStandard" && phys != "G4EmStandard_opt1" && phys != "G4EmStandard_opt2") {
739 msc->SetLateralDisplasmentFlag(true);
740 }
741 }
742 msc->SetLocked(true);
743}
@ fMinimal
@ fUseSafetyPlus
@ fUseDistanceToBoundary
void SetSkin(G4double)
Definition: G4VMscModel.hh:223
void SetRangeFactor(G4double)
Definition: G4VMscModel.hh:230
void SetLateralDisplasmentFlag(G4bool val)
Definition: G4VMscModel.hh:216
void SetStepLimitType(G4MscStepLimitType)
Definition: G4VMscModel.hh:258

References G4Electron::Electron(), fMinimal, fUseDistanceToBoundary, fUseSafetyPlus, G4Positron::Positron(), G4VMscModel::SetLateralDisplasmentFlag(), G4VEmModel::SetLocked(), G4VMscModel::SetRangeFactor(), G4VMscModel::SetSkin(), and G4VMscModel::SetStepLimitType().

Referenced by ActivateEmOptions(), and AddStandardScattering().

Field Documentation

◆ baseName

G4String G4EmModelActivator::baseName
private

Definition at line 86 of file G4EmModelActivator.hh.

Referenced by ActivateEmOptions().

◆ theParameters

G4EmParameters* G4EmModelActivator::theParameters
private

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