Geant4 11.1.1
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 "G4EmFluoDirectory.hh"
61#include "G4EmSaturation.hh"
62#include "G4ThreeVector.hh"
63#include <vector>
64
66{
67 fWVI = 0,
69 fDPWA
70};
71
73{
74 fDisabled = 0,
77};
78
80{
84};
85
91class G4VEmProcess;
92class G4StateManager;
93
95{
96public:
97
98 static G4EmParameters* Instance();
99
101
102 void SetDefaults();
103
104 // printing
105 void StreamInfo(std::ostream& os) const;
106 void Dump();
107 friend std::ostream& operator<< (std::ostream& os, const G4EmParameters&);
108
109 // boolean flags
110 void SetLossFluctuations(G4bool val);
111 G4bool LossFluctuation() const;
112
113 void SetBuildCSDARange(G4bool val);
114 G4bool BuildCSDARange() const;
115
116 void SetLPM(G4bool val);
117 G4bool LPM() const;
118
121
122 void SetApplyCuts(G4bool val);
123 G4bool ApplyCuts() const;
124
125 void SetFluo(G4bool val);
126 G4bool Fluo() const;
127
129
131 void SetBeardenFluoDir(G4bool val);
132 void SetANSTOFluoDir(G4bool val);
133 void SetXDB_EADLFluoDir(G4bool val);
134
137
138 void SetAuger(G4bool val);
139 void SetAugerCascade(G4bool val) { SetAuger(val); };
140 G4bool Auger() const;
141 G4bool AugerCascade() const { return Auger(); }
142
143 void SetPixe(G4bool val);
144 G4bool Pixe() const;
145
148
151
154
157
160
163
164 void SetIntegral(G4bool val);
165 G4bool Integral() const;
166
167 void SetBirksActive(G4bool val);
168 G4bool BirksActive() const;
169
170 void SetUseICRU90Data(G4bool val);
171 G4bool UseICRU90Data() const;
172
175
176 void SetDNAFast(G4bool val);
177 G4bool DNAFast() const;
178
179 void SetDNAStationary(G4bool val);
180 G4bool DNAStationary() const;
181
182 void SetDNAElectronMsc(G4bool val);
183 G4bool DNAElectronMsc() const;
184
185 // if general interaction is enabled then
186 // force interaction options should be disabled
189
192
195
198
201
204
207
210
211 // 5d
212 void SetOnIsolated(G4bool val);
213 G4bool OnIsolated() const;
214
215 void ActivateDNA();
216 void SetIsPrintedFlag(G4bool val);
217 G4bool IsPrintLocked() const;
218
219 // double parameters with values
220 void SetMinEnergy(G4double val);
221 G4double MinKinEnergy() const;
222
223 void SetMaxEnergy(G4double val);
224 G4double MaxKinEnergy() const;
225
228
231
234
237
240
245
246 void SetLambdaFactor(G4double val);
247 G4double LambdaFactor() const;
248
251
252 void SetMscThetaLimit(G4double val);
253 G4double MscThetaLimit() const;
254
255 void SetMscEnergyLimit(G4double val);
256 G4double MscEnergyLimit() const;
257
258 void SetMscRangeFactor(G4double val);
259 G4double MscRangeFactor() const;
260
263
264 void SetMscGeomFactor(G4double val);
265 G4double MscGeomFactor() const;
266
269
270 void SetMscLambdaLimit(G4double val);
271 G4double MscLambdaLimit() const;
272
273 void SetMscSkin(G4double val);
274 G4double MscSkin() const;
275
278
279 void SetMaxNIELEnergy(G4double val);
280 G4double MaxNIELEnergy() const;
281
284
285 void SetStepFunction(G4double v1, G4double v2);
290
293
296
297 // integer parameters
298
301 G4int NumberOfBins() const;
302
303 void SetVerbose(G4int val);
304 G4int Verbose() const;
305
306 void SetWorkerVerbose(G4int val);
307 G4int WorkerVerbose() const;
308
311
314
317
320
323
326
327 //5d
328 void SetConversionType(G4int val);
329 G4int GetConversionType() const;
330
331 // string parameters
334
337
338 void SetLivermoreDataDir(const G4String&);
339 const G4String& LivermoreDataDir();
340
341 // parameters per region or per process
342 void AddPAIModel(const G4String& particle,
343 const G4String& region,
344 const G4String& type);
345 const std::vector<G4String>& ParticlesPAI() const;
346 const std::vector<G4String>& RegionsPAI() const;
347 const std::vector<G4String>& TypesPAI() const;
348
349 void AddMicroElec(const G4String& region);
350 const std::vector<G4String>& RegionsMicroElec() const;
351
352 void AddDNA(const G4String& region, const G4String& type);
353 const std::vector<G4String>& RegionsDNA() const;
354 const std::vector<G4String>& TypesDNA() const;
355
356 void AddPhysics(const G4String& region, const G4String& type);
357 const std::vector<G4String>& RegionsPhysics() const;
358 const std::vector<G4String>& TypesPhysics() const;
359
360 void SetSubCutRegion(const G4String& region = "");
361
362 void SetDeexActiveRegion(const G4String& region, G4bool fdeex,
363 G4bool fauger, G4bool fpixe);
364
365 void SetProcessBiasingFactor(const G4String& procname,
366 G4double val, G4bool wflag);
367
368 void ActivateForcedInteraction(const G4String& procname,
369 const G4String& region,
370 G4double length,
371 G4bool wflag);
372
373 void ActivateSecondaryBiasing(const G4String& name,
374 const G4String& region,
375 G4double factor,
376 G4double energyLimit);
377
378 // define external saturation class
380 // create and access saturation class
382
383 // initialisation methods
387
389 G4EmParameters & operator=(const G4EmParameters &right) = delete;
390
391private:
392
394
395 void Initialise();
396
397 G4bool IsLocked() const;
398
399 void PrintWarning(G4ExceptionDescription& ed) const;
400
401 static G4EmParameters* theInstance;
402
403 G4EmParametersMessenger* theMessenger;
404 G4EmExtraParameters* fBParameters;
405 G4EmLowEParameters* fCParameters;
406 G4StateManager* fStateManager;
407 G4EmSaturation* emSaturation;
408
409 G4bool lossFluctuation;
410 G4bool buildCSDARange;
411 G4bool flagLPM;
412 G4bool cutAsFinalRange;
413 G4bool applyCuts;
414 G4bool lateralDisplacement;
415 G4bool lateralDisplacementAlg96;
416 G4bool muhadLateralDisplacement;
417 G4bool useAngGeneratorForIonisation;
418 G4bool useMottCorrection;
419 G4bool integral;
420 G4bool birks;
421 G4bool fICRU90;
422 G4bool gener;
423 G4bool fSamplingTable;
424 G4bool fPolarisation;
425 G4bool fMuDataFromFile;
426 G4bool fPEKShell;
427 G4bool fMscPosiCorr;
428 G4bool onIsolated; // 5d model conversion on free ions
429 G4bool fDNA;
430 G4bool fIsPrinted;
431
432 G4double minKinEnergy;
433 G4double maxKinEnergy;
434 G4double maxKinEnergyCSDA;
435 G4double max5DEnergyForMuPair;
436 G4double lowestElectronEnergy;
437 G4double lowestMuHadEnergy;
438 G4double lowestTripletEnergy;
439 G4double linLossLimit;
440 G4double bremsTh;
441 G4double bremsMuHadTh;
442 G4double lambdaFactor;
443 G4double factorForAngleLimit;
444 G4double thetaLimit;
445 G4double energyLimit;
446 G4double maxNIELEnergy;
447 G4double rangeFactor;
448 G4double rangeFactorMuHad;
449 G4double geomFactor;
450 G4double safetyFactor;
451 G4double lambdaLimit;
452 G4double skin;
453 G4double factorScreen;
454
455 G4int nbinsPerDecade;
456 G4int verbose;
457 G4int workerVerbose;
458 G4int tripletConv; // 5d model triplet generation type
459
460 G4TransportationWithMscType fTransportationWithMsc;
461 G4MscStepLimitType mscStepLimit;
462 G4MscStepLimitType mscStepLimitMuHad;
463 G4NuclearFormfactorType nucFormfactor;
465 G4EmFluctuationType fFluct;
466};
467
468//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
469
470#endif
471
G4DNAModelSubType
G4EmFluoDirectory
G4EmFluctuationType
@ fUniversalFluctuation
@ fUrbanFluctuation
@ fDummyFluctuation
G4eSingleScatteringType
@ fWVI
@ fMott
@ fDPWA
G4TransportationWithMscType
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4MscStepLimitType
G4NuclearFormfactorType
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)
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
G4EmParameters(G4EmParameters &)=delete
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
void SetDirectionalSplittingTarget(const G4ThreeVector &v)
G4bool DNAFast() const
G4DNAModelSubType DNAeSolvationSubType() const
G4bool RetrieveMuDataFromFile() const
G4bool PhotoeffectBelowKShell() 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
G4bool MscPositronCorrection() const
const std::vector< G4String > & RegionsPAI() 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 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)
G4bool ANSTOFluoDir()
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 &)
void SetFluoDirectory(G4EmFluoDirectory)
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
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)
G4EmParameters & operator=(const G4EmParameters &right)=delete
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
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