Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups 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 // $Id: G4EmModelActivator.cc 1651 2015-05-02 16:40:24Z vnivanch $
27 // GEANT4 tag $Name$
28 //
29 //---------------------------------------------------------------------------
30 //
31 // ClassName: G4EmModelActivator
32 //
33 // Author: V.Ivanchenko 01.06.2015
34 //
35 // Organisation: G4AI
36 // Contract: ESA contract 4000107387/12/NL/AT
37 // Customer: ESA/ESTEC
38 //
39 // Modified:
40 //
41 //----------------------------------------------------------------------------
42 //
43 
44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
46 
47 #include "G4EmModelActivator.hh"
48 #include "G4EmParameters.hh"
49 #include "G4ParticleDefinition.hh"
50 #include "G4ParticleTable.hh"
51 #include "G4RegionStore.hh"
52 #include "G4Region.hh"
53 #include "G4VEnergyLossProcess.hh"
54 #include "G4LossTableManager.hh"
55 #include "G4EmConfigurator.hh"
56 #include "G4PAIModel.hh"
57 #include "G4PAIPhotModel.hh"
58 #include "G4Gamma.hh"
59 #include "G4Electron.hh"
60 #include "G4Positron.hh"
61 #include "G4MuonPlus.hh"
62 #include "G4MuonMinus.hh"
63 #include "G4Proton.hh"
64 #include "G4GenericIon.hh"
65 #include "G4Alpha.hh"
66 #include "G4ProcessManager.hh"
67 #include "G4DummyModel.hh"
68 #include "G4EmProcessSubType.hh"
69 #include "G4PhysicsListHelper.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 {
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->RegionsMsc();
137  if(regnamesMSC.size() > 0)
138  {
139  ActivateEmOptions();
140  }
141 }
142 
143 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
144 
145 void G4EmModelActivator::ActivateEmOptions()
146 {
147  const std::vector<G4String> regnamesMSC = theParameters->RegionsMsc();
148  G4int nreg = regnamesMSC.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> typesMSC = theParameters->TypesMsc();
156 
157  // start configuration of models
160  const G4ParticleDefinition* phot = G4Gamma::Gamma();
162  G4EmConfigurator* em_config = man->EmConfigurator();
163  G4VAtomDeexcitation* adeexc = man->AtomDeexcitation();
164  G4VEmModel* mod;
165 
166  for(G4int i=0; i<nreg; ++i) {
167  G4String reg = regnamesMSC[i];
168  if(verbose > 0) {
169  G4cout << i << ". region <" << reg << ">; type <" << typesMSC[i] << "> "
170  << G4endl;
171  }
172 
173  if(baseName == typesMSC[i]) { continue; }
174 
175  if("G4EmStandard" == typesMSC[i]) {
176  G4UrbanMscModel* msc = new G4UrbanMscModel();
177  msc->SetLocked(true);
178  em_config->SetExtraEmModel("e-", "msc", msc, reg, 0.0, 100*MeV);
179  msc = new G4UrbanMscModel();
180  msc->SetLocked(true);
181  em_config->SetExtraEmModel("e+", "msc", msc, reg, 0.0, 100*MeV);
182 
183  } else if("G4EmStandard_opt1" == typesMSC[i] ||
184  "G4EmStandard_opt2" == typesMSC[i]) {
185  G4UrbanMscModel* msc = new G4UrbanMscModel();
187  msc->SetRangeFactor(0.2);
188  msc->SetLocked(true);
189  em_config->SetExtraEmModel("e-", "msc", msc, reg, 0.0, 100*MeV);
190  msc = new G4UrbanMscModel();
192  msc->SetRangeFactor(0.2);
193  msc->SetLocked(true);
194  em_config->SetExtraEmModel("e+", "msc", msc, reg, 0.0, 100*MeV);
195 
196  } else if("G4EmStandard_opt3" == typesMSC[i]) {
197  G4UrbanMscModel* msc = new G4UrbanMscModel();
199  msc->SetLocked(true);
200  em_config->SetExtraEmModel("e-", "msc", msc, reg);
201  msc = new G4UrbanMscModel();
203  msc->SetLocked(true);
204  em_config->SetExtraEmModel("e+", "msc", msc, reg);
205  msc = new G4UrbanMscModel();
207  msc->SetRangeFactor(0.2);
208  msc->SetLocked(true);
209  em_config->SetExtraEmModel("proton", "msc", msc, reg);
211  theParameters->SetDeexActiveRegion(reg, true, false, false);
212  }
213  theParameters->DefineRegParamForDeex(adeexc);
214  FindOrAddProcess(phot, "Rayl");
215  mod = new G4LivermoreRayleighModel();
216  em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
217  FindOrAddProcess(phot, "phot");
218  mod = new G4LivermorePhotoElectricModel();
219  em_config->SetExtraEmModel("gamma", "phot", mod, reg);
220  FindOrAddProcess(phot, "compt");
221  mod = new G4KleinNishinaModel();
222  em_config->SetExtraEmModel("gamma", "compt", mod, reg);
223 
224  } else if("G4EmStandard_opt4" == typesMSC[i]) {
225  G4UrbanMscModel* msc = new G4UrbanMscModel();
227  msc->SetRangeFactor(0.02);
228  msc->SetLocked(true);
229  em_config->SetExtraEmModel("e-", "msc", msc, reg, 0.0, 100*MeV);
230  msc = new G4UrbanMscModel();
232  msc->SetRangeFactor(0.02);
233  msc->SetLocked(true);
234  em_config->SetExtraEmModel("e+", "msc", msc, reg, 0.0, 100*MeV);
236  theParameters->SetDeexActiveRegion(reg, true, false, false);
237  }
238  theParameters->DefineRegParamForDeex(adeexc);
239  FindOrAddProcess(phot, "Rayl");
240  mod = new G4LivermoreRayleighModel();
241  em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
242  FindOrAddProcess(phot, "phot");
243  mod = new G4LivermorePhotoElectricModel();
244  em_config->SetExtraEmModel("gamma", "phot", mod, reg);
245  FindOrAddProcess(phot, "compt");
246  mod = new G4KleinNishinaModel();
247  em_config->SetExtraEmModel("gamma", "compt", mod, reg);
248  mod = new G4LowEPComptonModel();
249  mod->SetHighEnergyLimit(20*MeV);
250  em_config->SetExtraEmModel("gamma", "compt", mod, reg);
251 
252  } else if("G4EmStandardGS" == typesMSC[i]) {
255  msc->SetRangeFactor(0.1);
256  msc->SetLocked(true);
257  em_config->SetExtraEmModel("e-", "msc", msc, reg, 0.0, 100*MeV);
258  msc = new G4GoudsmitSaundersonMscModel();
260  msc->SetRangeFactor(0.1);
261  msc->SetLocked(true);
262  em_config->SetExtraEmModel("e+", "msc", msc, reg, 0.0, 100*MeV);
263 
264  } else if("G4EmStandardWVI" == typesMSC[i]) {
265  G4WentzelVIModel* msc = new G4WentzelVIModel();
266  msc->SetLocked(true);
267  em_config->SetExtraEmModel("e-", "msc", msc, reg, 0.0, 100*MeV);
268  msc = new G4WentzelVIModel();
269  msc->SetLocked(true);
270  em_config->SetExtraEmModel("e+", "msc", msc, reg, 0.0, 100*MeV);
271 
272  FindOrAddProcess(elec, "CoulombScat");
273  FindOrAddProcess(posi, "CoulombScat");
274  mod = new G4eCoulombScatteringModel();
275  em_config->SetExtraEmModel("e-", "CoulombScat", mod, reg, 0.0, 100*MeV);
276  mod = new G4eCoulombScatteringModel();
277  em_config->SetExtraEmModel("e+", "CoulombScat", mod, reg, 0.0, 100*MeV);
278 
280  theParameters->SetDeexActiveRegion(regnamesMSC[i], true, false, false);
281  }
282  theParameters->DefineRegParamForDeex(adeexc);
283 
284  } else if("G4EmStandardSS" == typesMSC[i] &&
285  baseName != "G4EmStandard_opt3") {
287  sc->SetPolarAngleLimit(0.0);
288  sc->SetLocked(true);
289  em_config->SetExtraEmModel("e-", "CoulombScat", sc, reg);
290 
291  sc = new G4eCoulombScatteringModel();
292  sc->SetPolarAngleLimit(0.0);
293  sc->SetLocked(true);
294  em_config->SetExtraEmModel("e+", "CoulombScat", sc, reg);
295 
296  sc = new G4eCoulombScatteringModel();
297  sc->SetPolarAngleLimit(0.0);
298  sc->SetLocked(true);
299  em_config->SetExtraEmModel("mu+", "CoulombScat", sc, reg);
300 
301  sc = new G4eCoulombScatteringModel();
302  sc->SetPolarAngleLimit(0.0);
303  sc->SetLocked(true);
304  em_config->SetExtraEmModel("pi+", "CoulombScat", sc, reg);
305 
306  sc = new G4eCoulombScatteringModel();
307  sc->SetPolarAngleLimit(0.0);
308  sc->SetLocked(true);
309  em_config->SetExtraEmModel("kaon+", "CoulombScat", sc, reg);
310 
311  sc = new G4eCoulombScatteringModel();
312  sc->SetPolarAngleLimit(0.0);
313  sc->SetLocked(true);
314  em_config->SetExtraEmModel("proton", "CoulombScat", sc, reg);
315 
317  theParameters->SetDeexActiveRegion(regnamesMSC[i], true, true, true);
318  }
319  theParameters->DefineRegParamForDeex(adeexc);
320 
321  } else if("G4EmLivermore" == typesMSC[i]) {
322 
323  mod = new G4LivermorePhotoElectricModel();
324  em_config->SetExtraEmModel("gamma", "phot", mod, reg);
325  mod = new G4LivermoreComptonModel();
326  em_config->SetExtraEmModel("gamma", "compt", mod, reg);
328  em_config->SetExtraEmModel("gamma", "conv", mod, reg);
329 
330  FindOrAddProcess(phot, "Rayl");
331  mod = new G4LivermoreRayleighModel();
332  em_config->SetExtraEmModel("gamma", "Rayl", mod, reg);
333 
334  mod = new G4LivermoreIonisationModel();
336  em_config->SetExtraEmModel("e-", "eIoni", mod, reg, 0.0, 0.1*MeV, uf);
337  mod = new G4LivermoreBremsstrahlungModel();
338  em_config->SetExtraEmModel("e-", "eBrem", mod, reg, 0.0, 1*GeV);
339 
341  theParameters->SetDeexActiveRegion(regnamesMSC[i], true, false, false);
342  }
343  theParameters->DefineRegParamForDeex(adeexc);
344 
345  } else if("G4EmPenelope" == typesMSC[i]) {
346  mod = new G4PenelopePhotoElectricModel();
347  em_config->SetExtraEmModel("gamma", "phot", mod, reg);
348  mod = new G4PenelopeComptonModel();
349  em_config->SetExtraEmModel("gamma", "compt", mod, reg);
350  mod = new G4PenelopeGammaConversionModel();
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, 1*GeV, uf);
360  mod = new G4PenelopeBremsstrahlungModel();
361  em_config->SetExtraEmModel("e-", "eBrem", mod, reg, 0.0, 1*GeV);
362 
363  mod = new G4PenelopeIonisationModel();
364  uf = new G4UniversalFluctuation();
365  em_config->SetExtraEmModel("e+", "eIoni", mod, reg, 0.0, 1*GeV, uf);
366  mod = new G4PenelopeBremsstrahlungModel();
367  em_config->SetExtraEmModel("e+", "eBrem", mod, reg, 0.0, 1*GeV);
368  mod = new G4PenelopeAnnihilationModel();
369  em_config->SetExtraEmModel("e+", "annihil", mod, reg, 0.0, 1*GeV);
370 
372  theParameters->SetDeexActiveRegion(regnamesMSC[i], true, false, false);
373  }
374  theParameters->DefineRegParamForDeex(adeexc);
375 
376  } else if("G4RadioactiveDecay" == typesMSC[i]) {
377 
379  theParameters->SetFluo(true);
380  theParameters->SetAuger(true);
381  theParameters->SetAugerCascade(true);
382  theParameters->SetDeexcitationIgnoreCut(true);
383  }
384  theParameters->DefineRegParamForDeex(adeexc);
385 
386  } else {
387  if(verbose > 0 && G4Threading::IsMasterThread()) {
388  G4cout << "### G4EmModelActivator::ActivateEmOptions WARNING: \n"
389  << " EM Physics configuration name <" << typesMSC[i]
390  << "> is not known - ignored" << G4endl;
391  }
392  }
393  }
394 }
395 
396 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
397 
398 void G4EmModelActivator::ActivatePAI()
399 {
400  const std::vector<G4String> regnamesPAI = theParameters->RegionsPAI();
401  G4int nreg = regnamesPAI.size();
402  if(0 == nreg) { return; }
403  G4int verbose = theParameters->Verbose() - 1;
404  if(verbose > 0) {
405  G4cout << "### G4EmModelActivator::ActivatePAI for " << nreg << " regions"
406  << G4endl;
407  }
408  const std::vector<G4String> particlesPAI = theParameters->ParticlesPAI();
409  const std::vector<G4String> typesPAI = theParameters->TypesPAI();
410 
411  const std::vector<G4VEnergyLossProcess*>& v = G4LossTableManager::Instance()
413  std::vector<G4VEnergyLossProcess*>::const_iterator itr;
414  G4RegionStore* regionStore = G4RegionStore::GetInstance();
415 
421 
422  for(G4int i = 0; i < nreg; ++i) {
423  const G4ParticleDefinition* p = nullptr;
424  if(particlesPAI[i] != "all") {
425  p = G4ParticleTable::GetParticleTable()->FindParticle(particlesPAI[i]);
426  if(!p) {
427  G4cout << "### WARNING: ActivatePAI::FindParticle fails to find "
428  << particlesPAI[i] << G4endl;
429  continue;
430  }
431  }
432  const G4Region* r = regionStore->GetRegion(regnamesPAI[i], false);
433  if(!r) {
434  G4cout << "### WARNING: ActivatePAI::GetRegion fails to find "
435  << regnamesPAI[i] << G4endl;
436  continue;
437  }
438 
439  G4String name = "hIoni";
440  if(p == elec || p == posi)
441  { name = "eIoni"; }
442  else if (p == mupl || p == mumi)
443  { name = "muIoni"; }
444  else if (p == gion)
445  { name = "ionIoni"; }
446 
447  for(itr = v.begin(); itr != v.end(); ++itr) {
448 
449  G4VEnergyLossProcess* proc = *itr;
450  if(!proc->IsIonisationProcess()) { continue; }
451 
452  G4String namep = proc->GetProcessName();
453  if(p) {
454  if(name != namep) { continue; }
455  } else {
456  if(namep != "hIoni" && namep != "muIoni" && namep != "eIoni")
457  { continue; }
458  }
459  G4VEmModel* em = nullptr;
460  G4VEmFluctuationModel* fm = nullptr;
461  if(typesPAI[i] == "PAIphoton") {
462  G4PAIPhotModel* mod = new G4PAIPhotModel(p,"PAIPhotModel");
463  em = mod;
464  fm = mod;
465  } else {
466  G4PAIModel* mod = new G4PAIModel(p,"PAIModel");
467  em = mod;
468  fm = mod;
469  }
470  proc->AddEmModel(0, em, fm, r);
471  if(verbose > 0) {
472  G4cout << "### G4EmModelActivator: add <" << typesPAI[i]
473  << "> model for " << particlesPAI[i]
474  << " in the " << regnamesPAI[i] << G4endl;
475  }
476  }
477  }
478 }
479 
480 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
481 
482 void G4EmModelActivator::ActivateMicroElec()
483 {
484  const std::vector<G4String> regnamesME = theParameters->RegionsMicroElec();
485  G4int nreg = regnamesME.size();
486  if(0 == nreg)
487  {
488  return;
489  }
490  G4int verbose = theParameters->Verbose() - 1;
491  if(verbose > 0)
492  {
493  G4cout << "### G4EmModelActivator::ActivateMicroElec for " << nreg
494  << " regions" << G4endl;
495  }
497 
499  const G4ParticleDefinition* prot = G4Proton::Proton();
501  G4ProcessManager* eman = elec->GetProcessManager();
502  G4ProcessManager* pman = prot->GetProcessManager();
503  G4ProcessManager* iman = gion->GetProcessManager();
504 
505  G4bool emsc = HasMsc(eman);
506 
507  // MicroElec elastic is not active in the world
508  G4MicroElecElastic* eElasProc =
509  new G4MicroElecElastic("e-G4MicroElecElastic");
510  eman->AddDiscreteProcess(eElasProc);
511 
512  G4MicroElecInelastic* eInelProc =
513  new G4MicroElecInelastic("e-G4MicroElecInelastic");
514  eman->AddDiscreteProcess(eInelProc);
515 
516  G4MicroElecInelastic* pInelProc =
517  new G4MicroElecInelastic("p_G4MicroElecInelastic");
518  pman->AddDiscreteProcess(pInelProc);
519 
520  G4MicroElecInelastic* iInelProc =
521  new G4MicroElecInelastic("ion_G4MicroElecInelastic");
522  iman->AddDiscreteProcess(iInelProc);
523 
524  // start configuration of models
525  G4EmConfigurator* em_config = man->EmConfigurator();
526  G4VEmModel* mod;
527 
528  // limits for MicroElec applicability
529  G4double elowest = 16.7 * eV;
530  G4double elimel = 100 * MeV;
531  G4double elimin = 9 * MeV;
532  G4double pmin = 50 * keV;
533  G4double pmax = 99.9 * MeV;
534 
535  G4LowECapture* ecap = new G4LowECapture(elowest);
536  eman->AddDiscreteProcess(ecap);
537 
538  for(G4int i = 0; i < nreg; ++i)
539  {
540 
541  G4String reg = regnamesME[i];
542  G4cout << "### MicroElec models are activated for G4Region " << reg
543  << G4endl
544  << " Energy limits for e- elastic: " << elowest/eV << " eV - "
545  << elimel/MeV << " MeV"
546  << G4endl
547  << " Energy limits for e- inelastic: " << elowest/eV << " eV - "
548  << elimin/MeV << " MeV"
549  << G4endl
550  << " Energy limits for hadrons/ions: " << pmin/MeV << " MeV - "
551  << pmax/MeV << " MeV"
552  << G4endl;
553 
554  // e-
555  if(emsc)
556  {
557  G4UrbanMscModel* msc = new G4UrbanMscModel();
558  msc->SetActivationLowEnergyLimit(elimel);
559  em_config->SetExtraEmModel("e-", "msc", msc, reg);
560  }
561  else
562  {
563  mod = new G4DummyModel();
564  em_config->SetExtraEmModel("e-", "CoulombScat", mod, reg, 0.0, elimel);
565  }
566 
567  mod = new G4MicroElecElasticModel();
568  em_config->SetExtraEmModel("e-",
569  "e-G4MicroElecElastic",
570  mod,
571  reg,
572  elowest,
573  elimin);
574 
575  mod = new G4MollerBhabhaModel();
576  mod->SetActivationLowEnergyLimit(elimin);
577  em_config->SetExtraEmModel("e-",
578  "eIoni",
579  mod,
580  reg,
581  0.0,
582  10 * TeV,
583  new G4UniversalFluctuation());
584 
585  mod = new G4MicroElecInelasticModel();
586  em_config->SetExtraEmModel("e-",
587  "e-G4MicroElecInelastic",
588  mod,
589  reg,
590  elowest,
591  elimin);
592 
593  // proton
594  mod = new G4BraggModel();
595  mod->SetActivationHighEnergyLimit(pmin);
596  em_config->SetExtraEmModel("proton",
597  "hIoni",
598  mod,
599  reg,
600  0.0,
601  2 * MeV,
602  new G4UniversalFluctuation());
603 
604  mod = new G4BetheBlochModel();
605  mod->SetActivationLowEnergyLimit(pmax);
606  em_config->SetExtraEmModel("proton",
607  "hIoni",
608  mod,
609  reg,
610  2 * MeV,
611  10 * TeV,
612  new G4UniversalFluctuation());
613 
614  mod = new G4MicroElecInelasticModel();
615  em_config->SetExtraEmModel("proton",
616  "p_G4MicroElecInelastic",
617  mod,
618  reg,
619  pmin,
620  pmax);
621 
622  // ions
623  mod = new G4BraggIonModel();
624  mod->SetActivationHighEnergyLimit(pmin);
625  em_config->SetExtraEmModel("GenericIon",
626  "ionIoni",
627  mod,
628  reg,
629  0.0,
630  2 * MeV,
631  new G4IonFluctuations());
632 
633  mod = new G4BetheBlochModel();
634  mod->SetActivationLowEnergyLimit(pmax);
635  em_config->SetExtraEmModel("GenericIon",
636  "ionIoni",
637  mod,
638  reg,
639  2 * MeV,
640  10 * TeV,
641  new G4IonFluctuations());
642 
643  mod = new G4MicroElecInelasticModel();
644  em_config->SetExtraEmModel("GenericIon",
645  "ion_G4MicroElecInelastic",
646  mod,
647  reg,
648  pmin,
649  pmax);
650  }
651 }
652 
653 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
654 
655 G4bool G4EmModelActivator::HasMsc(G4ProcessManager* pm) const
656 {
657  G4bool res = false;
658  G4ProcessVector* pv = pm->GetProcessList();
659  G4int nproc = pm->GetProcessListLength();
660  for(G4int i = 0; i < nproc; ++i)
661  {
662  if(((*pv)[i])->GetProcessSubType() == fMultipleScattering)
663  {
664  res = true;
665  break;
666  }
667  }
668  return res;
669 }
670 
671 void G4EmModelActivator::FindOrAddProcess(const G4ParticleDefinition* part,
672  const G4String& name)
673 {
674  G4ProcessManager* pm = part->GetProcessManager();
675  G4ProcessVector* pv = pm->GetProcessList();
676  G4int nproc = pm->GetProcessListLength();
677  for(G4int i = 0; i<nproc; ++i) {
678  if(((*pv)[i])->GetProcessName() == name) { return; }
679  }
680  if(name == "CoulombScat") {
682  cs->SetEmModel(new G4DummyModel());
683  pm->AddDiscreteProcess(cs);
684  } else if(name == "Rayl") {
686  rs->SetEmModel(new G4DummyModel());
687  pm->AddDiscreteProcess(rs);
688  }
689 }
690 
691 
692 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4EmConfigurator * EmConfigurator()
const XML_Char * name
Definition: expat.h:151
void SetActivationHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4LossTableManager * Instance()
void SetDeexcitationIgnoreCut(G4bool val)
const char * p
Definition: xmltok.h:285
void SetAuger(G4bool val)
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
void SetExtraEmModel(const G4String &particleName, const G4String &processName, G4VEmModel *, const G4String &regionName="", G4double emin=0.0, G4double emax=DBL_MAX, G4VEmFluctuationModel *fm=nullptr)
const std::vector< G4String > & ParticlesPAI() const
const std::vector< G4String > & RegionsMicroElec() const
const std::vector< G4String > & RegionsPAI() const
void SetStepLimitType(G4MscStepLimitType)
Definition: G4VMscModel.hh:223
int G4int
Definition: G4Types.hh:78
void SetDeexActiveRegion(const G4String &region, G4bool fdeex, G4bool fauger, G4bool fpixe)
const std::vector< G4String > & RegionsMsc() const
static const G4double reg
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
static constexpr double TeV
Definition: G4SIunits.hh:218
static G4RegionStore * GetInstance()
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=1)
G4GLOB_DLL std::ostream G4cout
G4int Verbose() const
void SetLocked(G4bool)
Definition: G4VEmModel.hh:847
bool G4bool
Definition: G4Types.hh:79
void SetRangeFactor(G4double)
Definition: G4VMscModel.hh:209
void SetAugerCascade(G4bool val)
static G4Proton * Proton()
Definition: G4Proton.cc:93
static constexpr double eV
Definition: G4SIunits.hh:215
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
tuple v
Definition: test.py:18
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:93
const std::vector< G4String > & TypesMsc() const
static G4ParticleTable * GetParticleTable()
G4ProcessManager * GetProcessManager() const
static G4EmParameters * Instance()
void DefineRegParamForDeex(G4VAtomDeexcitation *) const
static constexpr double GeV
Definition: G4SIunits.hh:217
G4EmModelActivator(const G4String &emphys="")
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4bool IsMasterThread()
Definition: G4Threading.cc:146
G4VAtomDeexcitation * AtomDeexcitation()
double G4double
Definition: G4Types.hh:76
G4int GetProcessListLength() const
static constexpr double keV
Definition: G4SIunits.hh:216
const std::vector< G4String > & TypesPAI() const
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:767
G4ProcessVector * GetProcessList() const
void SetFluo(G4bool val)
G4bool IsIonisationProcess() const