Geant4 9.6.0
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)
 
G4double GetKillBelowThreshold ()
 
- 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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
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
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

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 39 of file G4DNAScreenedRutherfordElasticModel.hh.

Constructor & Destructor Documentation

◆ G4DNAScreenedRutherfordElasticModel()

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

Definition at line 40 of file G4DNAScreenedRutherfordElasticModel.cc.

42 :G4VEmModel(nam),isInitialised(false)
43{
44 // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
45 fpWaterDensity = 0;
46
47 killBelowEnergy = 9*eV;
48 lowEnergyLimit = 0 * eV;
49 intermediateEnergyLimit = 200 * eV; // Switch between two final state models
50 highEnergyLimit = 1. * MeV;
51 SetLowEnergyLimit(lowEnergyLimit);
52 SetHighEnergyLimit(highEnergyLimit);
53
54 verboseLevel= 0;
55 // Verbosity scale:
56 // 0 = nothing
57 // 1 = warning for energy non-conservation
58 // 2 = details of energy budget
59 // 3 = calculation of cross sections, file openings, sampling of atoms
60 // 4 = entering in methods
61
62 if( verboseLevel>0 )
63 {
64 G4cout << "Screened Rutherford Elastic model is constructed " << G4endl
65 << "Energy range: "
66 << lowEnergyLimit / eV << " eV - "
67 << highEnergyLimit / MeV << " MeV"
68 << G4endl;
69 }
71}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592

◆ ~G4DNAScreenedRutherfordElasticModel()

G4DNAScreenedRutherfordElasticModel::~G4DNAScreenedRutherfordElasticModel ( )
virtual

Definition at line 75 of file G4DNAScreenedRutherfordElasticModel.cc.

76{}

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 156 of file G4DNAScreenedRutherfordElasticModel.cc.

161{
162 if (verboseLevel > 3)
163 G4cout << "Calling CrossSectionPerVolume() of G4DNAScreenedRutherfordElasticModel" << G4endl;
164
165 // Calculate total cross section for model
166
167 G4double sigma=0;
168
169 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
170
171 if(waterDensity!= 0.0)
172 // if (material == nistwater || material->GetBaseMaterial() == nistwater)
173 {
174
175 if (ekin < highEnergyLimit)
176 {
177
178 if (ekin < killBelowEnergy) return DBL_MAX;
179
180 G4double z = 10.;
181 G4double n = ScreeningFactor(ekin,z);
182 G4double crossSection = RutherfordCrossSection(ekin, z);
183 sigma = pi * crossSection / (n * (n + 1.));
184 }
185
186 if (verboseLevel > 2)
187 {
188 G4cout << "__________________________________" << G4endl;
189 G4cout << "°°° G4DNAScreenedRutherfordElasticModel - XS INFO START" << G4endl;
190 G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
191 G4cout << "°°° Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
192 G4cout << "°°° Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
193 // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
194 G4cout << "°°° G4DNAScreenedRutherfordElasticModel - XS INFO END" << G4endl;
195 }
196
197 }
198
199 return sigma*material->GetAtomicNumDensityVector()[1];
200}
double G4double
Definition: G4Types.hh:64
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:215
size_t GetIndex() const
Definition: G4Material.hh:261
const G4double pi
#define DBL_MAX
Definition: templates.hh:83

◆ GetKillBelowThreshold()

G4double G4DNAScreenedRutherfordElasticModel::GetKillBelowThreshold ( )
inline

Definition at line 64 of file G4DNAScreenedRutherfordElasticModel.hh.

64{ return killBelowEnergy; }

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 80 of file G4DNAScreenedRutherfordElasticModel.cc.

82{
83
84 if (verboseLevel > 3)
85 G4cout << "Calling G4DNAScreenedRutherfordElasticModel::Initialise()" << G4endl;
86
87 // Energy limits
88
89 if (LowEnergyLimit() < lowEnergyLimit)
90 {
91 G4cout << "G4DNAScreenedRutherfordElasticModel: low energy limit increased from " <<
92 LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
93 SetLowEnergyLimit(lowEnergyLimit);
94 }
95
96 if (HighEnergyLimit() > highEnergyLimit)
97 {
98 G4cout << "G4DNAScreenedRutherfordElasticModel: high energy limit decreased from " <<
99 HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl;
100 SetHighEnergyLimit(highEnergyLimit);
101 }
102
103 // Constants for final stae by Brenner & Zaider
104
105 betaCoeff.push_back(7.51525);
106 betaCoeff.push_back(-0.41912);
107 betaCoeff.push_back(7.2017E-3);
108 betaCoeff.push_back(-4.646E-5);
109 betaCoeff.push_back(1.02897E-7);
110
111 deltaCoeff.push_back(2.9612);
112 deltaCoeff.push_back(-0.26376);
113 deltaCoeff.push_back(4.307E-3);
114 deltaCoeff.push_back(-2.6895E-5);
115 deltaCoeff.push_back(5.83505E-8);
116
117 gamma035_10Coeff.push_back(-1.7013);
118 gamma035_10Coeff.push_back(-1.48284);
119 gamma035_10Coeff.push_back(0.6331);
120 gamma035_10Coeff.push_back(-0.10911);
121 gamma035_10Coeff.push_back(8.358E-3);
122 gamma035_10Coeff.push_back(-2.388E-4);
123
124 gamma10_100Coeff.push_back(-3.32517);
125 gamma10_100Coeff.push_back(0.10996);
126 gamma10_100Coeff.push_back(-4.5255E-3);
127 gamma10_100Coeff.push_back(5.8372E-5);
128 gamma10_100Coeff.push_back(-2.4659E-7);
129
130 gamma100_200Coeff.push_back(2.4775E-2);
131 gamma100_200Coeff.push_back(-2.96264E-5);
132 gamma100_200Coeff.push_back(-1.20655E-7);
133
134 //
135
136 if( verboseLevel>0 )
137 {
138 G4cout << "Screened Rutherford elastic model is initialized " << G4endl
139 << "Energy range: "
140 << LowEnergyLimit() / eV << " eV - "
141 << HighEnergyLimit() / MeV << " MeV"
142 << G4endl;
143 }
144
145 // Initialize water density pointer
147
148 if (isInitialised) { return; }
150 isInitialised = true;
151
152}
static G4DNAMolecularMaterial * Instance()
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:576
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 255 of file G4DNAScreenedRutherfordElasticModel.cc.

260{
261
262 if (verboseLevel > 3)
263 G4cout << "Calling SampleSecondaries() of G4DNAScreenedRutherfordElasticModel" << G4endl;
264
265 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
266
267 if (electronEnergy0 < killBelowEnergy)
268 {
272 return ;
273 }
274
275 G4double cosTheta = 0.;
276
277 if (electronEnergy0>= killBelowEnergy && electronEnergy0 < highEnergyLimit)
278 {
279 if (electronEnergy0<intermediateEnergyLimit)
280 {
281 if (verboseLevel > 3) G4cout << "---> Using Brenner & Zaider model" << G4endl;
282 cosTheta = BrennerZaiderRandomizeCosTheta(electronEnergy0);
283 }
284
285 if (electronEnergy0>=intermediateEnergyLimit)
286 {
287 if (verboseLevel > 3) G4cout << "---> Using Screened Rutherford model" << G4endl;
288 G4double z = 10.;
289 cosTheta = ScreenedRutherfordRandomizeCosTheta(electronEnergy0,z);
290 }
291
292 G4double phi = 2. * pi * G4UniformRand();
293
294 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
295 G4ThreeVector xVers = zVers.orthogonal();
296 G4ThreeVector yVers = zVers.cross(xVers);
297
298 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
299 G4double yDir = xDir;
300 xDir *= std::cos(phi);
301 yDir *= std::sin(phi);
302
303 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
304
306
308 }
309
310}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SetKillBelowThreshold()

void G4DNAScreenedRutherfordElasticModel::SetKillBelowThreshold ( G4double  threshold)
inline

Definition at line 108 of file G4DNAScreenedRutherfordElasticModel.hh.

109{
110 killBelowEnergy = threshold;
111 if (threshold < 9*CLHEP::eV)
112 G4Exception ("*** WARNING : the G4DNAScreenedRutherfordElasticModel class is not validated below 9 eV !","",JustWarning,"") ;
113}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAScreenedRutherfordElasticModel::fParticleChangeForGamma
protected

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