Geant4 11.1.1
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
95G4EmParameters::G4EmParameters()
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 fDNA = false;
140 fIsPrinted = false;
141
142 minKinEnergy = 0.1*CLHEP::keV;
143 maxKinEnergy = 100.0*CLHEP::TeV;
144 maxKinEnergyCSDA = 1.0*CLHEP::GeV;
145 max5DEnergyForMuPair = 0.0;
146 lowestElectronEnergy = 1.0*CLHEP::keV;
147 lowestMuHadEnergy = 1.0*CLHEP::keV;
148 lowestTripletEnergy = 1.0*CLHEP::MeV;
149 maxNIELEnergy = 0.0;
150 linLossLimit = 0.01;
151 bremsTh = bremsMuHadTh = maxKinEnergy;
152 lambdaFactor = 0.8;
153 factorForAngleLimit = 1.0;
154 thetaLimit = CLHEP::pi;
155 energyLimit = 100.0*CLHEP::MeV;
156 rangeFactor = 0.04;
157 rangeFactorMuHad = 0.2;
158 geomFactor = 2.5;
159 skin = 1.0;
160 safetyFactor = 0.6;
161 lambdaLimit = 1.0*CLHEP::mm;
162 factorScreen = 1.0;
163
164 nbinsPerDecade = 7;
165 verbose = 1;
166 workerVerbose = 0;
167 tripletConv = 0;
168
169 fTransportationWithMsc = G4TransportationWithMscType::fDisabled;
170 mscStepLimit = fUseSafety;
171 mscStepLimitMuHad = fMinimal;
172 nucFormfactor = fExponentialNF;
173 fSStype = fWVI;
174 fFluct = fUniversalFluctuation;
175}
176
178{
179 if(IsLocked()) { return; }
180 lossFluctuation = val;
181}
182
184{
185 return lossFluctuation;
186}
187
189{
190 if(IsLocked()) { return; }
191 buildCSDARange = val;
192}
193
195{
196 return buildCSDARange;
197}
198
200{
201 if(IsLocked()) { return; }
202 flagLPM = val;
203}
204
206{
207 return flagLPM;
208}
209
211{
212 if(IsLocked()) { return; }
213 cutAsFinalRange = val;
214}
215
217{
218 return cutAsFinalRange;
219}
220
222{
223 if(IsLocked()) { return; }
224 applyCuts = val;
225}
226
228{
229 return applyCuts;
230}
231
233{
234 if(IsLocked()) { return; }
235 fCParameters->SetFluo(val);
236}
237
239{
240 return fCParameters->Fluo();
241}
242
244{
245 return fCParameters->FluoDirectory();
246}
247
249{
250 if(IsLocked()) { return; }
251 fCParameters->SetFluoDirectory(val);
252}
253
255{
256 if(IsLocked()) { return; }
257 fCParameters->SetBeardenFluoDir(val);
258}
259
261{
262 if(IsLocked()) { return; }
263 fCParameters->SetANSTOFluoDir(val);
264}
265
267{
268 if(IsLocked()) { return; }
269 fCParameters->SetXDB_EADLFluoDir(val);
270}
271
273{
274 if(IsLocked()) { return; }
275 fCParameters->SetAuger(val);
276}
277
279{
280 auto dir = fCParameters->FluoDirectory();
281 return (dir == fluoBearden);
282}
283
285{
286 auto dir = fCParameters->FluoDirectory();
287 return (dir == fluoANSTO);
288}
289
291{
292 return fCParameters->Auger();
293}
294
296{
297 if(IsLocked()) { return; }
298 fCParameters->SetPixe(val);
299}
300
302{
303 return fCParameters->Pixe();
304}
305
307{
308 if(IsLocked()) { return; }
309 fCParameters->SetDeexcitationIgnoreCut(val);
310}
311
313{
314 return fCParameters->DeexcitationIgnoreCut();
315}
316
318{
319 if(IsLocked()) { return; }
320 lateralDisplacement = val;
321}
322
324{
325 return lateralDisplacement;
326}
327
329{
330 if(IsLocked()) { return; }
331 lateralDisplacementAlg96 = val;
332}
333
335{
336 return lateralDisplacementAlg96;
337}
338
340{
341 if(IsLocked()) { return; }
342 muhadLateralDisplacement = val;
343}
344
346{
347 return muhadLateralDisplacement;
348}
349
351{
352 if(IsLocked()) { return; }
353 useAngGeneratorForIonisation = val;
354}
355
357{
358 return useAngGeneratorForIonisation;
359}
360
362{
363 if(IsLocked()) { return; }
364 useMottCorrection = val;
365}
366
368{
369 return useMottCorrection;
370}
371
373{
374 if(IsLocked()) { return; }
375 integral = val;
376}
377
379{
380 return integral;
381}
382
384{
385 if(IsLocked()) { return; }
386 fPolarisation = val;
387}
388
390{
391 return fPolarisation;
392}
393
395{
396 if(IsLocked()) { return; }
397 birks = val;
398 if(birks && nullptr == emSaturation) { emSaturation = new G4EmSaturation(1); }
399}
400
402{
403 return birks;
404}
405
407{
408 if(IsLocked()) { return; }
409 fICRU90 = val;
410}
411
413{
414 return fICRU90;
415}
416
418{
419 if(IsLocked()) { return; }
420 fCParameters->SetDNAFast(val);
421 if(val) { ActivateDNA(); }
422}
423
425{
426 return fCParameters->DNAFast();
427}
428
430{
431 if(IsLocked()) { return; }
432 fCParameters->SetDNAStationary(val);
433 if(val) { ActivateDNA(); }
434}
435
437{
438 return fCParameters->DNAStationary();
439}
440
442{
443 if(IsLocked()) { return; }
444 fCParameters->SetDNAElectronMsc(val);
445 if(val) { ActivateDNA(); }
446}
447
449{
450 return fCParameters->DNAElectronMsc();
451}
452
454{
455 if(IsLocked()) { return; }
456 gener = val;
457}
458
460{
461 return gener;
462}
463
465{
466 if(IsLocked()) { return; }
467 birks = (nullptr != ptr);
468 if(emSaturation != ptr) {
469 delete emSaturation;
470 emSaturation = ptr;
471 }
472}
473
475{
476 return fMuDataFromFile;
477}
478
480{
481 fMuDataFromFile = v;
482}
483
485{
486 if(IsLocked()) { return; }
487 onIsolated = val;
488}
489
491{
492 return onIsolated;
493}
494
496{
497 if(IsLocked()) { return; }
498 fSamplingTable = val;
499}
500
502{
503 return fSamplingTable;
504}
505
507{
508 return fPEKShell;
509}
510
512{
513 if(IsLocked()) { return; }
514 fPEKShell = v;
515}
516
518{
519 return fMscPosiCorr;
520}
521
523{
524 if(IsLocked()) { return; }
525 fMscPosiCorr = v;
526}
527
529{
530 if(IsLocked()) { return; }
531 fDNA = true;
532}
533
535{
536 fIsPrinted = val;
537}
538
540{
541 return fIsPrinted;
542}
543
545{
546 if(nullptr == emSaturation) {
547#ifdef G4MULTITHREADED
548 G4MUTEXLOCK(&emParametersMutex);
549 if(nullptr == emSaturation) {
550#endif
551 emSaturation = new G4EmSaturation(1);
552#ifdef G4MULTITHREADED
553 }
554 G4MUTEXUNLOCK(&emParametersMutex);
555#endif
556 }
557 birks = true;
558 return emSaturation;
559}
560
562{
563 if(IsLocked()) { return; }
564 if(val > 1.e-3*CLHEP::eV && val < maxKinEnergy) {
565 minKinEnergy = val;
566 } else {
568 ed << "Value of MinKinEnergy - is out of range: " << val/CLHEP::MeV
569 << " MeV is ignored";
570 PrintWarning(ed);
571 }
572}
573
575{
576 return minKinEnergy;
577}
578
580{
581 if(IsLocked()) { return; }
582 if(val > std::max(minKinEnergy,9.99*CLHEP::MeV) && val < 1.e+7*CLHEP::TeV) {
583 maxKinEnergy = val;
584 } else {
586 ed << "Value of MaxKinEnergy is out of range: "
587 << val/CLHEP::GeV
588 << " GeV is ignored; allowed range 10 MeV - 1.e+7 TeV";
589 PrintWarning(ed);
590 }
591}
592
594{
595 return maxKinEnergy;
596}
597
599{
600 if(IsLocked()) { return; }
601 if(val > minKinEnergy && val <= 100*CLHEP::TeV) {
602 maxKinEnergyCSDA = val;
603 } else {
605 ed << "Value of MaxKinEnergyCSDA is out of range: "
606 << val/CLHEP::GeV << " GeV is ignored; allowed range "
607 << minKinEnergy << " MeV - 100 TeV";
608 PrintWarning(ed);
609 }
610}
611
613{
614 return maxKinEnergyCSDA;
615}
616
618{
619 if(IsLocked()) { return; }
620 if(val >= 0.0) { lowestElectronEnergy = val; }
621}
622
624{
625 return lowestElectronEnergy;
626}
627
629{
630 if(IsLocked()) { return; }
631 if(val >= 0.0) { lowestMuHadEnergy = val; }
632}
633
635{
636 return lowestMuHadEnergy;
637}
638
640{
641 if(IsLocked()) { return; }
642 if(val > 0.0) { lowestTripletEnergy = val; }
643}
644
646{
647 return lowestTripletEnergy;
648}
649
651{
652 if(IsLocked()) { return; }
653 if(val >= 0.0) { maxNIELEnergy = val; }
654}
655
657{
658 return maxNIELEnergy;
659}
660
662{
663 if(IsLocked()) { return; }
664 if(val > 0.0) { max5DEnergyForMuPair = val; }
665}
666
668{
669 return max5DEnergyForMuPair;
670}
671
673{
674 if(IsLocked()) { return; }
675 if(val > 0.0 && val < 0.5) {
676 linLossLimit = val;
677 } else {
679 ed << "Value of linLossLimit is out of range: " << val
680 << " is ignored";
681 PrintWarning(ed);
682 }
683}
684
686{
687 return linLossLimit;
688}
689
691{
692 if(IsLocked()) { return; }
693 if(val > 0.0) {
694 bremsTh = val;
695 } else {
697 ed << "Value of bremsstrahlung threshold is out of range: "
698 << val/GeV << " GeV is ignored";
699 PrintWarning(ed);
700 }
701}
702
704{
705 return bremsTh;
706}
707
709{
710 if(IsLocked()) { return; }
711 if(val > 0.0) {
712 bremsMuHadTh = val;
713 } else {
715 ed << "Value of bremsstrahlung threshold is out of range: "
716 << val/GeV << " GeV is ignored";
717 PrintWarning(ed);
718 }
719}
720
722{
723 return bremsMuHadTh;
724}
725
727{
728 if(IsLocked()) { return; }
729 if(val > 0.0 && val < 1.0) {
730 lambdaFactor = val;
731 } else {
733 ed << "Value of lambda factor is out of range: " << val
734 << " is ignored";
735 PrintWarning(ed);
736 }
737}
738
740{
741 return lambdaFactor;
742}
743
745{
746 if(IsLocked()) { return; }
747 if(val > 0.0) {
748 factorForAngleLimit = val;
749 } else {
751 ed << "Value of factor for enegry limit is out of range: "
752 << val << " is ignored";
753 PrintWarning(ed);
754 }
755}
756
758{
759 return factorForAngleLimit;
760}
761
763{
764 if(IsLocked()) { return; }
765 if(val >= 0.0 && val <= pi) {
766 thetaLimit = val;
767 } else {
769 ed << "Value of polar angle limit is out of range: "
770 << val << " is ignored";
771 PrintWarning(ed);
772 }
773}
774
776{
777 return thetaLimit;
778}
779
781{
782 if(IsLocked()) { return; }
783 if(val >= 0.0) {
784 energyLimit = val;
785 } else {
787 ed << "Value of msc energy limit is out of range: "
788 << val << " is ignored";
789 PrintWarning(ed);
790 }
791}
792
794{
795 return energyLimit;
796}
797
799{
800 if(IsLocked()) { return; }
801 if(val > 0.0 && val < 1.0) {
802 rangeFactor = val;
803 } else {
805 ed << "Value of rangeFactor is out of range: "
806 << val << " is ignored";
807 PrintWarning(ed);
808 }
809}
810
812{
813 return rangeFactor;
814}
815
817{
818 if(IsLocked()) { return; }
819 if(val > 0.0 && val < 1.0) {
820 rangeFactorMuHad = val;
821 } else {
823 ed << "Value of rangeFactorMuHad is out of range: "
824 << val << " is ignored";
825 PrintWarning(ed);
826 }
827}
828
830{
831 return rangeFactorMuHad;
832}
833
835{
836 if(IsLocked()) { return; }
837 if(val >= 1.0) {
838 geomFactor = val;
839 } else {
841 ed << "Value of geomFactor is out of range: "
842 << val << " is ignored";
843 PrintWarning(ed);
844 }
845}
846
848{
849 return geomFactor;
850}
851
853{
854 if(IsLocked()) { return; }
855 if(val >= 0.1) {
856 safetyFactor = val;
857 } else {
859 ed << "Value of safetyFactor is out of range: "
860 << val << " is ignored";
861 PrintWarning(ed);
862 }
863}
864
866{
867 return safetyFactor;
868}
869
871{
872 if(IsLocked()) { return; }
873 if(val >= 0.0) {
874 lambdaLimit = val;
875 } else {
877 ed << "Value of lambdaLimit is out of range: "
878 << val << " is ignored";
879 PrintWarning(ed);
880 }
881}
882
884{
885 return lambdaLimit;
886}
887
889{
890 if(IsLocked()) { return; }
891 if(val >= 1.0) {
892 skin = val;
893 } else {
895 ed << "Value of skin is out of range: "
896 << val << " is ignored";
897 PrintWarning(ed);
898 }
899}
900
902{
903 return skin;
904}
905
907{
908 if(IsLocked()) { return; }
909 if(val > 0.0) {
910 factorScreen = val;
911 } else {
913 ed << "Value of factorScreen is out of range: "
914 << val << " is ignored";
915 PrintWarning(ed);
916 }
917}
918
920{
921 return factorScreen;
922}
923
925{
926 if(IsLocked()) { return; }
927 fBParameters->SetStepFunction(v1, v2);
928}
929
931{
932 if(IsLocked()) { return; }
933 fBParameters->SetStepFunctionMuHad(v1, v2);
934}
935
937{
938 if(IsLocked()) { return; }
939 fBParameters->SetStepFunctionLightIons(v1, v2);
940}
941
943{
944 if(IsLocked()) { return; }
945 fBParameters->SetStepFunctionIons(v1, v2);
946}
947
949{
950 fBParameters->FillStepFunction(part, proc);
951}
952
954{
955 return nbinsPerDecade*G4lrint(std::log10(maxKinEnergy/minKinEnergy));
956}
957
959{
960 if(IsLocked()) { return; }
961 if(val >= 5 && val < 1000000) {
962 nbinsPerDecade = val;
963 } else {
965 ed << "Value of number of bins per decade is out of range: "
966 << val << " is ignored";
967 PrintWarning(ed);
968 }
969}
970
972{
973 return nbinsPerDecade;
974}
975
977{
978 if(IsLocked()) { return; }
979 verbose = val;
980 workerVerbose = std::min(workerVerbose, verbose);
981}
982
984{
985 return verbose;
986}
987
989{
990 if(IsLocked()) { return; }
991 workerVerbose = val;
992}
993
995{
996 return workerVerbose;
997}
998
1000{
1001 if(IsLocked()) { return; }
1002 fTransportationWithMsc = val;
1003}
1004
1006{
1007 return fTransportationWithMsc;
1008}
1009
1011{
1012 if(IsLocked()) { return; }
1013 fFluct = val;
1014}
1015
1017{
1018 return fFluct;
1019}
1020
1022{
1023 if(IsLocked()) { return; }
1024 mscStepLimit = val;
1025}
1026
1028{
1029 return mscStepLimit;
1030}
1031
1033{
1034 if(IsLocked()) { return; }
1035 mscStepLimitMuHad = val;
1036}
1037
1039{
1040 return mscStepLimitMuHad;
1041}
1042
1044{
1045 if(IsLocked()) { return; }
1046 fSStype = val;
1047}
1048
1050{
1051 return fSStype;
1052}
1053
1054void
1056{
1057 if(IsLocked()) { return; }
1058 nucFormfactor = val;
1059}
1060
1062{
1063 return nucFormfactor;
1064}
1065
1067{
1068 if(IsLocked()) { return; }
1069 fCParameters->SetDNAeSolvationSubType(val);
1070 ActivateDNA();
1071}
1072
1074{
1075 return fCParameters->DNAeSolvationSubType();
1076}
1077
1079{
1080 if(IsLocked()) { return; }
1081 tripletConv = val;
1082}
1083
1085{
1086 return tripletConv;
1087}
1088
1090{
1091 if(IsLocked()) { return; }
1092 fCParameters->SetPIXECrossSectionModel(sss);
1093}
1094
1096{
1097 return fCParameters->PIXECrossSectionModel();
1098}
1099
1101{
1102 if(IsLocked()) { return; }
1103 fCParameters->SetPIXEElectronCrossSectionModel(sss);
1104}
1105
1107{
1108 return fCParameters->PIXEElectronCrossSectionModel();
1109}
1110
1112{
1113 if(IsLocked()) { return; }
1114 fCParameters->SetLivermoreDataDir(sss);
1115}
1116
1118{
1119 return fCParameters->LivermoreDataDir();
1120}
1121
1122void G4EmParameters::PrintWarning(G4ExceptionDescription& ed) const
1123{
1124 G4Exception("G4EmParameters", "em0044", JustWarning, ed);
1125}
1126
1128 const G4String& region,
1129 const G4String& type)
1130{
1131 if(IsLocked()) { return; }
1132 fBParameters->AddPAIModel(particle, region, type);
1133}
1134
1135const std::vector<G4String>& G4EmParameters::ParticlesPAI() const
1136{
1137 return fBParameters->ParticlesPAI();
1138}
1139
1140const std::vector<G4String>& G4EmParameters::RegionsPAI() const
1141{
1142 return fBParameters->RegionsPAI();
1143}
1144
1145const std::vector<G4String>& G4EmParameters::TypesPAI() const
1146{
1147 return fBParameters->TypesPAI();
1148}
1149
1151{
1152 if(IsLocked()) { return; }
1153 fCParameters->AddMicroElec(region);
1154}
1155
1156const std::vector<G4String>& G4EmParameters::RegionsMicroElec() const
1157{
1158 return fCParameters->RegionsMicroElec();
1159}
1160
1161void G4EmParameters::AddDNA(const G4String& region, const G4String& type)
1162{
1163 if(IsLocked()) { return; }
1164 fCParameters->AddDNA(region, type);
1165 ActivateDNA();
1166}
1167
1168const std::vector<G4String>& G4EmParameters::RegionsDNA() const
1169{
1170 return fCParameters->RegionsDNA();
1171}
1172
1173const std::vector<G4String>& G4EmParameters::TypesDNA() const
1174{
1175 return fCParameters->TypesDNA();
1176}
1177
1178void G4EmParameters::AddPhysics(const G4String& region, const G4String& type)
1179{
1180 if(IsLocked()) { return; }
1181 fBParameters->AddPhysics(region, type);
1182}
1183
1184const std::vector<G4String>& G4EmParameters::RegionsPhysics() const
1185{
1186 return fBParameters->RegionsPhysics();
1187}
1188
1189const std::vector<G4String>& G4EmParameters::TypesPhysics() const
1190{
1191 return fBParameters->TypesPhysics();
1192}
1193
1195{
1196 if(IsLocked()) { return; }
1197 fBParameters->SetSubCutRegion(region);
1198}
1199
1200void
1202 G4bool aauger, G4bool apixe)
1203{
1204 if(IsLocked()) { return; }
1205 fCParameters->SetDeexActiveRegion(region, adeex, aauger, apixe);
1206}
1207
1208void
1210 G4double val, G4bool wflag)
1211{
1212 if(IsLocked()) { return; }
1213 fBParameters->SetProcessBiasingFactor(procname, val, wflag);
1214}
1215
1216void
1218 const G4String& region,
1219 G4double length,
1220 G4bool wflag)
1221{
1222 if(IsLocked() && !gener) { return; }
1223 fBParameters->ActivateForcedInteraction(procname, region, length, wflag);
1224}
1225
1226void
1228 const G4String& region,
1229 G4double factor,
1230 G4double energyLim)
1231{
1232 if(IsLocked()) { return; }
1233 fBParameters->ActivateSecondaryBiasing(procname, region, factor, energyLim);
1234}
1235
1237{
1238 fBParameters->DefineRegParamForLoss(ptr);
1239}
1240
1242{
1243 fBParameters->DefineRegParamForEM(ptr);
1244}
1245
1247{
1248 return fBParameters->QuantumEntanglement();
1249}
1250
1252{
1253 if(IsLocked()) { return; }
1254 fBParameters->SetQuantumEntanglement(v);
1255}
1256
1258 return fBParameters->GetDirectionalSplitting();
1259}
1260
1262{
1263 if(IsLocked()) { return; }
1264 fBParameters->SetDirectionalSplitting(v);
1265}
1266
1268{
1269 if(IsLocked()) { return; }
1270 fBParameters->SetDirectionalSplittingTarget(v);
1271}
1272
1274{
1275 return fBParameters->GetDirectionalSplittingTarget();
1276}
1277
1279{
1280 if(IsLocked()) { return; }
1281 fBParameters->SetDirectionalSplittingRadius(r);
1282}
1283
1285{
1286 return fBParameters->GetDirectionalSplittingRadius();
1287}
1288
1290{
1291 fCParameters->DefineRegParamForDeex(ptr);
1292}
1293
1294void G4EmParameters::StreamInfo(std::ostream& os) const
1295{
1296 G4long prec = os.precision(5);
1297 os << "=======================================================================" << "\n";
1298 os << "====== Electromagnetic Physics Parameters ========" << "\n";
1299 os << "=======================================================================" << "\n";
1300 os << "LPM effect enabled " <<flagLPM << "\n";
1301 os << "Enable creation and use of sampling tables " <<fSamplingTable << "\n";
1302 os << "Apply cuts on all EM processes " <<applyCuts << "\n";
1303 const char* transportationWithMsc = "Disabled";
1304 if(fTransportationWithMsc == G4TransportationWithMscType::fEnabled) {
1305 transportationWithMsc = "Enabled";
1306 } else if (fTransportationWithMsc == G4TransportationWithMscType::fMultipleSteps) {
1307 transportationWithMsc = "MultipleSteps";
1308 }
1309 os << "Use combined TransportationWithMsc " <<transportationWithMsc << "\n";
1310 os << "Use general process " <<gener << "\n";
1311 os << "Enable linear polarisation for gamma " <<fPolarisation << "\n";
1312 os << "Enable photoeffect sampling below K-shell " <<fPEKShell << "\n";
1313 os << "Enable sampling of quantum entanglement "
1314 <<fBParameters->QuantumEntanglement() << "\n";
1315 os << "X-section factor for integral approach " <<lambdaFactor << "\n";
1316 os << "Min kinetic energy for tables "
1317 <<G4BestUnit(minKinEnergy,"Energy") << "\n";
1318 os << "Max kinetic energy for tables "
1319 <<G4BestUnit(maxKinEnergy,"Energy") << "\n";
1320 os << "Number of bins per decade of a table " <<nbinsPerDecade << "\n";
1321 os << "Verbose level " <<verbose << "\n";
1322 os << "Verbose level for worker thread " <<workerVerbose << "\n";
1323 os << "Bremsstrahlung energy threshold above which \n"
1324 << " primary e+- is added to the list of secondary "
1325 <<G4BestUnit(bremsTh,"Energy") << "\n";
1326 os << "Bremsstrahlung energy threshold above which primary\n"
1327 << " muon/hadron is added to the list of secondary "
1328 <<G4BestUnit(bremsMuHadTh,"Energy") << "\n";
1329 os << "Lowest triplet kinetic energy "
1330 <<G4BestUnit(lowestTripletEnergy,"Energy") << "\n";
1331 os << "Enable sampling of gamma linear polarisation " <<fPolarisation << "\n";
1332 os << "5D gamma conversion model type " <<tripletConv << "\n";
1333 os << "5D gamma conversion model on isolated ion " <<onIsolated << "\n";
1334 if(max5DEnergyForMuPair>0.0) {
1335 os << "5D gamma conversion limit for muon pair "
1336 << max5DEnergyForMuPair/CLHEP::GeV << " GeV\n";
1337 }
1338 os << "Livermore data directory "
1339 << fCParameters->LivermoreDataDir() << "\n";
1340
1341 os << "=======================================================================" << "\n";
1342 os << "====== Ionisation Parameters ========" << "\n";
1343 os << "=======================================================================" << "\n";
1344 os << "Step function for e+- "
1345 <<"("<<fBParameters->GetStepFunctionP1() << ", "
1346 << fBParameters->GetStepFunctionP2()/CLHEP::mm << " mm)\n";
1347 os << "Step function for muons/hadrons "
1348 <<"("<<fBParameters->GetStepFunctionMuHadP1() << ", "
1349 << fBParameters->GetStepFunctionMuHadP2()/CLHEP::mm << " mm)\n";
1350 os << "Step function for light ions "
1351 <<"("<<fBParameters->GetStepFunctionLightIonsP1() << ", "
1352 << fBParameters->GetStepFunctionLightIonsP2()/CLHEP::mm << " mm)\n";
1353 os << "Step function for general ions "
1354 <<"("<<fBParameters->GetStepFunctionIonsP1() << ", "
1355 << fBParameters->GetStepFunctionIonsP2()/CLHEP::mm << " mm)\n";
1356 os << "Lowest e+e- kinetic energy "
1357 <<G4BestUnit(lowestElectronEnergy,"Energy") << "\n";
1358 os << "Lowest muon/hadron kinetic energy "
1359 <<G4BestUnit(lowestMuHadEnergy,"Energy") << "\n";
1360 os << "Use ICRU90 data " << fICRU90 << "\n";
1361 os << "Fluctuations of dE/dx are enabled " <<lossFluctuation << "\n";
1362 G4String namef = "Universal";
1363 if(fFluct == fUrbanFluctuation) { namef = "Urban"; }
1364 else if(fFluct == fDummyFluctuation) { namef = "Dummy"; }
1365 os << "Type of fluctuation model for leptons and hadrons " << namef << "\n";
1366 os << "Use built-in Birks satuaration " << birks << "\n";
1367 os << "Build CSDA range enabled " <<buildCSDARange << "\n";
1368 os << "Use cut as a final range enabled " <<cutAsFinalRange << "\n";
1369 os << "Enable angular generator interface "
1370 <<useAngGeneratorForIonisation << "\n";
1371 os << "Max kinetic energy for CSDA tables "
1372 <<G4BestUnit(maxKinEnergyCSDA,"Energy") << "\n";
1373 os << "Max kinetic energy for NIEL computation "
1374 <<G4BestUnit(maxNIELEnergy,"Energy") << "\n";
1375 os << "Linear loss limit " <<linLossLimit << "\n";
1376 os << "Read data from file for e+e- pair production by mu " <<fMuDataFromFile << "\n";
1377
1378 os << "=======================================================================" << "\n";
1379 os << "====== Multiple Scattering Parameters ========" << "\n";
1380 os << "=======================================================================" << "\n";
1381 os << "Type of msc step limit algorithm for e+- " <<mscStepLimit << "\n";
1382 os << "Type of msc step limit algorithm for muons/hadrons " <<mscStepLimitMuHad << "\n";
1383 os << "Msc lateral displacement for e+- enabled " <<lateralDisplacement << "\n";
1384 os << "Msc lateral displacement for muons and hadrons " <<muhadLateralDisplacement << "\n";
1385 os << "Urban msc model lateral displacement alg96 " <<lateralDisplacementAlg96 << "\n";
1386 os << "Range factor for msc step limit for e+- " <<rangeFactor << "\n";
1387 os << "Range factor for msc step limit for muons/hadrons " <<rangeFactorMuHad << "\n";
1388 os << "Geometry factor for msc step limitation of e+- " <<geomFactor << "\n";
1389 os << "Safety factor for msc step limit for e+- " <<safetyFactor << "\n";
1390 os << "Skin parameter for msc step limitation of e+- " <<skin << "\n";
1391 os << "Lambda limit for msc step limit for e+- " <<lambdaLimit/CLHEP::mm << " mm\n";
1392 os << "Use Mott correction for e- scattering " << useMottCorrection << "\n";
1393 os << "Factor used for dynamic computation of angular \n"
1394 << " limit between single and multiple scattering " << factorForAngleLimit << "\n";
1395 os << "Fixed angular limit between single \n"
1396 << " and multiple scattering "
1397 << thetaLimit/CLHEP::rad << " rad\n";
1398 os << "Upper energy limit for e+- multiple scattering "
1399 << energyLimit/CLHEP::MeV << " MeV\n";
1400 os << "Type of electron single scattering model " <<fSStype << "\n";
1401 os << "Type of nuclear form-factor " <<nucFormfactor << "\n";
1402 os << "Screening factor " <<factorScreen << "\n";
1403 os << "=======================================================================" << "\n";
1404
1405 if(fCParameters->Fluo()) {
1406 os << "====== Atomic Deexcitation Parameters ========" << "\n";
1407 os << "=======================================================================" << "\n";
1408 os << "Fluorescence enabled " <<fCParameters->Fluo() << "\n";
1409 G4String named = "fluor";
1411 if(fdir == fluoBearden) { named = "fluor_Bearden"; }
1412 else if(fdir == fluoANSTO) { named = "fluor_ANSTO"; }
1413 else if(fdir == fluoXDB_EADL) { named = "fluor_XDB_EADL"; }
1414 os << "Directory in G4LEDATA for fluorescence data files " << named << "\n";
1415 os << "Auger electron cascade enabled "
1416 <<fCParameters->Auger() << "\n";
1417 os << "PIXE atomic de-excitation enabled " <<fCParameters->Pixe() << "\n";
1418 os << "De-excitation module ignores cuts "
1419 <<fCParameters->DeexcitationIgnoreCut() << "\n";
1420 os << "Type of PIXE cross section for hadrons "
1421 <<fCParameters->PIXECrossSectionModel() << "\n";
1422 os << "Type of PIXE cross section for e+- "
1423 <<fCParameters->PIXEElectronCrossSectionModel() << "\n";
1424 os << "=======================================================================" << "\n";
1425 }
1426 if(fDNA) {
1427 os << "====== DNA Physics Parameters ========" << "\n";
1428 os << "=======================================================================" << "\n";
1429 os << "Use fast sampling in DNA models "
1430 << fCParameters->DNAFast() << "\n";
1431 os << "Use Stationary option in DNA models "
1432 << fCParameters->DNAStationary() << "\n";
1433 os << "Use DNA with multiple scattering of e- "
1434 << fCParameters->DNAElectronMsc() << "\n";
1435 os << "Use DNA e- solvation model type "
1436 << fCParameters->DNAeSolvationSubType() << "\n";
1437 os << "=======================================================================" << G4endl;
1438 }
1439 os.precision(prec);
1440}
1441
1443{
1444 if(fIsPrinted) return;
1445
1446#ifdef G4MULTITHREADED
1447 G4MUTEXLOCK(&emParametersMutex);
1448#endif
1450#ifdef G4MULTITHREADED
1451 G4MUTEXUNLOCK(&emParametersMutex);
1452#endif
1453}
1454
1455std::ostream& operator<< (std::ostream& os, const G4EmParameters& par)
1456{
1457 par.StreamInfo(os);
1458 return os;
1459}
1460
1461G4bool G4EmParameters::IsLocked() const
1462{
1463 return (!G4Threading::IsMasterThread() ||
1464 (fStateManager->GetCurrentState() != G4State_PreInit &&
1465 fStateManager->GetCurrentState() != G4State_Init &&
1466 fStateManager->GetCurrentState() != G4State_Idle));
1467}
1468
1469//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
@ G4State_Init
@ G4State_Idle
@ G4State_PreInit
G4DNAModelSubType
G4EmFluoDirectory
@ fluoBearden
@ fluoXDB_EADL
@ fluoANSTO
std::ostream & operator<<(std::ostream &os, const G4EmParameters &par)
G4EmFluctuationType
@ fUniversalFluctuation
@ fUrbanFluctuation
@ fDummyFluctuation
G4eSingleScatteringType
@ fWVI
G4TransportationWithMscType
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4MscStepLimitType
@ fMinimal
@ fUseSafety
G4NuclearFormfactorType
#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
long G4long
Definition: G4Types.hh:87
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)
void SetXDB_EADLFluoDir(G4bool val)
G4bool DNAStationary() const
void SetDeexcitationIgnoreCut(G4bool val)
const std::vector< G4String > & TypesDNA() const
void SetDNAElectronMsc(G4bool val)
void SetFluoDirectory(G4EmFluoDirectory val)
const G4String & LivermoreDataDir()
void SetFluo(G4bool val)
const std::vector< G4String > & RegionsMicroElec() const
G4bool DeexcitationIgnoreCut() const
const G4String & PIXECrossSectionModel()
G4EmFluoDirectory FluoDirectory() const
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
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)
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
void SetNumberOfBinsPerDecade(G4int val)
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 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
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 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)
G4bool ANSTOFluoDir()
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
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)
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)
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)
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
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
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)
G4bool GeneralProcessActive() const
static G4NistManager * Instance()
const G4ApplicationState & GetCurrentState() const
static G4StateManager * GetStateManager()
G4bool IsMasterThread()
Definition: G4Threading.cc:124
int G4lrint(double ad)
Definition: templates.hh:134