Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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() const;
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
105 void SetSpline(G4bool val);
106 G4bool Spline() const;
107
110
111 void SetApplyCuts(G4bool val);
112 G4bool ApplyCuts() const;
113
114 void SetFluo(G4bool val);
115 G4bool Fluo() const;
116
117 void SetBeardenFluoDir(G4bool val);
118 G4bool BeardenFluoDir() const;
119
120 void SetAuger(G4bool val);
121 G4bool Auger() const;
122
123 // obsolete methods
124 void SetAugerCascade(G4bool val);
125 G4bool AugerCascade() const;
126
127 void SetPixe(G4bool val);
128 G4bool Pixe() const;
129
132
135
138
141
144
147
148 void SetIntegral(G4bool val);
149 G4bool Integral() const;
150
151 void SetBirksActive(G4bool val);
152 G4bool BirksActive() const;
153
154 void SetUseICRU90Data(G4bool val);
155 G4bool UseICRU90Data() const;
156
157 void SetDNAFast(G4bool val);
158 G4bool DNAFast() const;
159
160 void SetDNAStationary(G4bool val);
161 G4bool DNAStationary() const;
162
163 void SetDNAElectronMsc(G4bool val);
164 G4bool DNAElectronMsc() const;
165
168
171
174
177
180
183
184 // 5d
185 void SetOnIsolated(G4bool val);
186 G4bool OnIsolated() const;
187
188 void ActivateDNA();
189
190 // double parameters with values
191 void SetMinSubRange(G4double val);
192 G4double MinSubRange() const;
193
194 void SetMinEnergy(G4double val);
195 G4double MinKinEnergy() const;
196
197 void SetMaxEnergy(G4double val);
198 G4double MaxKinEnergy() const;
199
202
205
208
211
214
219
220 void SetLambdaFactor(G4double val);
221 G4double LambdaFactor() const;
222
225
226 void SetMscThetaLimit(G4double val);
227 G4double MscThetaLimit() const;
228
229 void SetMscEnergyLimit(G4double val);
230 G4double MscEnergyLimit() const;
231
232 void SetMscRangeFactor(G4double val);
233 G4double MscRangeFactor() const;
234
237
238 void SetMscGeomFactor(G4double val);
239 G4double MscGeomFactor() const;
240
243
244 void SetMscLambdaLimit(G4double val);
245 G4double MscLambdaLimit() const;
246
247 void SetMscSkin(G4double val);
248 G4double MscSkin() const;
249
252
253 void SetMaxNIELEnergy(G4double val);
254 G4double MaxNIELEnergy() const;
255
258
259 void SetStepFunction(G4double v1, G4double v2);
264
267
270
271 // integer parameters
272 void SetNumberOfBins(G4int val);
273 G4int NumberOfBins() const;
274
277
278 void SetVerbose(G4int val);
279 G4int Verbose() const;
280
281 void SetWorkerVerbose(G4int val);
282 G4int WorkerVerbose() const;
283
286
289
292
295
298
299 //5d
300 void SetConversionType(G4int val);
301 G4int GetConversionType() const;
302
303 // string parameters
306
309
310 void SetLivermoreDataDir(const G4String&);
311 const G4String& LivermoreDataDir();
312
313 // parameters per region or per process
314 void AddPAIModel(const G4String& particle,
315 const G4String& region,
316 const G4String& type);
317 const std::vector<G4String>& ParticlesPAI() const;
318 const std::vector<G4String>& RegionsPAI() const;
319 const std::vector<G4String>& TypesPAI() const;
320
321 void AddMicroElec(const G4String& region);
322 const std::vector<G4String>& RegionsMicroElec() const;
323
324 void AddDNA(const G4String& region, const G4String& type);
325 const std::vector<G4String>& RegionsDNA() const;
326 const std::vector<G4String>& TypesDNA() const;
327
328 void AddPhysics(const G4String& region, const G4String& type);
329 const std::vector<G4String>& RegionsPhysics() const;
330 const std::vector<G4String>& TypesPhysics() const;
331
332 void SetSubCutoff(G4bool val, const G4String& region = "");
333
334 void SetDeexActiveRegion(const G4String& region, G4bool fdeex,
335 G4bool fauger, G4bool fpixe);
336
337 void SetProcessBiasingFactor(const G4String& procname,
338 G4double val, G4bool wflag);
339
340 void ActivateForcedInteraction(const G4String& procname,
341 const G4String& region,
342 G4double length,
343 G4bool wflag);
344
345 void ActivateSecondaryBiasing(const G4String& name,
346 const G4String& region,
347 G4double factor,
348 G4double energyLimit);
349
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
371 static G4EmParameters* theInstance;
372
373 G4EmParametersMessenger* theMessenger;
374 G4EmExtraParameters* fBParameters;
375 G4EmLowEParameters* fCParameters;
376 G4StateManager* fStateManager;
377 G4EmSaturation* emSaturation;
378
379 G4bool lossFluctuation;
380 G4bool buildCSDARange;
381 G4bool flagLPM;
382 G4bool spline;
383 G4bool cutAsFinalRange;
384 G4bool applyCuts;
385 G4bool lateralDisplacement;
386 G4bool lateralDisplacementAlg96;
387 G4bool muhadLateralDisplacement;
388 G4bool useAngGeneratorForIonisation;
389 G4bool useMottCorrection;
390 G4bool integral;
391 G4bool birks;
392 G4bool fICRU90;
393 G4bool gener;
394 G4bool fSamplingTable;
395 G4bool fPolarisation;
396 G4bool fMuDataFromFile;
397 G4bool onIsolated; // 5d model conversion on free ions
398 G4bool fDNA;
399
400 G4double minSubRange;
401 G4double minKinEnergy;
402 G4double maxKinEnergy;
403 G4double maxKinEnergyCSDA;
404 G4double max5DEnergyForMuPair;
405 G4double lowestElectronEnergy;
406 G4double lowestMuHadEnergy;
407 G4double lowestTripletEnergy;
408 G4double linLossLimit;
409 G4double bremsTh;
410 G4double bremsMuHadTh;
411 G4double lambdaFactor;
412 G4double factorForAngleLimit;
413 G4double thetaLimit;
414 G4double energyLimit;
415 G4double maxNIELEnergy;
416 G4double rangeFactor;
417 G4double rangeFactorMuHad;
418 G4double geomFactor;
419 G4double safetyFactor;
420 G4double lambdaLimit;
421 G4double skin;
422 G4double factorScreen;
423
424 G4int nbins;
425 G4int nbinsPerDecade;
426 G4int verbose;
427 G4int workerVerbose;
428 G4int tripletConv; // 5d model triplet generation type
429
430 G4MscStepLimitType mscStepLimit;
431 G4MscStepLimitType mscStepLimitMuHad;
432 G4NuclearFormfactorType nucFormfactor;
434
435#ifdef G4MULTITHREADED
436 static G4Mutex emParametersMutex;
437#endif
438};
439
440//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
441
442#endif
443
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)
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
G4EmParameters(G4EmParameters &)=delete
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
void SetDirectionalSplittingTarget(const G4ThreeVector &v)
G4bool BeardenFluoDir() const
G4bool DNAFast() const
G4DNAModelSubType DNAeSolvationSubType() const
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)
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)
G4EmParameters & operator=(const G4EmParameters &right)=delete
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