Geant4 11.3.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 "G4NistManager.hh"
54#include "G4RegionStore.hh"
55#include "G4Region.hh"
56#include "G4ApplicationState.hh"
57#include "G4StateManager.hh"
58#include "G4Threading.hh"
59#include "G4AutoLock.hh"
60
61G4EmParameters* G4EmParameters::theInstance = nullptr;
62
63namespace
64{
65 G4Mutex emParametersMutex = G4MUTEX_INITIALIZER;
66}
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
69
71{
72 if(nullptr == theInstance) {
73 G4AutoLock l(&emParametersMutex);
74 if(nullptr == theInstance) {
75 static G4EmParameters manager;
76 theInstance = &manager;
77 }
78 l.unlock();
79 }
80 return theInstance;
81}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
84
86{
87 delete theMessenger;
88 delete fBParameters;
89 delete fCParameters;
90 delete emSaturation;
91}
92
93//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
94
96{
98 theMessenger = new G4EmParametersMessenger(this);
99 Initialise();
100
101 fBParameters = new G4EmExtraParameters();
102 fCParameters = new G4EmLowEParameters();
103
104 fStateManager = G4StateManager::GetStateManager();
105 emSaturation = nullptr;
106}
107
109{
110 if(!IsLocked()) {
111 Initialise();
112 fBParameters->Initialise();
113 fCParameters->Initialise();
114 }
115}
116
117void G4EmParameters::Initialise()
118{
119 lossFluctuation = true;
120 buildCSDARange = false;
121 flagLPM = true;
122 cutAsFinalRange = false;
123 applyCuts = false;
124 lateralDisplacement = true;
125 lateralDisplacementAlg96 = true;
126 muhadLateralDisplacement = false;
127 useAngGeneratorForIonisation = false;
128 useMottCorrection = false;
129 integral = true;
130 birks = false;
131 fICRU90 = false;
132 gener = false;
133 onIsolated = false;
134 fSamplingTable = false;
135 fPolarisation = false;
136 fMuDataFromFile = false;
137 fPEKShell = true;
138 fMscPosiCorr = true;
139 fUseEPICS2017XS = false;
140 f3GammaAnnihilationOnFly = false;
141 fUseRiGePairProductionModel = false;
142 fDNA = false;
143 fIsPrinted = false;
144
145 minKinEnergy = 0.1*CLHEP::keV;
146 maxKinEnergy = 100.0*CLHEP::TeV;
147 maxKinEnergyCSDA = 1.0*CLHEP::GeV;
148 max5DEnergyForMuPair = 0.0;
149 lowestElectronEnergy = 1.0*CLHEP::keV;
150 lowestMuHadEnergy = 1.0*CLHEP::keV;
151 lowestTripletEnergy = 1.0*CLHEP::MeV;
152 maxNIELEnergy = 0.0;
153 linLossLimit = 0.01;
154 bremsTh = bremsMuHadTh = maxKinEnergy;
155 lambdaFactor = 0.8;
156 factorForAngleLimit = 1.0;
157 thetaLimit = CLHEP::pi;
158 energyLimit = 100.0*CLHEP::MeV;
159 rangeFactor = 0.04;
160 rangeFactorMuHad = 0.2;
161 geomFactor = 2.5;
162 skin = 1.0;
163 safetyFactor = 0.6;
164 lambdaLimit = 1.0*CLHEP::mm;
165 factorScreen = 1.0;
166
167 nbinsPerDecade = 7;
168 verbose = 1;
169 workerVerbose = 0;
170 nForFreeVector = 2;
171 tripletConv = 0;
172
173 fTransportationWithMsc = G4TransportationWithMscType::fDisabled;
174 mscStepLimit = fUseSafety;
175 mscStepLimitMuHad = fMinimal;
176 nucFormfactor = fExponentialNF;
177 fSStype = fWVI;
178 fFluct = fUniversalFluctuation;
179 fPositronium = fSimplePositronium;
180
181 const char* data_dir = G4FindDataDir("G4LEDATA");
182 if (nullptr != data_dir) {
183 fDirLEDATA = G4String(data_dir);
184 }
185 else {
186 G4Exception("G4EmParameters::Initialise()", "em0003", JustWarning,
187 "G4LEDATA data directory was not found.");
188 }
189}
190
192{
193 if(IsLocked()) { return; }
194 lossFluctuation = val;
195}
196
198{
199 return lossFluctuation;
200}
201
203{
204 if(IsLocked()) { return; }
205 buildCSDARange = val;
206}
207
209{
210 return buildCSDARange;
211}
212
214{
215 if(IsLocked()) { return; }
216 flagLPM = val;
217}
218
220{
221 return flagLPM;
222}
223
225{
226 if(IsLocked()) { return; }
227 cutAsFinalRange = val;
228}
229
231{
232 return cutAsFinalRange;
233}
234
236{
237 if(IsLocked()) { return; }
238 applyCuts = val;
239}
240
242{
243 return applyCuts;
244}
245
247{
248 if(IsLocked()) { return; }
249 fCParameters->SetFluo(val);
250}
251
253{
254 return fCParameters->Fluo();
255}
256
258{
259 return fCParameters->FluoDirectory();
260}
261
263{
264 if(IsLocked()) { return; }
265 fCParameters->SetFluoDirectory(val);
266}
267
269{
270 if(IsLocked()) { return; }
271 fCParameters->SetBeardenFluoDir(val);
272}
273
275{
276 if(IsLocked()) { return; }
277 fCParameters->SetANSTOFluoDir(val);
278}
279
281{
282 if(IsLocked()) { return; }
283 fCParameters->SetXDB_EADLFluoDir(val);
284}
285
287{
288 if(IsLocked()) { return; }
289 fCParameters->SetAuger(val);
290}
291
293{
294 auto dir = fCParameters->FluoDirectory();
295 return (dir == fluoBearden);
296}
297
299{
300 auto dir = fCParameters->FluoDirectory();
301 return (dir == fluoANSTO);
302}
303
305{
306 return fCParameters->Auger();
307}
308
310{
311 if(IsLocked()) { return; }
312 fCParameters->SetPixe(val);
313}
314
316{
317 return fCParameters->Pixe();
318}
319
321{
322 if(IsLocked()) { return; }
323 fCParameters->SetDeexcitationIgnoreCut(val);
324}
325
327{
328 return fCParameters->DeexcitationIgnoreCut();
329}
330
332{
333 if(IsLocked()) { return; }
334 lateralDisplacement = val;
335}
336
338{
339 return lateralDisplacement;
340}
341
343{
344 if(IsLocked()) { return; }
345 lateralDisplacementAlg96 = val;
346}
347
349{
350 return lateralDisplacementAlg96;
351}
352
354{
355 if(IsLocked()) { return; }
356 muhadLateralDisplacement = val;
357}
358
360{
361 return muhadLateralDisplacement;
362}
363
365{
366 if(IsLocked()) { return; }
367 useAngGeneratorForIonisation = val;
368}
369
371{
372 return useAngGeneratorForIonisation;
373}
374
376{
377 if(IsLocked()) { return; }
378 useMottCorrection = val;
379}
380
382{
383 return useMottCorrection;
384}
385
387{
388 if(IsLocked()) { return; }
389 integral = val;
390}
391
393{
394 return integral;
395}
396
398{
399 if(IsLocked()) { return; }
400 fPolarisation = val;
401}
402
404{
405 return fPolarisation;
406}
407
409{
410 if(IsLocked()) { return; }
411 birks = val;
412 if(birks && nullptr == emSaturation) { emSaturation = new G4EmSaturation(1); }
413}
414
416{
417 return birks;
418}
419
421{
422 if(IsLocked()) { return; }
423 fICRU90 = val;
424}
425
427{
428 return fICRU90;
429}
430
432{
433 if(IsLocked()) { return; }
434 fCParameters->SetDNAFast(val);
435 if(val) { ActivateDNA(); }
436}
437
439{
440 return fCParameters->DNAFast();
441}
442
444{
445 if(IsLocked()) { return; }
446 fCParameters->SetDNAStationary(val);
447 if(val) { ActivateDNA(); }
448}
449
451{
452 return fCParameters->DNAStationary();
453}
454
456{
457 if(IsLocked()) { return; }
458 fCParameters->SetDNAElectronMsc(val);
459 if(val) { ActivateDNA(); }
460}
461
463{
464 return fCParameters->DNAElectronMsc();
465}
466
468{
469 if(IsLocked()) { return; }
470 gener = val;
471}
472
474{
475 return gener;
476}
477
479{
480 if(IsLocked()) { return; }
481 birks = (nullptr != ptr);
482 if(emSaturation != ptr) {
483 delete emSaturation;
484 emSaturation = ptr;
485 }
486}
487
489{
490 return fMuDataFromFile;
491}
492
494{
495 fMuDataFromFile = v;
496}
497
499{
500 if(IsLocked()) { return; }
501 onIsolated = val;
502}
503
505{
506 return onIsolated;
507}
508
510{
511 if(IsLocked()) { return; }
512 fSamplingTable = val;
513}
514
516{
517 return fSamplingTable;
518}
519
521{
522 return fPEKShell;
523}
524
526{
527 if(IsLocked()) { return; }
528 fPEKShell = v;
529}
530
532{
533 return fMscPosiCorr;
534}
535
537{
538 if(IsLocked()) { return; }
539 fMscPosiCorr = v;
540}
541
543{
544 return fUseEPICS2017XS;
545}
546
548{
549 if(IsLocked()) { return; }
550 fUseEPICS2017XS = v;
551}
552
554{
555 return f3GammaAnnihilationOnFly;
556}
557
559{
560 if(IsLocked()) { return; }
561 f3GammaAnnihilationOnFly = v;
562}
563
565{
566 return fUseRiGePairProductionModel;
567}
568
570{
571 if (IsLocked()) { return; }
572 fUseRiGePairProductionModel = v;
573}
574
576{
577 if(IsLocked()) { return; }
578 fDNA = true;
579}
580
582{
583 fIsPrinted = val;
584}
585
587{
588 return fIsPrinted;
589}
590
592{
593 if(nullptr == emSaturation) {
594#ifdef G4MULTITHREADED
595 G4MUTEXLOCK(&emParametersMutex);
596 if(nullptr == emSaturation) {
597#endif
598 emSaturation = new G4EmSaturation(1);
599#ifdef G4MULTITHREADED
600 }
601 G4MUTEXUNLOCK(&emParametersMutex);
602#endif
603 }
604 birks = true;
605 return emSaturation;
606}
607
609{
610 if(IsLocked()) { return; }
611 if(val > 1.e-3*CLHEP::eV && val < maxKinEnergy) {
612 minKinEnergy = val;
613 } else {
615 ed << "Value of MinKinEnergy - is out of range: " << val/CLHEP::MeV
616 << " MeV is ignored";
617 PrintWarning(ed);
618 }
619}
620
622{
623 return minKinEnergy;
624}
625
627{
628 if(IsLocked()) { return; }
629 if(val > std::max(minKinEnergy,599.9*CLHEP::MeV) && val < 1.e+7*CLHEP::TeV) {
630 maxKinEnergy = val;
631 } else {
633 ed << "Value of MaxKinEnergy is out of range: "
634 << val/CLHEP::GeV
635 << " GeV is ignored; allowed range 600 MeV - 1.e+7 TeV";
636 PrintWarning(ed);
637 }
638}
639
641{
642 return maxKinEnergy;
643}
644
646{
647 if(IsLocked()) { return; }
648 if(val > minKinEnergy && val <= 100*CLHEP::TeV) {
649 maxKinEnergyCSDA = val;
650 } else {
652 ed << "Value of MaxKinEnergyCSDA is out of range: "
653 << val/CLHEP::GeV << " GeV is ignored; allowed range "
654 << minKinEnergy << " MeV - 100 TeV";
655 PrintWarning(ed);
656 }
657}
658
660{
661 return maxKinEnergyCSDA;
662}
663
665{
666 if(IsLocked()) { return; }
667 if(val >= 0.0) { lowestElectronEnergy = val; }
668}
669
671{
672 return lowestElectronEnergy;
673}
674
676{
677 if(IsLocked()) { return; }
678 if(val >= 0.0) { lowestMuHadEnergy = val; }
679}
680
682{
683 return lowestMuHadEnergy;
684}
685
687{
688 if(IsLocked()) { return; }
689 if(val > 0.0) { lowestTripletEnergy = val; }
690}
691
693{
694 return lowestTripletEnergy;
695}
696
698{
699 if(IsLocked()) { return; }
700 if(val >= 0.0) { maxNIELEnergy = val; }
701}
702
704{
705 return maxNIELEnergy;
706}
707
709{
710 if(IsLocked()) { return; }
711 if(val > 0.0) { max5DEnergyForMuPair = val; }
712}
713
715{
716 return max5DEnergyForMuPair;
717}
718
720{
721 if(IsLocked()) { return; }
722 if(val > 0.0 && val < 0.5) {
723 linLossLimit = val;
724 } else {
726 ed << "Value of linLossLimit is out of range: " << val
727 << " is ignored";
728 PrintWarning(ed);
729 }
730}
731
733{
734 return linLossLimit;
735}
736
738{
739 if(IsLocked()) { return; }
740 if(val > 0.0) {
741 bremsTh = val;
742 } else {
744 ed << "Value of bremsstrahlung threshold is out of range: "
745 << val/GeV << " GeV is ignored";
746 PrintWarning(ed);
747 }
748}
749
751{
752 return bremsTh;
753}
754
756{
757 if(IsLocked()) { return; }
758 if(val > 0.0) {
759 bremsMuHadTh = val;
760 } else {
762 ed << "Value of bremsstrahlung threshold is out of range: "
763 << val/GeV << " GeV is ignored";
764 PrintWarning(ed);
765 }
766}
767
769{
770 return bremsMuHadTh;
771}
772
774{
775 if(IsLocked()) { return; }
776 if(val > 0.0 && val < 1.0) {
777 lambdaFactor = val;
778 } else {
780 ed << "Value of lambda factor is out of range: " << val
781 << " is ignored";
782 PrintWarning(ed);
783 }
784}
785
787{
788 return lambdaFactor;
789}
790
792{
793 if(IsLocked()) { return; }
794 if(val > 0.0) {
795 factorForAngleLimit = val;
796 } else {
798 ed << "Value of factor for enegry limit is out of range: "
799 << val << " is ignored";
800 PrintWarning(ed);
801 }
802}
803
805{
806 return factorForAngleLimit;
807}
808
810{
811 if(IsLocked()) { return; }
812 if(val >= 0.0 && val <= pi) {
813 thetaLimit = val;
814 } else {
816 ed << "Value of polar angle limit is out of range: "
817 << val << " is ignored";
818 PrintWarning(ed);
819 }
820}
821
823{
824 return thetaLimit;
825}
826
828{
829 if(IsLocked()) { return; }
830 if(val >= 0.0) {
831 energyLimit = val;
832 } else {
834 ed << "Value of msc energy limit is out of range: "
835 << val << " is ignored";
836 PrintWarning(ed);
837 }
838}
839
841{
842 return energyLimit;
843}
844
846{
847 if(IsLocked()) { return; }
848 if(val > 0.0 && val < 1.0) {
849 rangeFactor = val;
850 } else {
852 ed << "Value of rangeFactor is out of range: "
853 << val << " is ignored";
854 PrintWarning(ed);
855 }
856}
857
859{
860 return rangeFactor;
861}
862
864{
865 if(IsLocked()) { return; }
866 if(val > 0.0 && val < 1.0) {
867 rangeFactorMuHad = val;
868 } else {
870 ed << "Value of rangeFactorMuHad is out of range: "
871 << val << " is ignored";
872 PrintWarning(ed);
873 }
874}
875
877{
878 return rangeFactorMuHad;
879}
880
882{
883 if(IsLocked()) { return; }
884 if(val >= 1.0) {
885 geomFactor = val;
886 } else {
888 ed << "Value of geomFactor is out of range: "
889 << val << " is ignored";
890 PrintWarning(ed);
891 }
892}
893
895{
896 return geomFactor;
897}
898
900{
901 if(IsLocked()) { return; }
902 if(val >= 0.1) {
903 safetyFactor = val;
904 } else {
906 ed << "Value of safetyFactor is out of range: "
907 << val << " is ignored";
908 PrintWarning(ed);
909 }
910}
911
913{
914 return safetyFactor;
915}
916
918{
919 if(IsLocked()) { return; }
920 if(val >= 0.0) {
921 lambdaLimit = val;
922 } else {
924 ed << "Value of lambdaLimit is out of range: "
925 << val << " is ignored";
926 PrintWarning(ed);
927 }
928}
929
931{
932 return lambdaLimit;
933}
934
936{
937 if(IsLocked()) { return; }
938 if(val >= 1.0) {
939 skin = val;
940 } else {
942 ed << "Value of skin is out of range: "
943 << val << " is ignored";
944 PrintWarning(ed);
945 }
946}
947
949{
950 return skin;
951}
952
954{
955 if(IsLocked()) { return; }
956 if(val > 0.0) {
957 factorScreen = val;
958 } else {
960 ed << "Value of factorScreen is out of range: "
961 << val << " is ignored";
962 PrintWarning(ed);
963 }
964}
965
967{
968 return factorScreen;
969}
970
972{
973 if(IsLocked()) { return; }
974 fBParameters->SetStepFunction(v1, v2);
975}
976
978{
979 if(IsLocked()) { return; }
980 fBParameters->SetStepFunctionMuHad(v1, v2);
981}
982
984{
985 if(IsLocked()) { return; }
986 fBParameters->SetStepFunctionLightIons(v1, v2);
987}
988
990{
991 if(IsLocked()) { return; }
992 fBParameters->SetStepFunctionIons(v1, v2);
993}
994
996{
997 fBParameters->FillStepFunction(part, proc);
998}
999
1001{
1002 return nbinsPerDecade*G4lrint(std::log10(maxKinEnergy/minKinEnergy));
1003}
1004
1006{
1007 if(IsLocked()) { return; }
1008 if(val >= 5 && val < 1000000) {
1009 nbinsPerDecade = val;
1010 } else {
1012 ed << "Value of number of bins per decade is out of range: "
1013 << val << " is ignored";
1014 PrintWarning(ed);
1015 }
1016}
1017
1019{
1020 return nbinsPerDecade;
1021}
1022
1024{
1025 if(IsLocked()) { return; }
1026 verbose = val;
1027 workerVerbose = std::min(workerVerbose, verbose);
1028}
1029
1031{
1032 return verbose;
1033}
1034
1036{
1037 if(IsLocked()) { return; }
1038 workerVerbose = val;
1039}
1040
1042{
1043 return workerVerbose;
1044}
1045
1047{
1048 if(IsLocked()) { return; }
1049 nForFreeVector = val;
1050}
1051
1053{
1054 return nForFreeVector;
1055}
1056
1058{
1059 if(IsLocked()) { return; }
1060 fTransportationWithMsc = val;
1061}
1062
1064{
1065 return fTransportationWithMsc;
1066}
1067
1069{
1070 if(IsLocked()) { return; }
1071 fFluct = val;
1072}
1073
1075{
1076 return fFluct;
1077}
1078
1080{
1081 if(IsLocked()) { return; }
1082 fPositronium = val;
1083}
1084
1086{
1087 return fPositronium;
1088}
1089
1091{
1092 if(IsLocked()) { return; }
1093 mscStepLimit = val;
1094}
1095
1097{
1098 return mscStepLimit;
1099}
1100
1102{
1103 if(IsLocked()) { return; }
1104 mscStepLimitMuHad = val;
1105}
1106
1108{
1109 return mscStepLimitMuHad;
1110}
1111
1113{
1114 if(IsLocked()) { return; }
1115 fSStype = val;
1116}
1117
1119{
1120 return fSStype;
1121}
1122
1123void
1125{
1126 if(IsLocked()) { return; }
1127 nucFormfactor = val;
1128}
1129
1131{
1132 return nucFormfactor;
1133}
1134
1136{
1137 if(IsLocked()) { return; }
1138 fCParameters->SetDNAeSolvationSubType(val);
1139 ActivateDNA();
1140}
1141
1143{
1144 return fCParameters->DNAeSolvationSubType();
1145}
1146
1148{
1149 if(IsLocked()) { return; }
1150 tripletConv = val;
1151}
1152
1154{
1155 return tripletConv;
1156}
1157
1159{
1160 if(IsLocked()) { return; }
1161 fCParameters->SetPIXECrossSectionModel(sss);
1162}
1163
1165{
1166 return fCParameters->PIXECrossSectionModel();
1167}
1168
1170{
1171 if(IsLocked()) { return; }
1172 fCParameters->SetPIXEElectronCrossSectionModel(sss);
1173}
1174
1176{
1177 return fCParameters->PIXEElectronCrossSectionModel();
1178}
1179
1181{
1182 if(IsLocked()) { return; }
1183 fCParameters->SetLivermoreDataDir(sss);
1184}
1185
1187{
1188 return fCParameters->LivermoreDataDir();
1189}
1190
1191void G4EmParameters::PrintWarning(G4ExceptionDescription& ed) const
1192{
1193 G4Exception("G4EmParameters", "em0044", JustWarning, ed);
1194}
1195
1197 const G4String& region,
1198 const G4String& type)
1199{
1200 if(IsLocked()) { return; }
1201 fBParameters->AddPAIModel(particle, region, type);
1202}
1203
1204const std::vector<G4String>& G4EmParameters::ParticlesPAI() const
1205{
1206 return fBParameters->ParticlesPAI();
1207}
1208
1209const std::vector<G4String>& G4EmParameters::RegionsPAI() const
1210{
1211 return fBParameters->RegionsPAI();
1212}
1213
1214const std::vector<G4String>& G4EmParameters::TypesPAI() const
1215{
1216 return fBParameters->TypesPAI();
1217}
1218
1220{
1221 if(IsLocked()) { return; }
1222 fCParameters->AddMicroElec(region);
1223}
1224
1225const std::vector<G4String>& G4EmParameters::RegionsMicroElec() const
1226{
1227 return fCParameters->RegionsMicroElec();
1228}
1229
1230void G4EmParameters::AddDNA(const G4String& region, const G4String& type)
1231{
1232 if(IsLocked()) { return; }
1233 fCParameters->AddDNA(region, type);
1234 ActivateDNA();
1235}
1236
1237const std::vector<G4String>& G4EmParameters::RegionsDNA() const
1238{
1239 return fCParameters->RegionsDNA();
1240}
1241
1242const std::vector<G4String>& G4EmParameters::TypesDNA() const
1243{
1244 return fCParameters->TypesDNA();
1245}
1246
1247void G4EmParameters::AddPhysics(const G4String& region, const G4String& type)
1248{
1249 if(IsLocked()) { return; }
1250 fBParameters->AddPhysics(region, type);
1251}
1252
1253const std::vector<G4String>& G4EmParameters::RegionsPhysics() const
1254{
1255 return fBParameters->RegionsPhysics();
1256}
1257
1258const std::vector<G4String>& G4EmParameters::TypesPhysics() const
1259{
1260 return fBParameters->TypesPhysics();
1261}
1262
1264{
1265 if(IsLocked()) { return; }
1266 fBParameters->SetSubCutRegion(region);
1267}
1268
1269void
1271 G4bool aauger, G4bool apixe)
1272{
1273 if(IsLocked()) { return; }
1274 fCParameters->SetDeexActiveRegion(region, adeex, aauger, apixe);
1275}
1276
1277void
1279 G4double val, G4bool wflag)
1280{
1281 if(IsLocked()) { return; }
1282 fBParameters->SetProcessBiasingFactor(procname, val, wflag);
1283}
1284
1285void
1287 const G4String& region,
1288 G4double length,
1289 G4bool wflag)
1290{
1291 if(IsLocked() && !gener) { return; }
1292 fBParameters->ActivateForcedInteraction(procname, region, length, wflag);
1293}
1294
1295void
1297 const G4String& region,
1298 G4double factor,
1299 G4double energyLim)
1300{
1301 if(IsLocked()) { return; }
1302 fBParameters->ActivateSecondaryBiasing(procname, region, factor, energyLim);
1303}
1304
1306{
1307 fBParameters->DefineRegParamForLoss(ptr);
1308}
1309
1311{
1312 fBParameters->DefineRegParamForEM(ptr);
1313}
1314
1316{
1317 return fBParameters->QuantumEntanglement();
1318}
1319
1321{
1322 if(IsLocked()) { return; }
1323 fBParameters->SetQuantumEntanglement(v);
1324}
1325
1327 return fBParameters->GetDirectionalSplitting();
1328}
1329
1331{
1332 if(IsLocked()) { return; }
1333 fBParameters->SetDirectionalSplitting(v);
1334}
1335
1337{
1338 if(IsLocked()) { return; }
1339 fBParameters->SetDirectionalSplittingTarget(v);
1340}
1341
1343{
1344 return fBParameters->GetDirectionalSplittingTarget();
1345}
1346
1348{
1349 if(IsLocked()) { return; }
1350 fBParameters->SetDirectionalSplittingRadius(r);
1351}
1352
1354{
1355 return fBParameters->GetDirectionalSplittingRadius();
1356}
1357
1359{
1360 fCParameters->DefineRegParamForDeex(ptr);
1361}
1362
1364{
1365 return fDirLEDATA;
1366}
1367
1368void G4EmParameters::StreamInfo(std::ostream& os) const
1369{
1370 G4long prec = os.precision(5);
1371 os << "=======================================================================" << "\n";
1372 os << "====== Electromagnetic Physics Parameters ========" << "\n";
1373 os << "=======================================================================" << "\n";
1374 os << "LPM effect enabled " <<flagLPM << "\n";
1375 os << "Enable creation and use of sampling tables " <<fSamplingTable << "\n";
1376 os << "Apply cuts on all EM processes " <<applyCuts << "\n";
1377 const char* transportationWithMsc = "Disabled";
1378 if(fTransportationWithMsc == G4TransportationWithMscType::fEnabled) {
1379 transportationWithMsc = "Enabled";
1380 } else if (fTransportationWithMsc == G4TransportationWithMscType::fMultipleSteps) {
1381 transportationWithMsc = "MultipleSteps";
1382 }
1383 os << "Use combined TransportationWithMsc " <<transportationWithMsc << "\n";
1384 os << "Use general process " <<gener << "\n";
1385 os << "Enable linear polarisation for gamma " <<fPolarisation << "\n";
1386 os << "Enable photoeffect sampling below K-shell " <<fPEKShell << "\n";
1387 os << "Enable sampling of quantum entanglement "
1388 <<fBParameters->QuantumEntanglement() << "\n";
1389 os << "X-section factor for integral approach " <<lambdaFactor << "\n";
1390 os << "Min kinetic energy for tables "
1391 <<G4BestUnit(minKinEnergy,"Energy") << "\n";
1392 os << "Max kinetic energy for tables "
1393 <<G4BestUnit(maxKinEnergy,"Energy") << "\n";
1394 os << "Number of bins per decade of a table " <<nbinsPerDecade << "\n";
1395 os << "Verbose level " <<verbose << "\n";
1396 os << "Verbose level for worker thread " <<workerVerbose << "\n";
1397 os << "Bremsstrahlung energy threshold above which \n"
1398 << " primary e+- is added to the list of secondary "
1399 <<G4BestUnit(bremsTh,"Energy") << "\n";
1400 os << "Bremsstrahlung energy threshold above which primary\n"
1401 << " muon/hadron is added to the list of secondary "
1402 <<G4BestUnit(bremsMuHadTh,"Energy") << "\n";
1403 G4String name3g = "SimplePositronium";
1404 if (fPositronium == fAllisonPositronium) { name3g = "AllisonPositronium"; }
1405 else if (fPositronium == fOrePowell) { name3g = "OrePowell"; }
1406 else if (fPositronium == fOrePowellPolar) { name3g = "OrePowellPolar"; }
1407 os << "Positron annihilation at rest model " << name3g << "\n";
1408
1409 os << "Enable 3 gamma annihilation on fly "
1410 << f3GammaAnnihilationOnFly << "\n";
1411 os << "Lowest triplet kinetic energy "
1412 <<G4BestUnit(lowestTripletEnergy,"Energy") << "\n";
1413 os << "Enable sampling of gamma linear polarisation " <<fPolarisation << "\n";
1414 os << "5D gamma conversion model type " <<tripletConv << "\n";
1415 os << "5D gamma conversion model on isolated ion " <<onIsolated << "\n";
1416 if(max5DEnergyForMuPair>0.0) {
1417 os << "5D gamma conversion limit for muon pair "
1418 << max5DEnergyForMuPair/CLHEP::GeV << " GeV\n";
1419 }
1420 os << "Use Ricardo-Gerardo pair production model "
1421 << fUseRiGePairProductionModel << "\n";
1422 os << "Livermore data directory "
1423 << fCParameters->LivermoreDataDir() << "\n";
1424
1425 os << "=======================================================================" << "\n";
1426 os << "====== Ionisation Parameters ========" << "\n";
1427 os << "=======================================================================" << "\n";
1428 os << "Step function for e+- "
1429 <<"("<<fBParameters->GetStepFunctionP1() << ", "
1430 << fBParameters->GetStepFunctionP2()/CLHEP::mm << " mm)\n";
1431 os << "Step function for muons/hadrons "
1432 <<"("<<fBParameters->GetStepFunctionMuHadP1() << ", "
1433 << fBParameters->GetStepFunctionMuHadP2()/CLHEP::mm << " mm)\n";
1434 os << "Step function for light ions "
1435 <<"("<<fBParameters->GetStepFunctionLightIonsP1() << ", "
1436 << fBParameters->GetStepFunctionLightIonsP2()/CLHEP::mm << " mm)\n";
1437 os << "Step function for general ions "
1438 <<"("<<fBParameters->GetStepFunctionIonsP1() << ", "
1439 << fBParameters->GetStepFunctionIonsP2()/CLHEP::mm << " mm)\n";
1440 os << "Lowest e+e- kinetic energy "
1441 <<G4BestUnit(lowestElectronEnergy,"Energy") << "\n";
1442 os << "Lowest muon/hadron kinetic energy "
1443 <<G4BestUnit(lowestMuHadEnergy,"Energy") << "\n";
1444 os << "Use ICRU90 data " << fICRU90 << "\n";
1445 os << "Fluctuations of dE/dx are enabled " <<lossFluctuation << "\n";
1446 G4String namef = "Universal";
1447 if(fFluct == fUrbanFluctuation) { namef = "Urban"; }
1448 else if(fFluct == fDummyFluctuation) { namef = "Dummy"; }
1449 os << "Type of fluctuation model for leptons and hadrons " << namef << "\n";
1450 os << "Use built-in Birks satuaration " << birks << "\n";
1451 os << "Build CSDA range enabled " <<buildCSDARange << "\n";
1452 os << "Use cut as a final range enabled " <<cutAsFinalRange << "\n";
1453 os << "Enable angular generator interface "
1454 <<useAngGeneratorForIonisation << "\n";
1455 os << "Max kinetic energy for CSDA tables "
1456 <<G4BestUnit(maxKinEnergyCSDA,"Energy") << "\n";
1457 os << "Max kinetic energy for NIEL computation "
1458 <<G4BestUnit(maxNIELEnergy,"Energy") << "\n";
1459 os << "Linear loss limit " <<linLossLimit << "\n";
1460 os << "Read data from file for e+e- pair production by mu " <<fMuDataFromFile << "\n";
1461
1462 os << "=======================================================================" << "\n";
1463 os << "====== Multiple Scattering Parameters ========" << "\n";
1464 os << "=======================================================================" << "\n";
1465 os << "Type of msc step limit algorithm for e+- " <<mscStepLimit << "\n";
1466 os << "Type of msc step limit algorithm for muons/hadrons " <<mscStepLimitMuHad << "\n";
1467 os << "Msc lateral displacement for e+- enabled " <<lateralDisplacement << "\n";
1468 os << "Msc lateral displacement for muons and hadrons " <<muhadLateralDisplacement << "\n";
1469 os << "Urban msc model lateral displacement alg96 " <<lateralDisplacementAlg96 << "\n";
1470 os << "Range factor for msc step limit for e+- " <<rangeFactor << "\n";
1471 os << "Range factor for msc step limit for muons/hadrons " <<rangeFactorMuHad << "\n";
1472 os << "Geometry factor for msc step limitation of e+- " <<geomFactor << "\n";
1473 os << "Safety factor for msc step limit for e+- " <<safetyFactor << "\n";
1474 os << "Skin parameter for msc step limitation of e+- " <<skin << "\n";
1475 os << "Lambda limit for msc step limit for e+- " <<lambdaLimit/CLHEP::mm << " mm\n";
1476 os << "Use Mott correction for e- scattering " << useMottCorrection << "\n";
1477 os << "Factor used for dynamic computation of angular \n"
1478 << " limit between single and multiple scattering " << factorForAngleLimit << "\n";
1479 os << "Fixed angular limit between single \n"
1480 << " and multiple scattering "
1481 << thetaLimit/CLHEP::rad << " rad\n";
1482 os << "Upper energy limit for e+- multiple scattering "
1483 << energyLimit/CLHEP::MeV << " MeV\n";
1484 os << "Type of electron single scattering model " <<fSStype << "\n";
1485 os << "Type of nuclear form-factor " <<nucFormfactor << "\n";
1486 os << "Screening factor " <<factorScreen << "\n";
1487 os << "=======================================================================" << "\n";
1488
1489 if(fCParameters->Fluo()) {
1490 os << "====== Atomic Deexcitation Parameters ========" << "\n";
1491 os << "=======================================================================" << "\n";
1492 os << "Fluorescence enabled " <<fCParameters->Fluo() << "\n";
1493 G4String named = "fluor";
1495 if(fdir == fluoBearden) { named = "fluor_Bearden"; }
1496 else if(fdir == fluoANSTO) { named = "fluor_ANSTO"; }
1497 else if(fdir == fluoXDB_EADL) { named = "fluor_XDB_EADL"; }
1498 os << "Directory in G4LEDATA for fluorescence data files " << named << "\n";
1499 os << "Auger electron cascade enabled "
1500 <<fCParameters->Auger() << "\n";
1501 os << "PIXE atomic de-excitation enabled " <<fCParameters->Pixe() << "\n";
1502 os << "De-excitation module ignores cuts "
1503 <<fCParameters->DeexcitationIgnoreCut() << "\n";
1504 os << "Type of PIXE cross section for hadrons "
1505 <<fCParameters->PIXECrossSectionModel() << "\n";
1506 os << "Type of PIXE cross section for e+- "
1507 <<fCParameters->PIXEElectronCrossSectionModel() << "\n";
1508 os << "=======================================================================" << "\n";
1509 }
1510 if(fDNA) {
1511 os << "====== DNA Physics Parameters ========" << "\n";
1512 os << "=======================================================================" << "\n";
1513 os << "Use fast sampling in DNA models "
1514 << fCParameters->DNAFast() << "\n";
1515 os << "Use Stationary option in DNA models "
1516 << fCParameters->DNAStationary() << "\n";
1517 os << "Use DNA with multiple scattering of e- "
1518 << fCParameters->DNAElectronMsc() << "\n";
1519 os << "Use DNA e- solvation model type "
1520 << fCParameters->DNAeSolvationSubType() << "\n";
1521 auto chemModel = fCParameters->GetChemTimeStepModel();
1522 if(fCParameters->GetChemTimeStepModel() != G4ChemTimeStepModel::Unknown)
1523 {
1524 std::vector<G4String> ChemModel{"Unknown","SBS","IRT","IRT_syn"};
1525 os << "Use DNA Chemistry model "
1526 << ChemModel.at((std::size_t)chemModel) << "\n";
1527 }
1528 os << "=======================================================================" << G4endl;
1529 }
1530 os.precision(prec);
1531}
1532
1534{
1535 if(fIsPrinted) return;
1536
1537#ifdef G4MULTITHREADED
1538 G4MUTEXLOCK(&emParametersMutex);
1539#endif
1541#ifdef G4MULTITHREADED
1542 G4MUTEXUNLOCK(&emParametersMutex);
1543#endif
1544}
1545
1546std::ostream& operator<< (std::ostream& os, const G4EmParameters& par)
1547{
1548 par.StreamInfo(os);
1549 return os;
1550}
1551
1552G4bool G4EmParameters::IsLocked() const
1553{
1554 return (!G4Threading::IsMasterThread() ||
1555 (fStateManager->GetCurrentState() != G4State_PreInit &&
1556 fStateManager->GetCurrentState() != G4State_Init &&
1557 fStateManager->GetCurrentState() != G4State_Idle));
1558}
1559
1560
1562{
1563 fCParameters-> SetChemTimeStepModel(model);
1564}
1565
1567{
1568 return fCParameters->GetChemTimeStepModel();
1569}
1570//....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 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)
void SetUseCutAsFinalRange(G4bool val)
void SetDeexcitationIgnoreCut(G4bool val)
void SetBirksActive(G4bool val)
G4double ScreeningFactor() const
G4bool UseMottCorrection() const
void SetFluo(G4bool val)
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
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)
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)
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 G4NistManager * Instance()
const G4ApplicationState & GetCurrentState() const
static G4StateManager * GetStateManager()
G4bool IsMasterThread()
int G4lrint(double ad)
Definition templates.hh:134