ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4EmModelActivator.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4EmModelActivator.cc
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"
47 #include "G4ParticleDefinition.hh"
48 #include "G4ParticleTable.hh"
49 #include "G4RegionStore.hh"
50 #include "G4Region.hh"
51 #include "G4VEnergyLossProcess.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"
68 #include "G4PhysicsListHelper.hh"
69 #include "G4EmParticleList.hh"
70 
71 #include "G4MicroElecElastic.hh"
73 
74 #include "G4MicroElecInelastic.hh"
76 
77 #include "G4BraggModel.hh"
78 #include "G4BraggIonModel.hh"
79 #include "G4BetheBlochModel.hh"
80 #include "G4UrbanMscModel.hh"
81 #include "G4MollerBhabhaModel.hh"
82 #include "G4IonFluctuations.hh"
84 #include "G4LowECapture.hh"
85 #include "G4hMultipleScattering.hh"
88 #include "G4ionIonisation.hh"
89 #include "G4KleinNishinaModel.hh"
90 
91 #include "G4CoulombScattering.hh"
93 #include "G4WentzelVIModel.hh"
95 #include "G4RayleighScattering.hh"
96 #include "G4UrbanMscModel.hh"
98 #include "G4LowEPComptonModel.hh"
99 
106 
108 #include "G4PenelopeComptonModel.hh"
114 
115 #include "G4SystemOfUnits.hh"
116 #include "G4Threading.hh"
117 #include <vector>
118 
119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
120 
122  : baseName(emphys)
123 {
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  {
135  }
136  const std::vector<G4String>& regnamesMSC = theParameters->RegionsPhysics();
137  if(regnamesMSC.size() > 0)
138  {
140  }
141 }
142 
143 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
144 
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();
161  const G4ParticleDefinition* prot = G4Proton::Proton();
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 = G4EmParameters::Instance()->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 = 100*TeV;
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  msc->SetRangeFactor(0.04);
189  msc->SetSkin(1);
191  AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy);
192 
193  msc = new G4UrbanMscModel();
194  msc->SetRangeFactor(0.04);
195  msc->SetSkin(1);
197  AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy);
198 
199  } else if("G4EmStandard_opt1" == typesPhys[i] ||
200  "G4EmStandard_opt2" == typesPhys[i]) {
201  G4UrbanMscModel* msc = new G4UrbanMscModel();
203  msc->SetRangeFactor(0.2);
204  AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy);
205 
206  msc = new G4UrbanMscModel();
208  msc->SetRangeFactor(0.2);
209  AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy);
210 
211  } else if("G4EmStandard_opt3" == typesPhys[i]) {
212 
213  G4DummyModel* dummy = new G4DummyModel();
214  G4UrbanMscModel* msc = new G4UrbanMscModel();
216  msc->SetLocked(true);
217  em_config->SetExtraEmModel("e-", "msc", msc, reg);
218  FindOrAddProcess(elec, "CoulombScat");
219  em_config->SetExtraEmModel("e-", "CoulombScat", dummy, reg);
220 
221  msc = new G4UrbanMscModel();
223  msc->SetLocked(true);
224  em_config->SetExtraEmModel("e+", "msc", msc, reg);
225  FindOrAddProcess(posi, "CoulombScat");
226  em_config->SetExtraEmModel("e+", "CoulombScat", dummy, reg);
227 
228  msc = new G4UrbanMscModel();
230  msc->SetRangeFactor(0.2);
231  msc->SetLocked(true);
232  em_config->SetExtraEmModel("proton", "msc", msc, reg);
233  FindOrAddProcess(prot, "CoulombScat");
234  em_config->SetExtraEmModel("proton", "CoulombScat", dummy, reg);
235 
238  theParameters->SetDeexActiveRegion(reg, true, false, false);
239  }
241  FindOrAddProcess(phot, "Rayl");
242  mod = new G4LivermoreRayleighModel();
243  em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
244  FindOrAddProcess(phot, "phot");
245  mod = new G4LivermorePhotoElectricModel();
246  em_config->SetExtraEmModel("gamma", "phot", mod, reg);
247  FindOrAddProcess(phot, "compt");
248  mod = new G4KleinNishinaModel();
249  em_config->SetExtraEmModel("gamma", "compt", mod, reg);
250 
251  } else if("G4EmStandard_opt4" == typesPhys[i]) {
252  G4VMscModel* msc = GetGSModel();
253  msc->SetRangeFactor(0.08);
254  msc->SetSkin(3);
256  AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy);
257 
258  msc = GetGSModel();
259  msc->SetRangeFactor(0.08);
260  msc->SetSkin(3);
262  AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy);
263 
267  theParameters->SetDeexActiveRegion(reg, true, false, false);
268  }
270 
271  FindOrAddProcess(phot, "Rayl");
272  mod = new G4LivermoreRayleighModel();
273  em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
274  FindOrAddProcess(phot, "compt");
275  mod = new G4KleinNishinaModel();
276  em_config->SetExtraEmModel("gamma", "compt", mod, reg);
277  mod = new G4LowEPComptonModel();
278  mod->SetHighEnergyLimit(20*MeV);
279  em_config->SetExtraEmModel("gamma", "compt", mod, reg);
280 
281  } else if("G4EmStandardGS" == typesPhys[i]) {
283  msc->SetRangeFactor(0.06);
284  AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy);
285 
286  msc = new G4GoudsmitSaundersonMscModel();
287  msc->SetRangeFactor(0.06);
288  AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy);
289 
290  } else if("G4EmStandardWVI" == typesPhys[i]) {
291  G4WentzelVIModel* msc = new G4WentzelVIModel();
292  AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy);
293 
294  msc = new G4WentzelVIModel();
295  AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy);
296 
298  theParameters->SetDeexActiveRegion(regnamesPhys[i], true, false, false);
299  }
301 
302  } else if("G4EmStandardSS" == typesPhys[i] &&
303  baseName != "G4EmStandard_opt3") {
304  G4EmParticleList emList;
305  for(const auto& particleName : emList.PartNames()) {
306  G4ParticleDefinition* particle = table->FindParticle(particleName);
307  if(particle && 0.0 != particle->GetPDGCharge()) {
308  FindOrAddProcess(particle, "CoulombScat");
310  sc->SetPolarAngleLimit(0.0);
311  sc->SetLocked(true);
312  em_config->SetExtraEmModel(particleName, "CoulombScat", sc, reg);
313  if(particleName == "mu+" || particleName == "mu-") {
314  em_config->SetExtraEmModel(particleName, "muMsc",
315  new G4DummyModel(), reg);
316  } else {
317  em_config->SetExtraEmModel(particleName, "msc",
318  new G4DummyModel(), reg);
319  }
320  }
321  }
323  theParameters->SetDeexActiveRegion(regnamesPhys[i], true, true, true);
324  }
326 
327  } else if("G4EmLivermore" == typesPhys[i]) {
328 
329  G4VMscModel* msc = GetGSModel();
330  msc->SetRangeFactor(0.08);
331  msc->SetSkin(3);
333  AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy);
334 
335  msc = GetGSModel();
336  msc->SetRangeFactor(0.08);
337  msc->SetSkin(3);
339  AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy);
340 
341  mod = new G4LivermorePhotoElectricModel();
342  em_config->SetExtraEmModel("gamma", "phot", mod, reg);
343  mod = new G4LivermoreComptonModel();
344  em_config->SetExtraEmModel("gamma", "compt", mod, reg);
346  em_config->SetExtraEmModel("gamma", "conv", mod, reg);
347 
348  FindOrAddProcess(phot, "Rayl");
349  mod = new G4LivermoreRayleighModel();
350  em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
351 
352  mod = new G4LivermoreIonisationModel();
354  em_config->SetExtraEmModel("e-", "eIoni", mod, reg, 0.0, 0.1*MeV, uf);
355  mod = new G4LivermoreBremsstrahlungModel();
356  em_config->SetExtraEmModel("e-", "eBrem", mod, reg, 0.0, highEnergyLimit);
357 
361  theParameters->SetDeexActiveRegion(regnamesPhys[i], true, false, false);
362  }
364 
365  } else if("G4EmPenelope" == typesPhys[i]) {
366 
367  G4VMscModel* msc = GetGSModel();
368  msc->SetRangeFactor(0.08);
369  msc->SetSkin(3);
371  AddStandardScattering(elec, em_config, msc, reg, mscEnergyLimit, highEnergy);
372 
373  msc = GetGSModel();
374  msc->SetRangeFactor(0.08);
375  msc->SetSkin(3);
377  AddStandardScattering(posi, em_config, msc, reg, mscEnergyLimit, highEnergy);
378 
379  mod = new G4PenelopePhotoElectricModel();
380  em_config->SetExtraEmModel("gamma", "phot", mod, reg);
381  mod = new G4PenelopeComptonModel();
382  em_config->SetExtraEmModel("gamma", "compt", mod, reg);
383  mod = new G4PenelopeGammaConversionModel();
384  em_config->SetExtraEmModel("gamma", "conv", mod, reg);
385 
386  FindOrAddProcess(phot, "Rayl");
387  mod = new G4PenelopeRayleighModel();
388  em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
389 
390  mod = new G4PenelopeIonisationModel();
392  em_config->SetExtraEmModel("e-", "eIoni", mod, reg, 0.0, highEnergyLimit, uf);
393  mod = new G4PenelopeBremsstrahlungModel();
394  em_config->SetExtraEmModel("e-", "eBrem", mod, reg, 0.0, highEnergyLimit);
395 
396  mod = new G4PenelopeIonisationModel();
397  uf = new G4UniversalFluctuation();
398  em_config->SetExtraEmModel("e+", "eIoni", mod, reg, 0.0, highEnergyLimit, uf);
399  mod = new G4PenelopeBremsstrahlungModel();
400  em_config->SetExtraEmModel("e+", "eBrem", mod, reg, 0.0, highEnergyLimit);
401  mod = new G4PenelopeAnnihilationModel();
402  em_config->SetExtraEmModel("e+", "annihil", mod, reg, 0.0, highEnergyLimit);
403 
407  theParameters->SetDeexActiveRegion(regnamesPhys[i], true, false, false);
408  }
410 
411  } else if("G4RadioactiveDecay" == typesPhys[i]) {
412 
416  }
418 
419  } else {
420  if(verbose > 0 && G4Threading::IsMasterThread()) {
421  G4cout << "### G4EmModelActivator::ActivateEmOptions WARNING: \n"
422  << " EM Physics configuration name <" << typesPhys[i]
423  << "> is not known - ignored" << G4endl;
424  }
425  }
426  }
427 }
428 
429 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
430 
432 {
433  const std::vector<G4String> regnamesPAI = theParameters->RegionsPAI();
434  G4int nreg = regnamesPAI.size();
435  if(0 == nreg) { return; }
436  G4int verbose = theParameters->Verbose() - 1;
437  if(verbose > 0) {
438  G4cout << "### G4EmModelActivator::ActivatePAI for " << nreg << " regions"
439  << G4endl;
440  }
441  const std::vector<G4String> particlesPAI = theParameters->ParticlesPAI();
442  const std::vector<G4String> typesPAI = theParameters->TypesPAI();
443 
444  const std::vector<G4VEnergyLossProcess*>& v = G4LossTableManager::Instance()
446 
447  G4RegionStore* regionStore = G4RegionStore::GetInstance();
448 
454 
455  for(G4int i = 0; i < nreg; ++i) {
456  const G4ParticleDefinition* p = nullptr;
457  if(particlesPAI[i] != "all") {
458  p = G4ParticleTable::GetParticleTable()->FindParticle(particlesPAI[i]);
459  if(!p) {
460  G4cout << "### WARNING: ActivatePAI::FindParticle fails to find "
461  << particlesPAI[i] << G4endl;
462  continue;
463  }
464  }
465  const G4Region* r = regionStore->GetRegion(regnamesPAI[i], false);
466  if(!r) {
467  G4cout << "### WARNING: ActivatePAI::GetRegion fails to find "
468  << regnamesPAI[i] << G4endl;
469  continue;
470  }
471 
472  G4String name = "hIoni";
473  if(p == elec || p == posi)
474  { name = "eIoni"; }
475  else if (p == mupl || p == mumi)
476  { name = "muIoni"; }
477  else if (p == gion)
478  { name = "ionIoni"; }
479 
480  for(auto proc : v) {
481 
482  if(!proc->IsIonisationProcess()) { continue; }
483 
484  G4String namep = proc->GetProcessName();
485  if(p) {
486  if(name != namep) { continue; }
487  } else {
488  if(namep != "hIoni" && namep != "muIoni" &&
489  namep != "eIoni" && namep != "ionIoni")
490  { continue; }
491  }
492  G4VEmModel* em = nullptr;
493  G4VEmFluctuationModel* fm = nullptr;
494  if(typesPAI[i] == "PAIphoton" || typesPAI[i] == "pai_photon") {
495  G4PAIPhotModel* mod = new G4PAIPhotModel(p,"PAIPhotModel");
496  em = mod;
497  fm = mod;
498  } else {
499  G4PAIModel* mod = new G4PAIModel(p,"PAIModel");
500  em = mod;
501  fm = mod;
502  }
503  proc->AddEmModel(0, em, fm, r);
504  if(verbose > 0) {
505  G4cout << "### G4EmModelActivator: add <" << typesPAI[i]
506  << "> model for " << particlesPAI[i]
507  << " in the " << regnamesPAI[i] << G4endl;
508  }
509  }
510  }
511 }
512 
513 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
514 
516 {
517  const std::vector<G4String> regnamesME = theParameters->RegionsMicroElec();
518  G4int nreg = regnamesME.size();
519  if(0 == nreg)
520  {
521  return;
522  }
523  G4int verbose = theParameters->Verbose() - 1;
524  if(verbose > 0)
525  {
526  G4cout << "### G4EmModelActivator::ActivateMicroElec for " << nreg
527  << " regions" << G4endl;
528  }
530 
532  const G4ParticleDefinition* prot = G4Proton::Proton();
534  G4ProcessManager* eman = elec->GetProcessManager();
535  G4ProcessManager* pman = prot->GetProcessManager();
536  G4ProcessManager* iman = gion->GetProcessManager();
537 
538  G4bool emsc = HasMsc(eman);
539 
540  // MicroElec elastic is not active in the world
541  G4MicroElecElastic* eElasProc =
542  new G4MicroElecElastic("e-G4MicroElecElastic");
543  eman->AddDiscreteProcess(eElasProc);
544 
545  G4MicroElecInelastic* eInelProc =
546  new G4MicroElecInelastic("e-G4MicroElecInelastic");
547  eman->AddDiscreteProcess(eInelProc);
548 
549  G4MicroElecInelastic* pInelProc =
550  new G4MicroElecInelastic("p_G4MicroElecInelastic");
551  pman->AddDiscreteProcess(pInelProc);
552 
553  G4MicroElecInelastic* iInelProc =
554  new G4MicroElecInelastic("ion_G4MicroElecInelastic");
555  iman->AddDiscreteProcess(iInelProc);
556 
557  // start configuration of models
558  G4EmConfigurator* em_config = man->EmConfigurator();
559  G4VEmModel* mod;
560 
561  // limits for MicroElec applicability
562  G4double elowest = 16.7 * eV;
563  G4double elimel = 100 * MeV;
564  G4double elimin = 9 * MeV;
565  G4double pmin = 50 * keV;
566  G4double pmax = 99.9 * MeV;
567 
568  G4LowECapture* ecap = new G4LowECapture(elowest);
569  eman->AddDiscreteProcess(ecap);
570 
571  for(G4int i = 0; i < nreg; ++i)
572  {
573 
574  G4String reg = regnamesME[i];
575  G4cout << "### MicroElec models are activated for G4Region " << reg
576  << G4endl
577  << " Energy limits for e- elastic: " << elowest/eV << " eV - "
578  << elimel/MeV << " MeV"
579  << G4endl
580  << " Energy limits for e- inelastic: " << elowest/eV << " eV - "
581  << elimin/MeV << " MeV"
582  << G4endl
583  << " Energy limits for hadrons/ions: " << pmin/MeV << " MeV - "
584  << pmax/MeV << " MeV"
585  << G4endl;
586 
587  // e-
588  if(emsc)
589  {
590  G4UrbanMscModel* msc = new G4UrbanMscModel();
591  msc->SetActivationLowEnergyLimit(elimel);
592  em_config->SetExtraEmModel("e-", "msc", msc, reg);
593  }
594  else
595  {
596  mod = new G4DummyModel();
597  em_config->SetExtraEmModel("e-", "CoulombScat", mod, reg, 0.0, elimel);
598  }
599 
600  mod = new G4MicroElecElasticModel();
601  em_config->SetExtraEmModel("e-",
602  "e-G4MicroElecElastic",
603  mod,
604  reg,
605  elowest,
606  elimin);
607 
608  mod = new G4MollerBhabhaModel();
609  mod->SetActivationLowEnergyLimit(elimin);
610  em_config->SetExtraEmModel("e-",
611  "eIoni",
612  mod,
613  reg,
614  0.0,
615  10 * TeV,
616  new G4UniversalFluctuation());
617 
618  mod = new G4MicroElecInelasticModel();
619  em_config->SetExtraEmModel("e-",
620  "e-G4MicroElecInelastic",
621  mod,
622  reg,
623  elowest,
624  elimin);
625 
626  // proton
627  mod = new G4BraggModel();
628  mod->SetActivationHighEnergyLimit(pmin);
629  em_config->SetExtraEmModel("proton",
630  "hIoni",
631  mod,
632  reg,
633  0.0,
634  2 * MeV,
635  new G4UniversalFluctuation());
636 
637  mod = new G4BetheBlochModel();
638  mod->SetActivationLowEnergyLimit(pmax);
639  em_config->SetExtraEmModel("proton",
640  "hIoni",
641  mod,
642  reg,
643  2 * MeV,
644  10 * TeV,
645  new G4UniversalFluctuation());
646 
647  mod = new G4MicroElecInelasticModel();
648  em_config->SetExtraEmModel("proton",
649  "p_G4MicroElecInelastic",
650  mod,
651  reg,
652  pmin,
653  pmax);
654 
655  // ions
656  mod = new G4BraggIonModel();
657  mod->SetActivationHighEnergyLimit(pmin);
658  em_config->SetExtraEmModel("GenericIon",
659  "ionIoni",
660  mod,
661  reg,
662  0.0,
663  2 * MeV,
664  new G4IonFluctuations());
665 
666  mod = new G4BetheBlochModel();
667  mod->SetActivationLowEnergyLimit(pmax);
668  em_config->SetExtraEmModel("GenericIon",
669  "ionIoni",
670  mod,
671  reg,
672  2 * MeV,
673  10 * TeV,
674  new G4IonFluctuations());
675 
676  mod = new G4MicroElecInelasticModel();
677  em_config->SetExtraEmModel("GenericIon",
678  "ion_G4MicroElecInelastic",
679  mod,
680  reg,
681  pmin,
682  pmax);
683  }
684 }
685 
686 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
687 
689 {
690  G4bool res = false;
691  G4ProcessVector* pv = pm->GetProcessList();
692  G4int nproc = pm->GetProcessListLength();
693  for(G4int i = 0; i < nproc; ++i)
694  {
695  if(((*pv)[i])->GetProcessSubType() == fMultipleScattering)
696  {
697  res = true;
698  break;
699  }
700  }
701  return res;
702 }
703 
704 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
705 
707 {
710  msc->SetRangeFactor(0.2);
711  msc->SetSkin(3);
712  msc->SetOptionMottCorrection(true);
713  msc->SetLocked(true);
714  return msc;
715 }
716 
717 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
718 
720  G4EmConfigurator* em_config,
721  G4VMscModel* mscmod,
722  const G4String& reg,
724 {
725  G4String pname = part->GetParticleName();
726 
727  // low-energy msc model
728  mscmod->SetLocked(true);
729  em_config->SetExtraEmModel(pname, "msc", mscmod, reg, 0.0, e1);
730 
731  // high energy msc model
732  G4WentzelVIModel* msc = new G4WentzelVIModel();
733  msc->SetLocked(true);
734  em_config->SetExtraEmModel(pname, "msc", msc, reg, e1, e2);
735 
736  // high energy single scattering model
737  FindOrAddProcess(part, "CoulombScat");
740  mod->SetLocked(true);
741  em_config->SetExtraEmModel(pname, "CoulombScat", mod, reg, 0.0, e2);
742 }
743 
744 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
745 
747  const G4String& name)
748 {
750  G4ProcessVector* pv = pm->GetProcessList();
751  G4int nproc = pm->GetProcessListLength();
752  for(G4int i = 0; i<nproc; ++i) {
753  if(((*pv)[i])->GetProcessName() == name) { return; }
754  }
755  if(name == "CoulombScat") {
757  cs->SetEmModel(new G4DummyModel());
758  pm->AddDiscreteProcess(cs);
759  } else if(name == "Rayl") {
761  rs->SetEmModel(new G4DummyModel());
762  pm->AddDiscreteProcess(rs);
763  }
764 }
765 
766 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......