Geant4 11.4.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 "G4EmFluoDirectory.hh"
61#include "G4EmSaturation.hh"
62#include "G4ThreeVector.hh"
64#include <vector>
65#include <map>
72
79
86
94
100class G4VEmProcess;
101class G4StateManager;
102
104{
105public:
106
107 static G4EmParameters* Instance();
108
110
111 void SetDefaults();
112
113 // printing
114 void StreamInfo(std::ostream& os) const;
115 void Dump();
116 friend std::ostream& operator<< (std::ostream& os, const G4EmParameters&);
117
118 // boolean flags
119 void SetLossFluctuations(G4bool val);
120 G4bool LossFluctuation() const;
121
122 void SetBuildCSDARange(G4bool val);
123 G4bool BuildCSDARange() const;
124
125 void SetLPM(G4bool val);
126 G4bool LPM() const;
127
130
131 void SetApplyCuts(G4bool val);
132 G4bool ApplyCuts() const;
133
134 void SetFluo(G4bool val);
135 G4bool Fluo() const;
136
138
140 void SetBeardenFluoDir(G4bool val);
141 void SetANSTOFluoDir(G4bool val);
142 void SetXDB_EADLFluoDir(G4bool val);
143
146
147 void SetAuger(G4bool val);
148 void SetAugerCascade(G4bool val) { SetAuger(val); };
149 G4bool Auger() const;
150 G4bool AugerCascade() const { return Auger(); }
151
152 void SetPixe(G4bool val);
153 G4bool Pixe() const;
154
157
160
163
166
169
172
173 void SetIntegral(G4bool val);
174 G4bool Integral() const;
175
176 void SetBirksActive(G4bool val);
177 G4bool BirksActive() const;
178
179 void SetUseICRU90Data(G4bool val);
180 G4bool UseICRU90Data() const;
181
184
187
188 void SetDNAFast(G4bool val);
189 G4bool DNAFast() const;
190
191 void SetDNAStationary(G4bool val);
192 G4bool DNAStationary() const;
193
194 void SetDNAElectronMsc(G4bool val);
195 G4bool DNAElectronMsc() const;
196
197 // if general interaction is enabled then
198 // force interaction options should be disabled
201
204
207
210
213
216
219
222
223 G4bool UseEPICS2017XS() const;
225
228
231
232 // 5d
233 void SetOnIsolated(G4bool val);
234 G4bool OnIsolated() const;
235
236 void ActivateDNA();
237 void SetIsPrintedFlag(G4bool val);
238 G4bool IsPrintLocked() const;
239
240 // double parameters with values
241 void SetMinEnergy(G4double val);
242 G4double MinKinEnergy() const;
243
244 void SetMaxEnergy(G4double val);
245 G4double MaxKinEnergy() const;
246
249
252
255
258
261
266
267 void SetLambdaFactor(G4double val);
268 G4double LambdaFactor() const;
269
272
273 void SetMscThetaLimit(G4double val);
274 G4double MscThetaLimit() const;
275
276 void SetMscEnergyLimit(G4double val);
277 G4double MscEnergyLimit() const;
278
279 void SetMscRangeFactor(G4double val);
280 G4double MscRangeFactor() const;
281
284
285 void SetMscGeomFactor(G4double val);
286 G4double MscGeomFactor() const;
287
290
291 void SetMscLambdaLimit(G4double val);
292 G4double MscLambdaLimit() const;
293
294 void SetMscSkin(G4double val);
295 G4double MscSkin() const;
296
299
300 void SetMaxNIELEnergy(G4double val);
301 G4double MaxNIELEnergy() const;
302
305
308
311
312 void SetStepFunction(G4double v1, G4double v2);
317
320
323
324 // integer parameters
325
328 G4int NumberOfBins() const;
329
330 void SetVerbose(G4int val);
331 G4int Verbose() const;
332
333 void SetWorkerVerbose(G4int val);
334 G4int WorkerVerbose() const;
335
338
341
344
347
350
353
356
357 //DNA chemistry model
358 void SetTimeStepModel(const G4ChemTimeStepModel& model);
360 //5d
361 void SetConversionType(G4int val);
362 G4int GetConversionType() const;
363
364 // string parameters
367
370
371 void SetLivermoreDataDir(const G4String&);
372 const G4String& LivermoreDataDir();
373
374 // parameters per region or per process
375 void AddPAIModel(const G4String& particle,
376 const G4String& region,
377 const G4String& type);
378 const std::vector<G4String>& ParticlesPAI() const;
379 const std::vector<G4String>& RegionsPAI() const;
380 const std::vector<G4String>& TypesPAI() const;
381
382 void AddMicroElec(const G4String& region);
383 const std::vector<G4String>& RegionsMicroElec() const;
384
385 void AddDNA(const G4String& region, const G4String& type);
386 const std::vector<G4String>& RegionsDNA() const;
387 const std::vector<G4String>& TypesDNA() const;
388
389 void AddPhysics(const G4String& region, const G4String& type);
390 const std::vector<G4String>& RegionsPhysics() const;
391 const std::vector<G4String>& TypesPhysics() const;
392
393 void SetSubCutRegion(const G4String& region = "");
394
395 void SetDeexActiveRegion(const G4String& region, G4bool fdeex,
396 G4bool fauger, G4bool fpixe);
397
398 void SetProcessBiasingFactor(const G4String& procname,
399 G4double val, G4bool wflag);
400
401 void ActivateForcedInteraction(const G4String& procname,
402 const G4String& region,
403 G4double length,
404 G4bool wflag);
405
406 void ActivateSecondaryBiasing(const G4String& name,
407 const G4String& region,
408 G4double factor,
409 G4double energyLimit);
410
411 // define external saturation class
413 // create and access saturation class
415
416 // defined fluctuations per G4Region
417 void SetFluctuationsForRegion(const G4String& regionName, G4bool flag);
418
419 // initialisation methods
423 void DefineFluctuationFlags(std::vector<G4bool>* theFluctFlags);
424
425 const G4String& GetDirLEDATA() const;
426
428 G4EmParameters & operator=(const G4EmParameters &right) = delete;
429
430private:
431
433
434 void Initialise();
435
436 G4bool IsLocked() const;
437
438 void PrintWarning(G4ExceptionDescription& ed) const;
439
440 static G4EmParameters* theInstance;
441
442 G4EmParametersMessenger* theMessenger;
443 G4EmExtraParameters* fBParameters;
444 G4EmLowEParameters* fCParameters;
445 G4StateManager* fStateManager;
446 G4EmSaturation* emSaturation;
447
448 G4bool lossFluctuation;
449 G4bool buildCSDARange;
450 G4bool flagLPM;
451 G4bool cutAsFinalRange;
452 G4bool applyCuts;
453 G4bool lateralDisplacement;
454 G4bool lateralDisplacementAlg96;
455 G4bool muhadLateralDisplacement;
456 G4bool useAngGeneratorForIonisation;
457 G4bool useMottCorrection;
458 G4bool integral;
459 G4bool birks;
460 G4bool fICRU90;
461 G4bool gener;
462 G4bool fSamplingTable;
463 G4bool fPolarisation;
464 G4bool fMuDataFromFile;
465 G4bool fPEKShell;
466 G4bool fMscPosiCorr;
467 G4bool fUseEPICS2017XS;
468 G4bool f3GammaAnnihilationOnFly;
469 G4bool fUseRiGePairProductionModel;
470 G4bool onIsolated; // 5d model conversion on free ions
471 G4bool fDNA;
472 G4bool fIsPrinted;
473
474 G4double minKinEnergy;
475 G4double maxKinEnergy;
476 G4double maxKinEnergyCSDA;
477 G4double max5DEnergyForMuPair;
478 G4double lowestElectronEnergy;
479 G4double lowestMuHadEnergy;
480 G4double lowestTripletEnergy;
481 G4double linLossLimit;
482 G4double bremsTh;
483 G4double bremsMuHadTh;
484 G4double lambdaFactor;
485 G4double factorForAngleLimit;
486 G4double thetaLimit;
487 G4double energyLimit;
488 G4double maxNIELEnergy;
489 G4double rangeFactor;
490 G4double rangeFactorMuHad;
491 G4double geomFactor;
492 G4double safetyFactor;
493 G4double lambdaLimit;
494 G4double skin;
495 G4double factorScreen;
496
497 G4int nbinsPerDecade;
498 G4int verbose;
499 G4int workerVerbose;
500 G4int nForFreeVector;
501 G4int tripletConv; // 5d model triplet generation type
502
503 G4TransportationWithMscType fTransportationWithMsc;
504 G4MscStepLimitType mscStepLimit;
505 G4MscStepLimitType mscStepLimitMuHad;
506 G4NuclearFormfactorType nucFormfactor;
508 G4EmFluctuationType fFluct;
509 G4PositronAtRestModelType fPositronium;
510
511 G4String fDirLEDATA;
512 std::vector<std::pair<G4String, G4bool> > fluctRegions;
513};
514
515//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
516
517#endif
518
G4ChemTimeStepModel
G4DNAModelSubType
G4EmFluoDirectory
G4EmFluctuationType
@ fUniversalFluctuation
@ fUrbanFluctuation
@ fDummyFluctuation
G4eSingleScatteringType
@ fWVI
@ fMott
@ fDPWA
G4TransportationWithMscType
G4PositronAtRestModelType
@ fSimplePositronium
@ fAllisonPositronium
@ fOrePowell
@ fOrePowellPolar
std::ostringstream G4ExceptionDescription
G4MscStepLimitType
CLHEP::Hep3Vector G4ThreeVector
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)
G4bool UseEPICS2017XS() const
static G4EmParameters * Instance()
void SetDirectionalSplittingTarget(const G4ThreeVector &v)
G4bool DNAFast() const
G4DNAModelSubType DNAeSolvationSubType() const
G4bool RetrieveMuDataFromFile() const
G4bool PhotoeffectBelowKShell() const
G4int NumberOfBins() const
G4double MscMuHadRangeFactor() const
void SetGeneralProcessActive(G4bool val)
void SetDNAFast(G4bool val)
void SetMscSafetyFactor(G4double val)
G4bool EnablePolarisation() const
void SetMaxDNAElectronEnergy(G4double val)
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 SetPositronAtRestModelType(G4PositronAtRestModelType val)
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
G4int NumberForFreeVector() 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 SetNumberForFreeVector(G4int val)
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)
G4double MaxDNAElectronEnergy() const
void SetUseCutAsFinalRange(G4bool val)
void SetDeexcitationIgnoreCut(G4bool val)
void SetBirksActive(G4bool val)
G4double ScreeningFactor() const
G4bool UseMottCorrection() const
void SetFluo(G4bool val)
G4double MaxDNAIonEnergy() const
void SetMuHadBremsstrahlungTh(G4double val)
const std::vector< G4String > & ParticlesPAI() const
G4bool Pixe() const
const std::vector< G4String > & RegionsPhysics() const
void SetFluctuationType(G4EmFluctuationType val)
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
void SetUseRiGePairProductionModel(G4bool v)
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)
G4PositronAtRestModelType PositronAtRestModelType() const
friend std::ostream & operator<<(std::ostream &os, const G4EmParameters &)
void SetMaxEnergyForCSDARange(G4double val)
G4bool AugerCascade() const
G4eSingleScatteringType SingleScatteringType() const
void SetMaxDNAIonEnergy(G4double val)
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)
void SetUseEPICS2017XS(G4bool v)
void DefineFluctuationFlags(std::vector< G4bool > *theFluctFlags)
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)
G4ChemTimeStepModel GetTimeStepModel() const
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
G4bool UseRiGePairProductionModel() const
void SetLowestMuHadEnergy(G4double val)
void SetFluctuationsForRegion(const G4String &regionName, G4bool flag)
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
G4bool Use3GammaAnnihilationOnFly() const
const G4String & GetDirLEDATA() 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)
void SetTimeStepModel(const G4ChemTimeStepModel &model)
void Set3GammaAnnihilationOnFly(G4bool v)
G4bool GeneralProcessActive() const