Geant4-11
G4EmParameters.hh
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 header file
29//
30// File name: G4EmParameters
31//
32// Author: Vladimir Ivanchenko for migration to MT
33//
34//
35// Creation date: 17.05.2013
36//
37// Modifications:
38//
39//
40// Class Description:
41//
42// A utility static class, responsable for keeping parameters
43// for all EM physics processes and models.
44//
45// It is initialized by the master thread but can be updated
46// at any moment. Parameters may be used in run time or at
47// initialisation
48//
49// -------------------------------------------------------------------
50//
51
52#ifndef G4EmParameters_h
53#define G4EmParameters_h 1
54
55#include "globals.hh"
56#include "G4ios.hh"
57#include "G4MscStepLimitType.hh"
59#include "G4DNAModelSubType.hh"
60#include "G4EmSaturation.hh"
61#include "G4ThreeVector.hh"
62#include "G4Threading.hh"
63#include <vector>
64
66 {
67 fWVI = 0,
69 fDPWA
70 };
71
77class G4VEmProcess;
78class G4StateManager;
79
81{
82public:
83
84 static G4EmParameters* Instance();
85
87
88 void SetDefaults();
89
90 // printing
91 void StreamInfo(std::ostream& os) const;
92 void Dump();
93 friend std::ostream& operator<< (std::ostream& os, const G4EmParameters&);
94
95 // boolean flags
97 G4bool LossFluctuation() const;
98
99 void SetBuildCSDARange(G4bool val);
100 G4bool BuildCSDARange() const;
101
102 void SetLPM(G4bool val);
103 G4bool LPM() const;
104
107
108 void SetApplyCuts(G4bool val);
109 G4bool ApplyCuts() const;
110
111 void SetFluo(G4bool val);
112 G4bool Fluo() const;
113
114 void SetBeardenFluoDir(G4bool val);
115 G4bool BeardenFluoDir() const;
116
117 void SetANSTOFluoDir(G4bool val);
118 G4bool ANSTOFluoDir() const;
119
120 void SetAuger(G4bool val);
121 void SetAugerCascade(G4bool val) { SetAuger(val); };
122 G4bool Auger() const;
123 G4bool AugerCascade() const { return Auger(); }
124
125 void SetPixe(G4bool val);
126 G4bool Pixe() const;
127
130
133
136
139
142
145
146 void SetIntegral(G4bool val);
147 G4bool Integral() const;
148
149 void SetBirksActive(G4bool val);
150 G4bool BirksActive() const;
151
152 void SetUseICRU90Data(G4bool val);
153 G4bool UseICRU90Data() const;
154
155 void SetDNAFast(G4bool val);
156 G4bool DNAFast() const;
157
158 void SetDNAStationary(G4bool val);
159 G4bool DNAStationary() const;
160
161 void SetDNAElectronMsc(G4bool val);
162 G4bool DNAElectronMsc() const;
163
164 // if general interaction is enabled then
165 // force interaction options should be disabled
168
171
174
177
180
183
184 // 5d
185 void SetOnIsolated(G4bool val);
186 G4bool OnIsolated() const;
187
188 void ActivateDNA();
189 void SetIsPrintedFlag(G4bool val);
190 G4bool IsPrintLocked() const;
191
192 // double parameters with values
193 void SetMinEnergy(G4double val);
194 G4double MinKinEnergy() const;
195
196 void SetMaxEnergy(G4double val);
197 G4double MaxKinEnergy() const;
198
201
204
207
210
213
218
219 void SetLambdaFactor(G4double val);
220 G4double LambdaFactor() const;
221
224
225 void SetMscThetaLimit(G4double val);
226 G4double MscThetaLimit() const;
227
228 void SetMscEnergyLimit(G4double val);
229 G4double MscEnergyLimit() const;
230
231 void SetMscRangeFactor(G4double val);
232 G4double MscRangeFactor() const;
233
236
237 void SetMscGeomFactor(G4double val);
238 G4double MscGeomFactor() const;
239
242
243 void SetMscLambdaLimit(G4double val);
244 G4double MscLambdaLimit() const;
245
246 void SetMscSkin(G4double val);
247 G4double MscSkin() const;
248
251
252 void SetMaxNIELEnergy(G4double val);
253 G4double MaxNIELEnergy() const;
254
257
258 void SetStepFunction(G4double v1, G4double v2);
263
266
269
270 // integer parameters
271
274 G4int NumberOfBins() const;
275
276 void SetVerbose(G4int val);
277 G4int Verbose() const;
278
279 void SetWorkerVerbose(G4int val);
280 G4int WorkerVerbose() const;
281
284
287
290
293
296
297 //5d
298 void SetConversionType(G4int val);
299 G4int GetConversionType() const;
300
301 // string parameters
304
307
308 void SetLivermoreDataDir(const G4String&);
309 const G4String& LivermoreDataDir();
310
311 // parameters per region or per process
312 void AddPAIModel(const G4String& particle,
313 const G4String& region,
314 const G4String& type);
315 const std::vector<G4String>& ParticlesPAI() const;
316 const std::vector<G4String>& RegionsPAI() const;
317 const std::vector<G4String>& TypesPAI() const;
318
319 void AddMicroElec(const G4String& region);
320 const std::vector<G4String>& RegionsMicroElec() const;
321
322 void AddDNA(const G4String& region, const G4String& type);
323 const std::vector<G4String>& RegionsDNA() const;
324 const std::vector<G4String>& TypesDNA() const;
325
326 void AddPhysics(const G4String& region, const G4String& type);
327 const std::vector<G4String>& RegionsPhysics() const;
328 const std::vector<G4String>& TypesPhysics() const;
329
330 void SetSubCutRegion(const G4String& region = "");
331
332 void SetDeexActiveRegion(const G4String& region, G4bool fdeex,
333 G4bool fauger, G4bool fpixe);
334
335 void SetProcessBiasingFactor(const G4String& procname,
336 G4double val, G4bool wflag);
337
338 void ActivateForcedInteraction(const G4String& procname,
339 const G4String& region,
340 G4double length,
341 G4bool wflag);
342
344 const G4String& region,
345 G4double factor,
347
348 // define external saturation class
350 // create and access saturation class
352
353 // initialisation methods
357
359 G4EmParameters & operator=(const G4EmParameters &right) = delete;
360
361private:
362
364
365 void Initialise();
366
367 G4bool IsLocked() const;
368
369 void PrintWarning(G4ExceptionDescription& ed) const;
370
372
378
396 G4bool onIsolated; // 5d model conversion on free ions
399
422
426 G4int tripletConv; // 5d model triplet generation type
427
432
433#ifdef G4MULTITHREADED
434 static G4Mutex emParametersMutex;
435#endif
436};
437
438//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
439
440#endif
441
G4DNAModelSubType
G4eSingleScatteringType
@ fWVI
@ fMott
@ fDPWA
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4MscStepLimitType
G4NuclearFormfactorType
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
void SetEmSaturation(G4EmSaturation *)
void SetBeardenFluoDir(G4bool val)
G4bool IsPrintLocked() const
void DefineRegParamForLoss(G4VEnergyLossProcess *) const
void SetLambdaFactor(G4double val)
static G4EmParameters * theInstance
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetBuildCSDARange(G4bool val)
G4double rangeFactorMuHad
void SetStepFunctionLightIons(G4double v1, G4double v2)
void PrintWarning(G4ExceptionDescription &ed) const
void SetEnablePolarisation(G4bool val)
G4StateManager * fStateManager
void AddDNA(const G4String &region, const G4String &type)
G4bool LateralDisplacementAlg96() const
void FillStepFunction(const G4ParticleDefinition *, G4VEnergyLossProcess *) const
G4EmParameters(G4EmParameters &)=delete
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
void SetDirectionalSplittingTarget(const G4ThreeVector &v)
G4bool BeardenFluoDir() const
G4bool DNAFast() const
G4DNAModelSubType DNAeSolvationSubType() const
G4EmSaturation * emSaturation
G4bool RetrieveMuDataFromFile() const
friend std::ostream & operator<<(std::ostream &os, const G4EmParameters &)
G4int NumberOfBins() const
G4double MscMuHadRangeFactor() const
void SetGeneralProcessActive(G4bool val)
void SetDNAFast(G4bool val)
void SetMscSafetyFactor(G4double val)
G4bool EnablePolarisation() const
void SetLateralDisplacementAlg96(G4bool val)
G4double factorScreen
void SetFactorForAngleLimit(G4double val)
const G4String & PIXECrossSectionModel()
const G4String & PIXEElectronCrossSectionModel()
G4double energyLimit
G4double MaxNIELEnergy() const
void SetRetrieveMuDataFromFile(G4bool v)
void SetDirectionalSplitting(G4bool v)
G4double lambdaLimit
const G4String & LivermoreDataDir()
G4bool OnIsolated() const
void SetMscMuHadRangeFactor(G4double val)
G4bool DNAElectronMsc() const
G4double lowestElectronEnergy
G4double MinKinEnergy() const
G4int NumberOfBinsPerDecade() const
G4MscStepLimitType mscStepLimit
G4bool BuildCSDARange() const
G4double MscThetaLimit() const
G4bool LossFluctuation() const
void SetLPM(G4bool val)
G4double MuHadBremsstrahlungTh() const
void AddPAIModel(const G4String &particle, const G4String &region, const G4String &type)
void SetDNAStationary(G4bool val)
G4bool ANSTOFluoDir() const
void SetDNAElectronMsc(G4bool val)
const std::vector< G4String > & TypesPhysics() const
void SetMaxEnergyFor5DMuPair(G4double val)
G4double GetDirectionalSplittingRadius()
G4bool muhadLateralDisplacement
void SetLinearLossLimit(G4double val)
void SetMscThetaLimit(G4double val)
G4int GetConversionType() const
G4MscStepLimitType MscMuHadStepLimitType() const
G4double MscEnergyLimit() const
G4bool lateralDisplacementAlg96
G4double factorForAngleLimit
void ActivateSecondaryBiasing(const G4String &name, const G4String &region, G4double factor, G4double energyLimit)
G4bool BirksActive() const
G4bool DNAStationary() const
G4bool IsLocked() const
G4double minKinEnergy
const std::vector< G4String > & RegionsPAI() const
G4double lowestTripletEnergy
G4double maxKinEnergyCSDA
void SetSubCutRegion(const G4String &region="")
void SetLossFluctuations(G4bool val)
void SetDeexActiveRegion(const G4String &region, G4bool fdeex, G4bool fauger, G4bool fpixe)
G4bool UseCutAsFinalRange() const
void SetPIXEElectronCrossSectionModel(const G4String &)
void SetLowestTripletEnergy(G4double val)
void SetMuHadLateralDisplacement(G4bool val)
G4bool Fluo() const
G4EmSaturation * GetEmSaturation()
void SetDNAeSolvationSubType(G4DNAModelSubType val)
G4bool useAngGeneratorForIonisation
void SetQuantumEntanglement(G4bool v)
void DefineRegParamForEM(G4VEmProcess *) const
void ActivateForcedInteraction(const G4String &procname, const G4String &region, G4double length, G4bool wflag)
G4double lowestMuHadEnergy
void ActivateAngularGeneratorForIonisation(G4bool val)
G4double maxKinEnergy
void SetScreeningFactor(G4double val)
void SetNuclearFormfactorType(G4NuclearFormfactorType val)
void SetStepFunction(G4double v1, G4double v2)
void SetLateralDisplacement(G4bool val)
G4int Verbose() const
G4bool QuantumEntanglement() const
void SetWorkerVerbose(G4int val)
void SetUseCutAsFinalRange(G4bool val)
void SetDeexcitationIgnoreCut(G4bool val)
void SetBirksActive(G4bool val)
G4double safetyFactor
G4double ScreeningFactor() const
G4bool UseMottCorrection() const
void SetFluo(G4bool val)
void SetMuHadBremsstrahlungTh(G4double val)
const std::vector< G4String > & ParticlesPAI() const
G4double max5DEnergyForMuPair
G4bool Pixe() const
const std::vector< G4String > & RegionsPhysics() const
G4double lambdaFactor
G4MscStepLimitType mscStepLimitMuHad
void DefineRegParamForDeex(G4VAtomDeexcitation *) const
void SetStepFunctionMuHad(G4double v1, G4double v2)
G4double MscSafetyFactor() const
G4NuclearFormfactorType nucFormfactor
void SetAugerCascade(G4bool val)
G4double rangeFactor
void SetVerbose(G4int val)
G4int WorkerVerbose() const
void SetMscGeomFactor(G4double val)
void SetMscLambdaLimit(G4double val)
void SetMscSkin(G4double val)
void SetApplyCuts(G4bool val)
const std::vector< G4String > & TypesDNA() const
G4MscStepLimitType MscStepLimitType() const
G4bool ApplyCuts() const
G4double linLossLimit
void SetEnableSamplingTable(G4bool val)
const std::vector< G4String > & RegionsDNA() const
void SetLivermoreDataDir(const G4String &)
G4double BremsstrahlungTh() const
void SetPixe(G4bool val)
void SetMaxNIELEnergy(G4double val)
void SetStepFunctionIons(G4double v1, G4double v2)
void SetMaxEnergyForCSDARange(G4double val)
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)
G4EmParameters & operator=(const G4EmParameters &right)=delete
void SetMscEnergyLimit(G4double val)
void SetBremsstrahlungTh(G4double val)
void SetAuger(G4bool val)
G4bool lateralDisplacement
G4bool GetDirectionalSplitting() const
G4EmLowEParameters * fCParameters
void SetIsPrintedFlag(G4bool val)
G4double MaxKinEnergy() const
G4bool UseICRU90Data() const
void SetDirectionalSplittingRadius(G4double r)
void SetConversionType(G4int val)
G4bool LateralDisplacement() const
void SetUseICRU90Data(G4bool val)
G4double bremsMuHadTh
G4bool MuHadLateralDisplacement() const
void SetOnIsolated(G4bool val)
G4bool EnableSamplingTable() const
G4EmParametersMessenger * theMessenger
void StreamInfo(std::ostream &os) const
G4double MscLambdaLimit() const
void SetPIXECrossSectionModel(const G4String &)
void SetIntegral(G4bool val)
G4ThreeVector GetDirectionalSplittingTarget() const
G4bool Auger() const
void SetUseMottCorrection(G4bool val)
void AddMicroElec(const G4String &region)
const std::vector< G4String > & RegionsMicroElec() const
G4double MaxEnergyFor5DMuPair() const
G4bool Integral() const
G4double MaxEnergyForCSDARange() const
G4EmExtraParameters * fBParameters
G4bool useMottCorrection
void SetLowestMuHadEnergy(G4double val)
const std::vector< G4String > & TypesPAI() const
G4bool LPM() const
G4bool UseAngularGeneratorForIonisation() const
G4double maxNIELEnergy
void SetMaxEnergy(G4double val)
G4double LinearLossLimit() const
G4NuclearFormfactorType NuclearFormfactorType() const
G4double LowestMuHadEnergy() const
G4double MscRangeFactor() const
G4double LambdaFactor() const
G4double FactorForAngleLimit() const
G4eSingleScatteringType fSStype
G4double MscSkin() const
void SetSingleScatteringType(G4eSingleScatteringType val)
void SetANSTOFluoDir(G4bool val)
G4double LowestTripletEnergy() const
G4double LowestElectronEnergy() const
void SetMscRangeFactor(G4double val)
G4bool GeneralProcessActive() const
const char * name(G4int ptype)