Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// -------------------------------------------------------------------
27//
28// GEANT4 Class file
29//
30// File name: G4EmParameters
31//
32// Author: Vladimir Ivanchenko
33//
34// Creation date: 18.05.2013
35//
36// Modifications:
37//
38// -------------------------------------------------------------------
39//
40//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
42
43#include "G4EmParameters.hh"
45#include "G4UnitsTable.hh"
46#include "G4SystemOfUnits.hh"
47#include "G4VEmProcess.hh"
51#include "G4EmLowEParameters.hh"
53#include "G4EmUtility.hh"
54#include "G4NistManager.hh"
55#include "G4RegionStore.hh"
56#include "G4Region.hh"
57#include "G4ApplicationState.hh"
58#include "G4StateManager.hh"
59#include "G4Threading.hh"
60#include "G4AutoLock.hh"
61
62G4EmParameters* G4EmParameters::theInstance = nullptr;
63
64namespace
65{
66 G4Mutex emParametersMutex = G4MUTEX_INITIALIZER;
67}
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
70
72{
73 if(nullptr == theInstance) {
74 G4AutoLock l(&emParametersMutex);
75 if(nullptr == theInstance) {
76 static G4EmParameters manager;
77 theInstance = &manager;
78 }
79 l.unlock();
80 }
81 return theInstance;
82}
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
85
87{
88 delete theMessenger;
89 delete fBParameters;
90 delete fCParameters;
91 delete emSaturation;
92}
93
94//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
95
97{
99 theMessenger = new G4EmParametersMessenger(this);
100 Initialise();
101
102 fBParameters = new G4EmExtraParameters();
103 fCParameters = new G4EmLowEParameters();
104
105 fStateManager = G4StateManager::GetStateManager();
106 emSaturation = nullptr;
107}
108
110{
111 if(!IsLocked()) {
112 Initialise();
113 fBParameters->Initialise();
114 fCParameters->Initialise();
115 }
116}
117
118void G4EmParameters::Initialise()
119{
120 lossFluctuation = true;
121 buildCSDARange = false;
122 flagLPM = true;
123 cutAsFinalRange = false;
124 applyCuts = false;
125 lateralDisplacement = true;
126 lateralDisplacementAlg96 = true;
127 muhadLateralDisplacement = false;
128 useAngGeneratorForIonisation = false;
129 useMottCorrection = false;
130 integral = true;
131 birks = false;
132 fICRU90 = false;
133 gener = false;
134 onIsolated = false;
135 fSamplingTable = false;
136 fPolarisation = false;
137 fMuDataFromFile = false;
138 fPEKShell = true;
139 fMscPosiCorr = true;
140 fUseEPICS2017XS = false;
141 f3GammaAnnihilationOnFly = false;
142 fUseRiGePairProductionModel = false;
143 fDNA = false;
144 fIsPrinted = false;
145
146 minKinEnergy = 0.1*CLHEP::keV;
147 maxKinEnergy = 100.0*CLHEP::TeV;
148 maxKinEnergyCSDA = 1.0*CLHEP::GeV;
149 max5DEnergyForMuPair = 0.0;
150 lowestElectronEnergy = 1.0*CLHEP::keV;
151 lowestMuHadEnergy = 1.0*CLHEP::keV;
152 lowestTripletEnergy = 1.0*CLHEP::MeV;
153 maxNIELEnergy = 0.0;
154 linLossLimit = 0.01;
155 bremsTh = bremsMuHadTh = maxKinEnergy;
156 lambdaFactor = 0.8;
157 factorForAngleLimit = 1.0;
158 thetaLimit = CLHEP::pi;
159 energyLimit = 100.0*CLHEP::MeV;
160 rangeFactor = 0.04;
161 rangeFactorMuHad = 0.2;
162 geomFactor = 2.5;
163 skin = 1.0;
164 safetyFactor = 0.6;
165 lambdaLimit = 1.0*CLHEP::mm;
166 factorScreen = 1.0;
167
168 nbinsPerDecade = 7;
169 verbose = 1;
170 workerVerbose = 0;
171 nForFreeVector = 2;
172 tripletConv = 0;
173
174 fTransportationWithMsc = G4TransportationWithMscType::fDisabled;
175 mscStepLimit = fUseSafety;
176 mscStepLimitMuHad = fMinimal;
177 nucFormfactor = fExponentialNF;
178 fSStype = fWVI;
179 fFluct = fUniversalFluctuation;
180 fPositronium = fSimplePositronium;
181
182 const char* data_dir = G4FindDataDir("G4LEDATA");
183 if (nullptr != data_dir) {
184 fDirLEDATA = G4String(data_dir);
185 }
186 else {
187 G4Exception("G4EmParameters::Initialise()", "em0003", JustWarning,
188 "G4LEDATA data directory was not found.");
189 }
190}
191
193{
194 if(IsLocked()) { return; }
195 lossFluctuation = val;
196}
197
199{
200 return lossFluctuation;
201}
202
204{
205 if(IsLocked()) { return; }
206 buildCSDARange = val;
207}
208
210{
211 return buildCSDARange;
212}
213
215{
216 if(IsLocked()) { return; }
217 flagLPM = val;
218}
219
221{
222 return flagLPM;
223}
224
226{
227 if(IsLocked()) { return; }
228 cutAsFinalRange = val;
229}
230
232{
233 return cutAsFinalRange;
234}
235
237{
238 if(IsLocked()) { return; }
239 applyCuts = val;
240}
241
243{
244 return applyCuts;
245}
246
248{
249 if(IsLocked()) { return; }
250 fCParameters->SetFluo(val);
251}
252
254{
255 return fCParameters->Fluo();
256}
257
259{
260 return fCParameters->FluoDirectory();
261}
262
264{
265 if(IsLocked()) { return; }
266 fCParameters->SetFluoDirectory(val);
267}
268
270{
271 if(IsLocked()) { return; }
272 fCParameters->SetBeardenFluoDir(val);
273}
274
276{
277 if(IsLocked()) { return; }
278 fCParameters->SetANSTOFluoDir(val);
279}
280
282{
283 if(IsLocked()) { return; }
284 fCParameters->SetXDB_EADLFluoDir(val);
285}
286
288{
289 if(IsLocked()) { return; }
290 fCParameters->SetAuger(val);
291}
292
294{
295 auto dir = fCParameters->FluoDirectory();
296 return (dir == fluoBearden);
297}
298
300{
301 auto dir = fCParameters->FluoDirectory();
302 return (dir == fluoANSTO);
303}
304
306{
307 return fCParameters->Auger();
308}
309
311{
312 if(IsLocked()) { return; }
313 fCParameters->SetPixe(val);
314}
315
317{
318 return fCParameters->Pixe();
319}
320
322{
323 if(IsLocked()) { return; }
324 fCParameters->SetDeexcitationIgnoreCut(val);
325}
326
328{
329 return fCParameters->DeexcitationIgnoreCut();
330}
331
333{
334 if(IsLocked()) { return; }
335 lateralDisplacement = val;
336}
337
339{
340 return lateralDisplacement;
341}
342
344{
345 if(IsLocked()) { return; }
346 lateralDisplacementAlg96 = val;
347}
348
350{
351 return lateralDisplacementAlg96;
352}
353
355{
356 if(IsLocked()) { return; }
357 muhadLateralDisplacement = val;
358}
359
361{
362 return muhadLateralDisplacement;
363}
364
366{
367 if(IsLocked()) { return; }
368 useAngGeneratorForIonisation = val;
369}
370
372{
373 return useAngGeneratorForIonisation;
374}
375
377{
378 if(IsLocked()) { return; }
379 useMottCorrection = val;
380}
381
383{
384 return useMottCorrection;
385}
386
388{
389 if(IsLocked()) { return; }
390 integral = val;
391}
392
394{
395 return integral;
396}
397
399{
400 if(IsLocked()) { return; }
401 fPolarisation = val;
402}
403
405{
406 return fPolarisation;
407}
408
410{
411 if(IsLocked()) { return; }
412 birks = val;
413 if(birks && nullptr == emSaturation) { emSaturation = new G4EmSaturation(1); }
414}
415
417{
418 return birks;
419}
420
422{
423 if(IsLocked()) { return; }
424 fICRU90 = val;
425}
426
428{
429 return fICRU90;
430}
431
433{
434 if(IsLocked()) { return; }
435 fCParameters->SetDNAFast(val);
436 if(val) { ActivateDNA(); }
437}
438
440{
441 return fCParameters->DNAFast();
442}
443
445{
446 if(IsLocked()) { return; }
447 fCParameters->SetDNAStationary(val);
448 if(val) { ActivateDNA(); }
449}
450
452{
453 return fCParameters->DNAStationary();
454}
455
457{
458 if(IsLocked()) { return; }
459 fCParameters->SetDNAElectronMsc(val);
460 if(val) { ActivateDNA(); }
461}
462
464{
465 return fCParameters->DNAElectronMsc();
466}
467
469{
470 if(IsLocked()) { return; }
471 gener = val;
472}
473
475{
476 return gener;
477}
478
480{
481 if(IsLocked()) { return; }
482 birks = (nullptr != ptr);
483 if(emSaturation != ptr) {
484 delete emSaturation;
485 emSaturation = ptr;
486 }
487}
488
490{
491 return fMuDataFromFile;
492}
493
495{
496 fMuDataFromFile = v;
497}
498
500{
501 if(IsLocked()) { return; }
502 onIsolated = val;
503}
504
506{
507 return onIsolated;
508}
509
511{
512 if(IsLocked()) { return; }
513 fSamplingTable = val;
514}
515
517{
518 return fSamplingTable;
519}
520
522{
523 return fPEKShell;
524}
525
527{
528 if(IsLocked()) { return; }
529 fPEKShell = v;
530}
531
533{
534 return fMscPosiCorr;
535}
536
538{
539 if(IsLocked()) { return; }
540 fMscPosiCorr = v;
541}
542
544{
545 return fUseEPICS2017XS;
546}
547
549{
550 if(IsLocked()) { return; }
551 fUseEPICS2017XS = v;
552}
553
555{
556 return f3GammaAnnihilationOnFly;
557}
558
560{
561 if(IsLocked()) { return; }
562 f3GammaAnnihilationOnFly = v;
563}
564
566{
567 return fUseRiGePairProductionModel;
568}
569
571{
572 if (IsLocked()) { return; }
573 fUseRiGePairProductionModel = v;
574}
575
577{
578 if(IsLocked()) { return; }
579 fDNA = true;
580}
581
583{
584 fIsPrinted = val;
585}
586
588{
589 return fIsPrinted;
590}
591
593{
594 if (nullptr == emSaturation) {
595 G4AutoLock l(&emParametersMutex);
596 if (nullptr == emSaturation) {
597 emSaturation = new G4EmSaturation(1);
598 }
599 l.unlock();
600 }
601 birks = true;
602 return emSaturation;
603}
604
606{
607 if(IsLocked()) { return; }
608 if(val > 1.e-3*CLHEP::eV && val < maxKinEnergy) {
609 minKinEnergy = val;
610 } else {
612 ed << "Value of MinKinEnergy - is out of range: " << val/CLHEP::MeV
613 << " MeV is ignored";
614 PrintWarning(ed);
615 }
616}
617
619{
620 return minKinEnergy;
621}
622
624{
625 if(IsLocked()) { return; }
626 if(val > std::max(minKinEnergy,599.9*CLHEP::MeV) && val < 1.e+7*CLHEP::TeV) {
627 maxKinEnergy = val;
628 } else {
630 ed << "Value of MaxKinEnergy is out of range: "
631 << val/CLHEP::GeV
632 << " GeV is ignored; allowed range 600 MeV - 1.e+7 TeV";
633 PrintWarning(ed);
634 }
635}
636
638{
639 return maxKinEnergy;
640}
641
643{
644 if(IsLocked()) { return; }
645 if(val > minKinEnergy && val <= 100*CLHEP::TeV) {
646 maxKinEnergyCSDA = val;
647 } else {
649 ed << "Value of MaxKinEnergyCSDA is out of range: "
650 << val/CLHEP::GeV << " GeV is ignored; allowed range "
651 << minKinEnergy << " MeV - 100 TeV";
652 PrintWarning(ed);
653 }
654}
655
657{
658 return maxKinEnergyCSDA;
659}
660
662{
663 if(IsLocked()) { return; }
664 if(val >= 0.0) { lowestElectronEnergy = val; }
665}
666
668{
669 return lowestElectronEnergy;
670}
671
673{
674 if(IsLocked()) { return; }
675 if(val >= 0.0) { lowestMuHadEnergy = val; }
676}
677
679{
680 return lowestMuHadEnergy;
681}
682
684{
685 if(IsLocked()) { return; }
686 if(val > 0.0) { lowestTripletEnergy = val; }
687}
688
690{
691 return lowestTripletEnergy;
692}
693
695{
696 if(IsLocked()) { return; }
697 if(val >= 0.0) { maxNIELEnergy = val; }
698}
699
701{
702 return maxNIELEnergy;
703}
704
706{
707 if(IsLocked()) { return; }
708 if(val > 0.0) { max5DEnergyForMuPair = val; }
709}
710
712{
713 return max5DEnergyForMuPair;
714}
715
717{
718 if(IsLocked()) { return; }
719 if(val > 0.0 && val < 1.0) {
720 linLossLimit = val;
721 } else {
723 ed << "Value of linLossLimit is out of range: " << val
724 << " is ignored";
725 PrintWarning(ed);
726 }
727}
728
730{
731 return linLossLimit;
732}
733
735{
736 if(IsLocked()) { return; }
737 if(val > 0.0) {
738 bremsTh = val;
739 } else {
741 ed << "Value of bremsstrahlung threshold is out of range: "
742 << val/GeV << " GeV is ignored";
743 PrintWarning(ed);
744 }
745}
746
748{
749 return bremsTh;
750}
751
753{
754 if(IsLocked()) { return; }
755 if(val > 0.0) {
756 bremsMuHadTh = val;
757 } else {
759 ed << "Value of bremsstrahlung threshold is out of range: "
760 << val/GeV << " GeV is ignored";
761 PrintWarning(ed);
762 }
763}
764
766{
767 return bremsMuHadTh;
768}
769
771{
772 if(IsLocked()) { return; }
773 if(val > 0.0 && val < 1.0) {
774 lambdaFactor = val;
775 } else {
777 ed << "Value of lambda factor is out of range: " << val
778 << " is ignored";
779 PrintWarning(ed);
780 }
781}
782
784{
785 return lambdaFactor;
786}
787
789{
790 if(IsLocked()) { return; }
791 if(val > 0.0) {
792 factorForAngleLimit = val;
793 } else {
795 ed << "Value of factor for enegry limit is out of range: "
796 << val << " is ignored";
797 PrintWarning(ed);
798 }
799}
800
802{
803 return factorForAngleLimit;
804}
805
807{
808 if(IsLocked()) { return; }
809 if(val >= 0.0 && val <= pi) {
810 thetaLimit = val;
811 } else {
813 ed << "Value of polar angle limit is out of range: "
814 << val << " is ignored";
815 PrintWarning(ed);
816 }
817}
818
820{
821 return thetaLimit;
822}
823
825{
826 if(IsLocked()) { return; }
827 if(val >= 0.0) {
828 energyLimit = val;
829 } else {
831 ed << "Value of msc energy limit is out of range: "
832 << val << " is ignored";
833 PrintWarning(ed);
834 }
835}
836
838{
839 return energyLimit;
840}
841
843{
844 if(IsLocked()) { return; }
845 if(val > 0.0 && val < 1.0) {
846 rangeFactor = val;
847 } else {
849 ed << "Value of rangeFactor is out of range: "
850 << val << " is ignored";
851 PrintWarning(ed);
852 }
853}
854
856{
857 return rangeFactor;
858}
859
861{
862 if(IsLocked()) { return; }
863 if(val > 0.0 && val < 1.0) {
864 rangeFactorMuHad = val;
865 } else {
867 ed << "Value of rangeFactorMuHad is out of range: "
868 << val << " is ignored";
869 PrintWarning(ed);
870 }
871}
872
874{
875 return rangeFactorMuHad;
876}
877
879{
880 if(IsLocked()) { return; }
881 if(val >= 1.0) {
882 geomFactor = val;
883 } else {
885 ed << "Value of geomFactor is out of range: "
886 << val << " is ignored";
887 PrintWarning(ed);
888 }
889}
890
892{
893 return geomFactor;
894}
895
897{
898 if(IsLocked()) { return; }
899 if(val >= 0.1) {
900 safetyFactor = val;
901 } else {
903 ed << "Value of safetyFactor is out of range: "
904 << val << " is ignored";
905 PrintWarning(ed);
906 }
907}
908
910{
911 return safetyFactor;
912}
913
915{
916 if(IsLocked()) { return; }
917 if(val >= 0.0) {
918 lambdaLimit = val;
919 } else {
921 ed << "Value of lambdaLimit is out of range: "
922 << val << " is ignored";
923 PrintWarning(ed);
924 }
925}
926
928{
929 return lambdaLimit;
930}
931
933{
934 if(IsLocked()) { return; }
935 if(val >= 1.0) {
936 skin = val;
937 } else {
939 ed << "Value of skin is out of range: "
940 << val << " is ignored";
941 PrintWarning(ed);
942 }
943}
944
946{
947 return skin;
948}
949
951{
952 if(IsLocked()) { return; }
953 if(val > 0.0) {
954 factorScreen = val;
955 } else {
957 ed << "Value of factorScreen is out of range: "
958 << val << " is ignored";
959 PrintWarning(ed);
960 }
961}
962
964{
965 return factorScreen;
966}
967
969{
970 if (IsLocked()) { return; }
971 fCParameters->SetMaxDNAElectronEnergy(val);
972}
973
975{
976 return fCParameters->MaxDNAElectronEnergy();
977}
978
980{
981 if (IsLocked()) { return; }
982 fCParameters->SetMaxDNAIonEnergy(val);
983}
984
986{
987 return fCParameters->MaxDNAIonEnergy();
988}
989
991{
992 if(IsLocked()) { return; }
993 fBParameters->SetStepFunction(v1, v2);
994}
995
997{
998 if(IsLocked()) { return; }
999 fBParameters->SetStepFunctionMuHad(v1, v2);
1000}
1001
1003{
1004 if(IsLocked()) { return; }
1005 fBParameters->SetStepFunctionLightIons(v1, v2);
1006}
1007
1009{
1010 if(IsLocked()) { return; }
1011 fBParameters->SetStepFunctionIons(v1, v2);
1012}
1013
1015{
1016 fBParameters->FillStepFunction(part, proc);
1017}
1018
1020{
1021 return nbinsPerDecade*G4lrint(std::log10(maxKinEnergy/minKinEnergy));
1022}
1023
1025{
1026 if(IsLocked()) { return; }
1027 if(val >= 5 && val < 1000000) {
1028 nbinsPerDecade = val;
1029 } else {
1031 ed << "Value of number of bins per decade is out of range: "
1032 << val << " is ignored";
1033 PrintWarning(ed);
1034 }
1035}
1036
1038{
1039 return nbinsPerDecade;
1040}
1041
1043{
1044 if(IsLocked()) { return; }
1045 verbose = val;
1046 workerVerbose = std::min(workerVerbose, verbose);
1047}
1048
1050{
1051 return verbose;
1052}
1053
1055{
1056 if(IsLocked()) { return; }
1057 workerVerbose = val;
1058}
1059
1061{
1062 return workerVerbose;
1063}
1064
1066{
1067 if(IsLocked()) { return; }
1068 nForFreeVector = val;
1069}
1070
1072{
1073 return nForFreeVector;
1074}
1075
1077{
1078 if(IsLocked()) { return; }
1079 fTransportationWithMsc = val;
1080}
1081
1083{
1084 return fTransportationWithMsc;
1085}
1086
1088{
1089 if(IsLocked()) { return; }
1090 fFluct = val;
1091}
1092
1094{
1095 return fFluct;
1096}
1097
1099{
1100 if(IsLocked()) { return; }
1101 fPositronium = val;
1102}
1103
1105{
1106 return fPositronium;
1107}
1108
1110{
1111 if(IsLocked()) { return; }
1112 mscStepLimit = val;
1113}
1114
1116{
1117 return mscStepLimit;
1118}
1119
1121{
1122 if(IsLocked()) { return; }
1123 mscStepLimitMuHad = val;
1124}
1125
1127{
1128 return mscStepLimitMuHad;
1129}
1130
1132{
1133 if(IsLocked()) { return; }
1134 fSStype = val;
1135}
1136
1138{
1139 return fSStype;
1140}
1141
1142void
1144{
1145 if(IsLocked()) { return; }
1146 nucFormfactor = val;
1147}
1148
1150{
1151 return nucFormfactor;
1152}
1153
1155{
1156 if(IsLocked()) { return; }
1157 fCParameters->SetDNAeSolvationSubType(val);
1158 ActivateDNA();
1159}
1160
1162{
1163 return fCParameters->DNAeSolvationSubType();
1164}
1165
1167{
1168 if(IsLocked()) { return; }
1169 tripletConv = val;
1170}
1171
1173{
1174 return tripletConv;
1175}
1176
1178{
1179 if(IsLocked()) { return; }
1180 fCParameters->SetPIXECrossSectionModel(sss);
1181}
1182
1184{
1185 return fCParameters->PIXECrossSectionModel();
1186}
1187
1189{
1190 if(IsLocked()) { return; }
1191 fCParameters->SetPIXEElectronCrossSectionModel(sss);
1192}
1193
1195{
1196 return fCParameters->PIXEElectronCrossSectionModel();
1197}
1198
1200{
1201 if(IsLocked()) { return; }
1202 fCParameters->SetLivermoreDataDir(sss);
1203}
1204
1206{
1207 return fCParameters->LivermoreDataDir();
1208}
1209
1210void G4EmParameters::PrintWarning(G4ExceptionDescription& ed) const
1211{
1212 G4Exception("G4EmParameters", "em0044", JustWarning, ed);
1213}
1214
1216 const G4String& region,
1217 const G4String& type)
1218{
1219 if(IsLocked()) { return; }
1220 fBParameters->AddPAIModel(particle, region, type);
1221}
1222
1223const std::vector<G4String>& G4EmParameters::ParticlesPAI() const
1224{
1225 return fBParameters->ParticlesPAI();
1226}
1227
1228const std::vector<G4String>& G4EmParameters::RegionsPAI() const
1229{
1230 return fBParameters->RegionsPAI();
1231}
1232
1233const std::vector<G4String>& G4EmParameters::TypesPAI() const
1234{
1235 return fBParameters->TypesPAI();
1236}
1237
1239{
1240 if(IsLocked()) { return; }
1241 fCParameters->AddMicroElec(region);
1242}
1243
1244const std::vector<G4String>& G4EmParameters::RegionsMicroElec() const
1245{
1246 return fCParameters->RegionsMicroElec();
1247}
1248
1249void G4EmParameters::AddDNA(const G4String& region, const G4String& type)
1250{
1251 if(IsLocked()) { return; }
1252 fCParameters->AddDNA(region, type);
1253 ActivateDNA();
1254}
1255
1256const std::vector<G4String>& G4EmParameters::RegionsDNA() const
1257{
1258 return fCParameters->RegionsDNA();
1259}
1260
1261const std::vector<G4String>& G4EmParameters::TypesDNA() const
1262{
1263 return fCParameters->TypesDNA();
1264}
1265
1266void G4EmParameters::AddPhysics(const G4String& region, const G4String& type)
1267{
1268 if(IsLocked()) { return; }
1269 fBParameters->AddPhysics(region, type);
1270}
1271
1272const std::vector<G4String>& G4EmParameters::RegionsPhysics() const
1273{
1274 return fBParameters->RegionsPhysics();
1275}
1276
1277const std::vector<G4String>& G4EmParameters::TypesPhysics() const
1278{
1279 return fBParameters->TypesPhysics();
1280}
1281
1283{
1284 if(IsLocked()) { return; }
1285 fBParameters->SetSubCutRegion(region);
1286}
1287
1288void
1290 G4bool aauger, G4bool apixe)
1291{
1292 if(IsLocked()) { return; }
1293 fCParameters->SetDeexActiveRegion(region, adeex, aauger, apixe);
1294}
1295
1296void
1298 G4double val, G4bool wflag)
1299{
1300 if(IsLocked()) { return; }
1301 fBParameters->SetProcessBiasingFactor(procname, val, wflag);
1302}
1303
1304void
1306 const G4String& region,
1307 G4double length,
1308 G4bool wflag)
1309{
1310 if(IsLocked() && !gener) { return; }
1311 fBParameters->ActivateForcedInteraction(procname, region, length, wflag);
1312}
1313
1314void
1316 const G4String& region,
1317 G4double factor,
1318 G4double energyLim)
1319{
1320 if(IsLocked()) { return; }
1321 fBParameters->ActivateSecondaryBiasing(procname, region, factor, energyLim);
1322}
1323
1325{
1326 fBParameters->DefineRegParamForLoss(ptr);
1327}
1328
1330{
1331 fBParameters->DefineRegParamForEM(ptr);
1332}
1333
1335{
1336 return fBParameters->QuantumEntanglement();
1337}
1338
1340{
1341 if(IsLocked()) { return; }
1342 fBParameters->SetQuantumEntanglement(v);
1343}
1344
1346 return fBParameters->GetDirectionalSplitting();
1347}
1348
1350{
1351 if(IsLocked()) { return; }
1352 fBParameters->SetDirectionalSplitting(v);
1353}
1354
1356{
1357 if(IsLocked()) { return; }
1358 fBParameters->SetDirectionalSplittingTarget(v);
1359}
1360
1362{
1363 return fBParameters->GetDirectionalSplittingTarget();
1364}
1365
1367{
1368 if(IsLocked()) { return; }
1369 fBParameters->SetDirectionalSplittingRadius(r);
1370}
1371
1373{
1374 return fBParameters->GetDirectionalSplittingRadius();
1375}
1376
1378{
1379 fCParameters->DefineRegParamForDeex(ptr);
1380}
1381
1383{
1384 return fDirLEDATA;
1385}
1386
1388{
1389 if (IsLocked()) { return; }
1390 G4String ss = fBParameters->CheckRegion(nam);
1391 if (!fluctRegions.empty()) {
1392 for (auto const& p : fluctRegions) { if (p.first == ss) { return; } }
1393 }
1394 fluctRegions.push_back(std::make_pair(ss, flag));
1395}
1396
1397void G4EmParameters::DefineFluctuationFlags(std::vector<G4bool>* theFlags)
1398{
1399 if (fluctRegions.empty()) { return; }
1400 G4EmUtility::FillFluctFlags(fluctRegions, theFlags);
1401}
1402
1403void G4EmParameters::StreamInfo(std::ostream& os) const
1404{
1405 G4long prec = os.precision(5);
1406 os << "=======================================================================" << "\n";
1407 os << "====== Electromagnetic Physics Parameters ========" << "\n";
1408 os << "=======================================================================" << "\n";
1409 os << "LPM effect enabled " <<flagLPM << "\n";
1410 os << "Enable creation and use of sampling tables " <<fSamplingTable << "\n";
1411 os << "Apply cuts on all EM processes " <<applyCuts << "\n";
1412 const char* transportationWithMsc = "Disabled";
1413 if(fTransportationWithMsc == G4TransportationWithMscType::fEnabled) {
1414 transportationWithMsc = "Enabled";
1415 } else if (fTransportationWithMsc == G4TransportationWithMscType::fMultipleSteps) {
1416 transportationWithMsc = "MultipleSteps";
1417 }
1418 os << "Use combined TransportationWithMsc " <<transportationWithMsc << "\n";
1419 os << "Use general process " <<gener << "\n";
1420 os << "Enable linear polarisation for gamma " <<fPolarisation << "\n";
1421 os << "Enable photoeffect sampling below K-shell " <<fPEKShell << "\n";
1422 os << "Enable sampling of quantum entanglement "
1423 <<fBParameters->QuantumEntanglement() << "\n";
1424 os << "X-section factor for integral approach " <<lambdaFactor << "\n";
1425 os << "Min kinetic energy for tables "
1426 <<G4BestUnit(minKinEnergy,"Energy") << "\n";
1427 os << "Max kinetic energy for tables "
1428 <<G4BestUnit(maxKinEnergy,"Energy") << "\n";
1429 os << "Number of bins per decade of a table " <<nbinsPerDecade << "\n";
1430 os << "Verbose level " <<verbose << "\n";
1431 os << "Verbose level for worker thread " <<workerVerbose << "\n";
1432 os << "Bremsstrahlung energy threshold above which \n"
1433 << " primary e+- is added to the list of secondary "
1434 <<G4BestUnit(bremsTh,"Energy") << "\n";
1435 os << "Bremsstrahlung energy threshold above which primary\n"
1436 << " muon/hadron is added to the list of secondary "
1437 <<G4BestUnit(bremsMuHadTh,"Energy") << "\n";
1438 G4String name3g = "SimplePositronium";
1439 if (fPositronium == fAllisonPositronium) { name3g = "AllisonPositronium"; }
1440 else if (fPositronium == fOrePowell) { name3g = "OrePowell"; }
1441 else if (fPositronium == fOrePowellPolar) { name3g = "OrePowellPolar"; }
1442 os << "Positron annihilation at rest model " << name3g << "\n";
1443
1444 os << "Enable 3 gamma annihilation on fly "
1445 << f3GammaAnnihilationOnFly << "\n";
1446 os << "Lowest triplet kinetic energy "
1447 <<G4BestUnit(lowestTripletEnergy,"Energy") << "\n";
1448 os << "Enable sampling of gamma linear polarisation " <<fPolarisation << "\n";
1449 os << "5D gamma conversion model type " <<tripletConv << "\n";
1450 os << "5D gamma conversion model on isolated ion " <<onIsolated << "\n";
1451 if(max5DEnergyForMuPair>0.0) {
1452 os << "5D gamma conversion limit for muon pair "
1453 << max5DEnergyForMuPair/CLHEP::GeV << " GeV\n";
1454 }
1455 os << "Use RiGe 5D e+e- pair production model by muons "
1456 << fUseRiGePairProductionModel << "\n";
1457 os << "Livermore data directory "
1458 << fCParameters->LivermoreDataDir() << "\n";
1459
1460 os << "=======================================================================" << "\n";
1461 os << "====== Ionisation Parameters ========" << "\n";
1462 os << "=======================================================================" << "\n";
1463 os << "Step function for e+- "
1464 <<"("<<fBParameters->GetStepFunctionP1() << ", "
1465 << fBParameters->GetStepFunctionP2()/CLHEP::mm << " mm)\n";
1466 os << "Step function for muons/hadrons "
1467 <<"("<<fBParameters->GetStepFunctionMuHadP1() << ", "
1468 << fBParameters->GetStepFunctionMuHadP2()/CLHEP::mm << " mm)\n";
1469 os << "Step function for light ions "
1470 <<"("<<fBParameters->GetStepFunctionLightIonsP1() << ", "
1471 << fBParameters->GetStepFunctionLightIonsP2()/CLHEP::mm << " mm)\n";
1472 os << "Step function for general ions "
1473 <<"("<<fBParameters->GetStepFunctionIonsP1() << ", "
1474 << fBParameters->GetStepFunctionIonsP2()/CLHEP::mm << " mm)\n";
1475 os << "Lowest e+e- kinetic energy "
1476 <<G4BestUnit(lowestElectronEnergy,"Energy") << "\n";
1477 os << "Lowest muon/hadron kinetic energy "
1478 <<G4BestUnit(lowestMuHadEnergy,"Energy") << "\n";
1479 os << "Use ICRU90 data " << fICRU90 << "\n";
1480 os << "Fluctuations of dE/dx are enabled " << lossFluctuation << "\n";
1481 G4String namef = "Universal";
1482 if(fFluct == fUrbanFluctuation) { namef = "Urban"; }
1483 else if(fFluct == fDummyFluctuation) { namef = "Dummy"; }
1484 os << "Type of fluctuation model for leptons and hadrons " << namef << "\n";
1485 if (!fluctRegions.empty()) {
1486 if (lossFluctuation) {
1487 os << "Fluctuations of dE/dx are disabled in G4Regions: ";
1488 } else {
1489 os << "Fluctuations of dE/dx are enabled in G4Regions: ";
1490 }
1491 G4int n = 0;
1492 for (auto const& p : fluctRegions) {
1493 if (p.second != lossFluctuation) {
1494 os << p.first << " ";
1495 if (n > 0 && (n/2)*2 == n) {
1496 os << "\n" << " ";
1497 }
1498 ++n;
1499 }
1500 }
1501 os << "\n";
1502 }
1503 os << "Use built-in Birks saturation " << birks << "\n";
1504 os << "Build CSDA range enabled " <<buildCSDARange << "\n";
1505 os << "Use cut as a final range enabled " <<cutAsFinalRange << "\n";
1506 os << "Enable angular generator interface "
1507 <<useAngGeneratorForIonisation << "\n";
1508 os << "Max kinetic energy for CSDA tables "
1509 <<G4BestUnit(maxKinEnergyCSDA,"Energy") << "\n";
1510 os << "Max kinetic energy for NIEL computation "
1511 <<G4BestUnit(maxNIELEnergy,"Energy") << "\n";
1512 os << "Linear loss limit " <<linLossLimit << "\n";
1513 os << "Read data from file for e+e- pair production by mu " <<fMuDataFromFile << "\n";
1514
1515 os << "=======================================================================" << "\n";
1516 os << "====== Multiple Scattering Parameters ========" << "\n";
1517 os << "=======================================================================" << "\n";
1518 os << "Type of msc step limit algorithm for e+- " <<mscStepLimit << "\n";
1519 os << "Type of msc step limit algorithm for muons/hadrons " <<mscStepLimitMuHad << "\n";
1520 os << "Msc lateral displacement for e+- enabled " <<lateralDisplacement << "\n";
1521 os << "Msc lateral displacement for muons and hadrons " <<muhadLateralDisplacement << "\n";
1522 os << "Urban msc model lateral displacement alg96 " <<lateralDisplacementAlg96 << "\n";
1523 os << "Range factor for msc step limit for e+- " <<rangeFactor << "\n";
1524 os << "Range factor for msc step limit for muons/hadrons " <<rangeFactorMuHad << "\n";
1525 os << "Geometry factor for msc step limitation of e+- " <<geomFactor << "\n";
1526 os << "Safety factor for msc step limit for e+- " <<safetyFactor << "\n";
1527 os << "Skin parameter for msc step limitation of e+- " <<skin << "\n";
1528 os << "Lambda limit for msc step limit for e+- " <<lambdaLimit/CLHEP::mm << " mm\n";
1529 os << "Use Mott correction for e- scattering " << useMottCorrection << "\n";
1530 os << "Factor used for dynamic computation of angular \n"
1531 << " limit between single and multiple scattering " << factorForAngleLimit << "\n";
1532 os << "Fixed angular limit between single \n"
1533 << " and multiple scattering "
1534 << thetaLimit/CLHEP::rad << " rad\n";
1535 os << "Upper energy limit for e+- multiple scattering "
1536 << energyLimit/CLHEP::MeV << " MeV\n";
1537 os << "Type of electron single scattering model " <<fSStype << "\n";
1538 os << "Type of nuclear form-factor " <<nucFormfactor << "\n";
1539 os << "Screening factor " <<factorScreen << "\n";
1540 os << "=======================================================================" << "\n";
1541
1542 if(fCParameters->Fluo()) {
1543 os << "====== Atomic Deexcitation Parameters ========" << "\n";
1544 os << "=======================================================================" << "\n";
1545 os << "Fluorescence enabled " <<fCParameters->Fluo() << "\n";
1546 G4String named = "fluor";
1548 if(fdir == fluoBearden) { named = "fluor_Bearden"; }
1549 else if(fdir == fluoANSTO) { named = "fluor_ANSTO"; }
1550 else if(fdir == fluoXDB_EADL) { named = "fluor_XDB_EADL"; }
1551 os << "Directory in G4LEDATA for fluorescence data files " << named << "\n";
1552 os << "Auger electron cascade enabled "
1553 <<fCParameters->Auger() << "\n";
1554 os << "PIXE atomic de-excitation enabled " <<fCParameters->Pixe() << "\n";
1555 os << "De-excitation module ignores cuts "
1556 <<fCParameters->DeexcitationIgnoreCut() << "\n";
1557 os << "Type of PIXE cross section for hadrons "
1558 <<fCParameters->PIXECrossSectionModel() << "\n";
1559 os << "Type of PIXE cross section for e+- "
1560 <<fCParameters->PIXEElectronCrossSectionModel() << "\n";
1561 os << "=======================================================================" << "\n";
1562 }
1563 if(fDNA) {
1564 os << "====== DNA Physics Parameters ========" << "\n";
1565 os << "=======================================================================" << "\n";
1566 os << "Max DNA electron energy "
1567 << fCParameters->MaxDNAElectronEnergy()/CLHEP::MeV << " MeV\n";
1568 os << "Max DNA ion energy "
1569 << fCParameters->MaxDNAIonEnergy()/CLHEP::MeV << " MeV\n";
1570 os << "Use fast sampling in DNA models "
1571 << fCParameters->DNAFast() << "\n";
1572 os << "Use Stationary option in DNA models "
1573 << fCParameters->DNAStationary() << "\n";
1574 os << "Use DNA with multiple scattering of e- "
1575 << fCParameters->DNAElectronMsc() << "\n";
1576 os << "Use DNA e- solvation model type "
1577 << fCParameters->DNAeSolvationSubType() << "\n";
1578 auto chemModel = fCParameters->GetChemTimeStepModel();
1579 if(fCParameters->GetChemTimeStepModel() != G4ChemTimeStepModel::Unknown)
1580 {
1581 std::vector<G4String> ChemModel{"Unknown","SBS","IRT","IRT_syn"};
1582 os << "Use DNA Chemistry model "
1583 << ChemModel.at((std::size_t)chemModel) << "\n";
1584 }
1585 os << "=======================================================================" << G4endl;
1586 }
1587 os.precision(prec);
1588}
1589
1591{
1592 if(fIsPrinted) return;
1593
1594#ifdef G4MULTITHREADED
1595 G4MUTEXLOCK(&emParametersMutex);
1596#endif
1598#ifdef G4MULTITHREADED
1599 G4MUTEXUNLOCK(&emParametersMutex);
1600#endif
1601}
1602
1603std::ostream& operator<< (std::ostream& os, const G4EmParameters& par)
1604{
1605 par.StreamInfo(os);
1606 return os;
1607}
1608
1609G4bool G4EmParameters::IsLocked() const
1610{
1611 return (!G4Threading::IsMasterThread() ||
1612 (fStateManager->GetCurrentState() != G4State_PreInit &&
1613 fStateManager->GetCurrentState() != G4State_Init &&
1614 fStateManager->GetCurrentState() != G4State_Idle));
1615}
1616
1617
1619{
1620 fCParameters-> SetChemTimeStepModel(model);
1621}
1622
1624{
1625 return fCParameters->GetChemTimeStepModel();
1626}
1627//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
@ G4State_Init
@ G4State_Idle
@ G4State_PreInit
G4TemplateAutoLock< G4Mutex > G4AutoLock
G4ChemTimeStepModel
G4DNAModelSubType
G4EmFluoDirectory
@ fluoBearden
@ fluoXDB_EADL
@ fluoANSTO
std::ostream & operator<<(std::ostream &os, const G4EmParameters &par)
G4EmFluctuationType
@ fUniversalFluctuation
@ fUrbanFluctuation
@ fDummyFluctuation
G4eSingleScatteringType
@ fWVI
G4TransportationWithMscType
G4PositronAtRestModelType
@ fSimplePositronium
@ fAllisonPositronium
@ fOrePowell
@ fOrePowellPolar
const char * G4FindDataDir(const char *)
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4MscStepLimitType
@ fUseSafety
#define G4BestUnit(a, b)
#define G4MUTEX_INITIALIZER
#define G4MUTEXLOCK(mutex)
#define G4MUTEXUNLOCK(mutex)
std::mutex G4Mutex
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
long G4long
Definition G4Types.hh:87
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void SetEmSaturation(G4EmSaturation *)
void SetBeardenFluoDir(G4bool val)
G4bool IsPrintLocked() const
void DefineRegParamForLoss(G4VEnergyLossProcess *) const
void SetLambdaFactor(G4double val)
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetBuildCSDARange(G4bool val)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void SetEnablePolarisation(G4bool val)
void AddDNA(const G4String &region, const G4String &type)
G4bool LateralDisplacementAlg96() const
void FillStepFunction(const G4ParticleDefinition *, G4VEnergyLossProcess *) const
G4EmParameters(G4EmParameters &)=delete
void SetNumberOfBinsPerDecade(G4int val)
G4bool UseEPICS2017XS() const
static G4EmParameters * Instance()
void SetDirectionalSplittingTarget(const G4ThreeVector &v)
G4bool DNAFast() const
G4DNAModelSubType DNAeSolvationSubType() const
G4bool RetrieveMuDataFromFile() const
G4bool PhotoeffectBelowKShell() const
G4int NumberOfBins() const
G4double MscMuHadRangeFactor() const
void SetGeneralProcessActive(G4bool val)
void SetDNAFast(G4bool val)
void SetMscSafetyFactor(G4double val)
G4bool EnablePolarisation() const
void SetMaxDNAElectronEnergy(G4double val)
void SetLateralDisplacementAlg96(G4bool val)
void SetFactorForAngleLimit(G4double val)
const G4String & PIXECrossSectionModel()
const G4String & PIXEElectronCrossSectionModel()
G4double MaxNIELEnergy() const
void SetRetrieveMuDataFromFile(G4bool v)
void SetDirectionalSplitting(G4bool v)
const G4String & LivermoreDataDir()
G4bool OnIsolated() const
void SetMscMuHadRangeFactor(G4double val)
G4bool DNAElectronMsc() const
G4double MinKinEnergy() const
G4int NumberOfBinsPerDecade() const
G4bool BuildCSDARange() const
G4double MscThetaLimit() const
G4bool LossFluctuation() const
void SetLPM(G4bool val)
G4double MuHadBremsstrahlungTh() const
void AddPAIModel(const G4String &particle, const G4String &region, const G4String &type)
void SetPositronAtRestModelType(G4PositronAtRestModelType val)
void SetDNAStationary(G4bool val)
void SetDNAElectronMsc(G4bool val)
const std::vector< G4String > & TypesPhysics() const
void SetMaxEnergyFor5DMuPair(G4double val)
G4double GetDirectionalSplittingRadius()
void SetLinearLossLimit(G4double val)
void SetMscThetaLimit(G4double val)
G4int GetConversionType() const
G4MscStepLimitType MscMuHadStepLimitType() const
G4double MscEnergyLimit() const
void ActivateSecondaryBiasing(const G4String &name, const G4String &region, G4double factor, G4double energyLimit)
G4bool BirksActive() const
G4bool DNAStationary() const
G4bool MscPositronCorrection() const
const std::vector< G4String > & RegionsPAI() const
G4int NumberForFreeVector() const
void SetSubCutRegion(const G4String &region="")
void SetLossFluctuations(G4bool val)
void SetDeexActiveRegion(const G4String &region, G4bool fdeex, G4bool fauger, G4bool fpixe)
G4bool UseCutAsFinalRange() const
void SetPIXEElectronCrossSectionModel(const G4String &)
void SetLowestTripletEnergy(G4double val)
void SetMuHadLateralDisplacement(G4bool val)
G4bool Fluo() const
G4EmSaturation * GetEmSaturation()
void SetDNAeSolvationSubType(G4DNAModelSubType val)
void SetQuantumEntanglement(G4bool v)
void SetNumberForFreeVector(G4int val)
void DefineRegParamForEM(G4VEmProcess *) const
void ActivateForcedInteraction(const G4String &procname, const G4String &region, G4double length, G4bool wflag)
void SetXDB_EADLFluoDir(G4bool val)
void ActivateAngularGeneratorForIonisation(G4bool val)
void SetScreeningFactor(G4double val)
void SetNuclearFormfactorType(G4NuclearFormfactorType val)
void SetStepFunction(G4double v1, G4double v2)
void SetLateralDisplacement(G4bool val)
G4int Verbose() const
G4bool QuantumEntanglement() const
void SetWorkerVerbose(G4int val)
G4double MaxDNAElectronEnergy() const
void SetUseCutAsFinalRange(G4bool val)
void SetDeexcitationIgnoreCut(G4bool val)
void SetBirksActive(G4bool val)
G4double ScreeningFactor() const
G4bool UseMottCorrection() const
void SetFluo(G4bool val)
G4double MaxDNAIonEnergy() const
void SetMuHadBremsstrahlungTh(G4double val)
const std::vector< G4String > & ParticlesPAI() const
G4bool Pixe() const
const std::vector< G4String > & RegionsPhysics() const
void SetFluctuationType(G4EmFluctuationType val)
void DefineRegParamForDeex(G4VAtomDeexcitation *) const
void SetStepFunctionMuHad(G4double v1, G4double v2)
G4double MscSafetyFactor() const
void SetVerbose(G4int val)
G4int WorkerVerbose() const
void SetMscGeomFactor(G4double val)
void SetMscLambdaLimit(G4double val)
void SetMscSkin(G4double val)
void SetApplyCuts(G4bool val)
const std::vector< G4String > & TypesDNA() const
void SetUseRiGePairProductionModel(G4bool v)
G4MscStepLimitType MscStepLimitType() const
G4bool ApplyCuts() const
void SetEnableSamplingTable(G4bool val)
const std::vector< G4String > & RegionsDNA() const
void SetLivermoreDataDir(const G4String &)
void SetFluoDirectory(G4EmFluoDirectory)
G4double BremsstrahlungTh() const
void SetPixe(G4bool val)
void SetMaxNIELEnergy(G4double val)
void SetStepFunctionIons(G4double v1, G4double v2)
G4PositronAtRestModelType PositronAtRestModelType() const
void SetMaxEnergyForCSDARange(G4double val)
G4eSingleScatteringType SingleScatteringType() const
void SetMaxDNAIonEnergy(G4double val)
G4bool DeexcitationIgnoreCut() const
G4TransportationWithMscType TransportationWithMsc() const
void SetProcessBiasingFactor(const G4String &procname, G4double val, G4bool wflag)
G4double MscGeomFactor() const
void SetMscMuHadStepLimitType(G4MscStepLimitType val)
G4EmFluctuationType FluctuationType() const
void SetMscStepLimitType(G4MscStepLimitType val)
void AddPhysics(const G4String &region, const G4String &type)
void SetMscEnergyLimit(G4double val)
void SetBremsstrahlungTh(G4double val)
void SetAuger(G4bool val)
void SetUseEPICS2017XS(G4bool v)
void DefineFluctuationFlags(std::vector< G4bool > *theFluctFlags)
G4bool GetDirectionalSplitting() const
void SetIsPrintedFlag(G4bool val)
G4double MaxKinEnergy() const
G4bool UseICRU90Data() const
void SetDirectionalSplittingRadius(G4double r)
void SetConversionType(G4int val)
G4bool LateralDisplacement() const
void SetUseICRU90Data(G4bool val)
G4ChemTimeStepModel GetTimeStepModel() const
G4bool MuHadLateralDisplacement() const
void SetOnIsolated(G4bool val)
void SetTransportationWithMsc(G4TransportationWithMscType val)
G4bool EnableSamplingTable() const
G4EmFluoDirectory FluoDirectory() const
void StreamInfo(std::ostream &os) const
G4double MscLambdaLimit() const
void SetPIXECrossSectionModel(const G4String &)
void SetIntegral(G4bool val)
G4ThreeVector GetDirectionalSplittingTarget() const
G4bool Auger() const
void SetUseMottCorrection(G4bool val)
void AddMicroElec(const G4String &region)
const std::vector< G4String > & RegionsMicroElec() const
G4double MaxEnergyFor5DMuPair() const
G4bool Integral() const
G4bool BeardenFluoDir()
void SetMscPositronCorrection(G4bool v)
G4double MaxEnergyForCSDARange() const
G4bool UseRiGePairProductionModel() const
void SetLowestMuHadEnergy(G4double val)
void SetFluctuationsForRegion(const G4String &regionName, G4bool flag)
const std::vector< G4String > & TypesPAI() const
G4bool LPM() const
G4bool UseAngularGeneratorForIonisation() const
void SetMaxEnergy(G4double val)
G4double LinearLossLimit() const
G4NuclearFormfactorType NuclearFormfactorType() const
G4double LowestMuHadEnergy() const
G4double MscRangeFactor() const
G4double LambdaFactor() const
G4double FactorForAngleLimit() const
G4bool Use3GammaAnnihilationOnFly() const
const G4String & GetDirLEDATA() const
G4double MscSkin() const
void SetSingleScatteringType(G4eSingleScatteringType val)
void SetANSTOFluoDir(G4bool val)
G4double LowestTripletEnergy() const
void SetPhotoeffectBelowKShell(G4bool v)
G4double LowestElectronEnergy() const
void SetMscRangeFactor(G4double val)
void SetTimeStepModel(const G4ChemTimeStepModel &model)
void Set3GammaAnnihilationOnFly(G4bool v)
G4bool GeneralProcessActive() const
static void FillFluctFlags(std::vector< std::pair< G4String, G4bool > > &reg, std::vector< G4bool > *flags)
static G4NistManager * Instance()
const G4ApplicationState & GetCurrentState() const
static G4StateManager * GetStateManager()
G4bool IsMasterThread()
int G4lrint(double ad)
Definition templates.hh:134