Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAUeharaScreenedRutherfordElasticModel Class Reference

#include <G4DNAUeharaScreenedRutherfordElasticModel.hh>

+ Inheritance diagram for G4DNAUeharaScreenedRutherfordElasticModel:

Public Member Functions

 G4DNAUeharaScreenedRutherfordElasticModel (const G4ParticleDefinition *p=0, const G4String &nam="DNAUeharaScreenedRutherfordElasticModel")
 
 ~G4DNAUeharaScreenedRutherfordElasticModel () override=default
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SelectFasterComputation (G4bool input)
 
void SetKillBelowThreshold (G4double threshold)
 
G4double GetKillBelowThreshold ()
 
void SelectHighEnergyLimit (G4double threshold)
 
G4DNAUeharaScreenedRutherfordElasticModeloperator= (const G4DNAUeharaScreenedRutherfordElasticModel &right)=delete
 
 G4DNAUeharaScreenedRutherfordElasticModel (const G4DNAUeharaScreenedRutherfordElasticModel &)=delete
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma = nullptr
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Detailed Description

Definition at line 36 of file G4DNAUeharaScreenedRutherfordElasticModel.hh.

Constructor & Destructor Documentation

◆ G4DNAUeharaScreenedRutherfordElasticModel() [1/2]

G4DNAUeharaScreenedRutherfordElasticModel::G4DNAUeharaScreenedRutherfordElasticModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNAUeharaScreenedRutherfordElasticModel" 
)

Definition at line 40 of file G4DNAUeharaScreenedRutherfordElasticModel.cc.

42 : G4VEmModel(nam)
43{
44 // Energy limits of the models
45 intermediateEnergyLimit = 200. * eV;
46 iLowEnergyLimit = 9.*eV;
47 iHighEnergyLimit = 1.*MeV;
48
49 verboseLevel = 0;
50 // Verbosity scale:
51 // 0 = nothing
52 // 1 = warning for energy non-conservation
53 // 2 = details of energy budget
54 // 3 = calculation of cross sections, file openings, sampling of atoms
55 // 4 = entering in methods
56
57 fasterCode = false;
58}

◆ ~G4DNAUeharaScreenedRutherfordElasticModel()

G4DNAUeharaScreenedRutherfordElasticModel::~G4DNAUeharaScreenedRutherfordElasticModel ( )
overridedefault

◆ G4DNAUeharaScreenedRutherfordElasticModel() [2/2]

G4DNAUeharaScreenedRutherfordElasticModel::G4DNAUeharaScreenedRutherfordElasticModel ( const G4DNAUeharaScreenedRutherfordElasticModel )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4DNAUeharaScreenedRutherfordElasticModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  ekin,
G4double  emin,
G4double  emax 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 143 of file G4DNAUeharaScreenedRutherfordElasticModel.cc.

149{
150#ifdef UEHARA_VERBOSE
151 if (verboseLevel > 3)
152 {
153 G4cout
154 << "Calling CrossSectionPerVolume() of G4DNAUeharaScreenedRutherfordElasticModel"
155 << G4endl;
156 }
157#endif
158
159 // Calculate total cross section for model
160
161 G4double sigma = 0.;
162 if(ekin < iLowEnergyLimit || ekin > iHighEnergyLimit) { return sigma; }
163
164 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
165
166 G4double z = 7.42; // FROM PMB 37 (1992) 1841-1858 p1842
167 G4double n = ScreeningFactor(ekin,z);
168 G4double crossSection = RutherfordCrossSection(ekin, z);
169 sigma = pi * crossSection / (n * (n + 1.));
170
171#ifdef UEHARA_VERBOSE
172 if (verboseLevel > 2)
173 {
174 G4cout << "__________________________________" << G4endl;
175 G4cout << "=== G4DNAUeharaScreenedRutherfordElasticModel - XS INFO START"
176 << G4endl;
177 G4cout << "=== Kinetic energy(eV)=" << ekin/eV
178 << " particle : " << particleDefinition->GetParticleName() << G4endl;
179 G4cout << "=== Cross section per water molecule (cm^2)=" << sigma/cm/cm
180 << G4endl;
181 G4cout << "=== Cross section per water molecule (cm^-1)="
182 << sigma*waterDensity/(1./cm) << G4endl;
183 G4cout << "=== G4DNAUeharaScreenedRutherfordElasticModel - XS INFO END"
184 << G4endl;
185 }
186#endif
187
188 return sigma*waterDensity;
189}
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
size_t GetIndex() const
Definition: G4Material.hh:255
const G4double pi

◆ GetKillBelowThreshold()

G4double G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold ( )
inline

Definition at line 152 of file G4DNAUeharaScreenedRutherfordElasticModel.hh.

153{
155 errMsg << "*** WARNING : "
156 << "G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold"
157 << "is deprecated, the returned value is nonsense";
158
160 "G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold",
161 "DEPRECATED", JustWarning, errMsg);
162
163 return -1;
164}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40

◆ Initialise()

void G4DNAUeharaScreenedRutherfordElasticModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector  
)
overridevirtual

Implements G4VEmModel.

Definition at line 63 of file G4DNAUeharaScreenedRutherfordElasticModel.cc.

66{
67 if(isInitialised) { return; }
68 if (verboseLevel > 3)
69 {
70 }
71
72 if(particle->GetParticleName() != "e-")
73 {
74 G4Exception("*** WARNING: the G4DNAUeharaScreenedRutherfordElasticModel is "
75 "not intented to be used with another particle than the electron",
76 "",FatalException,"") ;
77 }
78
79
80 if( verboseLevel>1 )
81 {
82 G4cout << "G4DNAUeharaScreenedRutherfordElasticModel::Initialise()"
83 << G4endl;
84 G4cout << "Energy range: "
85 << LowEnergyLimit() / eV << " eV - "
86 << HighEnergyLimit() / MeV << " MeV"
87 << G4endl;
88 }
89 // Constants for final state by Brenner & Zaider
90 // Note: the instantiation must be placed after if (isInitialised)
91
92 betaCoeff=
93 {
94 7.51525,
95 -0.41912,
96 7.2017E-3,
97 -4.646E-5,
98 1.02897E-7};
99
100 deltaCoeff=
101 {
102 2.9612,
103 -0.26376,
104 4.307E-3,
105 -2.6895E-5,
106 5.83505E-8};
107
108 gamma035_10Coeff=
109 {
110 -1.7013,
111 -1.48284,
112 0.6331,
113 -0.10911,
114 8.358E-3,
115 -2.388E-4};
116
117 gamma10_100Coeff =
118 {
119 -3.32517,
120 0.10996,
121 -4.5255E-3,
122 5.8372E-5,
123 -2.4659E-7};
124
125 gamma100_200Coeff=
126 {
127 2.4775E-2,
128 -2.96264E-5,
129 -1.20655E-7};
130
131 // Initialize water density pointer
132 fpWaterDensity =
134 GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
135
137 isInitialised = true;
138}
@ FatalException
static G4DNAMolecularMaterial * Instance()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:691
const G4String & GetParticleName() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634

◆ operator=()

G4DNAUeharaScreenedRutherfordElasticModel & G4DNAUeharaScreenedRutherfordElasticModel::operator= ( const G4DNAUeharaScreenedRutherfordElasticModel right)
delete

◆ SampleSecondaries()

void G4DNAUeharaScreenedRutherfordElasticModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicElectron,
G4double  tmin,
G4double  maxEnergy 
)
overridevirtual

Implements G4VEmModel.

Definition at line 248 of file G4DNAUeharaScreenedRutherfordElasticModel.cc.

254{
255#ifdef UEHARA_VERBOSE
256 if (verboseLevel > 3)
257 {
258 G4cout
259 << "Calling SampleSecondaries() of G4DNAUeharaScreenedRutherfordElasticModel"
260 << G4endl;
261 }
262#endif
263
264 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
265 if(electronEnergy0 < iLowEnergyLimit || electronEnergy0 > iHighEnergyLimit)
266 return;
267
268 G4double cosTheta = 0.;
269
270 if (electronEnergy0<intermediateEnergyLimit)
271 {
272#ifdef UEHARA_VERBOSE
273 if (verboseLevel > 3)
274 G4cout << "---> Using Brenner & Zaider model" << G4endl;
275#endif
276 cosTheta = BrennerZaiderRandomizeCosTheta(electronEnergy0);
277 }
278 else //if (electronEnergy0>=intermediateEnergyLimit)
279 {
280#ifdef UEHARA_VERBOSE
281 if (verboseLevel > 3)
282 G4cout << "---> Using Screened Rutherford model" << G4endl;
283#endif
284 G4double z = 7.42; // FROM PMB 37 (1992) 1841-1858 p1842
285 cosTheta = ScreenedRutherfordRandomizeCosTheta(electronEnergy0,z);
286 }
287
288 G4double phi = 2. * pi * G4UniformRand();
289
290 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
291 G4ThreeVector xVers = zVers.orthogonal();
292 G4ThreeVector yVers = zVers.cross(xVers);
293
294 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
295 G4double yDir = xDir;
296 xDir *= std::cos(phi);
297 yDir *= std::sin(phi);
298
299 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
300
302
304}
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)

◆ SelectFasterComputation()

void G4DNAUeharaScreenedRutherfordElasticModel::SelectFasterComputation ( G4bool  input)
inline

Definition at line 114 of file G4DNAUeharaScreenedRutherfordElasticModel.hh.

116{
117 fasterCode = input;
118}

◆ SelectHighEnergyLimit()

void G4DNAUeharaScreenedRutherfordElasticModel::SelectHighEnergyLimit ( G4double  threshold)
inline

Definition at line 124 of file G4DNAUeharaScreenedRutherfordElasticModel.hh.

126{
127 if(threshold > 10. * CLHEP::keV)
128 {
130 "*** WARNING : the G4DNAUeharaScreenedRutherfordElasticModel class is "
131 "used above 10 keV !",
132 "", JustWarning, "");
133 }
134
135 SetHighEnergyLimit(threshold);
136}
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746

◆ SetKillBelowThreshold()

void G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold ( G4double  threshold)
inline

Definition at line 139 of file G4DNAUeharaScreenedRutherfordElasticModel.hh.

140{
142 errMsg << "*** WARNING : "
143 << "G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold"
144 << "is deprecated, the kill threshold won't be taken into account";
145
147 "G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold",
148 "DEPRECATED", JustWarning, errMsg);
149}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAUeharaScreenedRutherfordElasticModel::fParticleChangeForGamma = nullptr
protected

Definition at line 89 of file G4DNAUeharaScreenedRutherfordElasticModel.hh.

Referenced by Initialise(), and SampleSecondaries().


The documentation for this class was generated from the following files: