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

#include <G4DNAScreenedRutherfordElasticModel.hh>

+ Inheritance diagram for G4DNAScreenedRutherfordElasticModel:

Public Member Functions

 G4DNAScreenedRutherfordElasticModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="DNAScreenedRutherfordElasticModel")
 
 ~G4DNAScreenedRutherfordElasticModel () override
 
G4DNAScreenedRutherfordElasticModeloperator= (const G4DNAScreenedRutherfordElasticModel &right)=delete
 
 G4DNAScreenedRutherfordElasticModel (const G4DNAScreenedRutherfordElasticModel &)=delete
 
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 SetKillBelowThreshold (G4double threshold)
 
void SelectFasterComputation (G4bool input)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
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 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 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 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)
 
void SetLPMFlag (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- 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
 
std::size_t currentCoupleIndex = 0
 
std::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 38 of file G4DNAScreenedRutherfordElasticModel.hh.

Constructor & Destructor Documentation

◆ G4DNAScreenedRutherfordElasticModel() [1/2]

G4DNAScreenedRutherfordElasticModel::G4DNAScreenedRutherfordElasticModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNAScreenedRutherfordElasticModel" )

Definition at line 41 of file G4DNAScreenedRutherfordElasticModel.cc.

43 :
44G4VEmModel(nam)
45{
46 fpWaterDensity = nullptr;
47
48 lowEnergyLimit = 0 * eV;
49 intermediateEnergyLimit = 200 * eV; // Switch between two final state models
50 highEnergyLimit = 1. * MeV;
51
52 SetLowEnergyLimit(lowEnergyLimit);
53 SetHighEnergyLimit(highEnergyLimit);
54
55 verboseLevel = 0;
56 // Verbosity scale:
57 // 0 = nothing
58 // 1 = warning for energy non-conservation
59 // 2 = details of energy budget
60 // 3 = calculation of cross sections, file openings, sampling of atoms
61 // 4 = entering in methods
62
63#ifdef SR_VERBOSE
64 if (verboseLevel > 0)
65 {
66 G4cout << "Screened Rutherford Elastic model is constructed "
67 << G4endl
68 << "Energy range: "
69 << lowEnergyLimit / eV << " eV - "
70 << highEnergyLimit / MeV << " MeV"
71 << G4endl;
72 }
73#endif
75
76 // Selection of computation method
77 // We do not recommend "true" usage with the current cumul. proba. settings
78 fasterCode = false;
79}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void SetHighEnergyLimit(G4double)
void SetLowEnergyLimit(G4double)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4DNAScreenedRutherfordElasticModel()

G4DNAScreenedRutherfordElasticModel::~G4DNAScreenedRutherfordElasticModel ( )
overridedefault

◆ G4DNAScreenedRutherfordElasticModel() [2/2]

G4DNAScreenedRutherfordElasticModel::G4DNAScreenedRutherfordElasticModel ( const G4DNAScreenedRutherfordElasticModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 189 of file G4DNAScreenedRutherfordElasticModel.cc.

199{
200#ifdef SR_VERBOSE
201 if (verboseLevel > 3)
202 {
203 G4cout << "Calling CrossSectionPerVolume() of "
204 "G4DNAScreenedRutherfordElasticModel"
205 << G4endl;
206 }
207#endif
208
209 // Calculate total cross section for model
210
211 G4double sigma=0.;
212 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
213
214 if(ekin <= HighEnergyLimit() && ekin >= LowEnergyLimit())
215 {
216 G4double z = 10.;
217 G4double n = ScreeningFactor(ekin,z);
218 G4double crossSection = RutherfordCrossSection(ekin, z);
219 sigma = pi * crossSection / (n * (n + 1.));
220 }
221
222#ifdef SR_VERBOSE
223 if (verboseLevel > 2)
224 {
225 G4cout << "__________________________________" << G4endl;
226 G4cout << "=== G4DNAScreenedRutherfordElasticModel - XS INFO START"
227 << G4endl;
228 G4cout << "=== Kinetic energy(eV)=" << ekin/eV
229 << " particle : " << particleDefinition->GetParticleName()
230 << G4endl;
231 G4cout << "=== Cross section per water molecule (cm^2)=" << sigma/cm/cm
232 << G4endl;
233 G4cout << "=== Cross section per water molecule (cm^-1)="
234 << sigma*waterDensity/(1./cm) << G4endl;
235 G4cout << "=== G4DNAScreenedRutherfordElasticModel - XS INFO END"
236 << G4endl;
237 }
238#endif
239
240 return sigma*waterDensity;
241}
double G4double
Definition G4Types.hh:83
std::size_t GetIndex() const
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const

◆ Initialise()

void G4DNAScreenedRutherfordElasticModel::Initialise ( const G4ParticleDefinition * particle,
const G4DataVector &  )
overridevirtual

Implements G4VEmModel.

Definition at line 88 of file G4DNAScreenedRutherfordElasticModel.cc.

91{
92#ifdef SR_VERBOSE
93 if (verboseLevel > 3)
94 {
95 G4cout << "Calling G4DNAScreenedRutherfordElasticModel::Initialise()"
96 << G4endl;
97 }
98#endif
99
100 if(particle->GetParticleName() != "e-")
101 {
102 G4Exception ("*** WARNING: the G4DNAScreenedRutherfordElasticModel is not "
103 "intented to be used with another particle than the electron",
104 "",FatalException,"") ;
105 }
106
107 // Energy limits
108 if (LowEnergyLimit() < 9*eV)
109 {
110 G4Exception("*** WARNING: the G4DNAScreenedRutherfordElasticModel class is "
111 "not validated below 9 eV",
112 "",JustWarning,"") ;
113 }
114
115 if (HighEnergyLimit() > 1*MeV)
116 {
117 G4Exception("*** WARNING: the G4DNAScreenedRutherfordElasticModel class is "
118 "not validated above 1 MeV",
119 "",JustWarning,"") ;
120 }
121
122 //
123#ifdef SR_VERBOSE
124 if( verboseLevel>0 )
125 {
126 G4cout << "Screened Rutherford elastic model is initialized " << G4endl
127 << "Energy range: "
128 << LowEnergyLimit() / eV << " eV - "
129 << HighEnergyLimit() / MeV << " MeV"
130 << G4endl;
131 }
132#endif
133
134 if (isInitialised) { return; } // return here, prevent reinit consts + pointer
135
136 // Initialize water density pointer
137 fpWaterDensity = G4DNAMolecularMaterial::Instance()->
138 GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
139
141 isInitialised = true;
142
143 // Constants for final state by Brenner & Zaider
144 // note: if called after if(isInitialised) no need for clear and resetting
145 // the values at every call
146
147 betaCoeff=
148 {
149 7.51525,
150 -0.41912,
151 7.2017E-3,
152 -4.646E-5,
153 1.02897E-7};
154
155 deltaCoeff=
156 {
157 2.9612,
158 -0.26376,
159 4.307E-3,
160 -2.6895E-5,
161 5.83505E-8};
162
163 gamma035_10Coeff =
164 {
165 -1.7013,
166 -1.48284,
167 0.6331,
168 -0.10911,
169 8.358E-3,
170 -2.388E-4};
171
172 gamma10_100Coeff =
173 {
174 -3.32517,
175 0.10996,
176 -4.5255E-3,
177 5.8372E-5,
178 -2.4659E-7};
179
180 gamma100_200Coeff =
181 {
182 2.4775E-2,
183 -2.96264E-5,
184 -1.20655E-7};
185}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
static G4DNAMolecularMaterial * Instance()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
const G4String & GetParticleName() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 300 of file G4DNAScreenedRutherfordElasticModel.cc.

306{
307#ifdef SR_VERBOSE
308 if (verboseLevel > 3)
309 {
310 G4cout << "Calling SampleSecondaries() of "
311 "G4DNAScreenedRutherfordElasticModel"
312 << G4endl;
313 }
314#endif
315
316 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
317 G4double cosTheta = 0.;
318
319 if (electronEnergy0<intermediateEnergyLimit)
320 {
321#ifdef SR_VERBOSE
322 if (verboseLevel > 3)
323 {G4cout << "---> Using Brenner & Zaider model" << G4endl;}
324#endif
325 cosTheta = BrennerZaiderRandomizeCosTheta(electronEnergy0);
326 }
327
328 if (electronEnergy0>=intermediateEnergyLimit)
329 {
330#ifdef SR_VERBOSE
331 if (verboseLevel > 3)
332 {G4cout << "---> Using Screened Rutherford model" << G4endl;}
333#endif
334 G4double z = 10.;
335 cosTheta = ScreenedRutherfordRandomizeCosTheta(electronEnergy0,z);
336 }
337
338 G4double phi = 2. * pi * G4UniformRand();
339
340 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
341 G4ThreeVector xVers = zVers.orthogonal();
342 G4ThreeVector yVers = zVers.cross(xVers);
343
344 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
345 G4double yDir = xDir;
346 xDir *= std::cos(phi);
347 yDir *= std::sin(phi);
348
349 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
350
352
354}
#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 G4DNAScreenedRutherfordElasticModel::SelectFasterComputation ( G4bool input)
inline

Definition at line 119 of file G4DNAScreenedRutherfordElasticModel.hh.

120{
121 fasterCode = input;
122}

◆ SetKillBelowThreshold()

void G4DNAScreenedRutherfordElasticModel::SetKillBelowThreshold ( G4double threshold)
inline

Definition at line 106 of file G4DNAScreenedRutherfordElasticModel.hh.

107{
109 errMsg << "The method G4DNAScreenedRutherfordElasticModel::"
110 "SetKillBelowThreshold is deprecated";
111
112 G4Exception("G4DNAScreenedRutherfordElasticModel::SetKillBelowThreshold",
113 "deprecated",
115 errMsg);
116}
std::ostringstream G4ExceptionDescription

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAScreenedRutherfordElasticModel::fParticleChangeForGamma
protected

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