Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4DNAScreenedRutherfordElasticModel Class Reference

#include <G4DNAScreenedRutherfordElasticModel.hh>

Inheritance diagram for G4DNAScreenedRutherfordElasticModel:
Collaboration 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 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 SetFluctuationFlag (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
 
G4bool lossFlucFlag
 

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

Constructor & Destructor Documentation

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

Definition at line 44 of file G4DNAScreenedRutherfordElasticModel.cc.

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

Here is the call graph for this function:

G4DNAScreenedRutherfordElasticModel::~G4DNAScreenedRutherfordElasticModel ( )
virtual

Definition at line 84 of file G4DNAScreenedRutherfordElasticModel.cc.

85 {
86 }

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 192 of file G4DNAScreenedRutherfordElasticModel.cc.

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

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

Implements G4VEmModel.

Definition at line 91 of file G4DNAScreenedRutherfordElasticModel.cc.

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

Implements G4VEmModel.

Definition at line 306 of file G4DNAScreenedRutherfordElasticModel.cc.

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

Definition at line 120 of file G4DNAScreenedRutherfordElasticModel.hh.

121 {
122  fasterCode = input;
123 }
void G4DNAScreenedRutherfordElasticModel::SetKillBelowThreshold ( G4double  threshold)
inline

Definition at line 107 of file G4DNAScreenedRutherfordElasticModel.hh.

108 {
109  G4ExceptionDescription errMsg;
110  errMsg << "The method G4DNAScreenedRutherfordElasticModel::"
111  "SetKillBelowThreshold is deprecated";
112 
113  G4Exception("G4DNAScreenedRutherfordElasticModel::SetKillBelowThreshold",
114  "deprecated",
115  JustWarning,
116  errMsg);
117 }
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* G4DNAScreenedRutherfordElasticModel::fParticleChangeForGamma
protected

Definition at line 66 of file G4DNAScreenedRutherfordElasticModel.hh.


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