Geant4 11.2.2
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 nForFreeVector = 2;
168 tripletConv = 0;
169
170 fTransportationWithMsc = G4TransportationWithMscType::fDisabled;
171 mscStepLimit = fUseSafety;
172 mscStepLimitMuHad = fMinimal;
173 nucFormfactor = fExponentialNF;
174 fSStype = fWVI;
175 fFluct = fUniversalFluctuation;
176
177 fDirLEDATA = G4String(G4FindDataDir("G4LEDATA"));
178}
179
181{
182 if(IsLocked()) { return; }
183 lossFluctuation = val;
184}
185
187{
188 return lossFluctuation;
189}
190
192{
193 if(IsLocked()) { return; }
194 buildCSDARange = val;
195}
196
198{
199 return buildCSDARange;
200}
201
203{
204 if(IsLocked()) { return; }
205 flagLPM = val;
206}
207
209{
210 return flagLPM;
211}
212
214{
215 if(IsLocked()) { return; }
216 cutAsFinalRange = val;
217}
218
220{
221 return cutAsFinalRange;
222}
223
225{
226 if(IsLocked()) { return; }
227 applyCuts = val;
228}
229
231{
232 return applyCuts;
233}
234
236{
237 if(IsLocked()) { return; }
238 fCParameters->SetFluo(val);
239}
240
242{
243 return fCParameters->Fluo();
244}
245
247{
248 return fCParameters->FluoDirectory();
249}
250
252{
253 if(IsLocked()) { return; }
254 fCParameters->SetFluoDirectory(val);
255}
256
258{
259 if(IsLocked()) { return; }
260 fCParameters->SetBeardenFluoDir(val);
261}
262
264{
265 if(IsLocked()) { return; }
266 fCParameters->SetANSTOFluoDir(val);
267}
268
270{
271 if(IsLocked()) { return; }
272 fCParameters->SetXDB_EADLFluoDir(val);
273}
274
276{
277 if(IsLocked()) { return; }
278 fCParameters->SetAuger(val);
279}
280
282{
283 auto dir = fCParameters->FluoDirectory();
284 return (dir == fluoBearden);
285}
286
288{
289 auto dir = fCParameters->FluoDirectory();
290 return (dir == fluoANSTO);
291}
292
294{
295 return fCParameters->Auger();
296}
297
299{
300 if(IsLocked()) { return; }
301 fCParameters->SetPixe(val);
302}
303
305{
306 return fCParameters->Pixe();
307}
308
310{
311 if(IsLocked()) { return; }
312 fCParameters->SetDeexcitationIgnoreCut(val);
313}
314
316{
317 return fCParameters->DeexcitationIgnoreCut();
318}
319
321{
322 if(IsLocked()) { return; }
323 lateralDisplacement = val;
324}
325
327{
328 return lateralDisplacement;
329}
330
332{
333 if(IsLocked()) { return; }
334 lateralDisplacementAlg96 = val;
335}
336
338{
339 return lateralDisplacementAlg96;
340}
341
343{
344 if(IsLocked()) { return; }
345 muhadLateralDisplacement = val;
346}
347
349{
350 return muhadLateralDisplacement;
351}
352
354{
355 if(IsLocked()) { return; }
356 useAngGeneratorForIonisation = val;
357}
358
360{
361 return useAngGeneratorForIonisation;
362}
363
365{
366 if(IsLocked()) { return; }
367 useMottCorrection = val;
368}
369
371{
372 return useMottCorrection;
373}
374
376{
377 if(IsLocked()) { return; }
378 integral = val;
379}
380
382{
383 return integral;
384}
385
387{
388 if(IsLocked()) { return; }
389 fPolarisation = val;
390}
391
393{
394 return fPolarisation;
395}
396
398{
399 if(IsLocked()) { return; }
400 birks = val;
401 if(birks && nullptr == emSaturation) { emSaturation = new G4EmSaturation(1); }
402}
403
405{
406 return birks;
407}
408
410{
411 if(IsLocked()) { return; }
412 fICRU90 = val;
413}
414
416{
417 return fICRU90;
418}
419
421{
422 if(IsLocked()) { return; }
423 fCParameters->SetDNAFast(val);
424 if(val) { ActivateDNA(); }
425}
426
428{
429 return fCParameters->DNAFast();
430}
431
433{
434 if(IsLocked()) { return; }
435 fCParameters->SetDNAStationary(val);
436 if(val) { ActivateDNA(); }
437}
438
440{
441 return fCParameters->DNAStationary();
442}
443
445{
446 if(IsLocked()) { return; }
447 fCParameters->SetDNAElectronMsc(val);
448 if(val) { ActivateDNA(); }
449}
450
452{
453 return fCParameters->DNAElectronMsc();
454}
455
457{
458 if(IsLocked()) { return; }
459 gener = val;
460}
461
463{
464 return gener;
465}
466
468{
469 if(IsLocked()) { return; }
470 birks = (nullptr != ptr);
471 if(emSaturation != ptr) {
472 delete emSaturation;
473 emSaturation = ptr;
474 }
475}
476
478{
479 return fMuDataFromFile;
480}
481
483{
484 fMuDataFromFile = v;
485}
486
488{
489 if(IsLocked()) { return; }
490 onIsolated = val;
491}
492
494{
495 return onIsolated;
496}
497
499{
500 if(IsLocked()) { return; }
501 fSamplingTable = val;
502}
503
505{
506 return fSamplingTable;
507}
508
510{
511 return fPEKShell;
512}
513
515{
516 if(IsLocked()) { return; }
517 fPEKShell = v;
518}
519
521{
522 return fMscPosiCorr;
523}
524
526{
527 if(IsLocked()) { return; }
528 fMscPosiCorr = v;
529}
530
532{
533 if(IsLocked()) { return; }
534 fDNA = true;
535}
536
538{
539 fIsPrinted = val;
540}
541
543{
544 return fIsPrinted;
545}
546
548{
549 if(nullptr == emSaturation) {
550#ifdef G4MULTITHREADED
551 G4MUTEXLOCK(&emParametersMutex);
552 if(nullptr == emSaturation) {
553#endif
554 emSaturation = new G4EmSaturation(1);
555#ifdef G4MULTITHREADED
556 }
557 G4MUTEXUNLOCK(&emParametersMutex);
558#endif
559 }
560 birks = true;
561 return emSaturation;
562}
563
565{
566 if(IsLocked()) { return; }
567 if(val > 1.e-3*CLHEP::eV && val < maxKinEnergy) {
568 minKinEnergy = val;
569 } else {
571 ed << "Value of MinKinEnergy - is out of range: " << val/CLHEP::MeV
572 << " MeV is ignored";
573 PrintWarning(ed);
574 }
575}
576
578{
579 return minKinEnergy;
580}
581
583{
584 if(IsLocked()) { return; }
585 if(val > std::max(minKinEnergy,9.99*CLHEP::MeV) && val < 1.e+7*CLHEP::TeV) {
586 maxKinEnergy = val;
587 } else {
589 ed << "Value of MaxKinEnergy is out of range: "
590 << val/CLHEP::GeV
591 << " GeV is ignored; allowed range 10 MeV - 1.e+7 TeV";
592 PrintWarning(ed);
593 }
594}
595
597{
598 return maxKinEnergy;
599}
600
602{
603 if(IsLocked()) { return; }
604 if(val > minKinEnergy && val <= 100*CLHEP::TeV) {
605 maxKinEnergyCSDA = val;
606 } else {
608 ed << "Value of MaxKinEnergyCSDA is out of range: "
609 << val/CLHEP::GeV << " GeV is ignored; allowed range "
610 << minKinEnergy << " MeV - 100 TeV";
611 PrintWarning(ed);
612 }
613}
614
616{
617 return maxKinEnergyCSDA;
618}
619
621{
622 if(IsLocked()) { return; }
623 if(val >= 0.0) { lowestElectronEnergy = val; }
624}
625
627{
628 return lowestElectronEnergy;
629}
630
632{
633 if(IsLocked()) { return; }
634 if(val >= 0.0) { lowestMuHadEnergy = val; }
635}
636
638{
639 return lowestMuHadEnergy;
640}
641
643{
644 if(IsLocked()) { return; }
645 if(val > 0.0) { lowestTripletEnergy = val; }
646}
647
649{
650 return lowestTripletEnergy;
651}
652
654{
655 if(IsLocked()) { return; }
656 if(val >= 0.0) { maxNIELEnergy = val; }
657}
658
660{
661 return maxNIELEnergy;
662}
663
665{
666 if(IsLocked()) { return; }
667 if(val > 0.0) { max5DEnergyForMuPair = val; }
668}
669
671{
672 return max5DEnergyForMuPair;
673}
674
676{
677 if(IsLocked()) { return; }
678 if(val > 0.0 && val < 0.5) {
679 linLossLimit = val;
680 } else {
682 ed << "Value of linLossLimit is out of range: " << val
683 << " is ignored";
684 PrintWarning(ed);
685 }
686}
687
689{
690 return linLossLimit;
691}
692
694{
695 if(IsLocked()) { return; }
696 if(val > 0.0) {
697 bremsTh = val;
698 } else {
700 ed << "Value of bremsstrahlung threshold is out of range: "
701 << val/GeV << " GeV is ignored";
702 PrintWarning(ed);
703 }
704}
705
707{
708 return bremsTh;
709}
710
712{
713 if(IsLocked()) { return; }
714 if(val > 0.0) {
715 bremsMuHadTh = val;
716 } else {
718 ed << "Value of bremsstrahlung threshold is out of range: "
719 << val/GeV << " GeV is ignored";
720 PrintWarning(ed);
721 }
722}
723
725{
726 return bremsMuHadTh;
727}
728
730{
731 if(IsLocked()) { return; }
732 if(val > 0.0 && val < 1.0) {
733 lambdaFactor = val;
734 } else {
736 ed << "Value of lambda factor is out of range: " << val
737 << " is ignored";
738 PrintWarning(ed);
739 }
740}
741
743{
744 return lambdaFactor;
745}
746
748{
749 if(IsLocked()) { return; }
750 if(val > 0.0) {
751 factorForAngleLimit = val;
752 } else {
754 ed << "Value of factor for enegry limit is out of range: "
755 << val << " is ignored";
756 PrintWarning(ed);
757 }
758}
759
761{
762 return factorForAngleLimit;
763}
764
766{
767 if(IsLocked()) { return; }
768 if(val >= 0.0 && val <= pi) {
769 thetaLimit = val;
770 } else {
772 ed << "Value of polar angle limit is out of range: "
773 << val << " is ignored";
774 PrintWarning(ed);
775 }
776}
777
779{
780 return thetaLimit;
781}
782
784{
785 if(IsLocked()) { return; }
786 if(val >= 0.0) {
787 energyLimit = val;
788 } else {
790 ed << "Value of msc energy limit is out of range: "
791 << val << " is ignored";
792 PrintWarning(ed);
793 }
794}
795
797{
798 return energyLimit;
799}
800
802{
803 if(IsLocked()) { return; }
804 if(val > 0.0 && val < 1.0) {
805 rangeFactor = val;
806 } else {
808 ed << "Value of rangeFactor is out of range: "
809 << val << " is ignored";
810 PrintWarning(ed);
811 }
812}
813
815{
816 return rangeFactor;
817}
818
820{
821 if(IsLocked()) { return; }
822 if(val > 0.0 && val < 1.0) {
823 rangeFactorMuHad = val;
824 } else {
826 ed << "Value of rangeFactorMuHad is out of range: "
827 << val << " is ignored";
828 PrintWarning(ed);
829 }
830}
831
833{
834 return rangeFactorMuHad;
835}
836
838{
839 if(IsLocked()) { return; }
840 if(val >= 1.0) {
841 geomFactor = val;
842 } else {
844 ed << "Value of geomFactor is out of range: "
845 << val << " is ignored";
846 PrintWarning(ed);
847 }
848}
849
851{
852 return geomFactor;
853}
854
856{
857 if(IsLocked()) { return; }
858 if(val >= 0.1) {
859 safetyFactor = val;
860 } else {
862 ed << "Value of safetyFactor is out of range: "
863 << val << " is ignored";
864 PrintWarning(ed);
865 }
866}
867
869{
870 return safetyFactor;
871}
872
874{
875 if(IsLocked()) { return; }
876 if(val >= 0.0) {
877 lambdaLimit = val;
878 } else {
880 ed << "Value of lambdaLimit is out of range: "
881 << val << " is ignored";
882 PrintWarning(ed);
883 }
884}
885
887{
888 return lambdaLimit;
889}
890
892{
893 if(IsLocked()) { return; }
894 if(val >= 1.0) {
895 skin = val;
896 } else {
898 ed << "Value of skin is out of range: "
899 << val << " is ignored";
900 PrintWarning(ed);
901 }
902}
903
905{
906 return skin;
907}
908
910{
911 if(IsLocked()) { return; }
912 if(val > 0.0) {
913 factorScreen = val;
914 } else {
916 ed << "Value of factorScreen is out of range: "
917 << val << " is ignored";
918 PrintWarning(ed);
919 }
920}
921
923{
924 return factorScreen;
925}
926
928{
929 if(IsLocked()) { return; }
930 fBParameters->SetStepFunction(v1, v2);
931}
932
934{
935 if(IsLocked()) { return; }
936 fBParameters->SetStepFunctionMuHad(v1, v2);
937}
938
940{
941 if(IsLocked()) { return; }
942 fBParameters->SetStepFunctionLightIons(v1, v2);
943}
944
946{
947 if(IsLocked()) { return; }
948 fBParameters->SetStepFunctionIons(v1, v2);
949}
950
952{
953 fBParameters->FillStepFunction(part, proc);
954}
955
957{
958 return nbinsPerDecade*G4lrint(std::log10(maxKinEnergy/minKinEnergy));
959}
960
962{
963 if(IsLocked()) { return; }
964 if(val >= 5 && val < 1000000) {
965 nbinsPerDecade = val;
966 } else {
968 ed << "Value of number of bins per decade is out of range: "
969 << val << " is ignored";
970 PrintWarning(ed);
971 }
972}
973
975{
976 return nbinsPerDecade;
977}
978
980{
981 if(IsLocked()) { return; }
982 verbose = val;
983 workerVerbose = std::min(workerVerbose, verbose);
984}
985
987{
988 return verbose;
989}
990
992{
993 if(IsLocked()) { return; }
994 workerVerbose = val;
995}
996
998{
999 return workerVerbose;
1000}
1001
1003{
1004 if(IsLocked()) { return; }
1005 nForFreeVector = val;
1006}
1007
1009{
1010 return nForFreeVector;
1011}
1012
1014{
1015 if(IsLocked()) { return; }
1016 fTransportationWithMsc = val;
1017}
1018
1020{
1021 return fTransportationWithMsc;
1022}
1023
1025{
1026 if(IsLocked()) { return; }
1027 fFluct = val;
1028}
1029
1031{
1032 return fFluct;
1033}
1034
1036{
1037 if(IsLocked()) { return; }
1038 mscStepLimit = val;
1039}
1040
1042{
1043 return mscStepLimit;
1044}
1045
1047{
1048 if(IsLocked()) { return; }
1049 mscStepLimitMuHad = val;
1050}
1051
1053{
1054 return mscStepLimitMuHad;
1055}
1056
1058{
1059 if(IsLocked()) { return; }
1060 fSStype = val;
1061}
1062
1064{
1065 return fSStype;
1066}
1067
1068void
1070{
1071 if(IsLocked()) { return; }
1072 nucFormfactor = val;
1073}
1074
1076{
1077 return nucFormfactor;
1078}
1079
1081{
1082 if(IsLocked()) { return; }
1083 fCParameters->SetDNAeSolvationSubType(val);
1084 ActivateDNA();
1085}
1086
1088{
1089 return fCParameters->DNAeSolvationSubType();
1090}
1091
1093{
1094 if(IsLocked()) { return; }
1095 tripletConv = val;
1096}
1097
1099{
1100 return tripletConv;
1101}
1102
1104{
1105 if(IsLocked()) { return; }
1106 fCParameters->SetPIXECrossSectionModel(sss);
1107}
1108
1110{
1111 return fCParameters->PIXECrossSectionModel();
1112}
1113
1115{
1116 if(IsLocked()) { return; }
1117 fCParameters->SetPIXEElectronCrossSectionModel(sss);
1118}
1119
1124
1126{
1127 if(IsLocked()) { return; }
1128 fCParameters->SetLivermoreDataDir(sss);
1129}
1130
1132{
1133 return fCParameters->LivermoreDataDir();
1134}
1135
1136void G4EmParameters::PrintWarning(G4ExceptionDescription& ed) const
1137{
1138 G4Exception("G4EmParameters", "em0044", JustWarning, ed);
1139}
1140
1142 const G4String& region,
1143 const G4String& type)
1144{
1145 if(IsLocked()) { return; }
1146 fBParameters->AddPAIModel(particle, region, type);
1147}
1148
1149const std::vector<G4String>& G4EmParameters::ParticlesPAI() const
1150{
1151 return fBParameters->ParticlesPAI();
1152}
1153
1154const std::vector<G4String>& G4EmParameters::RegionsPAI() const
1155{
1156 return fBParameters->RegionsPAI();
1157}
1158
1159const std::vector<G4String>& G4EmParameters::TypesPAI() const
1160{
1161 return fBParameters->TypesPAI();
1162}
1163
1165{
1166 if(IsLocked()) { return; }
1167 fCParameters->AddMicroElec(region);
1168}
1169
1170const std::vector<G4String>& G4EmParameters::RegionsMicroElec() const
1171{
1172 return fCParameters->RegionsMicroElec();
1173}
1174
1175void G4EmParameters::AddDNA(const G4String& region, const G4String& type)
1176{
1177 if(IsLocked()) { return; }
1178 fCParameters->AddDNA(region, type);
1179 ActivateDNA();
1180}
1181
1182const std::vector<G4String>& G4EmParameters::RegionsDNA() const
1183{
1184 return fCParameters->RegionsDNA();
1185}
1186
1187const std::vector<G4String>& G4EmParameters::TypesDNA() const
1188{
1189 return fCParameters->TypesDNA();
1190}
1191
1192void G4EmParameters::AddPhysics(const G4String& region, const G4String& type)
1193{
1194 if(IsLocked()) { return; }
1195 fBParameters->AddPhysics(region, type);
1196}
1197
1198const std::vector<G4String>& G4EmParameters::RegionsPhysics() const
1199{
1200 return fBParameters->RegionsPhysics();
1201}
1202
1203const std::vector<G4String>& G4EmParameters::TypesPhysics() const
1204{
1205 return fBParameters->TypesPhysics();
1206}
1207
1209{
1210 if(IsLocked()) { return; }
1211 fBParameters->SetSubCutRegion(region);
1212}
1213
1214void
1216 G4bool aauger, G4bool apixe)
1217{
1218 if(IsLocked()) { return; }
1219 fCParameters->SetDeexActiveRegion(region, adeex, aauger, apixe);
1220}
1221
1222void
1224 G4double val, G4bool wflag)
1225{
1226 if(IsLocked()) { return; }
1227 fBParameters->SetProcessBiasingFactor(procname, val, wflag);
1228}
1229
1230void
1232 const G4String& region,
1233 G4double length,
1234 G4bool wflag)
1235{
1236 if(IsLocked() && !gener) { return; }
1237 fBParameters->ActivateForcedInteraction(procname, region, length, wflag);
1238}
1239
1240void
1242 const G4String& region,
1243 G4double factor,
1244 G4double energyLim)
1245{
1246 if(IsLocked()) { return; }
1247 fBParameters->ActivateSecondaryBiasing(procname, region, factor, energyLim);
1248}
1249
1251{
1252 fBParameters->DefineRegParamForLoss(ptr);
1253}
1254
1256{
1257 fBParameters->DefineRegParamForEM(ptr);
1258}
1259
1261{
1262 return fBParameters->QuantumEntanglement();
1263}
1264
1266{
1267 if(IsLocked()) { return; }
1268 fBParameters->SetQuantumEntanglement(v);
1269}
1270
1274
1276{
1277 if(IsLocked()) { return; }
1278 fBParameters->SetDirectionalSplitting(v);
1279}
1280
1282{
1283 if(IsLocked()) { return; }
1284 fBParameters->SetDirectionalSplittingTarget(v);
1285}
1286
1291
1293{
1294 if(IsLocked()) { return; }
1295 fBParameters->SetDirectionalSplittingRadius(r);
1296}
1297
1302
1304{
1305 fCParameters->DefineRegParamForDeex(ptr);
1306}
1307
1309{
1310 return fDirLEDATA;
1311}
1312
1313void G4EmParameters::StreamInfo(std::ostream& os) const
1314{
1315 G4long prec = os.precision(5);
1316 os << "=======================================================================" << "\n";
1317 os << "====== Electromagnetic Physics Parameters ========" << "\n";
1318 os << "=======================================================================" << "\n";
1319 os << "LPM effect enabled " <<flagLPM << "\n";
1320 os << "Enable creation and use of sampling tables " <<fSamplingTable << "\n";
1321 os << "Apply cuts on all EM processes " <<applyCuts << "\n";
1322 const char* transportationWithMsc = "Disabled";
1323 if(fTransportationWithMsc == G4TransportationWithMscType::fEnabled) {
1324 transportationWithMsc = "Enabled";
1325 } else if (fTransportationWithMsc == G4TransportationWithMscType::fMultipleSteps) {
1326 transportationWithMsc = "MultipleSteps";
1327 }
1328 os << "Use combined TransportationWithMsc " <<transportationWithMsc << "\n";
1329 os << "Use general process " <<gener << "\n";
1330 os << "Enable linear polarisation for gamma " <<fPolarisation << "\n";
1331 os << "Enable photoeffect sampling below K-shell " <<fPEKShell << "\n";
1332 os << "Enable sampling of quantum entanglement "
1333 <<fBParameters->QuantumEntanglement() << "\n";
1334 os << "X-section factor for integral approach " <<lambdaFactor << "\n";
1335 os << "Min kinetic energy for tables "
1336 <<G4BestUnit(minKinEnergy,"Energy") << "\n";
1337 os << "Max kinetic energy for tables "
1338 <<G4BestUnit(maxKinEnergy,"Energy") << "\n";
1339 os << "Number of bins per decade of a table " <<nbinsPerDecade << "\n";
1340 os << "Verbose level " <<verbose << "\n";
1341 os << "Verbose level for worker thread " <<workerVerbose << "\n";
1342 os << "Bremsstrahlung energy threshold above which \n"
1343 << " primary e+- is added to the list of secondary "
1344 <<G4BestUnit(bremsTh,"Energy") << "\n";
1345 os << "Bremsstrahlung energy threshold above which primary\n"
1346 << " muon/hadron is added to the list of secondary "
1347 <<G4BestUnit(bremsMuHadTh,"Energy") << "\n";
1348 os << "Lowest triplet kinetic energy "
1349 <<G4BestUnit(lowestTripletEnergy,"Energy") << "\n";
1350 os << "Enable sampling of gamma linear polarisation " <<fPolarisation << "\n";
1351 os << "5D gamma conversion model type " <<tripletConv << "\n";
1352 os << "5D gamma conversion model on isolated ion " <<onIsolated << "\n";
1353 if(max5DEnergyForMuPair>0.0) {
1354 os << "5D gamma conversion limit for muon pair "
1355 << max5DEnergyForMuPair/CLHEP::GeV << " GeV\n";
1356 }
1357 os << "Livermore data directory "
1358 << fCParameters->LivermoreDataDir() << "\n";
1359
1360 os << "=======================================================================" << "\n";
1361 os << "====== Ionisation Parameters ========" << "\n";
1362 os << "=======================================================================" << "\n";
1363 os << "Step function for e+- "
1364 <<"("<<fBParameters->GetStepFunctionP1() << ", "
1365 << fBParameters->GetStepFunctionP2()/CLHEP::mm << " mm)\n";
1366 os << "Step function for muons/hadrons "
1367 <<"("<<fBParameters->GetStepFunctionMuHadP1() << ", "
1368 << fBParameters->GetStepFunctionMuHadP2()/CLHEP::mm << " mm)\n";
1369 os << "Step function for light ions "
1370 <<"("<<fBParameters->GetStepFunctionLightIonsP1() << ", "
1371 << fBParameters->GetStepFunctionLightIonsP2()/CLHEP::mm << " mm)\n";
1372 os << "Step function for general ions "
1373 <<"("<<fBParameters->GetStepFunctionIonsP1() << ", "
1374 << fBParameters->GetStepFunctionIonsP2()/CLHEP::mm << " mm)\n";
1375 os << "Lowest e+e- kinetic energy "
1376 <<G4BestUnit(lowestElectronEnergy,"Energy") << "\n";
1377 os << "Lowest muon/hadron kinetic energy "
1378 <<G4BestUnit(lowestMuHadEnergy,"Energy") << "\n";
1379 os << "Use ICRU90 data " << fICRU90 << "\n";
1380 os << "Fluctuations of dE/dx are enabled " <<lossFluctuation << "\n";
1381 G4String namef = "Universal";
1382 if(fFluct == fUrbanFluctuation) { namef = "Urban"; }
1383 else if(fFluct == fDummyFluctuation) { namef = "Dummy"; }
1384 os << "Type of fluctuation model for leptons and hadrons " << namef << "\n";
1385 os << "Use built-in Birks satuaration " << birks << "\n";
1386 os << "Build CSDA range enabled " <<buildCSDARange << "\n";
1387 os << "Use cut as a final range enabled " <<cutAsFinalRange << "\n";
1388 os << "Enable angular generator interface "
1389 <<useAngGeneratorForIonisation << "\n";
1390 os << "Max kinetic energy for CSDA tables "
1391 <<G4BestUnit(maxKinEnergyCSDA,"Energy") << "\n";
1392 os << "Max kinetic energy for NIEL computation "
1393 <<G4BestUnit(maxNIELEnergy,"Energy") << "\n";
1394 os << "Linear loss limit " <<linLossLimit << "\n";
1395 os << "Read data from file for e+e- pair production by mu " <<fMuDataFromFile << "\n";
1396
1397 os << "=======================================================================" << "\n";
1398 os << "====== Multiple Scattering Parameters ========" << "\n";
1399 os << "=======================================================================" << "\n";
1400 os << "Type of msc step limit algorithm for e+- " <<mscStepLimit << "\n";
1401 os << "Type of msc step limit algorithm for muons/hadrons " <<mscStepLimitMuHad << "\n";
1402 os << "Msc lateral displacement for e+- enabled " <<lateralDisplacement << "\n";
1403 os << "Msc lateral displacement for muons and hadrons " <<muhadLateralDisplacement << "\n";
1404 os << "Urban msc model lateral displacement alg96 " <<lateralDisplacementAlg96 << "\n";
1405 os << "Range factor for msc step limit for e+- " <<rangeFactor << "\n";
1406 os << "Range factor for msc step limit for muons/hadrons " <<rangeFactorMuHad << "\n";
1407 os << "Geometry factor for msc step limitation of e+- " <<geomFactor << "\n";
1408 os << "Safety factor for msc step limit for e+- " <<safetyFactor << "\n";
1409 os << "Skin parameter for msc step limitation of e+- " <<skin << "\n";
1410 os << "Lambda limit for msc step limit for e+- " <<lambdaLimit/CLHEP::mm << " mm\n";
1411 os << "Use Mott correction for e- scattering " << useMottCorrection << "\n";
1412 os << "Factor used for dynamic computation of angular \n"
1413 << " limit between single and multiple scattering " << factorForAngleLimit << "\n";
1414 os << "Fixed angular limit between single \n"
1415 << " and multiple scattering "
1416 << thetaLimit/CLHEP::rad << " rad\n";
1417 os << "Upper energy limit for e+- multiple scattering "
1418 << energyLimit/CLHEP::MeV << " MeV\n";
1419 os << "Type of electron single scattering model " <<fSStype << "\n";
1420 os << "Type of nuclear form-factor " <<nucFormfactor << "\n";
1421 os << "Screening factor " <<factorScreen << "\n";
1422 os << "=======================================================================" << "\n";
1423
1424 if(fCParameters->Fluo()) {
1425 os << "====== Atomic Deexcitation Parameters ========" << "\n";
1426 os << "=======================================================================" << "\n";
1427 os << "Fluorescence enabled " <<fCParameters->Fluo() << "\n";
1428 G4String named = "fluor";
1430 if(fdir == fluoBearden) { named = "fluor_Bearden"; }
1431 else if(fdir == fluoANSTO) { named = "fluor_ANSTO"; }
1432 else if(fdir == fluoXDB_EADL) { named = "fluor_XDB_EADL"; }
1433 os << "Directory in G4LEDATA for fluorescence data files " << named << "\n";
1434 os << "Auger electron cascade enabled "
1435 <<fCParameters->Auger() << "\n";
1436 os << "PIXE atomic de-excitation enabled " <<fCParameters->Pixe() << "\n";
1437 os << "De-excitation module ignores cuts "
1438 <<fCParameters->DeexcitationIgnoreCut() << "\n";
1439 os << "Type of PIXE cross section for hadrons "
1440 <<fCParameters->PIXECrossSectionModel() << "\n";
1441 os << "Type of PIXE cross section for e+- "
1442 <<fCParameters->PIXEElectronCrossSectionModel() << "\n";
1443 os << "=======================================================================" << "\n";
1444 }
1445 if(fDNA) {
1446 os << "====== DNA Physics Parameters ========" << "\n";
1447 os << "=======================================================================" << "\n";
1448 os << "Use fast sampling in DNA models "
1449 << fCParameters->DNAFast() << "\n";
1450 os << "Use Stationary option in DNA models "
1451 << fCParameters->DNAStationary() << "\n";
1452 os << "Use DNA with multiple scattering of e- "
1453 << fCParameters->DNAElectronMsc() << "\n";
1454 os << "Use DNA e- solvation model type "
1455 << fCParameters->DNAeSolvationSubType() << "\n";
1456 os << "=======================================================================" << G4endl;
1457 }
1458 os.precision(prec);
1459}
1460
1462{
1463 if(fIsPrinted) return;
1464
1465#ifdef G4MULTITHREADED
1466 G4MUTEXLOCK(&emParametersMutex);
1467#endif
1469#ifdef G4MULTITHREADED
1470 G4MUTEXUNLOCK(&emParametersMutex);
1471#endif
1472}
1473
1474std::ostream& operator<< (std::ostream& os, const G4EmParameters& par)
1475{
1476 par.StreamInfo(os);
1477 return os;
1478}
1479
1480G4bool G4EmParameters::IsLocked() const
1481{
1482 return (!G4Threading::IsMasterThread() ||
1483 (fStateManager->GetCurrentState() != G4State_PreInit &&
1484 fStateManager->GetCurrentState() != G4State_Init &&
1485 fStateManager->GetCurrentState() != G4State_Idle));
1486}
1487
1488//....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
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
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
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
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
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
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)
G4bool GeneralProcessActive() const
static G4NistManager * Instance()
const G4ApplicationState & GetCurrentState() const
static G4StateManager * GetStateManager()
G4bool IsMasterThread()
int G4lrint(double ad)
Definition templates.hh:134