Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4EmParameters.hh
이 파일의 문서화 페이지로 가기
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 // $Id: G4EmParameters.hh 66885 2013-01-16 17:37:13Z gunter $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4EmParameters
34 //
35 // Author: Vladimir Ivanchenko for migration to MT
36 //
37 //
38 // Creation date: 17.05.2013
39 //
40 // Modifications:
41 //
42 //
43 // Class Description:
44 //
45 // A utility static class, responsable for keeping parameters
46 // for all EM physics processes and models.
47 //
48 // It is initialized by the master thread but can be updated
49 // at any moment. Parameters may be used in run time or at
50 // initialisation
51 //
52 // -------------------------------------------------------------------
53 //
54 
55 #ifndef G4EmParameters_h
56 #define G4EmParameters_h 1
57 
58 #include "globals.hh"
59 #include "G4ios.hh"
60 #include "G4MscStepLimitType.hh"
62 #include "G4DNAModelSubType.hh"
63 #include "G4EmSaturation.hh"
64 #include "G4Threading.hh"
65 #include <vector>
66 
69 class G4VEmProcess;
71 class G4StateManager;
72 
74 {
75 public:
76 
77  static G4EmParameters* Instance();
78 
80 
81  void SetDefaults();
82 
83  // printing
84  std::ostream& StreamInfo(std::ostream& os) const;
85  void Dump() const;
86  friend std::ostream& operator<< (std::ostream& os, const G4EmParameters&);
87 
88  // boolean flags
89  void SetLossFluctuations(G4bool val);
90  G4bool LossFluctuation() const;
91 
92  void SetBuildCSDARange(G4bool val);
93  G4bool BuildCSDARange() const;
94 
95  void SetLPM(G4bool val);
96  G4bool LPM() const;
97 
98  void SetSpline(G4bool val);
99  G4bool Spline() const;
100 
101  void SetUseCutAsFinalRange(G4bool val);
102  G4bool UseCutAsFinalRange() const;
103 
104  void SetApplyCuts(G4bool val);
105  G4bool ApplyCuts() const;
106 
107  void SetFluo(G4bool val);
108  G4bool Fluo() const;
109 
110  void SetBeardenFluoDir(G4bool val);
111  G4bool BeardenFluoDir() const;
112 
113  void SetAuger(G4bool val);
114  G4bool Auger() const;
115 
116  void SetAugerCascade(G4bool val);
117  G4bool AugerCascade() const;
118 
119  void SetPixe(G4bool val);
120  G4bool Pixe() const;
121 
124 
125  void SetLateralDisplacement(G4bool val);
126  G4bool LateralDisplacement() const;
127 
130 
133 
136 
139 
140  void SetUseMottCorrection(G4bool val);
141  G4bool UseMottCorrection() const;
142 
143  void SetIntegral(G4bool val);
144  G4bool Integral() const;
145 
146  void SetBirksActive(G4bool val);
147  G4bool BirksActive() const;
148 
149  void SetDNAFast(G4bool val);
150  G4bool DNAFast() const;
151 
152  void SetDNAStationary(G4bool val);
153  G4bool DNAStationary() const;
154 
155  void SetDNAElectronMsc(G4bool val);
156  G4bool DNAElectronMsc() const;
157 
158  void SetGammaSharkActive(G4bool val);
159  G4bool GammaSharkActive() const;
160 
163 
164  // 5d
165  void SetOnIsolated(G4bool val);
166  bool OnIsolated() const;
167 
168  // double parameters with values
169  void SetMinSubRange(G4double val);
170  G4double MinSubRange() const;
171 
172  void SetMinEnergy(G4double val);
173  G4double MinKinEnergy() const;
174 
175  void SetMaxEnergy(G4double val);
176  G4double MaxKinEnergy() const;
177 
180 
183 
184  void SetLowestMuHadEnergy(G4double val);
185  G4double LowestMuHadEnergy() const;
186 
189 
190  void SetLinearLossLimit(G4double val);
191  G4double LinearLossLimit() const;
192 
193  void SetBremsstrahlungTh(G4double val);
194  G4double BremsstrahlungTh() const;
195 
196  void SetLambdaFactor(G4double val);
197  G4double LambdaFactor() const;
198 
201 
202  void SetMscThetaLimit(G4double val);
203  G4double MscThetaLimit() const;
204 
205  void SetMscEnergyLimit(G4double val);
206  G4double MscEnergyLimit() const;
207 
208  void SetMscRangeFactor(G4double val);
209  G4double MscRangeFactor() const;
210 
213 
214  void SetMscGeomFactor(G4double val);
215  G4double MscGeomFactor() const;
216 
217  void SetMscSkin(G4double val);
218  G4double MscSkin() const;
219 
220  void SetScreeningFactor(G4double val);
221  G4double ScreeningFactor() const;
222 
223  void SetStepFunction(G4double v1, G4double v2);
224 
226 
227  // integer parameters
228  void SetNumberOfBins(G4int val);
229  G4int NumberOfBins() const;
230 
233 
234  void SetVerbose(G4int val);
235  G4int Verbose() const;
236 
237  void SetWorkerVerbose(G4int val);
238  G4int WorkerVerbose() const;
239 
242 
245 
248 
251 
252  //5d
253  void SetConversionType(G4int val);
254  G4int GetConversionType() const;
255 
256  // string parameters
257  void SetPIXECrossSectionModel(const G4String&);
259 
262 
263  // parameters per region or per process
264  void AddPAIModel(const G4String& particle,
265  const G4String& region,
266  const G4String& type);
267  const std::vector<G4String>& ParticlesPAI() const;
268  const std::vector<G4String>& RegionsPAI() const;
269  const std::vector<G4String>& TypesPAI() const;
270 
271  void AddMicroElec(const G4String& region);
272  const std::vector<G4String>& RegionsMicroElec() const;
273 
274  void AddDNA(const G4String& region, const G4String& type);
275  const std::vector<G4String>& RegionsDNA() const;
276  const std::vector<G4String>& TypesDNA() const;
277 
278  void AddMsc(const G4String& region, const G4String& type);
279  const std::vector<G4String>& RegionsMsc() const;
280  const std::vector<G4String>& TypesMsc() const;
281 
282  void AddPhysics(const G4String& region, const G4String& type);
283  const std::vector<G4String>& RegionsPhysics() const;
284  const std::vector<G4String>& TypesPhysics() const;
285 
286  void SetSubCutoff(G4bool val, const G4String& region = "");
287 
288  void SetDeexActiveRegion(const G4String& region, G4bool fdeex,
289  G4bool fauger, G4bool fpixe);
290 
291  void SetProcessBiasingFactor(const G4String& procname,
292  G4double val, G4bool wflag);
293 
294  void ActivateForcedInteraction(const G4String& procname,
295  const G4String& region,
296  G4double length,
297  G4bool wflag);
298 
300  const G4String& region,
301  G4double factor,
303 
304  // initialisation methods
306  G4bool isElectron) const;
307  void DefineRegParamForEM(G4VEmProcess*) const;
309 
310  G4EmParameters(G4EmParameters &) = delete;
311  G4EmParameters & operator=(const G4EmParameters &right) = delete;
312 
313 private:
314 
315  G4EmParameters();
316 
317  void Initialise();
318 
319  G4bool IsLocked() const;
320 
321  G4String CheckRegion(const G4String&) const;
322 
323  void PrintWarning(G4ExceptionDescription& ed) const;
324 
326 
328 
330 
332 
357  G4bool onIsolated; // 5d model conversion on free ions
358 
381 
386  G4int tripletConv; // 5d model triplet generation type
387 
392 
395 
396  std::vector<G4String> m_particlesPAI;
397  std::vector<G4String> m_regnamesPAI;
398  std::vector<G4String> m_typesPAI;
399 
400  std::vector<G4String> m_regnamesME;
401 
402  std::vector<G4String> m_regnamesDNA;
403  std::vector<G4String> m_typesDNA;
404 
405  std::vector<G4String> m_regnamesMsc;
406  std::vector<G4String> m_typesMsc;
407 
408  std::vector<G4String> m_regnamesSubCut;
409  std::vector<G4bool> m_subCuts;
410 
411  std::vector<G4String> m_regnamesDeex;
412  std::vector<G4bool> m_fluo;
413  std::vector<G4bool> m_auger;
414  std::vector<G4bool> m_pixe;
415 
416  std::vector<G4String> m_procBiasedXS;
417  std::vector<G4double> m_factBiasedXS;
418  std::vector<G4bool> m_weightBiasedXS;
419 
420  std::vector<G4String> m_procForced;
421  std::vector<G4String> m_regnamesForced;
422  std::vector<G4double> m_lengthForced;
423  std::vector<G4bool> m_weightForced;
424 
425  std::vector<G4String> m_procBiasedSec;
426  std::vector<G4String> m_regnamesBiasedSec;
427  std::vector<G4double> m_factBiasedSec;
428  std::vector<G4double> m_elimBiasedSec;
429 
430 #ifdef G4MULTITHREADED
431  static G4Mutex emParametersMutex;
432 #endif
433 };
434 
435 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
436 
437 #endif
438 
std::vector< G4double > m_factBiasedXS
std::vector< G4bool > m_weightForced
G4String CheckRegion(const G4String &) const
std::vector< G4String > m_procBiasedXS
void SetMaxEnergy(G4double val)
void SetMscEnergyLimit(G4double val)
void SetDeexcitationIgnoreCut(G4bool val)
G4double BremsstrahlungTh() const
const XML_Char * name
Definition: expat.h:151
G4MscStepLimitType MscMuHadStepLimitType() const
G4DNAModelSubType DNAeSolvationSubType() const
const std::vector< G4String > & RegionsPhysics() const
void SetBremsstrahlungTh(G4double val)
friend std::ostream & operator<<(std::ostream &os, const G4EmParameters &)
std::vector< G4bool > m_weightBiasedXS
void SetMinSubRange(G4double val)
G4int NumberOfBinsPerDecade() const
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:45
void SetMscMuHadStepLimitType(G4MscStepLimitType val)
G4MscStepLimitType
G4bool LateralDisplacementAlg96() const
G4MscStepLimitType MscStepLimitType() const
void SetLambdaFactor(G4double val)
void SetSpline(G4bool val)
std::vector< G4String > m_regnamesDeex
G4bool UseAngularGeneratorForIonisation() const
G4double lambdaFactor
G4EmSaturation * GetEmSaturation()
G4bool lateralDisplacement
void SetUseCutAsFinalRange(G4bool val)
G4double linLossLimit
void SetPIXEElectronCrossSectionModel(const G4String &)
G4double rangeFactor
void DefineRegParamForEM(G4VEmProcess *) const
void PrintWarning(G4ExceptionDescription &ed) const
G4bool GammaSharkActive() const
G4bool LateralDisplacement() const
G4bool ApplyCuts() const
void DefineRegParamForLoss(G4VEnergyLossProcess *, G4bool isElectron) const
G4bool BeardenFluoDir() const
void SetDNAElectronMsc(G4bool val)
void SetMuHadLateralDisplacement(G4bool val)
G4double MaxKinEnergy() const
G4bool BuildCSDARange() const
void SetBeardenFluoDir(G4bool val)
G4int GetConversionType() const
G4bool DNAStationary() const
void AddPhysics(const G4String &region, const G4String &type)
G4bool LPM() const
G4bool muhadLateralDisplacement
void SetVerbose(G4int val)
std::vector< G4String > m_regnamesForced
void SetStepFunctionMuHad(G4double v1, G4double v2)
void SetLPM(G4bool val)
G4bool UseCutAsFinalRange() const
const std::vector< G4String > & TypesDNA() const
G4bool DNAFast() const
std::vector< G4bool > m_fluo
G4bool MuHadLateralDisplacement() const
void SetOnIsolated(G4bool val)
void AddPAIModel(const G4String &particle, const G4String &region, const G4String &type)
void SetEmSaturation(G4EmSaturation *)
void ActivateSecondaryBiasing(const G4String &name, const G4String &region, G4double factor, G4double energyLimit)
G4double MscMuHadRangeFactor() const
G4bool lateralDisplacementAlg96
G4StateManager * fStateManager
void SetBuildCSDARange(G4bool val)
void SetLossFluctuations(G4bool val)
void SetLowestTripletEnergy(G4double val)
G4double factorScreen
G4bool BirksActive() const
void ActivateAngularGeneratorForIonisation(G4bool val)
const std::vector< G4String > & RegionsPAI() const
void SetLateralDisplacement(G4bool val)
const std::vector< G4String > & RegionsMicroElec() const
G4EmSaturation * emSaturation
std::vector< G4bool > m_auger
void SetMinEnergy(G4double val)
G4DNAModelSubType dnaElectronSolvation
void SetMscSkin(G4double val)
G4double LinearLossLimit() const
G4bool latDisplacementBeyondSafety
double G4double
Definition: G4Types.hh:76
G4double MaxEnergyForCSDARange() const
bool G4bool
Definition: G4Types.hh:79
void SetDNAStationary(G4bool val)
G4double rangeFactorMuHad
void SetMscThetaLimit(G4double val)
G4EmParameters & operator=(const G4EmParameters &right)=delete
G4NuclearFormfactorType nucFormfactor
G4bool LatDisplacementBeyondSafety() const
void SetMscMuHadRangeFactor(G4double val)
void SetConversionType(G4int val)
G4double MinSubRange() const
void SetStepFunction(G4double v1, G4double v2)
void SetLowestElectronEnergy(G4double val)
void SetMscRangeFactor(G4double val)
void SetPIXECrossSectionModel(const G4String &)
const std::vector< G4String > & TypesPAI() const
G4double LowestMuHadEnergy() const
G4double MscSkin() const
void SetSubCutoff(G4bool val, const G4String &region="")
void SetFactorForAngleLimit(G4double val)
std::vector< G4String > m_regnamesDNA
G4bool UseMottCorrection() const
void SetLatDisplacementBeyondSafety(G4bool val)
G4int WorkerVerbose() const
void SetFluo(G4bool val)
G4double finalRangeMuHad
G4double lowestTripletEnergy
G4double LambdaFactor() const
void SetPixe(G4bool val)
G4MscStepLimitType mscStepLimit
void AddMicroElec(const G4String &region)
G4double energyLimit
void SetAuger(G4bool val)
void SetNumberOfBinsPerDecade(G4int val)
G4int Verbose() const
void AddDNA(const G4String &region, const G4String &type)
G4MscStepLimitType mscStepLimitMuHad
G4bool Pixe() const
void SetIntegral(G4bool val)
G4double MinKinEnergy() const
G4double MscRangeFactor() const
G4double ScreeningFactor() const
const std::vector< G4String > & TypesMsc() const
G4double lowestMuHadEnergy
G4bool AugerCascade() const
G4bool useAngGeneratorForIonisation
const G4String & PIXECrossSectionModel()
void SetBirksActive(G4bool val)
void ActivateForcedInteraction(const G4String &procname, const G4String &region, G4double length, G4bool wflag)
void SetLowestMuHadEnergy(G4double val)
G4double LowestTripletEnergy() const
void SetMscStepLimitType(G4MscStepLimitType val)
static G4EmParameters * theInstance
void SetAugerCascade(G4bool val)
G4double dRoverRange
void SetDNAFast(G4bool val)
int G4int
Definition: G4Types.hh:78
G4double MscThetaLimit() const
void AddMsc(const G4String &region, const G4String &type)
G4NuclearFormfactorType NuclearFormfactorType() const
void SetNuclearFormfactorType(G4NuclearFormfactorType val)
std::vector< G4double > m_lengthForced
std::vector< G4String > m_particlesPAI
G4double maxKinEnergy
std::vector< G4bool > m_subCuts
G4NuclearFormfactorType
G4bool DeexcitationIgnoreCut() const
std::vector< G4String > m_procForced
const G4String & PIXEElectronCrossSectionModel()
std::vector< G4String > m_regnamesBiasedSec
G4double minSubRange
std::vector< G4double > m_elimBiasedSec
void DefineRegParamForDeex(G4VAtomDeexcitation *) const
void SetLinearLossLimit(G4double val)
std::vector< G4String > m_typesMsc
G4double MscGeomFactor() const
void SetWorkerVerbose(G4int val)
std::vector< G4String > m_regnamesMsc
G4double factorForAngleLimit
std::vector< G4String > m_procBiasedSec
G4bool Auger() const
void SetMscGeomFactor(G4double val)
G4double maxKinEnergyCSDA
G4double dRoverRangeMuHad
std::vector< G4String > m_typesDNA
void SetMaxEnergyForCSDARange(G4double val)
const std::vector< G4String > & TypesPhysics() const
G4EmParametersMessenger * theMessenger
G4double FactorForAngleLimit() const
void SetGammaSharkActive(G4bool val)
G4double minKinEnergy
void SetDeexActiveRegion(const G4String &region, G4bool fdeex, G4bool fauger, G4bool fpixe)
void SetLateralDisplacementAlg96(G4bool val)
std::vector< G4double > m_factBiasedSec
void SetScreeningFactor(G4double val)
G4bool isElectron(G4int ityp)
std::vector< G4String > m_regnamesSubCut
G4bool LossFluctuation() const
void SetProcessBiasingFactor(const G4String &procname, G4double val, G4bool wflag)
G4int NumberOfBins() const
void SetNumberOfBins(G4int val)
std::vector< G4String > m_typesPAI
G4double MscEnergyLimit() const
void SetUseMottCorrection(G4bool val)
const std::vector< G4String > & RegionsDNA() const
const std::vector< G4String > & RegionsMsc() const
std::vector< G4bool > m_pixe
G4bool DNAElectronMsc() const
void Dump() const
G4bool IsLocked() const
G4DNAModelSubType
G4double LowestElectronEnergy() const
G4bool Spline() const
G4String nameElectronPIXE
const std::vector< G4String > & ParticlesPAI() const
bool OnIsolated() const
std::vector< G4String > m_regnamesME
std::ostream & StreamInfo(std::ostream &os) const
static G4EmParameters * Instance()
void SetDNAeSolvationSubType(G4DNAModelSubType val)
G4bool Integral() const
std::vector< G4String > m_regnamesPAI
void SetApplyCuts(G4bool val)
std::mutex G4Mutex
Definition: G4Threading.hh:84
G4bool useMottCorrection
G4bool Fluo() const
G4double lowestElectronEnergy