Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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"
101
108
116
117#include "G4SystemOfUnits.hh"
118#include "G4Threading.hh"
119#include <vector>
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122
124 : baseName(emphys)
125{
126 theParameters = G4EmParameters::Instance();
127
128 const std::vector<G4String>& regnamesPAI = theParameters->RegionsPAI();
129 if(regnamesPAI.size() > 0)
130 {
131 ActivatePAI();
132 }
133 const std::vector<G4String>& regnamesME = theParameters->RegionsMicroElec();
134 if(regnamesME.size() > 0)
135 {
136 ActivateMicroElec();
137 }
138 const std::vector<G4String>& regnamesMSC = theParameters->RegionsPhysics();
139 if(regnamesMSC.size() > 0)
140 {
141 ActivateEmOptions();
142 }
143}
144
145//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
146
147void G4EmModelActivator::ActivateEmOptions()
148{
149 const std::vector<G4String>& regnamesPhys = theParameters->RegionsPhysics();
150 std::size_t nreg = regnamesPhys.size();
151 if(0 == nreg) { return; }
152 G4int verbose = theParameters->Verbose() - 1;
153 if(verbose > 0) {
154 G4cout << "### G4EmModelActivator::ActivateEmOptions for " << nreg << " regions"
155 << G4endl;
156 }
157 const std::vector<G4String>& typesPhys = theParameters->TypesPhysics();
158
159 // start configuration of models
162 const G4ParticleDefinition* phot = G4Gamma::Gamma();
165 G4EmConfigurator* em_config = man->EmConfigurator();
166 G4VAtomDeexcitation* adeexc = man->AtomDeexcitation();
168 G4VEmModel* mod;
169
170 // high energy limit for low-energy e+- model of msc
171 G4double mscEnergyLimit = theParameters->MscEnergyLimit();
172
173 // high energy limit for Livermore and Penelope models
174 G4double highEnergyLimit = 1*GeV;
175
176 // general high energy limit
177 G4double highEnergy = theParameters->MaxKinEnergy();
178
179 for(std::size_t i=0; i<nreg; ++i) {
180 G4String reg = regnamesPhys[i];
181 if(verbose > 0) {
182 G4cout << i << ". region <" << reg << ">; type <" << typesPhys[i] << "> "
183 << G4endl;
184 }
185
186 if(baseName == typesPhys[i]) { continue; }
187
188 if("G4EmStandard" == typesPhys[i]) {
189 G4UrbanMscModel* msc = new G4UrbanMscModel();
190 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
191
192 msc = new G4UrbanMscModel();
193 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
194
195 } else if("G4EmStandard_opt1" == typesPhys[i] || "G4EmStandard_opt2" == typesPhys[i]) {
196 G4UrbanMscModel* msc = new G4UrbanMscModel();
197 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
198
199 msc = new G4UrbanMscModel();
200 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
201
202 } else if("G4EmStandard_opt3" == typesPhys[i]) {
203
204 G4DummyModel* dummy = new G4DummyModel();
205 G4UrbanMscModel* msc = new G4UrbanMscModel();
206 SetMscParameters(elec, msc, typesPhys[i]);
207 em_config->SetExtraEmModel("e-", "msc", msc, reg);
208 FindOrAddProcess(elec, "CoulombScat");
209 em_config->SetExtraEmModel("e-", "CoulombScat", dummy, reg);
210
211 msc = new G4UrbanMscModel();
212 SetMscParameters(posi, msc, typesPhys[i]);
213 em_config->SetExtraEmModel("e+", "msc", msc, reg);
214 FindOrAddProcess(posi, "CoulombScat");
215 em_config->SetExtraEmModel("e+", "CoulombScat", dummy, reg);
216
217 msc = new G4UrbanMscModel();
218 SetMscParameters(prot, msc, typesPhys[i]);
219 em_config->SetExtraEmModel("proton", "msc", msc, reg);
220 FindOrAddProcess(prot, "CoulombScat");
221 em_config->SetExtraEmModel("proton", "CoulombScat", dummy, reg);
222
223 theParameters->SetNumberOfBinsPerDecade(20);
225 theParameters->SetDeexActiveRegion(reg, true, false, false);
226 }
227 theParameters->DefineRegParamForDeex(adeexc);
228 FindOrAddProcess(phot, "Rayl");
229 mod = new G4LivermoreRayleighModel();
230 em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
231 FindOrAddProcess(phot, "phot");
233 em_config->SetExtraEmModel("gamma", "phot", mod, reg);
234 FindOrAddProcess(phot, "compt");
235 mod = new G4KleinNishinaModel();
236 em_config->SetExtraEmModel("gamma", "compt", mod, reg);
237
238 } else if("G4EmStandard_opt4" == typesPhys[i]) {
240 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
241
243 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
244
245 theParameters->SetNumberOfBinsPerDecade(20);
246 theParameters->SetUseMottCorrection(true);
248 theParameters->SetDeexActiveRegion(reg, true, false, false);
249 }
250 theParameters->DefineRegParamForDeex(adeexc);
251
252 FindOrAddProcess(phot, "Rayl");
253 mod = new G4LivermoreRayleighModel();
254 em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
255 FindOrAddProcess(phot, "phot");
257 FindOrAddProcess(phot, "compt");
258 mod = new G4KleinNishinaModel();
259 em_config->SetExtraEmModel("gamma", "compt", mod, reg);
260 mod = new G4LowEPComptonModel();
261 mod->SetHighEnergyLimit(20*MeV);
262 em_config->SetExtraEmModel("gamma", "compt", mod, reg);
263 FindOrAddProcess(phot, "conv");
264 mod = new G4BetheHeitler5DModel();
265 em_config->SetExtraEmModel("gamma", "conv", mod, reg);
266
267 } else if("G4EmStandardGS" == typesPhys[i]) {
269 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
270
272 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
273
274 } else if("G4EmStandardWVI" == typesPhys[i]) {
276 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
277
278 msc = new G4WentzelVIModel();
279 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
280
281 theParameters->SetMscThetaLimit(0.15);
282
284 theParameters->SetDeexActiveRegion(regnamesPhys[i], true, false, false);
285 }
286 theParameters->DefineRegParamForDeex(adeexc);
287
288 } else if("G4EmStandardSS" == typesPhys[i] &&
289 baseName != "G4EmStandard_opt3") {
290 G4EmParticleList emList;
291 for(const auto& particleName : emList.PartNames()) {
292 G4ParticleDefinition* particle = table->FindParticle(particleName);
293 if(particle && 0.0 != particle->GetPDGCharge()) {
294 FindOrAddProcess(particle, "CoulombScat");
296 sc->SetPolarAngleLimit(0.0);
297 sc->SetLocked(true);
298 em_config->SetExtraEmModel(particleName, "CoulombScat", sc, reg);
299 if(particleName == "mu+" || particleName == "mu-") {
300 em_config->SetExtraEmModel(particleName, "muMsc",
301 new G4DummyModel(), reg);
302 } else {
303 em_config->SetExtraEmModel(particleName, "msc",
304 new G4DummyModel(), reg);
305 }
306 }
307 }
309 theParameters->SetDeexActiveRegion(regnamesPhys[i], true, true, true);
310 }
311 theParameters->DefineRegParamForDeex(adeexc);
312
313 } else if("G4EmLivermore" == typesPhys[i]) {
314
316 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
317
319 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
320
322 em_config->SetExtraEmModel("gamma", "phot", mod, reg);
323 mod = new G4LivermoreComptonModel();
324 em_config->SetExtraEmModel("gamma", "compt", mod, reg);
326 em_config->SetExtraEmModel("gamma", "conv", mod, reg);
327
328 FindOrAddProcess(phot, "Rayl");
329 mod = new G4LivermoreRayleighModel();
330 em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
331
332 mod = new G4LivermoreIonisationModel();
334 em_config->SetExtraEmModel("e-", "eIoni", mod, reg, 0.0, 0.1*MeV, uf);
336 em_config->SetExtraEmModel("e-", "eBrem", mod, reg, 0.0, highEnergyLimit);
337
338 theParameters->SetNumberOfBinsPerDecade(20);
339 theParameters->SetUseMottCorrection(true);
341 theParameters->SetDeexActiveRegion(regnamesPhys[i], true, false, false);
342 }
343 theParameters->DefineRegParamForDeex(adeexc);
344
345 } else if("G4EmPenelope" == typesPhys[i]) {
346
348 AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
349
351 AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy, typesPhys[i]);
352
354 em_config->SetExtraEmModel("gamma", "phot", mod, reg);
355 mod = new G4PenelopeComptonModel();
356 em_config->SetExtraEmModel("gamma", "compt", mod, reg);
358 em_config->SetExtraEmModel("gamma", "conv", mod, reg);
359
360 FindOrAddProcess(phot, "Rayl");
361 mod = new G4PenelopeRayleighModel();
362 em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
363
364 mod = new G4PenelopeIonisationModel();
366 em_config->SetExtraEmModel("e-", "eIoni", mod, reg, 0.0, highEnergyLimit, uf);
368 em_config->SetExtraEmModel("e-", "eBrem", mod, reg, 0.0, highEnergyLimit);
369
370 mod = new G4PenelopeIonisationModel();
371 uf = new G4UniversalFluctuation();
372 em_config->SetExtraEmModel("e+", "eIoni", mod, reg, 0.0, highEnergyLimit, uf);
374 em_config->SetExtraEmModel("e+", "eBrem", mod, reg, 0.0, highEnergyLimit);
375 mod = new G4PenelopeAnnihilationModel();
376 em_config->SetExtraEmModel("e+", "annihil", mod, reg, 0.0, highEnergyLimit);
377
378 theParameters->SetNumberOfBinsPerDecade(20);
379 theParameters->SetUseMottCorrection(true);
381 theParameters->SetDeexActiveRegion(regnamesPhys[i], true, false, false);
382 }
383 theParameters->DefineRegParamForDeex(adeexc);
384
385 } else {
386 if(verbose > 0 && G4Threading::IsMasterThread()) {
387 G4cout << "### G4EmModelActivator::ActivateEmOptions WARNING: \n"
388 << " EM Physics configuration name <" << typesPhys[i]
389 << "> is not known - ignored" << G4endl;
390 }
391 }
392 }
393}
394
395//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
396
397void G4EmModelActivator::ActivatePAI()
398{
399 const std::vector<G4String> regnamesPAI = theParameters->RegionsPAI();
400 std::size_t nreg = regnamesPAI.size();
401 if(0 == nreg) { return; }
402 G4int verbose = theParameters->Verbose() - 1;
403 if(verbose > 0) {
404 G4cout << "### G4EmModelActivator::ActivatePAI for " << nreg << " regions"
405 << G4endl;
406 }
407 const std::vector<G4String> particlesPAI = theParameters->ParticlesPAI();
408 const std::vector<G4String> typesPAI = theParameters->TypesPAI();
409
410 const std::vector<G4VEnergyLossProcess*>& v = G4LossTableManager::Instance()
412
414
420
421 for(std::size_t i = 0; i < nreg; ++i) {
422 const G4ParticleDefinition* p = nullptr;
423 if(particlesPAI[i] != "all") {
424 p = G4ParticleTable::GetParticleTable()->FindParticle(particlesPAI[i]);
425 if(!p) {
426 G4cout << "### WARNING: ActivatePAI::FindParticle fails to find "
427 << particlesPAI[i] << G4endl;
428 continue;
429 }
430 }
431 const G4Region* r = regionStore->GetRegion(regnamesPAI[i], false);
432 if(!r) {
433 G4cout << "### WARNING: ActivatePAI::GetRegion fails to find "
434 << regnamesPAI[i] << G4endl;
435 continue;
436 }
437
438 G4String name = "hIoni";
439 if(p == elec || p == posi)
440 { name = "eIoni"; }
441 else if (p == mupl || p == mumi)
442 { name = "muIoni"; }
443 else if (p == gion)
444 { name = "ionIoni"; }
445
446 for(auto proc : v) {
447
448 if(!proc->IsIonisationProcess()) { continue; }
449
450 G4String namep = proc->GetProcessName();
451 if(p) {
452 if(name != namep) { continue; }
453 } else {
454 if(namep != "hIoni" && namep != "muIoni" &&
455 namep != "eIoni" && namep != "ionIoni")
456 { continue; }
457 }
458 G4double emin = 50*CLHEP::keV;
459 if(namep == "eIoni") emin = 110*CLHEP::eV;
460 else if(namep == "muIoni") emin = 5*CLHEP::keV;
461
462 G4VEmModel* em = nullptr;
463 G4VEmFluctuationModel* fm = nullptr;
464 if(typesPAI[i] == "PAIphoton" || typesPAI[i] == "pai_photon") {
465 G4PAIPhotModel* mod = new G4PAIPhotModel(p,"PAIPhotModel");
466 em = mod;
467 fm = mod;
468 } else {
469 G4PAIModel* mod = new G4PAIModel(p,"PAIModel");
470 em = mod;
471 fm = mod;
472 }
473 // first added the default model for world
474 // second added low energy model below PAI threshold in the region
475 // finally added PAI for the region
476 G4VEmModel* em0 = nullptr;
477 G4VEmFluctuationModel* fm0 = nullptr;
478 if(namep == "eIoni") {
479 fm0 = new G4UniversalFluctuation();
480 proc->SetEmModel(new G4MollerBhabhaModel());
481 proc->SetFluctModel(fm0);
482 em0 = new G4MollerBhabhaModel();
483 } else if(namep == "ionIoni") {
484 fm0 = new G4IonFluctuations();
485 proc->SetEmModel(new G4LindhardSorensenIonModel());
486 proc->SetFluctModel(fm0);
487 em0 = new G4LindhardSorensenIonModel();
488 } else {
489 fm0 = new G4UniversalFluctuation();
490 proc->SetEmModel(new G4BraggModel());
491 proc->SetEmModel(new G4BetheBlochModel());
492 proc->SetFluctModel(fm0);
493 em0 = new G4BraggModel();
494 }
495 em0->SetHighEnergyLimit(emin);
496 proc->AddEmModel(-1, em0, fm0, r);
497 em->SetLowEnergyLimit(emin);
498 proc->AddEmModel(-1, em, fm, r);
499 if(verbose > 0) {
500 G4cout << "### G4EmModelActivator: add <" << typesPAI[i]
501 << "> model for " << particlesPAI[i]
502 << " in the " << regnamesPAI[i]
503 << " Emin(keV)= " << emin/CLHEP::keV << G4endl;
504 }
505 }
506 }
507}
508
509//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
510
511void G4EmModelActivator::ActivateMicroElec()
512{
513 const std::vector<G4String> regnamesME = theParameters->RegionsMicroElec();
514 std::size_t nreg = regnamesME.size();
515 if(0 == nreg)
516 {
517 return;
518 }
519 G4int verbose = theParameters->Verbose() - 1;
520 if(verbose > 0)
521 {
522 G4cout << "### G4EmModelActivator::ActivateMicroElec for " << nreg
523 << " regions" << G4endl;
524 }
526
530 G4ProcessManager* eman = elec->GetProcessManager();
531 G4ProcessManager* pman = prot->GetProcessManager();
532 G4ProcessManager* iman = gion->GetProcessManager();
533
534 G4bool emsc = HasMsc(eman);
535
536 // MicroElec elastic is not active in the world
537 G4MicroElecElastic* eElasProc =
538 new G4MicroElecElastic("e-G4MicroElecElastic");
539 eman->AddDiscreteProcess(eElasProc);
540
541 G4MicroElecInelastic* eInelProc =
542 new G4MicroElecInelastic("e-G4MicroElecInelastic");
543 eman->AddDiscreteProcess(eInelProc);
544
545 G4MicroElecInelastic* pInelProc =
546 new G4MicroElecInelastic("p_G4MicroElecInelastic");
547 pman->AddDiscreteProcess(pInelProc);
548
549 G4MicroElecInelastic* iInelProc =
550 new G4MicroElecInelastic("ion_G4MicroElecInelastic");
551 iman->AddDiscreteProcess(iInelProc);
552
553 // start configuration of models
554 G4EmConfigurator* em_config = man->EmConfigurator();
555 G4VEmModel* mod;
556
557 // limits for MicroElec applicability
558 G4double elowest = 16.7 * eV;
559 G4double elimel = 100 * MeV;
560 G4double elimin = 9 * MeV;
561 G4double pmin = 50 * keV;
562 G4double pmax = 99.9 * MeV;
563
564 G4LowECapture* ecap = new G4LowECapture(elowest);
565 eman->AddDiscreteProcess(ecap);
566
567 for(std::size_t i = 0; i < nreg; ++i)
568 {
569
570 G4String reg = regnamesME[i];
571 G4cout << "### MicroElec models are activated for G4Region " << reg
572 << G4endl
573 << " Energy limits for e- elastic: " << elowest/eV << " eV - "
574 << elimel/MeV << " MeV"
575 << G4endl
576 << " Energy limits for e- inelastic: " << elowest/eV << " eV - "
577 << elimin/MeV << " MeV"
578 << G4endl
579 << " Energy limits for hadrons/ions: " << pmin/MeV << " MeV - "
580 << pmax/MeV << " MeV"
581 << G4endl;
582
583 // e-
584 if(emsc)
585 {
586 G4UrbanMscModel* msc = new G4UrbanMscModel();
587 msc->SetActivationLowEnergyLimit(elimel);
588 em_config->SetExtraEmModel("e-", "msc", msc, reg);
589 }
590 else
591 {
592 mod = new G4DummyModel();
593 em_config->SetExtraEmModel("e-", "CoulombScat", mod, reg, 0.0, elimel);
594 }
595
596 mod = new G4MicroElecElasticModel();
597 em_config->SetExtraEmModel("e-",
598 "e-G4MicroElecElastic",
599 mod,
600 reg,
601 elowest,
602 elimin);
603
604 mod = new G4MollerBhabhaModel();
605 mod->SetActivationLowEnergyLimit(elimin);
606 em_config->SetExtraEmModel("e-",
607 "eIoni",
608 mod,
609 reg,
610 0.0,
611 10 * TeV,
613
614 mod = new G4MicroElecInelasticModel();
615 em_config->SetExtraEmModel("e-",
616 "e-G4MicroElecInelastic",
617 mod,
618 reg,
619 elowest,
620 elimin);
621
622 // proton
623 mod = new G4BraggModel();
625 em_config->SetExtraEmModel("proton",
626 "hIoni",
627 mod,
628 reg,
629 0.0,
630 2 * MeV,
632
633 mod = new G4BetheBlochModel();
635 em_config->SetExtraEmModel("proton",
636 "hIoni",
637 mod,
638 reg,
639 2 * MeV,
640 10 * TeV,
642
643 mod = new G4MicroElecInelasticModel();
644 em_config->SetExtraEmModel("proton",
645 "p_G4MicroElecInelastic",
646 mod,
647 reg,
648 pmin,
649 pmax);
650
651 // ions
652 mod = new G4BraggIonModel();
654 em_config->SetExtraEmModel("GenericIon",
655 "ionIoni",
656 mod,
657 reg,
658 0.0,
659 2 * MeV,
660 new G4IonFluctuations());
661
662 mod = new G4BetheBlochModel();
664 em_config->SetExtraEmModel("GenericIon",
665 "ionIoni",
666 mod,
667 reg,
668 2 * MeV,
669 10 * TeV,
670 new G4IonFluctuations());
671
672 mod = new G4MicroElecInelasticModel();
673 em_config->SetExtraEmModel("GenericIon",
674 "ion_G4MicroElecInelastic",
675 mod,
676 reg,
677 pmin,
678 pmax);
679 }
680}
681
682//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
683
684G4bool G4EmModelActivator::HasMsc(G4ProcessManager* pm) const
685{
686 G4bool res = false;
688 G4int nproc = pm->GetProcessListLength();
689 for(G4int i = 0; i < nproc; ++i)
690 {
691 if(((*pv)[i])->GetProcessSubType() == fMultipleScattering)
692 {
693 res = true;
694 break;
695 }
696 }
697 return res;
698}
699
700//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
701
702void G4EmModelActivator::AddStandardScattering(const G4ParticleDefinition* part,
703 G4EmConfigurator* em_config,
704 G4VMscModel* mscmod,
705 const G4String& reg,
706 G4double e1, G4double e2,
707 const G4String& type)
708{
709 G4String pname = part->GetParticleName();
710
711 // low-energy msc model
712 SetMscParameters(part, mscmod, type);
713 em_config->SetExtraEmModel(pname, "msc", mscmod, reg, 0.0, e1);
714
715 // high energy msc model
717 SetMscParameters(part, msc, type);
718 em_config->SetExtraEmModel(pname, "msc", msc, reg, e1, e2);
719
720 // high energy single scattering model
721 FindOrAddProcess(part, "CoulombScat");
724 mod->SetLocked(true);
725 em_config->SetExtraEmModel(pname, "CoulombScat", mod, reg, 0.0, e2);
726}
727
728//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
729
730void G4EmModelActivator::SetMscParameters(const G4ParticleDefinition* part,
731 G4VMscModel* msc, const G4String& phys)
732{
733 if(part == G4Electron::Electron() || part == G4Positron::Positron()) {
734 if(phys == "G4EmStandard_opt1" || phys == "G4EmStandard_opt2") {
735 msc->SetRangeFactor(0.2);
737 } else if(phys == "G4EmStandard_opt3") {
739 } else if(phys == "G4EmStandard_opt4" || phys == "G4EmLivermore" || phys == "G4EmPenelope") {
740 msc->SetRangeFactor(0.08);
742 msc->SetSkin(3);
743 } else if(phys == "G4EmStandardGS") {
744 msc->SetRangeFactor(0.06);
745 }
746 } else {
747 if(phys != "G4EmStandard" && phys != "G4EmStandard_opt1" && phys != "G4EmStandard_opt2") {
748 msc->SetLateralDisplasmentFlag(true);
749 }
750 }
751 msc->SetLocked(true);
752}
753
754//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
755
756void G4EmModelActivator::FindOrAddProcess(const G4ParticleDefinition* part,
757 const G4String& name)
758{
761 G4int nproc = pm->GetProcessListLength();
762 for(G4int i = 0; i<nproc; ++i) {
763 if(((*pv)[i])->GetProcessName() == name) { return; }
764 }
765 if(name == "CoulombScat") {
767 cs->SetEmModel(new G4DummyModel());
768 pm->AddDiscreteProcess(cs);
769 } else if(name == "Rayl") {
771 rs->SetEmModel(new G4DummyModel());
772 pm->AddDiscreteProcess(rs);
773 }
774}
775
776//....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:781
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
void SetActivationHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:760
void SetLocked(G4bool)
Definition: G4VEmModel.hh:863
void SetEmModel(G4VEmModel *, G4int index=0)
void SetSkin(G4double)
Definition: G4VMscModel.hh:229
void SetRangeFactor(G4double)
Definition: G4VMscModel.hh:236
void SetLateralDisplasmentFlag(G4bool val)
Definition: G4VMscModel.hh:222
void SetStepLimitType(G4MscStepLimitType)
Definition: G4VMscModel.hh:264
const char * name(G4int ptype)
G4bool IsMasterThread()
Definition: G4Threading.cc:124