Geant4 11.1.1
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=0, const G4String &nam="DNAScreenedRutherfordElasticModel")
 
virtual ~G4DNAScreenedRutherfordElasticModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SetKillBelowThreshold (G4double threshold)
 
void SelectFasterComputation (G4bool input)
 
- 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
 
- 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 38 of file G4DNAScreenedRutherfordElasticModel.hh.

Constructor & Destructor Documentation

◆ G4DNAScreenedRutherfordElasticModel()

G4DNAScreenedRutherfordElasticModel::G4DNAScreenedRutherfordElasticModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNAScreenedRutherfordElasticModel" 
)

Definition at line 41 of file G4DNAScreenedRutherfordElasticModel.cc.

43 :
44G4VEmModel(nam), isInitialised(false)
45{
46 fpWaterDensity = 0;
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:57
G4GLOB_DLL std::ostream G4cout
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753

◆ ~G4DNAScreenedRutherfordElasticModel()

G4DNAScreenedRutherfordElasticModel::~G4DNAScreenedRutherfordElasticModel ( )
virtual

Definition at line 83 of file G4DNAScreenedRutherfordElasticModel.cc.

84{
85}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 190 of file G4DNAScreenedRutherfordElasticModel.cc.

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

◆ Initialise()

void G4DNAScreenedRutherfordElasticModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 89 of file G4DNAScreenedRutherfordElasticModel.cc.

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 301 of file G4DNAScreenedRutherfordElasticModel.cc.

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

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAScreenedRutherfordElasticModel::fParticleChangeForGamma
protected

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