Geant4-11
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
60
61#ifdef G4MULTITHREADED
62 G4Mutex G4EmParameters::emParametersMutex = G4MUTEX_INITIALIZER;
63#endif
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
66
68{
69 if(nullptr == theInstance) {
70#ifdef G4MULTITHREADED
71 G4MUTEXLOCK(&emParametersMutex);
72 if(nullptr == theInstance) {
73#endif
74 static G4EmParameters manager;
75 theInstance = &manager;
76#ifdef G4MULTITHREADED
77 }
78 G4MUTEXUNLOCK(&emParametersMutex);
79#endif
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{
100 Initialise();
101
104
106 emSaturation = nullptr;
107}
108
110{
111 if(!IsLocked()) {
112 Initialise();
115 }
116}
117
119{
120 lossFluctuation = true;
121 buildCSDARange = false;
122 flagLPM = true;
123 cutAsFinalRange = false;
124 applyCuts = false;
125 lateralDisplacement = true;
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 fDNA = false;
139 fIsPrinted = false;
140
142 maxKinEnergy = 100.0*CLHEP::TeV;
148 maxNIELEnergy = 0.0;
149 linLossLimit = 0.01;
151 lambdaFactor = 0.8;
154 energyLimit = 100.0*CLHEP::MeV;
155 rangeFactor = 0.04;
156 rangeFactorMuHad = 0.2;
157 geomFactor = 2.5;
158 skin = 1.0;
159 safetyFactor = 0.6;
161 factorScreen = 1.0;
162
163 nbinsPerDecade = 7;
164 verbose = 1;
165 workerVerbose = 0;
166 tripletConv = 0;
167
171 fSStype = fWVI;
172}
173
175{
176 if(IsLocked()) { return; }
177 lossFluctuation = val;
178}
179
181{
182 return lossFluctuation;
183}
184
186{
187 if(IsLocked()) { return; }
188 buildCSDARange = val;
189}
190
192{
193 return buildCSDARange;
194}
195
197{
198 if(IsLocked()) { return; }
199 flagLPM = val;
200}
201
203{
204 return flagLPM;
205}
206
208{
209 if(IsLocked()) { return; }
210 cutAsFinalRange = val;
211}
212
214{
215 return cutAsFinalRange;
216}
217
219{
220 if(IsLocked()) { return; }
221 applyCuts = val;
222}
223
225{
226 return applyCuts;
227}
228
230{
231 if(IsLocked()) { return; }
232 fCParameters->SetFluo(val);
233}
234
236{
237 return fCParameters->Fluo();
238}
239
241{
242 if(IsLocked()) { return; }
244}
245
247{
249}
250
252{
253 if(IsLocked()) { return; }
255}
256
258{
259 return fCParameters->ANSTOFluoDir();
260}
261
263{
264 if(IsLocked()) { return; }
266}
267
269{
270 return fCParameters->Auger();
271}
272
274{
275 if(IsLocked()) { return; }
276 fCParameters->SetPixe(val);
277}
278
280{
281 return fCParameters->Pixe();
282}
283
285{
286 if(IsLocked()) { return; }
288}
289
291{
293}
294
296{
297 if(IsLocked()) { return; }
299}
300
302{
303 return lateralDisplacement;
304}
305
307{
308 if(IsLocked()) { return; }
310}
311
313{
315}
316
318{
319 if(IsLocked()) { return; }
321}
322
324{
326}
327
329{
330 if(IsLocked()) { return; }
332}
333
335{
337}
338
340{
341 if(IsLocked()) { return; }
342 useMottCorrection = val;
343}
344
346{
347 return useMottCorrection;
348}
349
351{
352 if(IsLocked()) { return; }
353 integral = val;
354}
355
357{
358 return integral;
359}
360
362{
363 if(IsLocked()) { return; }
364 fPolarisation = val;
365}
366
368{
369 return fPolarisation;
370}
371
373{
374 if(IsLocked()) { return; }
375 birks = val;
376 if(birks && nullptr == emSaturation) { emSaturation = new G4EmSaturation(1); }
377}
378
380{
381 return birks;
382}
383
385{
386 if(IsLocked()) { return; }
387 fICRU90 = val;
388}
389
391{
392 return fICRU90;
393}
394
396{
397 if(IsLocked()) { return; }
399 if(val) { ActivateDNA(); }
400}
401
403{
404 return fCParameters->DNAFast();
405}
406
408{
409 if(IsLocked()) { return; }
411 if(val) { ActivateDNA(); }
412}
413
415{
416 return fCParameters->DNAStationary();
417}
418
420{
421 if(IsLocked()) { return; }
423 if(val) { ActivateDNA(); }
424}
425
427{
429}
430
432{
433 if(IsLocked()) { return; }
434 gener = val;
435}
436
438{
439 return gener;
440}
441
443{
444 if(IsLocked()) { return; }
445 birks = (nullptr != ptr);
446 if(emSaturation != ptr) {
447 delete emSaturation;
448 emSaturation = ptr;
449 }
450}
451
453{
454 return fMuDataFromFile;
455}
456
458{
459 fMuDataFromFile = v;
460}
461
463{
464 if(IsLocked()) { return; }
465 onIsolated = val;
466}
467
469{
470 return onIsolated;
471}
472
474{
475 if(IsLocked()) { return; }
476 fSamplingTable = val;
477}
478
480{
481 return fSamplingTable;
482}
483
485{
486 if(IsLocked()) { return; }
487 fDNA = true;
488}
489
491{
492 fIsPrinted = val;
493}
494
496{
497 return fIsPrinted;
498}
499
501{
502 if(nullptr == emSaturation) {
503#ifdef G4MULTITHREADED
504 G4MUTEXLOCK(&emParametersMutex);
505 if(nullptr == emSaturation) {
506#endif
508#ifdef G4MULTITHREADED
509 }
510 G4MUTEXUNLOCK(&emParametersMutex);
511#endif
512 }
513 birks = true;
514 return emSaturation;
515}
516
518{
519 if(IsLocked()) { return; }
520 if(val > 1.e-3*CLHEP::eV && val < maxKinEnergy) {
521 minKinEnergy = val;
522 } else {
524 ed << "Value of MinKinEnergy - is out of range: " << val/CLHEP::MeV
525 << " MeV is ignored";
526 PrintWarning(ed);
527 }
528}
529
531{
532 return minKinEnergy;
533}
534
536{
537 if(IsLocked()) { return; }
538 if(val > std::max(minKinEnergy,9.99*CLHEP::MeV) && val < 1.e+7*CLHEP::TeV) {
539 maxKinEnergy = val;
540 } else {
542 ed << "Value of MaxKinEnergy is out of range: "
543 << val/CLHEP::GeV
544 << " GeV is ignored; allowed range 10 MeV - 1.e+7 TeV";
545 PrintWarning(ed);
546 }
547}
548
550{
551 return maxKinEnergy;
552}
553
555{
556 if(IsLocked()) { return; }
557 if(val > minKinEnergy && val <= 100*CLHEP::TeV) {
558 maxKinEnergyCSDA = val;
559 } else {
561 ed << "Value of MaxKinEnergyCSDA is out of range: "
562 << val/CLHEP::GeV << " GeV is ignored; allowed range "
563 << minKinEnergy << " MeV - 100 TeV";
564 PrintWarning(ed);
565 }
566}
567
569{
570 return maxKinEnergyCSDA;
571}
572
574{
575 if(IsLocked()) { return; }
576 if(val >= 0.0) { lowestElectronEnergy = val; }
577}
578
580{
581 return lowestElectronEnergy;
582}
583
585{
586 if(IsLocked()) { return; }
587 if(val >= 0.0) { lowestMuHadEnergy = val; }
588}
589
591{
592 return lowestMuHadEnergy;
593}
594
596{
597 if(IsLocked()) { return; }
598 if(val > 0.0) { lowestTripletEnergy = val; }
599}
600
602{
603 return lowestTripletEnergy;
604}
605
607{
608 if(IsLocked()) { return; }
609 if(val >= 0.0) { maxNIELEnergy = val; }
610}
611
613{
614 return maxNIELEnergy;
615}
616
618{
619 if(IsLocked()) { return; }
620 if(val > 0.0) { max5DEnergyForMuPair = val; }
621}
622
624{
626}
627
629{
630 if(IsLocked()) { return; }
631 if(val > 0.0 && val < 0.5) {
632 linLossLimit = val;
633 } else {
635 ed << "Value of linLossLimit is out of range: " << val
636 << " is ignored";
637 PrintWarning(ed);
638 }
639}
640
642{
643 return linLossLimit;
644}
645
647{
648 if(IsLocked()) { return; }
649 if(val > 0.0) {
650 bremsTh = val;
651 } else {
653 ed << "Value of bremsstrahlung threshold is out of range: "
654 << val/GeV << " GeV is ignored";
655 PrintWarning(ed);
656 }
657}
658
660{
661 return bremsTh;
662}
663
665{
666 if(IsLocked()) { return; }
667 if(val > 0.0) {
668 bremsMuHadTh = val;
669 } else {
671 ed << "Value of bremsstrahlung threshold is out of range: "
672 << val/GeV << " GeV is ignored";
673 PrintWarning(ed);
674 }
675}
676
678{
679 return bremsMuHadTh;
680}
681
683{
684 if(IsLocked()) { return; }
685 if(val > 0.0 && val < 1.0) {
686 lambdaFactor = val;
687 } else {
689 ed << "Value of lambda factor is out of range: " << val
690 << " is ignored";
691 PrintWarning(ed);
692 }
693}
694
696{
697 return lambdaFactor;
698}
699
701{
702 if(IsLocked()) { return; }
703 if(val > 0.0) {
705 } else {
707 ed << "Value of factor for enegry limit is out of range: "
708 << val << " is ignored";
709 PrintWarning(ed);
710 }
711}
712
714{
715 return factorForAngleLimit;
716}
717
719{
720 if(IsLocked()) { return; }
721 if(val >= 0.0 && val <= pi) {
722 thetaLimit = val;
723 } else {
725 ed << "Value of polar angle limit is out of range: "
726 << val << " is ignored";
727 PrintWarning(ed);
728 }
729}
730
732{
733 return thetaLimit;
734}
735
737{
738 if(IsLocked()) { return; }
739 if(val >= 0.0) {
740 energyLimit = val;
741 } else {
743 ed << "Value of msc energy limit is out of range: "
744 << val << " is ignored";
745 PrintWarning(ed);
746 }
747}
748
750{
751 return energyLimit;
752}
753
755{
756 if(IsLocked()) { return; }
757 if(val > 0.0 && val < 1.0) {
758 rangeFactor = val;
759 } else {
761 ed << "Value of rangeFactor is out of range: "
762 << val << " is ignored";
763 PrintWarning(ed);
764 }
765}
766
768{
769 return rangeFactor;
770}
771
773{
774 if(IsLocked()) { return; }
775 if(val > 0.0 && val < 1.0) {
776 rangeFactorMuHad = val;
777 } else {
779 ed << "Value of rangeFactorMuHad is out of range: "
780 << val << " is ignored";
781 PrintWarning(ed);
782 }
783}
784
786{
787 return rangeFactorMuHad;
788}
789
791{
792 if(IsLocked()) { return; }
793 if(val >= 1.0) {
794 geomFactor = val;
795 } else {
797 ed << "Value of geomFactor is out of range: "
798 << val << " is ignored";
799 PrintWarning(ed);
800 }
801}
802
804{
805 return geomFactor;
806}
807
809{
810 if(IsLocked()) { return; }
811 if(val >= 0.1) {
812 safetyFactor = val;
813 } else {
815 ed << "Value of safetyFactor is out of range: "
816 << val << " is ignored";
817 PrintWarning(ed);
818 }
819}
820
822{
823 return safetyFactor;
824}
825
827{
828 if(IsLocked()) { return; }
829 if(val >= 0.0) {
830 lambdaLimit = val;
831 } else {
833 ed << "Value of lambdaLimit is out of range: "
834 << val << " is ignored";
835 PrintWarning(ed);
836 }
837}
838
840{
841 return lambdaLimit;
842}
843
845{
846 if(IsLocked()) { return; }
847 if(val >= 1.0) {
848 skin = val;
849 } else {
851 ed << "Value of skin is out of range: "
852 << val << " is ignored";
853 PrintWarning(ed);
854 }
855}
856
858{
859 return skin;
860}
861
863{
864 if(IsLocked()) { return; }
865 if(val > 0.0) {
866 factorScreen = val;
867 } else {
869 ed << "Value of factorScreen is out of range: "
870 << val << " is ignored";
871 PrintWarning(ed);
872 }
873}
874
876{
877 return factorScreen;
878}
879
881{
882 if(IsLocked()) { return; }
884}
885
887{
888 if(IsLocked()) { return; }
890}
891
893{
894 if(IsLocked()) { return; }
896}
897
899{
900 if(IsLocked()) { return; }
902}
903
905{
906 fBParameters->FillStepFunction(part, proc);
907}
908
910{
912}
913
915{
916 if(IsLocked()) { return; }
917 if(val >= 5 && val < 1000000) {
918 nbinsPerDecade = val;
919 } else {
921 ed << "Value of number of bins per decade is out of range: "
922 << val << " is ignored";
923 PrintWarning(ed);
924 }
925}
926
928{
929 return nbinsPerDecade;
930}
931
933{
934 if(IsLocked()) { return; }
935 verbose = val;
937}
938
940{
941 return verbose;
942}
943
945{
946 if(IsLocked()) { return; }
947 workerVerbose = val;
948}
949
951{
952 return workerVerbose;
953}
954
956{
957 if(IsLocked()) { return; }
958 mscStepLimit = val;
959}
960
962{
963 return mscStepLimit;
964}
965
967{
968 if(IsLocked()) { return; }
969 mscStepLimitMuHad = val;
970}
971
973{
974 return mscStepLimitMuHad;
975}
976
978{
979 if(IsLocked()) { return; }
980 fSStype = val;
981}
982
984{
985 return fSStype;
986}
987
988void
990{
991 if(IsLocked()) { return; }
992 nucFormfactor = val;
993}
994
996{
997 return nucFormfactor;
998}
999
1001{
1002 if(IsLocked()) { return; }
1004 ActivateDNA();
1005}
1006
1008{
1010}
1011
1013{
1014 if(IsLocked()) { return; }
1015 tripletConv = val;
1016}
1017
1019{
1020 return tripletConv;
1021}
1022
1024{
1025 if(IsLocked()) { return; }
1027}
1028
1030{
1032}
1033
1035{
1036 if(IsLocked()) { return; }
1038}
1039
1041{
1043}
1044
1046{
1047 if(IsLocked()) { return; }
1049}
1050
1052{
1054}
1055
1057{
1058 G4Exception("G4EmParameters", "em0044", JustWarning, ed);
1059}
1060
1062 const G4String& region,
1063 const G4String& type)
1064{
1065 if(IsLocked()) { return; }
1066 fBParameters->AddPAIModel(particle, region, type);
1067}
1068
1069const std::vector<G4String>& G4EmParameters::ParticlesPAI() const
1070{
1071 return fBParameters->ParticlesPAI();
1072}
1073
1074const std::vector<G4String>& G4EmParameters::RegionsPAI() const
1075{
1076 return fBParameters->RegionsPAI();
1077}
1078
1079const std::vector<G4String>& G4EmParameters::TypesPAI() const
1080{
1081 return fBParameters->TypesPAI();
1082}
1083
1085{
1086 if(IsLocked()) { return; }
1087 fCParameters->AddMicroElec(region);
1088}
1089
1090const std::vector<G4String>& G4EmParameters::RegionsMicroElec() const
1091{
1093}
1094
1095void G4EmParameters::AddDNA(const G4String& region, const G4String& type)
1096{
1097 if(IsLocked()) { return; }
1098 fCParameters->AddDNA(region, type);
1099 ActivateDNA();
1100}
1101
1102const std::vector<G4String>& G4EmParameters::RegionsDNA() const
1103{
1104 return fCParameters->RegionsDNA();
1105}
1106
1107const std::vector<G4String>& G4EmParameters::TypesDNA() const
1108{
1109 return fCParameters->TypesDNA();
1110}
1111
1112void G4EmParameters::AddPhysics(const G4String& region, const G4String& type)
1113{
1114 if(IsLocked()) { return; }
1115 fBParameters->AddPhysics(region, type);
1116}
1117
1118const std::vector<G4String>& G4EmParameters::RegionsPhysics() const
1119{
1120 return fBParameters->RegionsPhysics();
1121}
1122
1123const std::vector<G4String>& G4EmParameters::TypesPhysics() const
1124{
1125 return fBParameters->TypesPhysics();
1126}
1127
1129{
1130 if(IsLocked()) { return; }
1132}
1133
1134void
1136 G4bool aauger, G4bool apixe)
1137{
1138 if(IsLocked()) { return; }
1139 fCParameters->SetDeexActiveRegion(region, adeex, aauger, apixe);
1140}
1141
1142void
1144 G4double val, G4bool wflag)
1145{
1146 if(IsLocked()) { return; }
1147 fBParameters->SetProcessBiasingFactor(procname, val, wflag);
1148}
1149
1150void
1152 const G4String& region,
1153 G4double length,
1154 G4bool wflag)
1155{
1156 if(IsLocked() && !gener) { return; }
1157 fBParameters->ActivateForcedInteraction(procname, region, length, wflag);
1158}
1159
1160void
1162 const G4String& region,
1163 G4double factor,
1164 G4double energyLim)
1165{
1166 if(IsLocked()) { return; }
1167 fBParameters->ActivateSecondaryBiasing(procname, region, factor, energyLim);
1168}
1169
1171{
1173}
1174
1176{
1178}
1179
1181{
1183}
1184
1186{
1187 if(IsLocked()) { return; }
1189}
1190
1193}
1194
1196{
1197 if(IsLocked()) { return; }
1199}
1200
1202{
1203 if(IsLocked()) { return; }
1205}
1206
1208{
1210}
1211
1213{
1214 if(IsLocked()) { return; }
1216}
1217
1219{
1221}
1222
1224{
1226}
1227
1228void G4EmParameters::StreamInfo(std::ostream& os) const
1229{
1230 G4int prec = os.precision(5);
1231 os << "=======================================================================" << "\n";
1232 os << "====== Electromagnetic Physics Parameters ========" << "\n";
1233 os << "=======================================================================" << "\n";
1234 os << "LPM effect enabled " <<flagLPM << "\n";
1235 os << "Enable creation and use of sampling tables " <<fSamplingTable << "\n";
1236 os << "Apply cuts on all EM processes " <<applyCuts << "\n";
1237 os << "Use general process " <<gener << "\n";
1238 os << "Enable linear polarisation for gamma " <<fPolarisation << "\n";
1239 os << "Enable sampling of quantum entanglement "
1240 <<fBParameters->QuantumEntanglement() << "\n";
1241 os << "X-section factor for integral approach " <<lambdaFactor << "\n";
1242 os << "Min kinetic energy for tables "
1243 <<G4BestUnit(minKinEnergy,"Energy") << "\n";
1244 os << "Max kinetic energy for tables "
1245 <<G4BestUnit(maxKinEnergy,"Energy") << "\n";
1246 os << "Number of bins per decade of a table " <<nbinsPerDecade << "\n";
1247 os << "Verbose level " <<verbose << "\n";
1248 os << "Verbose level for worker thread " <<workerVerbose << "\n";
1249 os << "Bremsstrahlung energy threshold above which \n"
1250 << " primary e+- is added to the list of secondary "
1251 <<G4BestUnit(bremsTh,"Energy") << "\n";
1252 os << "Bremsstrahlung energy threshold above which primary\n"
1253 << " muon/hadron is added to the list of secondary "
1254 <<G4BestUnit(bremsMuHadTh,"Energy") << "\n";
1255 os << "Lowest triplet kinetic energy "
1256 <<G4BestUnit(lowestTripletEnergy,"Energy") << "\n";
1257 os << "Enable sampling of gamma linear polarisation " <<fPolarisation << "\n";
1258 os << "5D gamma conversion model type " <<tripletConv << "\n";
1259 os << "5D gamma conversion model on isolated ion " <<onIsolated << "\n";
1260 if(max5DEnergyForMuPair>0.0) {
1261 os << "5D gamma conversion limit for muon pair "
1262 << max5DEnergyForMuPair/CLHEP::GeV << " GeV\n";
1263 }
1264 os << "Livermore data directory "
1265 << fCParameters->LivermoreDataDir() << "\n";
1266
1267 os << "=======================================================================" << "\n";
1268 os << "====== Ionisation Parameters ========" << "\n";
1269 os << "=======================================================================" << "\n";
1270 os << "Step function for e+- "
1271 <<"("<<fBParameters->GetStepFunctionP1() << ", "
1272 << fBParameters->GetStepFunctionP2()/CLHEP::mm << " mm)\n";
1273 os << "Step function for muons/hadrons "
1274 <<"("<<fBParameters->GetStepFunctionMuHadP1() << ", "
1276 os << "Step function for light ions "
1277 <<"("<<fBParameters->GetStepFunctionLightIonsP1() << ", "
1279 os << "Step function for general ions "
1280 <<"("<<fBParameters->GetStepFunctionIonsP1() << ", "
1282 os << "Lowest e+e- kinetic energy "
1283 <<G4BestUnit(lowestElectronEnergy,"Energy") << "\n";
1284 os << "Lowest muon/hadron kinetic energy "
1285 <<G4BestUnit(lowestMuHadEnergy,"Energy") << "\n";
1286 os << "Fluctuations of dE/dx are enabled " <<lossFluctuation << "\n";
1287 os << "Use ICRU90 data " << fICRU90 << "\n";
1288 os << "Use built-in Birks satuaration " << birks << "\n";
1289 os << "Build CSDA range enabled " <<buildCSDARange << "\n";
1290 os << "Use cut as a final range enabled " <<cutAsFinalRange << "\n";
1291 os << "Enable angular generator interface "
1293 os << "Max kinetic energy for CSDA tables "
1294 <<G4BestUnit(maxKinEnergyCSDA,"Energy") << "\n";
1295 os << "Max kinetic energy for NIEL computation "
1296 <<G4BestUnit(maxNIELEnergy,"Energy") << "\n";
1297 os << "Linear loss limit " <<linLossLimit << "\n";
1298 os << "Read data from file for e+e- pair production by mu " <<fMuDataFromFile << "\n";
1299
1300 os << "=======================================================================" << "\n";
1301 os << "====== Multiple Scattering Parameters ========" << "\n";
1302 os << "=======================================================================" << "\n";
1303 os << "Type of msc step limit algorithm for e+- " <<mscStepLimit << "\n";
1304 os << "Type of msc step limit algorithm for muons/hadrons " <<mscStepLimitMuHad << "\n";
1305 os << "Msc lateral displacement for e+- enabled " <<lateralDisplacement << "\n";
1306 os << "Msc lateral displacement for muons and hadrons " <<muhadLateralDisplacement << "\n";
1307 os << "Urban msc model lateral displacement alg96 " <<lateralDisplacementAlg96 << "\n";
1308 os << "Range factor for msc step limit for e+- " <<rangeFactor << "\n";
1309 os << "Range factor for msc step limit for muons/hadrons " <<rangeFactorMuHad << "\n";
1310 os << "Geometry factor for msc step limitation of e+- " <<geomFactor << "\n";
1311 os << "Safety factor for msc step limit for e+- " <<safetyFactor << "\n";
1312 os << "Skin parameter for msc step limitation of e+- " <<skin << "\n";
1313 os << "Lambda limit for msc step limit for e+- " <<lambdaLimit/CLHEP::mm << " mm\n";
1314 os << "Use Mott correction for e- scattering " << useMottCorrection << "\n";
1315 os << "Factor used for dynamic computation of angular \n"
1316 << " limit between single and multiple scattering " << factorForAngleLimit << "\n";
1317 os << "Fixed angular limit between single \n"
1318 << " and multiple scattering "
1319 << thetaLimit/CLHEP::rad << " rad\n";
1320 os << "Upper energy limit for e+- multiple scattering "
1321 << energyLimit/CLHEP::MeV << " MeV\n";
1322 os << "Type of electron single scattering model " <<fSStype << "\n";
1323 os << "Type of nuclear form-factor " <<nucFormfactor << "\n";
1324 os << "Screening factor " <<factorScreen << "\n";
1325 os << "=======================================================================" << "\n";
1326
1327 if(fCParameters->Fluo()) {
1328 os << "====== Atomic Deexcitation Parameters ========" << "\n";
1329 os << "=======================================================================" << "\n";
1330 os << "Fluorescence enabled " <<fCParameters->Fluo() << "\n";
1331 os << "Fluorescence Bearden data files enabled "
1332 <<fCParameters->BeardenFluoDir() << "\n";
1333 os << "Fluorescence ANSTO data files enabled "
1334 <<fCParameters->ANSTOFluoDir() << "\n";
1335 os << "Auger electron cascade enabled "
1336 <<fCParameters->Auger() << "\n";
1337 os << "PIXE atomic de-excitation enabled " <<fCParameters->Pixe() << "\n";
1338 os << "De-excitation module ignores cuts "
1340 os << "Type of PIXE cross section for hadrons "
1342 os << "Type of PIXE cross section for e+- "
1344 os << "=======================================================================" << "\n";
1345 }
1346 if(fDNA) {
1347 os << "====== DNA Physics Parameters ========" << "\n";
1348 os << "=======================================================================" << "\n";
1349 os << "Use fast sampling in DNA models "
1350 << fCParameters->DNAFast() << "\n";
1351 os << "Use Stationary option in DNA models "
1352 << fCParameters->DNAStationary() << "\n";
1353 os << "Use DNA with multiple scattering of e- "
1354 << fCParameters->DNAElectronMsc() << "\n";
1355 os << "Use DNA e- solvation model type "
1356 << fCParameters->DNAeSolvationSubType() << "\n";
1357 os << "=======================================================================" << G4endl;
1358 }
1359 os.precision(prec);
1360}
1361
1363{
1364 if(fIsPrinted) return;
1365
1366#ifdef G4MULTITHREADED
1367 G4MUTEXLOCK(&emParametersMutex);
1368#endif
1370#ifdef G4MULTITHREADED
1371 G4MUTEXUNLOCK(&emParametersMutex);
1372#endif
1373}
1374
1375std::ostream& operator<< (std::ostream& os, const G4EmParameters& par)
1376{
1377 par.StreamInfo(os);
1378 return os;
1379}
1380
1382{
1383 return (!G4Threading::IsMasterThread() ||
1387}
1388
1389//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
static const char sss[MAX_N_PAR+2]
Definition: Evaluator.cc:63
@ G4State_Init
@ G4State_Idle
@ G4State_PreInit
G4DNAModelSubType
G4eSingleScatteringType
@ fWVI
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4MscStepLimitType
@ fMinimal
@ fUseSafety
G4NuclearFormfactorType
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double pi
Definition: G4SIunits.hh:55
#define G4BestUnit(a, b)
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
std::mutex G4Mutex
Definition: G4Threading.hh:81
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetStepFunctionP1() const
G4double GetStepFunctionP2() const
void SetQuantumEntanglement(G4bool v)
void SetStepFunctionMuHad(G4double v1, G4double v2)
void SetDirectionalSplittingRadius(G4double r)
G4double GetStepFunctionIonsP2() const
void SetStepFunctionIons(G4double v1, G4double v2)
const std::vector< G4String > & ParticlesPAI() const
G4double GetStepFunctionMuHadP1() const
void ActivateForcedInteraction(const G4String &procname, const G4String &region, G4double length, G4bool wflag)
void ActivateSecondaryBiasing(const G4String &name, const G4String &region, G4double factor, G4double energyLimit)
void SetDirectionalSplitting(G4bool v)
const std::vector< G4String > & TypesPhysics() const
const std::vector< G4String > & TypesPAI() const
void DefineRegParamForEM(G4VEmProcess *) const
void SetDirectionalSplittingTarget(const G4ThreeVector &v)
const std::vector< G4String > & RegionsPAI() const
void SetSubCutRegion(const G4String &region)
G4double GetStepFunctionLightIonsP1() const
void FillStepFunction(const G4ParticleDefinition *, G4VEnergyLossProcess *) const
void SetProcessBiasingFactor(const G4String &procname, G4double val, G4bool wflag)
G4ThreeVector GetDirectionalSplittingTarget() const
void DefineRegParamForLoss(G4VEnergyLossProcess *) const
const std::vector< G4String > & RegionsPhysics() const
G4double GetStepFunctionIonsP1() const
void SetStepFunction(G4double v1, G4double v2)
G4double GetStepFunctionLightIonsP2() const
G4double GetStepFunctionMuHadP2() const
G4double GetDirectionalSplittingRadius()
void SetStepFunctionLightIons(G4double v1, G4double v2)
void AddPAIModel(const G4String &particle, const G4String &region, const G4String &type)
void AddPhysics(const G4String &region, const G4String &type)
void SetAuger(G4bool val)
void SetDeexActiveRegion(const G4String &region, G4bool fdeex, G4bool fauger, G4bool fpixe)
void SetLivermoreDataDir(const G4String &)
void SetDNAFast(G4bool val)
G4bool DNAStationary() const
void SetDeexcitationIgnoreCut(G4bool val)
const std::vector< G4String > & TypesDNA() const
G4bool BeardenFluoDir() const
void SetDNAElectronMsc(G4bool val)
const G4String & LivermoreDataDir()
void SetFluo(G4bool val)
const std::vector< G4String > & RegionsMicroElec() const
G4bool DeexcitationIgnoreCut() const
const G4String & PIXECrossSectionModel()
void DefineRegParamForDeex(G4VAtomDeexcitation *) const
G4DNAModelSubType DNAeSolvationSubType() const
void AddDNA(const G4String &region, const G4String &type)
void SetDNAStationary(G4bool val)
void SetDNAeSolvationSubType(G4DNAModelSubType val)
const G4String & PIXEElectronCrossSectionModel()
void SetBeardenFluoDir(G4bool val)
void SetANSTOFluoDir(G4bool val)
void SetPIXECrossSectionModel(const G4String &)
G4bool DNAElectronMsc() const
void SetPixe(G4bool val)
const std::vector< G4String > & RegionsDNA() const
G4bool ANSTOFluoDir() const
void SetPIXEElectronCrossSectionModel(const G4String &)
void AddMicroElec(const G4String &region)
void SetEmSaturation(G4EmSaturation *)
void SetBeardenFluoDir(G4bool val)
G4bool IsPrintLocked() const
void DefineRegParamForLoss(G4VEnergyLossProcess *) const
void SetLambdaFactor(G4double val)
static G4EmParameters * theInstance
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetBuildCSDARange(G4bool val)
G4double rangeFactorMuHad
void SetStepFunctionLightIons(G4double v1, G4double v2)
void PrintWarning(G4ExceptionDescription &ed) const
void SetEnablePolarisation(G4bool val)
G4StateManager * fStateManager
void AddDNA(const G4String &region, const G4String &type)
G4bool LateralDisplacementAlg96() const
void FillStepFunction(const G4ParticleDefinition *, G4VEnergyLossProcess *) const
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
void SetDirectionalSplittingTarget(const G4ThreeVector &v)
G4bool BeardenFluoDir() const
G4bool DNAFast() const
G4DNAModelSubType DNAeSolvationSubType() const
G4EmSaturation * emSaturation
G4bool RetrieveMuDataFromFile() 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)
G4double factorScreen
void SetFactorForAngleLimit(G4double val)
const G4String & PIXECrossSectionModel()
const G4String & PIXEElectronCrossSectionModel()
G4double energyLimit
G4double MaxNIELEnergy() const
void SetRetrieveMuDataFromFile(G4bool v)
void SetDirectionalSplitting(G4bool v)
G4double lambdaLimit
const G4String & LivermoreDataDir()
G4bool OnIsolated() const
void SetMscMuHadRangeFactor(G4double val)
G4bool DNAElectronMsc() const
G4double lowestElectronEnergy
G4double MinKinEnergy() const
G4int NumberOfBinsPerDecade() const
G4MscStepLimitType mscStepLimit
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 SetDNAStationary(G4bool val)
G4bool ANSTOFluoDir() const
void SetDNAElectronMsc(G4bool val)
const std::vector< G4String > & TypesPhysics() const
void SetMaxEnergyFor5DMuPair(G4double val)
G4double GetDirectionalSplittingRadius()
G4bool muhadLateralDisplacement
void SetLinearLossLimit(G4double val)
void SetMscThetaLimit(G4double val)
G4int GetConversionType() const
G4MscStepLimitType MscMuHadStepLimitType() const
G4double MscEnergyLimit() const
G4bool lateralDisplacementAlg96
G4double factorForAngleLimit
void ActivateSecondaryBiasing(const G4String &name, const G4String &region, G4double factor, G4double energyLimit)
G4bool BirksActive() const
G4bool DNAStationary() const
G4bool IsLocked() const
G4double minKinEnergy
const std::vector< G4String > & RegionsPAI() const
G4double lowestTripletEnergy
G4double maxKinEnergyCSDA
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)
G4bool useAngGeneratorForIonisation
void SetQuantumEntanglement(G4bool v)
void DefineRegParamForEM(G4VEmProcess *) const
void ActivateForcedInteraction(const G4String &procname, const G4String &region, G4double length, G4bool wflag)
G4double lowestMuHadEnergy
void ActivateAngularGeneratorForIonisation(G4bool val)
G4double maxKinEnergy
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 safetyFactor
G4double ScreeningFactor() const
G4bool UseMottCorrection() const
void SetFluo(G4bool val)
void SetMuHadBremsstrahlungTh(G4double val)
const std::vector< G4String > & ParticlesPAI() const
G4double max5DEnergyForMuPair
G4bool Pixe() const
const std::vector< G4String > & RegionsPhysics() const
G4double lambdaFactor
G4MscStepLimitType mscStepLimitMuHad
void DefineRegParamForDeex(G4VAtomDeexcitation *) const
void SetStepFunctionMuHad(G4double v1, G4double v2)
G4double MscSafetyFactor() const
G4NuclearFormfactorType nucFormfactor
G4double rangeFactor
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
G4MscStepLimitType MscStepLimitType() const
G4bool ApplyCuts() const
G4double linLossLimit
void SetEnableSamplingTable(G4bool val)
const std::vector< G4String > & RegionsDNA() const
void SetLivermoreDataDir(const G4String &)
G4double BremsstrahlungTh() const
void SetPixe(G4bool val)
void SetMaxNIELEnergy(G4double val)
void SetStepFunctionIons(G4double v1, G4double v2)
void SetMaxEnergyForCSDARange(G4double val)
G4eSingleScatteringType SingleScatteringType() const
G4bool DeexcitationIgnoreCut() const
void SetProcessBiasingFactor(const G4String &procname, G4double val, G4bool wflag)
G4double MscGeomFactor() const
void SetMscMuHadStepLimitType(G4MscStepLimitType val)
void SetMscStepLimitType(G4MscStepLimitType val)
void AddPhysics(const G4String &region, const G4String &type)
void SetMscEnergyLimit(G4double val)
void SetBremsstrahlungTh(G4double val)
void SetAuger(G4bool val)
G4bool lateralDisplacement
G4bool GetDirectionalSplitting() const
G4EmLowEParameters * fCParameters
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)
G4double bremsMuHadTh
G4bool MuHadLateralDisplacement() const
void SetOnIsolated(G4bool val)
G4bool EnableSamplingTable() const
G4EmParametersMessenger * theMessenger
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
G4double MaxEnergyForCSDARange() const
G4EmExtraParameters * fBParameters
G4bool useMottCorrection
void SetLowestMuHadEnergy(G4double val)
const std::vector< G4String > & TypesPAI() const
G4bool LPM() const
G4bool UseAngularGeneratorForIonisation() const
G4double maxNIELEnergy
void SetMaxEnergy(G4double val)
G4double LinearLossLimit() const
G4NuclearFormfactorType NuclearFormfactorType() const
G4double LowestMuHadEnergy() const
G4double MscRangeFactor() const
G4double LambdaFactor() const
G4double FactorForAngleLimit() const
G4eSingleScatteringType fSStype
G4double MscSkin() const
void SetSingleScatteringType(G4eSingleScatteringType val)
void SetANSTOFluoDir(G4bool val)
G4double LowestTripletEnergy() const
G4double LowestElectronEnergy() const
void SetMscRangeFactor(G4double val)
G4bool GeneralProcessActive() const
static G4NistManager * Instance()
const G4ApplicationState & GetCurrentState() const
static G4StateManager * GetStateManager()
std::ostream & operator<<(std::ostream &, const BasicVector3D< float > &)
static constexpr double mm
Definition: SystemOfUnits.h:96
static const double prec
Definition: RanecuEngine.cc:61
static constexpr double TeV
static constexpr double GeV
static constexpr double keV
static constexpr double MeV
static constexpr double eV
static constexpr double rad
static constexpr double pi
Definition: SystemOfUnits.h:55
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4bool IsMasterThread()
Definition: G4Threading.cc:124
int G4lrint(double ad)
Definition: templates.hh:134