Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmParameters.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// -------------------------------------------------------------------
27//
28// GEANT4 Class file
29//
30// File name: G4EmParameters
31//
32// Author: Vladimir Ivanchenko
33//
34// Creation date: 18.05.2013
35//
36// Modifications:
37//
38// -------------------------------------------------------------------
39//
40//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
42
43#include "G4EmParameters.hh"
45#include "G4UnitsTable.hh"
46#include "G4SystemOfUnits.hh"
47#include "G4VEmProcess.hh"
51#include "G4EmLowEParameters.hh"
53#include "G4NistManager.hh"
54#include "G4RegionStore.hh"
55#include "G4Region.hh"
56#include "G4ApplicationState.hh"
57#include "G4StateManager.hh"
58
59G4EmParameters* G4EmParameters::theInstance = nullptr;
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
96G4EmParameters::G4EmParameters()
97{
99 theMessenger = new G4EmParametersMessenger(this);
100 Initialise();
101
102 fBParameters = new G4EmExtraParameters();
103 fCParameters = new G4EmLowEParameters();
104
105 fStateManager = G4StateManager::GetStateManager();
106 emSaturation = nullptr;
107}
108
110{
111 if(!IsLocked()) {
112 Initialise();
113 fBParameters->Initialise();
114 fCParameters->Initialise();
115 }
116}
117
118void G4EmParameters::Initialise()
119{
120 lossFluctuation = true;
121 buildCSDARange = false;
122 flagLPM = true;
123 spline = true;
124 cutAsFinalRange = false;
125 applyCuts = false;
126 lateralDisplacement = true;
127 lateralDisplacementAlg96 = true;
128 muhadLateralDisplacement = false;
129 useAngGeneratorForIonisation = false;
130 useMottCorrection = false;
131 integral = true;
132 birks = false;
133 fICRU90 = false;
134 gener = false;
135 onIsolated = false;
136 fSamplingTable = false;
137 fPolarisation = false;
138 fMuDataFromFile = false;
139 fDNA = false;
140
141 minSubRange = 1.0;
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 nbins = 84;
165 nbinsPerDecade = 7;
166 verbose = 1;
167 workerVerbose = 0;
168 tripletConv = 0;
169
170 mscStepLimit = fUseSafety;
171 mscStepLimitMuHad = fMinimal;
172 nucFormfactor = fExponentialNF;
173 fSStype = fWVI;
174}
175
177{
178 if(IsLocked()) { return; }
179 lossFluctuation = val;
180}
181
183{
184 return lossFluctuation;
185}
186
188{
189 if(IsLocked()) { return; }
190 buildCSDARange = val;
191}
192
194{
195 return buildCSDARange;
196}
197
199{
200 if(IsLocked()) { return; }
201 flagLPM = val;
202}
203
205{
206 return flagLPM;
207}
208
210{
211 if(IsLocked()) { return; }
212 spline = val;
213}
214
216{
217 return spline;
218}
219
221{
222 if(IsLocked()) { return; }
223 cutAsFinalRange = val;
224}
225
227{
228 return cutAsFinalRange;
229}
230
232{
233 if(IsLocked()) { return; }
234 applyCuts = val;
235}
236
238{
239 return applyCuts;
240}
241
243{
244 if(IsLocked()) { return; }
245 fCParameters->SetFluo(val);
246}
247
249{
250 return fCParameters->Fluo();
251}
252
254{
255 if(IsLocked()) { return; }
256 fCParameters->SetBeardenFluoDir(val);
257}
258
260{
261 return fCParameters->BeardenFluoDir();
262}
263
265{
266 if(IsLocked()) { return; }
267 fCParameters->SetAuger(val);
268}
269
271{
272 return fCParameters->Auger();
273}
274
276{
277 if(IsLocked()) { return; }
278 fCParameters->SetAuger(val);
279}
280
282{
283 return fCParameters->Auger();
284}
285
287{
288 if(IsLocked()) { return; }
289 fCParameters->SetPixe(val);
290}
291
293{
294 return fCParameters->Pixe();
295}
296
298{
299 if(IsLocked()) { return; }
300 fCParameters->SetDeexcitationIgnoreCut(val);
301}
302
304{
305 return fCParameters->DeexcitationIgnoreCut();
306}
307
309{
310 if(IsLocked()) { return; }
311 lateralDisplacement = val;
312}
313
315{
316 return lateralDisplacement;
317}
318
320{
321 if(IsLocked()) { return; }
322 lateralDisplacementAlg96 = val;
323}
324
326{
327 return lateralDisplacementAlg96;
328}
329
331{
332 if(IsLocked()) { return; }
333 muhadLateralDisplacement = val;
334}
335
337{
338 return muhadLateralDisplacement;
339}
340
342{
343 if(IsLocked()) { return; }
344 useAngGeneratorForIonisation = val;
345}
346
348{
349 return useAngGeneratorForIonisation;
350}
351
353{
354 if(IsLocked()) { return; }
355 useMottCorrection = val;
356}
357
359{
360 return useMottCorrection;
361}
362
364{
365 if(IsLocked()) { return; }
366 integral = val;
367}
368
370{
371 return integral;
372}
373
375{
376 if(IsLocked()) { return; }
377 fPolarisation = val;
378}
379
381{
382 return fPolarisation;
383}
384
386{
387 birks = val;
388#ifdef G4MULTITHREADED
389 G4MUTEXLOCK(&G4EmParameters::emParametersMutex);
390#endif
391 if(birks) {
392 if(!emSaturation) { emSaturation = new G4EmSaturation(1); }
393 emSaturation->InitialiseG4Saturation();
394 }
395#ifdef G4MULTITHREADED
396 G4MUTEXUNLOCK(&G4EmParameters::emParametersMutex);
397#endif
398}
399
401{
402 return birks;
403}
404
406{
407 if(IsLocked()) { return; }
408 fICRU90 = val;
409}
410
412{
413 return fICRU90;
414}
415
417{
418 if(IsLocked()) { return; }
419 fCParameters->SetDNAFast(val);
420 if(val) { ActivateDNA(); }
421}
422
424{
425 return fCParameters->DNAFast();
426}
427
429{
430 if(IsLocked()) { return; }
431 fCParameters->SetDNAStationary(val);
432 if(val) { ActivateDNA(); }
433}
434
436{
437 return fCParameters->DNAStationary();
438}
439
441{
442 if(IsLocked()) { return; }
443 fCParameters->SetDNAElectronMsc(val);
444 if(val) { ActivateDNA(); }
445}
446
448{
449 return fCParameters->DNAElectronMsc();
450}
451
453{
454 if(IsLocked()) { return; }
455 gener = val;
456 // if general interaction is enabled then sub-cutoff and
457 // force interaction options should be disabled
458 if(gener) { fBParameters->Initialise(); }
459}
460
462{
463 return gener;
464}
465
467{
468 if(emSaturation != ptr) {
469 delete emSaturation;
470 emSaturation = ptr;
471 SetBirksActive(true);
472 }
473}
474
476{
477 return fMuDataFromFile;
478}
479
481{
482 fMuDataFromFile = v;
483}
484
486{
487 if(IsLocked()) { return; }
488 onIsolated = val;
489}
490
492{
493 return onIsolated;
494}
495
497{
498 if(IsLocked()) { return; }
499 fSamplingTable = val;
500}
501
503{
504 return fSamplingTable;
505}
506
508{
509 if(IsLocked()) { return; }
510 fDNA = true;
511}
512
514{
515 if(!emSaturation) { SetBirksActive(true); }
516 return emSaturation;
517}
518
520{
521 if(IsLocked()) { return; }
522 if(val > 0.0 && val < 1.0) {
523 minSubRange = val;
524 } else {
526 ed << "Value of MinSubRange is out of range (0 - 1): " << val
527 << " is ignored";
528 PrintWarning(ed);
529 }
530}
531
533{
534 return minSubRange;
535}
536
538{
539 if(IsLocked()) { return; }
540 if(val > 1.e-3*eV && val < maxKinEnergy) {
541 minKinEnergy = val;
542 nbins = nbinsPerDecade*G4lrint(std::log10(maxKinEnergy/minKinEnergy));
543 } else {
545 ed << "Value of MinKinEnergy - is out of range: " << val/MeV
546 << " MeV is ignored";
547 PrintWarning(ed);
548 }
549}
550
552{
553 return minKinEnergy;
554}
555
557{
558 if(IsLocked()) { return; }
559 if(val > std::max(minKinEnergy,9.99*MeV) && val < 1.e+7*TeV) {
560 maxKinEnergy = val;
561 nbins = nbinsPerDecade*G4lrint(std::log10(maxKinEnergy/minKinEnergy));
562 } else {
564 ed << "Value of MaxKinEnergy is out of range: "
565 << val/GeV << " GeV is ignored; allowed range 10 MeV - 1.e+7 TeV";
566 PrintWarning(ed);
567 }
568}
569
571{
572 return maxKinEnergy;
573}
574
576{
577 if(IsLocked()) { return; }
578 if(val > minKinEnergy && val <= 100*TeV) {
579 maxKinEnergyCSDA = val;
580 } else {
582 ed << "Value of MaxKinEnergyCSDA is out of range: "
583 << val/GeV << " GeV is ignored; allowed range "
584 << minKinEnergy << " MeV - 100 TeV";
585 PrintWarning(ed);
586 }
587}
588
590{
591 return maxKinEnergyCSDA;
592}
593
595{
596 if(IsLocked()) { return; }
597 if(val >= 0.0) { lowestElectronEnergy = val; }
598}
599
601{
602 return lowestElectronEnergy;
603}
604
606{
607 if(IsLocked()) { return; }
608 if(val >= 0.0) { lowestMuHadEnergy = val; }
609}
610
612{
613 return lowestMuHadEnergy;
614}
615
617{
618 if(IsLocked()) { return; }
619 if(val > 0.0) { lowestTripletEnergy = val; }
620}
621
623{
624 return lowestTripletEnergy;
625}
626
628{
629 if(IsLocked()) { return; }
630 if(val >= 0.0) { maxNIELEnergy = val; }
631}
632
634{
635 return maxNIELEnergy;
636}
637
639{
640 if(IsLocked()) { return; }
641 if(val > 0.0) { max5DEnergyForMuPair = val; }
642}
643
645{
646 return max5DEnergyForMuPair;
647}
648
650{
651 if(IsLocked()) { return; }
652 if(val > 0.0 && val < 0.5) {
653 linLossLimit = val;
654 } else {
656 ed << "Value of linLossLimit is out of range: " << val
657 << " is ignored";
658 PrintWarning(ed);
659 }
660}
661
663{
664 return linLossLimit;
665}
666
668{
669 if(IsLocked()) { return; }
670 if(val > 0.0) {
671 bremsTh = val;
672 } else {
674 ed << "Value of bremsstrahlung threshold is out of range: "
675 << val/GeV << " GeV is ignored";
676 PrintWarning(ed);
677 }
678}
679
681{
682 return bremsTh;
683}
684
686{
687 if(IsLocked()) { return; }
688 if(val > 0.0) {
689 bremsMuHadTh = val;
690 } else {
692 ed << "Value of bremsstrahlung threshold is out of range: "
693 << val/GeV << " GeV is ignored";
694 PrintWarning(ed);
695 }
696}
697
699{
700 return bremsMuHadTh;
701}
702
704{
705 if(IsLocked()) { return; }
706 if(val > 0.0 && val < 1.0) {
707 lambdaFactor = val;
708 } else {
710 ed << "Value of lambda factor is out of range: " << val
711 << " is ignored";
712 PrintWarning(ed);
713 }
714}
715
717{
718 return lambdaFactor;
719}
720
722{
723 if(IsLocked()) { return; }
724 if(val > 0.0) {
725 factorForAngleLimit = val;
726 } else {
728 ed << "Value of factor for enegry limit is out of range: "
729 << val << " is ignored";
730 PrintWarning(ed);
731 }
732}
733
735{
736 return factorForAngleLimit;
737}
738
740{
741 if(IsLocked()) { return; }
742 if(val >= 0.0 && val <= pi) {
743 thetaLimit = val;
744 } else {
746 ed << "Value of polar angle limit is out of range: "
747 << val << " is ignored";
748 PrintWarning(ed);
749 }
750}
751
753{
754 return thetaLimit;
755}
756
758{
759 if(IsLocked()) { return; }
760 if(val >= 0.0) {
761 energyLimit = val;
762 } else {
764 ed << "Value of msc energy limit is out of range: "
765 << val << " is ignored";
766 PrintWarning(ed);
767 }
768}
769
771{
772 return energyLimit;
773}
774
776{
777 if(IsLocked()) { return; }
778 if(val > 0.0 && val < 1.0) {
779 rangeFactor = val;
780 } else {
782 ed << "Value of rangeFactor is out of range: "
783 << val << " is ignored";
784 PrintWarning(ed);
785 }
786}
787
789{
790 return rangeFactor;
791}
792
794{
795 if(IsLocked()) { return; }
796 if(val > 0.0 && val < 1.0) {
797 rangeFactorMuHad = val;
798 } else {
800 ed << "Value of rangeFactorMuHad is out of range: "
801 << val << " is ignored";
802 PrintWarning(ed);
803 }
804}
805
807{
808 return rangeFactorMuHad;
809}
810
812{
813 if(IsLocked()) { return; }
814 if(val >= 1.0) {
815 geomFactor = val;
816 } else {
818 ed << "Value of geomFactor is out of range: "
819 << val << " is ignored";
820 PrintWarning(ed);
821 }
822}
823
825{
826 return geomFactor;
827}
828
830{
831 if(IsLocked()) { return; }
832 if(val >= 0.1) {
833 safetyFactor = val;
834 } else {
836 ed << "Value of safetyFactor is out of range: "
837 << val << " is ignored";
838 PrintWarning(ed);
839 }
840}
841
843{
844 return safetyFactor;
845}
846
848{
849 if(IsLocked()) { return; }
850 if(val >= 0.0) {
851 lambdaLimit = val;
852 } else {
854 ed << "Value of lambdaLimit is out of range: "
855 << val << " is ignored";
856 PrintWarning(ed);
857 }
858}
859
861{
862 return lambdaLimit;
863}
864
866{
867 if(IsLocked()) { return; }
868 if(val >= 1.0) {
869 skin = val;
870 } else {
872 ed << "Value of skin is out of range: "
873 << val << " is ignored";
874 PrintWarning(ed);
875 }
876}
877
879{
880 return skin;
881}
882
884{
885 if(IsLocked()) { return; }
886 if(val > 0.0) {
887 factorScreen = val;
888 } else {
890 ed << "Value of factorScreen is out of range: "
891 << val << " is ignored";
892 PrintWarning(ed);
893 }
894}
895
897{
898 return factorScreen;
899}
900
902{
903 if(IsLocked()) { return; }
904 fBParameters->SetStepFunction(v1, v2);
905}
906
908{
909 if(IsLocked()) { return; }
910 fBParameters->SetStepFunctionMuHad(v1, v2);
911}
912
914{
915 if(IsLocked()) { return; }
916 fBParameters->SetStepFunctionLightIons(v1, v2);
917}
918
920{
921 if(IsLocked()) { return; }
922 fBParameters->SetStepFunctionIons(v1, v2);
923}
924
926{
927 fBParameters->FillStepFunction(part, proc);
928}
929
931{
932 if(IsLocked()) { return; }
933 if(val >= 5 && val < 10000000) {
934 nbins = val;
935 nbinsPerDecade = G4lrint(nbins/std::log10(maxKinEnergy/minKinEnergy));
936 } else {
938 ed << "Value of number of bins is out of range: "
939 << val << " is ignored";
940 PrintWarning(ed);
941 }
942}
943
945{
946 return nbins;
947}
948
950{
951 if(IsLocked()) { return; }
952 if(val >= 5 && val < 1000000) {
953 nbinsPerDecade = val;
954 nbins = nbinsPerDecade*G4lrint(std::log10(maxKinEnergy/minKinEnergy));
955 } else {
957 ed << "Value of number of bins per decade is out of range: "
958 << val << " is ignored";
959 PrintWarning(ed);
960 }
961}
962
964{
965 return nbinsPerDecade;
966}
967
969{
970 if(IsLocked()) { return; }
971 verbose = val;
972 workerVerbose = std::min(workerVerbose, verbose);
973}
974
976{
977 return verbose;
978}
979
981{
982 if(IsLocked()) { return; }
983 workerVerbose = val;
984}
985
987{
988 return workerVerbose;
989}
990
992{
993 if(IsLocked()) { return; }
994 mscStepLimit = val;
995}
996
998{
999 return mscStepLimit;
1000}
1001
1003{
1004 if(IsLocked()) { return; }
1005 mscStepLimitMuHad = val;
1006}
1007
1009{
1010 return mscStepLimitMuHad;
1011}
1012
1014{
1015 if(IsLocked()) { return; }
1016 fSStype = val;
1017}
1018
1020{
1021 return fSStype;
1022}
1023
1024void
1026{
1027 if(IsLocked()) { return; }
1028 nucFormfactor = val;
1029}
1030
1032{
1033 return nucFormfactor;
1034}
1035
1037{
1038 if(IsLocked()) { return; }
1039 fCParameters->SetDNAeSolvationSubType(val);
1040 ActivateDNA();
1041}
1042
1044{
1045 return fCParameters->DNAeSolvationSubType();
1046}
1047
1049{
1050 if(IsLocked()) { return; }
1051 tripletConv = val;
1052}
1053
1055{
1056 return tripletConv;
1057}
1058
1060{
1061 if(IsLocked()) { return; }
1062 fCParameters->SetPIXECrossSectionModel(sss);
1063}
1064
1066{
1067 return fCParameters->PIXECrossSectionModel();
1068}
1069
1071{
1072 if(IsLocked()) { return; }
1073 fCParameters->SetPIXEElectronCrossSectionModel(sss);
1074}
1075
1077{
1078 return fCParameters->PIXEElectronCrossSectionModel();
1079}
1080
1082{
1083 if(IsLocked()) { return; }
1084 fCParameters->SetLivermoreDataDir(sss);
1085}
1086
1088{
1089 return fCParameters->LivermoreDataDir();
1090}
1091
1092void G4EmParameters::PrintWarning(G4ExceptionDescription& ed) const
1093{
1094 G4Exception("G4EmParameters", "em0044", JustWarning, ed);
1095}
1096
1098 const G4String& region,
1099 const G4String& type)
1100{
1101 if(IsLocked()) { return; }
1102 fBParameters->AddPAIModel(particle, region, type);
1103}
1104
1105const std::vector<G4String>& G4EmParameters::ParticlesPAI() const
1106{
1107 return fBParameters->ParticlesPAI();
1108}
1109
1110const std::vector<G4String>& G4EmParameters::RegionsPAI() const
1111{
1112 return fBParameters->RegionsPAI();
1113}
1114
1115const std::vector<G4String>& G4EmParameters::TypesPAI() const
1116{
1117 return fBParameters->TypesPAI();
1118}
1119
1121{
1122 if(IsLocked()) { return; }
1123 fCParameters->AddMicroElec(region);
1124}
1125
1126const std::vector<G4String>& G4EmParameters::RegionsMicroElec() const
1127{
1128 return fCParameters->RegionsMicroElec();
1129}
1130
1131void G4EmParameters::AddDNA(const G4String& region, const G4String& type)
1132{
1133 if(IsLocked()) { return; }
1134 fCParameters->AddDNA(region, type);
1135 ActivateDNA();
1136}
1137
1138const std::vector<G4String>& G4EmParameters::RegionsDNA() const
1139{
1140 return fCParameters->RegionsDNA();
1141}
1142
1143const std::vector<G4String>& G4EmParameters::TypesDNA() const
1144{
1145 return fCParameters->TypesDNA();
1146}
1147
1148void G4EmParameters::AddPhysics(const G4String& region, const G4String& type)
1149{
1150 if(IsLocked()) { return; }
1151 fBParameters->AddPhysics(region, type);
1152}
1153
1154const std::vector<G4String>& G4EmParameters::RegionsPhysics() const
1155{
1156 return fBParameters->RegionsPhysics();
1157}
1158
1159const std::vector<G4String>& G4EmParameters::TypesPhysics() const
1160{
1161 return fBParameters->TypesPhysics();
1162}
1163
1165{
1166 if(IsLocked() && !gener) { return; }
1167 fBParameters->SetSubCutoff(val, region);
1168}
1169
1170void
1172 G4bool aauger, G4bool apixe)
1173{
1174 if(IsLocked()) { return; }
1175 fCParameters->SetDeexActiveRegion(region, adeex, aauger, apixe);
1176}
1177
1178void
1180 G4double val, G4bool wflag)
1181{
1182 if(IsLocked()) { return; }
1183 fBParameters->SetProcessBiasingFactor(procname, val, wflag);
1184}
1185
1186void
1188 const G4String& region,
1189 G4double length,
1190 G4bool wflag)
1191{
1192 if(IsLocked() && !gener) { return; }
1193 fBParameters->ActivateForcedInteraction(procname, region, length, wflag);
1194}
1195
1196void
1198 const G4String& region,
1199 G4double factor,
1200 G4double energyLim)
1201{
1202 if(IsLocked()) { return; }
1203 fBParameters->ActivateSecondaryBiasing(procname, region, factor, energyLim);
1204}
1205
1207{
1208 fBParameters->DefineRegParamForLoss(ptr);
1209}
1210
1212{
1213 fBParameters->DefineRegParamForEM(ptr);
1214}
1215
1217{
1218 return fBParameters->QuantumEntanglement();
1219}
1220
1222{
1223 if(IsLocked()) { return; }
1224 fBParameters->SetQuantumEntanglement(v);
1225}
1226
1228 return fBParameters->GetDirectionalSplitting();
1229}
1230
1232{
1233 if(IsLocked()) { return; }
1234 fBParameters->SetDirectionalSplitting(v);
1235}
1236
1238{
1239 if(IsLocked()) { return; }
1240 fBParameters->SetDirectionalSplittingTarget(v);
1241}
1242
1244{
1245 return fBParameters->GetDirectionalSplittingTarget();
1246}
1247
1249{
1250 if(IsLocked()) { return; }
1251 fBParameters->SetDirectionalSplittingRadius(r);
1252}
1253
1255{
1256 return fBParameters->GetDirectionalSplittingRadius();
1257}
1258
1260{
1261 fCParameters->DefineRegParamForDeex(ptr);
1262}
1263
1264void G4EmParameters::StreamInfo(std::ostream& os) const
1265{
1266 G4int prec = os.precision(5);
1267 os << "=======================================================================" << "\n";
1268 os << "====== Electromagnetic Physics Parameters ========" << "\n";
1269 os << "=======================================================================" << "\n";
1270 os << "LPM effect enabled " <<flagLPM << "\n";
1271 os << "Spline of EM tables enabled " <<spline << "\n";
1272 os << "Enable creation and use of sampling tables " <<fSamplingTable << "\n";
1273 os << "Apply cuts on all EM processes " <<applyCuts << "\n";
1274 os << "Use integral approach for tracking " <<integral << "\n";
1275 os << "Use general process " <<gener << "\n";
1276 os << "Enable linear polarisation for gamma " <<fPolarisation << "\n";
1277 os << "Enable sampling of quantum entanglement "
1278 <<fBParameters->QuantumEntanglement() << "\n";
1279 os << "X-section factor for integral approach " <<lambdaFactor << "\n";
1280 os << "Min kinetic energy for tables "
1281 <<G4BestUnit(minKinEnergy,"Energy") << "\n";
1282 os << "Max kinetic energy for tables "
1283 <<G4BestUnit(maxKinEnergy,"Energy") << "\n";
1284 os << "Number of bins in tables " <<nbins << "\n";
1285 os << "Number of bins per decade of a table " <<nbinsPerDecade << "\n";
1286 os << "Verbose level " <<verbose << "\n";
1287 os << "Verbose level for worker thread " <<workerVerbose << "\n";
1288 os << "Bremsstrahlung energy threshold above which \n"
1289 << " primary e+- is added to the list of secondary "
1290 <<G4BestUnit(bremsTh,"Energy") << "\n";
1291 os << "Bremsstrahlung energy threshold above which primary\n"
1292 << " muon/hadron is added to the list of secondary "
1293 <<G4BestUnit(bremsMuHadTh,"Energy") << "\n";
1294 os << "Lowest triplet kinetic energy "
1295 <<G4BestUnit(lowestTripletEnergy,"Energy") << "\n";
1296 os << "Enable sampling of gamma linear polarisation " <<fPolarisation << "\n";
1297 os << "5D gamma conversion model type " <<tripletConv << "\n";
1298 os << "5D gamma conversion model on isolated ion " <<onIsolated << "\n";
1299 if(max5DEnergyForMuPair>0.0) {
1300 os << "5D gamma conversion limit for muon pair "
1301 << max5DEnergyForMuPair/CLHEP::GeV << " GeV\n";
1302 }
1303 os << "Livermore data directory "
1304 << fCParameters->LivermoreDataDir() << "\n";
1305
1306 os << "=======================================================================" << "\n";
1307 os << "====== Ionisation Parameters ========" << "\n";
1308 os << "=======================================================================" << "\n";
1309 os << "Step function for e+- "
1310 <<"("<<fBParameters->GetStepFunctionP1() << ", "
1311 << fBParameters->GetStepFunctionP2()/CLHEP::mm << " mm)\n";
1312 os << "Step function for muons/hadrons "
1313 <<"("<<fBParameters->GetStepFunctionMuHadP1() << ", "
1314 << fBParameters->GetStepFunctionMuHadP2()/CLHEP::mm << " mm)\n";
1315 os << "Step function for light ions "
1316 <<"("<<fBParameters->GetStepFunctionLightIonsP1() << ", "
1317 << fBParameters->GetStepFunctionLightIonsP2()/CLHEP::mm << " mm)\n";
1318 os << "Step function for general ions "
1319 <<"("<<fBParameters->GetStepFunctionIonsP1() << ", "
1320 << fBParameters->GetStepFunctionIonsP2()/CLHEP::mm << " mm)\n";
1321 os << "Lowest e+e- kinetic energy "
1322 <<G4BestUnit(lowestElectronEnergy,"Energy") << "\n";
1323 os << "Lowest muon/hadron kinetic energy "
1324 <<G4BestUnit(lowestMuHadEnergy,"Energy") << "\n";
1325 os << "Fluctuations of dE/dx are enabled " <<lossFluctuation << "\n";
1326 os << "Use ICRU90 data " << fICRU90 << "\n";
1327 os << "Use built-in Birks satuaration " << birks << "\n";
1328 os << "Build CSDA range enabled " <<buildCSDARange << "\n";
1329 os << "Use cut as a final range enabled " <<cutAsFinalRange << "\n";
1330 os << "Enable angular generator interface "
1331 <<useAngGeneratorForIonisation << "\n";
1332 os << "Factor of cut reduction for sub-cutoff method " << minSubRange << "\n";
1333 os << "Max kinetic energy for CSDA tables "
1334 <<G4BestUnit(maxKinEnergyCSDA,"Energy") << "\n";
1335 os << "Max kinetic energy for NIEL computation "
1336 <<G4BestUnit(maxNIELEnergy,"Energy") << "\n";
1337 os << "Linear loss limit " <<linLossLimit << "\n";
1338 os << "Read data from file for e+e- pair production by mu " <<fMuDataFromFile << "\n";
1339
1340 os << "=======================================================================" << "\n";
1341 os << "====== Multiple Scattering Parameters ========" << "\n";
1342 os << "=======================================================================" << "\n";
1343 os << "Type of msc step limit algorithm for e+- " <<mscStepLimit << "\n";
1344 os << "Type of msc step limit algorithm for muons/hadrons " <<mscStepLimitMuHad << "\n";
1345 os << "Msc lateral displacement for e+- enabled " <<lateralDisplacement << "\n";
1346 os << "Msc lateral displacement for muons and hadrons " <<muhadLateralDisplacement << "\n";
1347 os << "Urban msc model lateral displacement alg96 " <<lateralDisplacementAlg96 << "\n";
1348 os << "Range factor for msc step limit for e+- " <<rangeFactor << "\n";
1349 os << "Range factor for msc step limit for muons/hadrons " <<rangeFactorMuHad << "\n";
1350 os << "Geometry factor for msc step limitation of e+- " <<geomFactor << "\n";
1351 os << "Safety factor for msc step limit for e+- " <<safetyFactor << "\n";
1352 os << "Skin parameter for msc step limitation of e+- " <<skin << "\n";
1353 os << "Lambda limit for msc step limit for e+- " <<lambdaLimit/CLHEP::mm << " mm\n";
1354 os << "Use Mott correction for e- scattering " << useMottCorrection << "\n";
1355 os << "Factor used for dynamic computation of angular \n"
1356 << " limit between single and multiple scattering " << factorForAngleLimit << "\n";
1357 os << "Fixed angular limit between single \n"
1358 << " and multiple scattering "
1359 << thetaLimit/CLHEP::rad << " rad\n";
1360 os << "Upper energy limit for e+- multiple scattering "
1361 << energyLimit/CLHEP::MeV << " MeV\n";
1362 os << "Type of electron single scattering model " <<fSStype << "\n";
1363 os << "Type of nuclear form-factor " <<nucFormfactor << "\n";
1364 os << "Screening factor " <<factorScreen << "\n";
1365 os << "=======================================================================" << "\n";
1366
1367 if(fCParameters->Fluo()) {
1368 os << "====== Atomic Deexcitation Parameters ========" << "\n";
1369 os << "=======================================================================" << "\n";
1370 os << "Fluorescence enabled " <<fCParameters->Fluo() << "\n";
1371 os << "Fluorescence Bearden data files enabled "
1372 <<fCParameters->BeardenFluoDir() << "\n";
1373 os << "Auger electron cascade enabled "
1374 <<fCParameters->Auger() << "\n";
1375 os << "PIXE atomic de-excitation enabled " <<fCParameters->Pixe() << "\n";
1376 os << "De-excitation module ignores cuts "
1377 <<fCParameters->DeexcitationIgnoreCut() << "\n";
1378 os << "Type of PIXE cross section for hadrons "
1379 <<fCParameters->PIXECrossSectionModel() << "\n";
1380 os << "Type of PIXE cross section for e+- "
1381 <<fCParameters->PIXEElectronCrossSectionModel() << "\n";
1382 os << "=======================================================================" << "\n";
1383 }
1384 if(fDNA) {
1385 os << "====== DNA Physics Parameters ========" << "\n";
1386 os << "=======================================================================" << "\n";
1387 os << "Use fast sampling in DNA models "
1388 << fCParameters->DNAFast() << "\n";
1389 os << "Use Stationary option in DNA models "
1390 << fCParameters->DNAStationary() << "\n";
1391 os << "Use DNA with multiple scattering of e- "
1392 << fCParameters->DNAElectronMsc() << "\n";
1393 os << "Use DNA e- solvation model type "
1394 << fCParameters->DNAeSolvationSubType() << "\n";
1395 os << "=======================================================================" << "\n";
1396 }
1397 os.precision(prec);
1398}
1399
1401{
1402#ifdef G4MULTITHREADED
1403 G4MUTEXLOCK(&emParametersMutex);
1404#endif
1406#ifdef G4MULTITHREADED
1407 G4MUTEXUNLOCK(&emParametersMutex);
1408#endif
1409}
1410
1411std::ostream& operator<< (std::ostream& os, const G4EmParameters& par)
1412{
1413 par.StreamInfo(os);
1414 return os;
1415}
1416
1417G4bool G4EmParameters::IsLocked() const
1418{
1419 return (!G4Threading::IsMasterThread() ||
1420 (fStateManager->GetCurrentState() != G4State_PreInit &&
1421 fStateManager->GetCurrentState() != G4State_Init &&
1422 fStateManager->GetCurrentState() != G4State_Idle));
1423}
1424
1425//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
@ G4State_Init
@ G4State_Idle
@ G4State_PreInit
G4DNAModelSubType
std::ostream & operator<<(std::ostream &os, const G4EmParameters &par)
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
#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
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 SetSubCutoff(G4bool val, const G4String &region="")
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
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 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)
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 Dump() 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
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)
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
void SetSubCutoff(G4bool val, const G4String &region="")
const std::vector< G4String > & RegionsPAI() const
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 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)
G4double MinSubRange() const
const std::vector< G4String > & ParticlesPAI() const
G4bool Pixe() const
const std::vector< G4String > & RegionsPhysics() const
void DefineRegParamForDeex(G4VAtomDeexcitation *) const
void SetStepFunctionMuHad(G4double v1, G4double v2)
G4double MscSafetyFactor() const
void SetAugerCascade(G4bool val)
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 &)
G4double BremsstrahlungTh() const
void SetPixe(G4bool val)
void SetMaxNIELEnergy(G4double val)
void SetStepFunctionIons(G4double v1, G4double v2)
void SetMaxEnergyForCSDARange(G4double val)
G4bool AugerCascade() const
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 GetDirectionalSplitting() const
G4double MaxKinEnergy() const
G4bool UseICRU90Data() const
void SetDirectionalSplittingRadius(G4double r)
void SetConversionType(G4int val)
G4bool LateralDisplacement() const
G4bool Spline() const
void SetUseICRU90Data(G4bool val)
G4bool MuHadLateralDisplacement() const
void SetOnIsolated(G4bool val)
G4bool EnableSamplingTable() 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
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
void SetNumberOfBins(G4int val)
G4NuclearFormfactorType NuclearFormfactorType() const
G4double LowestMuHadEnergy() const
G4double MscRangeFactor() const
G4double LambdaFactor() const
G4double FactorForAngleLimit() const
void SetSpline(G4bool val)
G4double MscSkin() const
void SetMinSubRange(G4double val)
void SetSingleScatteringType(G4eSingleScatteringType val)
G4double LowestTripletEnergy() const
G4double LowestElectronEnergy() const
void SetMscRangeFactor(G4double val)
G4bool GeneralProcessActive() const
void InitialiseG4Saturation()
static G4NistManager * Instance()
const G4ApplicationState & GetCurrentState() const
static G4StateManager * GetStateManager()
G4bool IsMasterThread()
Definition: G4Threading.cc:124
int G4lrint(double ad)
Definition: templates.hh:134