Geant4  10.02
G4EmParameters.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: G4EmParameters.cc 69320 2013-04-30 15:59:36Z vnivanch $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4EmParameters
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 18.05.2013
38 //
39 // Modifications:
40 //
41 //
42 //
43 // -------------------------------------------------------------------
44 //
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47 
48 #include "G4EmParameters.hh"
49 #include "G4PhysicalConstants.hh"
50 #include "G4UnitsTable.hh"
51 #include "G4SystemOfUnits.hh"
53 #include "G4NistManager.hh"
54 
56 
57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
58 
60 {
61  if(0 == theInstance) {
62  static G4EmParameters manager;
63  theInstance = &manager;
64  }
65  return theInstance;
66 }
67 
68 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
69 
71 {
72  delete theMessenger;
73 }
74 
75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
76 
78 {
81 
82  SetDefaults();
83 }
84 
85 #include "G4AutoLock.hh"
86 namespace { G4Mutex EmParametersMutex = G4MUTEX_INITIALIZER; }
87 
89 {
90  G4AutoLock l(&EmParametersMutex);
91 
92  lossFluctuation = true;
93  buildCSDARange = false;
94  flagLPM = true;
95  spline = true;
96  finalRange = false;
97  applyCuts = false;
98  fluo = false;
99  beardenFluoDir = false;
100  auger = false;
101  augerCascade = false;
102  pixe = false;
103  deexIgnoreCut = false;
104  lateralDisplacement = true;
105  muhadLateralDisplacement = false;
108  useMottCorrection = false;
109 
110  minSubRange = 1.0;
111  minKinEnergy = 0.1*keV;
112  maxKinEnergy = 10.0*TeV;
113  maxKinEnergyCSDA = 1.0*GeV;
115  lowestMuHadEnergy = 1.0*keV;
116  linLossLimit = 0.01;
118  lambdaFactor = 0.8;
119  factorForAngleLimit = 1.0;
121  rangeFactor = 0.04;
122  rangeFactorMuHad = 0.2;
123  geomFactor = 2.5;
124  skin = 1.0;
125 
126  nbins = 77;
127  nbinsPerDecade = 7;
128  verbose = 1;
129  workerVerbose = 0;
130 
133 
134  namePIXE = "Empirical";
135  nameElectronPIXE = "Livermore";
136 }
137 
139 {
140  G4AutoLock l(&EmParametersMutex);
141  lossFluctuation = val;
142 }
143 
145 {
146  return lossFluctuation;
147 }
148 
150 {
151  G4AutoLock l(&EmParametersMutex);
152  buildCSDARange = val;
153 }
154 
156 {
157  return buildCSDARange;
158 }
159 
161 {
162  G4AutoLock l(&EmParametersMutex);
163  flagLPM = val;
164 }
165 
167 {
168  return flagLPM;
169 }
170 
172 {
173  G4AutoLock l(&EmParametersMutex);
174  spline = val;
175 }
176 
178 {
179  return spline;
180 }
181 
183 {
184  G4AutoLock l(&EmParametersMutex);
185  finalRange = val;
186 }
187 
189 {
190  return finalRange;
191 }
192 
194 {
195  G4AutoLock l(&EmParametersMutex);
196  applyCuts = val;
197 }
198 
200 {
201  return applyCuts;
202 }
203 
205 {
206  G4AutoLock l(&EmParametersMutex);
207  fluo = val;
208 }
209 
211 {
212  return fluo;
213 }
214 
216 {
217  G4AutoLock l(&EmParametersMutex);
218  beardenFluoDir = val;
219 }
220 
222 {
223  return beardenFluoDir;
224 }
225 
227 {
228  G4AutoLock l(&EmParametersMutex);
229  auger = val;
230  if(val) { fluo = true; }
231 }
232 
234 {
235  return auger;
236 }
237 
239 {
240  G4AutoLock l(&EmParametersMutex);
241  augerCascade = val;
242  if(val) { fluo = true; auger = true; }
243 }
244 
246 {
247  return augerCascade;
248 }
249 
251 {
252  G4AutoLock l(&EmParametersMutex);
253  pixe = val;
254  if(val) { fluo = true; }
255 }
256 
258 {
259  return pixe;
260 }
261 
263 {
264  G4AutoLock l(&EmParametersMutex);
265  deexIgnoreCut = val;
266 }
267 
269 {
270  return deexIgnoreCut;
271 }
272 
274 {
275  G4AutoLock l(&EmParametersMutex);
276  lateralDisplacement = val;
277 }
278 
280 {
281  return lateralDisplacement;
282 }
283 
285 {
286  G4AutoLock l(&EmParametersMutex);
288 }
289 
291 {
293 }
294 
296 {
297  G4AutoLock l(&EmParametersMutex);
299 }
300 
302 {
304 }
305 
307 {
309 }
310 
312 {
314 }
315 
317 {
318  useMottCorrection = val;
319 }
320 
322 {
323  return useMottCorrection;
324 }
325 
327 {
328  G4AutoLock l(&EmParametersMutex);
329  if(val > 0.0 && val < 1.0) {
330  minSubRange = val;
331  } else {
333  ed << "Value of MinSubRange is out of range (0 - 1): " << val
334  << " is ignored";
335  PrintWarning(ed);
336  }
337 }
338 
340 {
341  return minSubRange;
342 }
343 
345 {
346  G4AutoLock l(&EmParametersMutex);
347  if(val > 1.e-3*eV && val < maxKinEnergy) {
348  minKinEnergy = val;
350  } else {
352  ed << "Value of MinKinEnergy is out of range: " << val/MeV
353  << " MeV is ignored";
354  PrintWarning(ed);
355  }
356 }
357 
359 {
360  return minKinEnergy;
361 }
362 
364 {
365  G4AutoLock l(&EmParametersMutex);
366  if(val > minKinEnergy && val < 1.e+7*TeV) {
367  maxKinEnergy = val;
369  } else {
371  ed << "Value of MaxKinEnergy is out of range: "
372  << val/GeV << " GeV is ignored";
373  PrintWarning(ed);
374  }
375 }
376 
378 {
379  return maxKinEnergy;
380 }
381 
383 {
384  G4AutoLock l(&EmParametersMutex);
385  if(val > minKinEnergy && val <= 100*TeV) {
386  maxKinEnergyCSDA = val;
387  } else {
389  ed << "Value of MaxKinEnergyCSDA is out of range: "
390  << val/GeV << " GeV is ignored";
391  PrintWarning(ed);
392  }
393 }
394 
396 {
397  return maxKinEnergyCSDA;
398 }
399 
401 {
402  G4AutoLock l(&EmParametersMutex);
403  if(val >= 0.0) {
404  lowestElectronEnergy = val;
405  } else {
407  ed << "Value of lowestElectronEnergy is out of range: "
408  << val/MeV << " MeV is ignored";
409  PrintWarning(ed);
410  }
411 }
412 
414 {
415  return lowestElectronEnergy;
416 }
417 
419 {
420  G4AutoLock l(&EmParametersMutex);
421  if(val >= 0.0) {
422  lowestMuHadEnergy = val;
423  } else {
425  ed << "Value of lowestMuHadEnergy is out of range: "
426  << val/MeV << " MeV is ignored";
427  PrintWarning(ed);
428  }
429 }
430 
432 {
433  return lowestMuHadEnergy;
434 }
435 
437 {
438  G4AutoLock l(&EmParametersMutex);
439  if(val > 0.0 && val < 0.5) {
440  linLossLimit = val;
441  } else {
443  ed << "Value of linLossLimit is out of range: " << val
444  << " is ignored";
445  PrintWarning(ed);
446  }
447 }
448 
450 {
451  return linLossLimit;
452 }
453 
455 {
456  G4AutoLock l(&EmParametersMutex);
457  if(val > 0.0) {
458  bremsTh = val;
459  } else {
461  ed << "Value of bremsstrahlung threshold is out of range: "
462  << val/GeV << " GeV is ignored";
463  PrintWarning(ed);
464  }
465 }
466 
468 {
469  return bremsTh;
470 }
471 
473 {
474  G4AutoLock l(&EmParametersMutex);
475  if(val > 0.0 && val < 1.0) {
476  lambdaFactor = val;
477  } else {
479  ed << "Value of lambda factor is out of range: " << val
480  << " is ignored";
481  PrintWarning(ed);
482  }
483 }
484 
486 {
487  return lambdaFactor;
488 }
489 
491 {
492  G4AutoLock l(&EmParametersMutex);
493  if(val > 0.0) {
494  factorForAngleLimit = val;
495  } else {
497  ed << "Value of factor for enegry limit is out of range: "
498  << val << " is ignored";
499  PrintWarning(ed);
500  }
501 }
502 
504 {
505  return factorForAngleLimit;
506 }
507 
509 {
510  G4AutoLock l(&EmParametersMutex);
511  if(val >= 0.0 && val <= pi) {
512  thetaLimit = val;
513  } else {
515  ed << "Value of polar angle limit is out of range: "
516  << val << " is ignored";
517  PrintWarning(ed);
518  }
519 }
520 
522 {
523  return thetaLimit;
524 }
525 
527 {
528  G4AutoLock l(&EmParametersMutex);
529  if(val > 0.0 && val < 1.0) {
530  //G4cout << " G4EmParameters::SetMscRangeFactor: " << val << G4endl;
531  rangeFactor = val;
532  } else {
534  ed << "Value of rangeFactor is out of range: "
535  << val << " is ignored";
536  PrintWarning(ed);
537  }
538 }
539 
541 {
542  return rangeFactor;
543 }
544 
546 {
547  G4AutoLock l(&EmParametersMutex);
548  if(val > 0.0 && val < 1.0) {
549  rangeFactorMuHad = val;
550  } else {
552  ed << "Value of rangeFactorMuHad is out of range: "
553  << val << " is ignored";
554  PrintWarning(ed);
555  }
556 }
557 
559 {
560  return rangeFactorMuHad;
561 }
562 
564 {
565  G4AutoLock l(&EmParametersMutex);
566  if(val >= 1.0) {
567  geomFactor = val;
568  } else {
570  ed << "Value of geomFactor is out of range: "
571  << val << " is ignored";
572  PrintWarning(ed);
573  }
574 }
575 
577 {
578  return geomFactor;
579 }
580 
582 {
583  G4AutoLock l(&EmParametersMutex);
584  if(val >= 1.0) {
585  skin = val;
586  } else {
588  ed << "Value of skin is out of range: "
589  << val << " is ignored";
590  PrintWarning(ed);
591  }
592 }
593 
595 {
596  return skin;
597 }
598 
600 {
601  G4AutoLock l(&EmParametersMutex);
602  if(val >= 5 && val < 10000000) {
603  nbins = val;
605  } else {
607  ed << "Value of number of bins is out of range: "
608  << val << " is ignored";
609  PrintWarning(ed);
610  }
611 }
612 
614 {
615  return nbins;
616 }
617 
619 {
620  G4AutoLock l(&EmParametersMutex);
621  if(val >= 5 && val < 1000000) {
622  nbinsPerDecade = val;
624  } else {
626  ed << "Value of number of bins per decade is out of range: "
627  << val << " is ignored";
628  PrintWarning(ed);
629  }
630 }
631 
633 {
634  return nbinsPerDecade;
635 }
636 
638 {
639  G4AutoLock l(&EmParametersMutex);
640  verbose = val;
642 }
643 
645 {
646  return verbose;
647 }
648 
650 {
651  G4AutoLock l(&EmParametersMutex);
652  workerVerbose = val;
653 }
654 
656 {
657  return workerVerbose;
658 }
659 
661 {
662  G4AutoLock l(&EmParametersMutex);
663  mscStepLimit = val;
664 }
665 
667 {
668  return mscStepLimit;
669 }
670 
672 {
673  G4AutoLock l(&EmParametersMutex);
674  mscStepLimitMuHad = val;
675 }
676 
678 {
679  return mscStepLimitMuHad;
680 }
681 
683 {
684  G4cout << "G4EmParameters::SetPIXECrossSectionModel " << sss << G4endl;
685  G4AutoLock l(&EmParametersMutex);
686  namePIXE = sss;
687 }
688 
690 {
691  return namePIXE;
692 }
693 
695 {
696  G4AutoLock l(&EmParametersMutex);
698 }
699 
701 {
702  return nameElectronPIXE;
703 }
704 
706 {
707  G4Exception("G4EmParameters", "em0044", JustWarning, ed);
708 }
709 
711  const G4String& region,
712  const G4String& type)
713 {
714  G4String r = region;
715  if(r == "" || r == "world" || r == "World") r = "DefaultRegionForTheWorld";
716  G4int nreg = m_regnamesPAI.size();
717  for(G4int i=0; i<nreg; ++i) {
718  if((m_particlesPAI[i] == particle ||
719  m_particlesPAI[i] == "all" ||
720  particle == "all") &&
721  (m_regnamesPAI[i] == r ||
722  m_regnamesPAI[i] == "DefaultRegionForTheWorld" ||
723  r == "DefaultRegionForTheWorld") ) {
724 
725  m_typesPAI[i] = type;
726  if(particle == "all") m_particlesPAI[i] = particle;
727  if(r == "DefaultRegionForTheWorld") { m_regnamesPAI[i] = r; }
728  return;
729  }
730  }
731  m_particlesPAI.push_back(particle);
732  m_regnamesPAI.push_back(r);
733  m_typesPAI.push_back(type);
734 }
735 
736 const std::vector<G4String>& G4EmParameters::ParticlesPAI() const
737 {
738  return m_particlesPAI;
739 }
740 
741 const std::vector<G4String>& G4EmParameters::RegionsPAI() const
742 {
743  return m_regnamesPAI;
744 }
745 
746 const std::vector<G4String>& G4EmParameters::TypesPAI() const
747 {
748  return m_typesPAI;
749 }
750 
752 {
753  G4String r = region;
754  if(r == "" || r == "world" || r == "World") r = "DefaultRegionForTheWorld";
755  G4int nreg = m_regnamesME.size();
756  for(G4int i=0; i<nreg; ++i) {
757  if(r == m_regnamesME[i]) { return; }
758  }
759  m_regnamesME.push_back(r);
760 }
761 
762 const std::vector<G4String>& G4EmParameters::RegionsMicroElec() const
763 {
764  return m_regnamesME;
765 }
766 
767 void G4EmParameters::AddDNA(const G4String& region, const G4String& type)
768 {
769  G4String r = region;
770  if(r == "" || r == "world" || r == "World") r = "DefaultRegionForTheWorld";
771  G4int nreg = m_regnamesDNA.size();
772  for(G4int i=0; i<nreg; ++i) {
773  if(r == m_regnamesDNA[i]) { return; }
774  }
775  m_regnamesDNA.push_back(r);
776  m_typesDNA.push_back(type);
777 }
778 
779 const std::vector<G4String>& G4EmParameters::RegionsDNA() const
780 {
781  return m_regnamesDNA;
782 }
783 
784 const std::vector<G4String>& G4EmParameters::TypesDNA() const
785 {
786  return m_typesDNA;
787 }
788 
789 std::ostream& G4EmParameters::StreamInfo(std::ostream& os) const
790 {
791  G4int prec = os.precision(5);
792  os << "=======================================================================" << "\n";
793  os << "====== Electromagnetic Physics Parameters ========" << "\n";
794  os << "=======================================================================" << "\n";
795  os << "Fluctuations of dE/dx are enabled " <<lossFluctuation << "\n";
796  os << "Build CSDA range enabled " <<buildCSDARange << "\n";
797  os << "LPM effect enabled " <<flagLPM << "\n";
798  os << "Spline of EM tables enabled " <<spline << "\n";
799  os << "Use cut as a final range enabled " <<finalRange << "\n";
800  os << "Apply cuts on all EM processes " <<applyCuts << "\n";
801  os << "Fluorescence enabled " <<fluo << "\n";
802  os << "Fluorescence Bearden data files enabled " <<beardenFluoDir << "\n";
803  os << "Auger electron production enabled " <<auger << "\n";
804  os << "Auger cascade enabled " <<augerCascade << "\n";
805  os << "PIXE atomic de-excitation enabled " <<pixe << "\n";
806  os << "De-excitation module ignores cuts " <<deexIgnoreCut << "\n";
807  os << "Msc lateraral displacement for e+- enabled " <<lateralDisplacement << "\n";
808  os << "Msc lateraral displacement for muons and hadrons " <<muhadLateralDisplacement << "\n";
809  os << "Msc lateraral displacement beyond geometry safety " <<latDisplacementBeyondSafety << "\n";
810  os << "Enable angular generator interface "
812  os << "Use Mott correction for e- scattering "
813  <<useMottCorrection << "\n";
814 
815  os << "Factor of cut reduction for sub-cutoff method " <<minSubRange << "\n";
816  os << "Min kinetic energy for tables "
817  <<G4BestUnit(minKinEnergy,"Energy") << "\n";
818  os << "Max kinetic energy for tables "
819  <<G4BestUnit(maxKinEnergy,"Energy") << "\n";
820  os << "Max kinetic energy for CSDA tables "
821  <<G4BestUnit(maxKinEnergyCSDA,"Energy") << "\n";
822  os << "Lowest e+e- kinetic energy "
823  <<G4BestUnit(lowestElectronEnergy,"Energy") << "\n";
824  os << "Lowest muon/hadron kinetic energy "
825  <<G4BestUnit(lowestMuHadEnergy,"Energy") << "\n";
826  os << "Linear loss limit " <<linLossLimit << "\n";
827  os << "Bremsstrahlung energy threshold above which \n"
828  << " primary is added to the list of secondary "
829  <<G4BestUnit(bremsTh,"Energy") << "\n";
830  os << "X-section factor for integral approach " <<lambdaFactor << "\n";
831  os << "Factor used for dynamic computation of angular \n"
832  << " limit between single and multiple scattering " << factorForAngleLimit << "\n";
833  os << "Fixed angular limit between single \n"
834  << " and multiple scattering "
835  <<thetaLimit/rad << " rad" << "\n";
836  os << "Range factor for msc step limit for e+- " <<rangeFactor << "\n";
837  os << "Range factor for msc step limit for muons/hadrons " <<rangeFactorMuHad << "\n";
838  os << "Geometry factor for msc step limitation of e+- " <<geomFactor << "\n";
839  os << "Skin parameter for msc step limitation of e+- " <<skin << "\n";
840 
841  os << "Number of bins in tables " <<nbins << "\n";
842  os << "Number of bins per decade of a table " <<nbinsPerDecade << "\n";
843  os << "Verbose level " <<verbose << "\n";
844  os << "Verbose level for worker thread " <<workerVerbose << "\n";
845 
846  os << "Type of msc step limit algorithm for e+- " <<mscStepLimit << "\n";
847  os << "Type of msc step limit algorithm for muons/hadrons " <<mscStepLimitMuHad << "\n";
848 
849  os << "Type of PIXE cross section for hadrons " <<namePIXE << "\n";
850  os << "Type of PIXE cross section for e+- " <<nameElectronPIXE << "\n";
851  os << "=======================================================================" << "\n";
852  os.precision(prec);
853  return os;
854 }
855 
857 {
859 }
860 
861 std::ostream& operator<< (std::ostream& os, const G4EmParameters& par)
862 {
863  return par.StreamInfo(os);
864 }
865 
866 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
867 
G4bool UseCutAsFinalRange() const
static G4EmParameters * theInstance
void SetLossFluctuations(G4bool val)
G4int NumberOfBinsPerDecade() const
void SetApplyCuts(G4bool val)
G4MscStepLimitType mscStepLimit
static const double MeV
Definition: G4SIunits.hh:211
G4int NumberOfBins() const
void SetVerbose(G4int val)
G4bool Spline() const
G4int WorkerVerbose() const
G4double MaxKinEnergy() const
std::vector< G4String > m_typesDNA
void SetDeexcitationIgnoreCut(G4bool val)
void SetUseMottCorrection(G4bool val)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
std::vector< G4String > m_typesPAI
void SetLowestElectronEnergy(G4double val)
void SetLatDisplacementBeyondSafety(G4bool val)
G4MscStepLimitType MscMuHadStepLimitType() const
void SetMscStepLimitType(G4MscStepLimitType val)
G4double linLossLimit
G4bool AugerCascade() const
G4double LowestElectronEnergy() const
void SetBeardenFluoDir(G4bool val)
void SetLinearLossLimit(G4double val)
G4double MscGeomFactor() const
G4double MscMuHadRangeFactor() const
G4double MscThetaLimit() const
std::vector< G4String > m_regnamesME
void SetAuger(G4bool val)
const std::vector< G4String > & ParticlesPAI() const
G4String nameElectronPIXE
const std::vector< G4String > & RegionsMicroElec() const
G4double maxKinEnergy
G4bool muhadLateralDisplacement
std::ostream & StreamInfo(std::ostream &os) const
void SetNumberOfBins(G4int val)
void SetMinSubRange(G4double val)
G4bool LPM() const
void SetMaxEnergyForCSDARange(G4double val)
const std::vector< G4String > & RegionsPAI() const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4bool ApplyCuts() const
const std::vector< G4String > & TypesDNA() const
G4bool Fluo() const
void SetPIXEElectronCrossSectionModel(const G4String &)
int G4int
Definition: G4Types.hh:78
G4bool useMottCorrection
void SetMaxEnergy(G4double val)
void SetBremsstrahlungTh(G4double val)
static G4NistManager * Instance()
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:175
G4double MinSubRange() const
G4bool DeexcitationIgnoreCut() const
G4bool LatDisplacementBeyondSafety() const
void AddPAIModel(const G4String &particle, const G4String &region, const G4String &type)
void SetLateralDisplacement(G4bool val)
G4bool BuildCSDARange() const
void SetWorkerVerbose(G4int val)
G4bool lateralDisplacement
G4double LinearLossLimit() const
void SetPIXECrossSectionModel(const G4String &)
G4double LambdaFactor() const
static const double prec
Definition: RanecuEngine.cc:58
G4GLOB_DLL std::ostream G4cout
G4int Verbose() const
G4bool MuHadLateralDisplacement() const
const G4String & PIXECrossSectionModel()
G4double LowestMuHadEnergy() const
G4double lowestMuHadEnergy
bool G4bool
Definition: G4Types.hh:79
void SetMscRangeFactor(G4double val)
void SetAugerCascade(G4bool val)
void Dump() const
void SetNumberOfBinsPerDecade(G4int val)
G4bool LateralDisplacement() const
G4double minSubRange
static const double GeV
Definition: G4SIunits.hh:214
void SetLowestMuHadEnergy(G4double val)
void SetMscGeomFactor(G4double val)
std::vector< G4String > m_regnamesDNA
G4double MinKinEnergy() const
void SetMscMuHadStepLimitType(G4MscStepLimitType val)
std::vector< G4String > m_regnamesPAI
void AddDNA(const G4String &region, const G4String &type)
std::vector< G4String > m_particlesPAI
G4double lowestElectronEnergy
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const double rad
Definition: G4SIunits.hh:148
G4MscStepLimitType mscStepLimitMuHad
G4double lambdaFactor
static const double eV
Definition: G4SIunits.hh:212
void SetBuildCSDARange(G4bool val)
G4int G4Mutex
Definition: G4Threading.hh:173
void AddMicroElec(const G4String &region)
G4bool Auger() const
G4bool LossFluctuation() const
void SetMscMuHadRangeFactor(G4double val)
void SetPixe(G4bool val)
void SetSpline(G4bool val)
static const char sss[MAX_N_PAR+2]
Definition: Evaluator.cc:64
static const double pi
Definition: G4SIunits.hh:74
void SetMuHadLateralDisplacement(G4bool val)
int G4lrint(double ad)
Definition: templates.hh:163
G4double maxKinEnergyCSDA
void SetMinEnergy(G4double val)
void SetLPM(G4bool val)
G4bool UseAngularGeneratorForIonisation() const
G4MscStepLimitType
G4double BremsstrahlungTh() const
static G4EmParameters * Instance()
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double MscRangeFactor() const
void ActivateAngularGeneratorForIonisation(G4bool val)
void SetUseCutAsFinalRange(G4bool val)
G4bool latDisplacementBeyondSafety
G4bool UseMottCorrection() const
const std::vector< G4String > & RegionsDNA() const
#define G4endl
Definition: G4ios.hh:61
static const double TeV
Definition: G4SIunits.hh:215
G4bool Pixe() const
G4double minKinEnergy
static const double keV
Definition: G4SIunits.hh:213
G4MscStepLimitType MscStepLimitType() const
void SetMscThetaLimit(G4double val)
void SetLambdaFactor(G4double val)
void SetFactorForAngleLimit(G4double val)
G4double rangeFactorMuHad
double G4double
Definition: G4Types.hh:76
G4bool useAngGeneratorForIonisation
void PrintWarning(G4ExceptionDescription &ed)
G4EmParametersMessenger * theMessenger
G4double MaxEnergyForCSDARange() const
G4double MscSkin() const
const G4String & PIXEElectronCrossSectionModel()
const std::vector< G4String > & TypesPAI() const
void SetFluo(G4bool val)
G4double FactorForAngleLimit() const
G4bool BeardenFluoDir() const
std::ostream & operator<<(std::ostream &os, const G4EmParameters &par)
void SetMscSkin(G4double val)
G4double rangeFactor
G4double factorForAngleLimit