Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNAUeharaScreenedRutherfordElasticModel Class Reference

#include <G4DNAUeharaScreenedRutherfordElasticModel.hh>

Inheritance diagram for G4DNAUeharaScreenedRutherfordElasticModel:
Collaboration diagram for G4DNAUeharaScreenedRutherfordElasticModel:

Public Member Functions

 G4DNAUeharaScreenedRutherfordElasticModel (const G4ParticleDefinition *p=0, const G4String &nam="DNAUeharaScreenedRutherfordElasticModel")
 
virtual ~G4DNAUeharaScreenedRutherfordElasticModel ()
 
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 SelectFasterComputation (G4bool input)
 
void SetKillBelowThreshold (G4double threshold)
 
G4double GetKillBelowThreshold ()
 
void SelectHighEnergyLimit (G4double threshold)
 
- 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 *, G4double &eloss, G4double &niel, G4double length)
 
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 ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector
< G4EmElementSelector * > * 
GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
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
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
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 SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 

Protected Attributes

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

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 *)
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 36 of file G4DNAUeharaScreenedRutherfordElasticModel.hh.

Constructor & Destructor Documentation

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

Definition at line 41 of file G4DNAUeharaScreenedRutherfordElasticModel.cc.

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

Here is the call graph for this function:

G4DNAUeharaScreenedRutherfordElasticModel::~G4DNAUeharaScreenedRutherfordElasticModel ( )
virtual

Definition at line 81 of file G4DNAUeharaScreenedRutherfordElasticModel.cc.

82 {
83 }

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 190 of file G4DNAUeharaScreenedRutherfordElasticModel.cc.

195 {
196 #ifdef UEHARA_VERBOSE
197  if (verboseLevel > 3)
198  {
199  G4cout
200  << "Calling CrossSectionPerVolume() of G4DNAUeharaScreenedRutherfordElasticModel"
201  << G4endl;
202  }
203 #endif
204 
205  // Calculate total cross section for model
206 
207  G4double sigma = 0.;
208  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
209 
210  // Use activation boundaries instead
211  //if(ekin > HighEnergyLimit() || ekin < LowEnergyLimit()) return 0.;
212 
213  if(waterDensity!= 0.0)
214  {
215  G4double z = 7.42; // FROM PMB 37 (1992) 1841-1858 p1842
216  G4double n = ScreeningFactor(ekin,z);
217  G4double crossSection = RutherfordCrossSection(ekin, z);
218  sigma = pi * crossSection / (n * (n + 1.));
219 
220 #ifdef UEHARA_VERBOSE
221  if (verboseLevel > 2)
222  {
223  G4cout << "__________________________________" << G4endl;
224  G4cout << "=== G4DNAUeharaScreenedRutherfordElasticModel - XS INFO START"
225  << G4endl;
226  G4cout << "=== Kinetic energy(eV)=" << ekin/eV
227  << " particle : " << particleDefinition->GetParticleName() << G4endl;
228  G4cout << "=== Cross section per water molecule (cm^2)=" << sigma/cm/cm
229  << G4endl;
230  G4cout << "=== Cross section per water molecule (cm^-1)="
231  << sigma*waterDensity/(1./cm) << G4endl;
232  G4cout << "=== G4DNAUeharaScreenedRutherfordElasticModel - XS INFO END"
233  << G4endl;
234  }
235 #endif
236  }
237 
238  return sigma*waterDensity;
239 }
size_t GetIndex() const
Definition: G4Material.hh:262
G4GLOB_DLL std::ostream G4cout
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
#define G4endl
Definition: G4ios.hh:61
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold ( )
inline

Definition at line 147 of file G4DNAUeharaScreenedRutherfordElasticModel.hh.

148 {
149  G4ExceptionDescription errMsg;
150  errMsg << "*** WARNING : "
151  << "G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold"
152  << "is deprecated, the returned value is nonsense";
153 
154  G4Exception (
155  "G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold",
156  "DEPRECATED", JustWarning, errMsg);
157 
158  return -1;
159 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 89 of file G4DNAUeharaScreenedRutherfordElasticModel.cc.

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

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 299 of file G4DNAUeharaScreenedRutherfordElasticModel.cc.

304 {
305 #ifdef UEHARA_VERBOSE
306  if (verboseLevel > 3)
307  {
308  G4cout
309  << "Calling SampleSecondaries() of G4DNAUeharaScreenedRutherfordElasticModel"
310  << G4endl;
311  }
312 #endif
313 
314  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
315 
316  G4double cosTheta = 0.;
317 
318  if (electronEnergy0<intermediateEnergyLimit)
319  {
320 #ifdef UEHARA_VERBOSE
321  if (verboseLevel > 3)
322  G4cout << "---> Using Brenner & Zaider model" << G4endl;
323 #endif
324  cosTheta = BrennerZaiderRandomizeCosTheta(electronEnergy0);
325  }
326  else //if (electronEnergy0>=intermediateEnergyLimit)
327  {
328 #ifdef UEHARA_VERBOSE
329  if (verboseLevel > 3)
330  G4cout << "---> Using Screened Rutherford model" << G4endl;
331 #endif
332  G4double z = 7.42; // FROM PMB 37 (1992) 1841-1858 p1842
333  cosTheta = ScreenedRutherfordRandomizeCosTheta(electronEnergy0,z);
334  }
335 
336  G4double phi = 2. * pi * G4UniformRand();
337 
338  G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
339  G4ThreeVector xVers = zVers.orthogonal();
340  G4ThreeVector yVers = zVers.cross(xVers);
341 
342  G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
343  G4double yDir = xDir;
344  xDir *= std::cos(phi);
345  yDir *= std::sin(phi);
346 
347  G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
348 
350 
352 }
G4double GetKineticEnergy() const
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector orthogonal() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
static constexpr double pi
Definition: G4SIunits.hh:75
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4DNAUeharaScreenedRutherfordElasticModel::SelectFasterComputation ( G4bool  input)
inline

Definition at line 110 of file G4DNAUeharaScreenedRutherfordElasticModel.hh.

111 {
112  fasterCode = input;
113 }
void G4DNAUeharaScreenedRutherfordElasticModel::SelectHighEnergyLimit ( G4double  threshold)
inline

Definition at line 119 of file G4DNAUeharaScreenedRutherfordElasticModel.hh.

121 {
122  if(threshold > 10. * CLHEP::keV)
123  {
124  G4Exception (
125  "*** WARNING : the G4DNAUeharaScreenedRutherfordElasticModel class is "
126  "used above 10 keV !",
127  "", JustWarning, "");
128  }
129 
130  SetHighEnergyLimit(threshold);
131 }
static constexpr double keV
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

void G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold ( G4double  threshold)
inline

Definition at line 134 of file G4DNAUeharaScreenedRutherfordElasticModel.hh.

135 {
136  G4ExceptionDescription errMsg;
137  errMsg << "*** WARNING : "
138  << "G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold"
139  << "is deprecated, the kill threshold won't be taken into account";
140 
141  G4Exception (
142  "G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold",
143  "DEPRECATED", JustWarning, errMsg);
144 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

Member Data Documentation

G4ParticleChangeForGamma* G4DNAUeharaScreenedRutherfordElasticModel::fParticleChangeForGamma
protected

Definition at line 81 of file G4DNAUeharaScreenedRutherfordElasticModel.hh.


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