Geant4  10.02.p03
G4ICRU73QOModel Class Reference

#include <G4ICRU73QOModel.hh>

Inheritance diagram for G4ICRU73QOModel:
Collaboration diagram for G4ICRU73QOModel:

Public Member Functions

 G4ICRU73QOModel (const G4ParticleDefinition *p=0, const G4String &nam="ICRU73QO")
 
virtual ~G4ICRU73QOModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle *> *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
- 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 GetPartialCrossSection (const G4Material *, G4int, const G4ParticleDefinition *, G4double)
 
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 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 *> *)
 
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=0)
 
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 Member Functions

virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLoss * GetParticleChangeForLoss ()
 
G4ParticleChangeForGamma * GetParticleChangeForGamma ()
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Private Member Functions

void SetParticle (const G4ParticleDefinition *p)
 
void SetLowestKinEnergy (const G4double val)
 
G4double DEDX (const G4Material *material, G4double kineticEnergy)
 
G4double DEDXPerElement (G4int Z, G4double kineticEnergy)
 
G4ICRU73QOModeloperator= (const G4ICRU73QOModel &right)
 
 G4ICRU73QOModel (const G4ICRU73QOModel &)
 
G4int GetNumberOfShells (G4int Z) const
 
G4double GetShellEnergy (G4int Z, G4int nbOfTheShell) const
 
G4double GetOscillatorEnergy (G4int Z, G4int nbOfTheShell) const
 
G4double GetShellStrength (G4int Z, G4int nbOfTheShell) const
 
G4double GetL0 (G4double normEnergy) const
 
G4double GetL1 (G4double normEnergy) const
 
G4double GetL2 (G4double normEnergy) const
 

Private Attributes

const G4ParticleDefinitionparticle
 
G4ParticleDefinitiontheElectron
 
G4ParticleChangeForLoss * fParticleChange
 
G4DensityEffectDatadenEffData
 
G4double mass
 
G4double charge
 
G4double chargeSquare
 
G4double massRate
 
G4double ratio
 
G4double lowestKinEnergy
 
G4bool isInitialised
 
G4int indexZ [100]
 
G4int sizeL0
 
G4int sizeL1
 
G4int sizeL2
 

Static Private Attributes

static const G4int NQOELEM = 26
 
static const G4int NQODATA = 130
 
static const G4int ZElementAvailable [NQOELEM]
 
static const G4int startElemIndex [NQOELEM]
 
static const G4int nbofShellsForElement [NQOELEM]
 
static const G4double ShellEnergy [NQODATA]
 
static const G4double SubShellOccupation [NQODATA]
 
static const G4double L0 [67][2]
 
static const G4double L1 [22][2]
 
static const G4double L2 [14][2]
 
static const G4double factorBethe [99]
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChange * pParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 71 of file G4ICRU73QOModel.hh.

Constructor & Destructor Documentation

◆ G4ICRU73QOModel() [1/2]

G4ICRU73QOModel::G4ICRU73QOModel ( const G4ParticleDefinition p = 0,
const G4String nam = "ICRU73QO" 
)

Definition at line 67 of file G4ICRU73QOModel.cc.

68  : G4VEmModel(nam),
69  particle(0),
70  isInitialised(false)
71 {
72  mass = charge = chargeSquare = massRate = ratio = 0.0;
73  if(p) { SetParticle(p); }
74  SetHighEnergyLimit(10.0*MeV);
75 
76  lowestKinEnergy = 5.0*keV;
77 
78  sizeL0 = 67;
79  sizeL1 = 22;
80  sizeL2 = 14;
81 
83 
84  for (G4int i = 0; i < 100; ++i)
85  {
86  indexZ[i] = -1;
87  }
88  for(G4int i = 0; i < NQOELEM; ++i)
89  {
90  if(ZElementAvailable[i] > 0) {
91  indexZ[ZElementAvailable[i]] = i;
92  }
93  }
94  fParticleChange = 0;
95  denEffData = 0;
96 }
static const double MeV
Definition: G4SIunits.hh:211
G4ParticleChangeForLoss * fParticleChange
int G4int
Definition: G4Types.hh:78
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
G4ParticleDefinition * theElectron
G4DensityEffectData * denEffData
void SetParticle(const G4ParticleDefinition *p)
static const G4int NQOELEM
G4double lowestKinEnergy
static G4Electron * Electron()
Definition: G4Electron.cc:94
static const double keV
Definition: G4SIunits.hh:213
static const G4int ZElementAvailable[NQOELEM]
const G4ParticleDefinition * particle
Here is the call graph for this function:

◆ ~G4ICRU73QOModel()

G4ICRU73QOModel::~G4ICRU73QOModel ( )
virtual

Definition at line 100 of file G4ICRU73QOModel.cc.

101 {}

◆ G4ICRU73QOModel() [2/2]

G4ICRU73QOModel::G4ICRU73QOModel ( const G4ICRU73QOModel )
private

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4ICRU73QOModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  Z,
G4double  A,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 154 of file G4ICRU73QOModel.cc.

160 {
162  (p,kineticEnergy,cutEnergy,maxEnergy);
163  return cross;
164 }
Float_t Z
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ ComputeCrossSectionPerElectron()

G4double G4ICRU73QOModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Definition at line 129 of file G4ICRU73QOModel.cc.

134 {
135  G4double cross = 0.0;
136  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
137  G4double maxEnergy = std::min(tmax,maxKinEnergy);
138  if(cutEnergy < maxEnergy) {
139 
140  G4double energy = kineticEnergy + mass;
141  G4double energy2 = energy*energy;
142  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
143  cross = (maxEnergy - cutEnergy)/(cutEnergy*maxEnergy)
144  - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
145 
146  cross *= CLHEP::twopi_mc2_rcl2*chargeSquare/beta2;
147  }
148 
149  return cross;
150 }
double energy
Definition: plottest35.C:25
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static const double twopi_mc2_rcl2
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeDEDXPerVolume()

G4double G4ICRU73QOModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy 
)
virtual

Reimplemented from G4VEmModel.

Reimplemented in G4ICRU73NoDeltaModel.

Definition at line 183 of file G4ICRU73QOModel.cc.

187 {
188  SetParticle(p);
189  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
190  G4double tkin = kineticEnergy/massRate;
191  G4double dedx = 0.0;
192  if(tkin > lowestKinEnergy) { dedx = DEDX(material, tkin); }
193  else { dedx = DEDX(material, lowestKinEnergy)*sqrt(tkin/lowestKinEnergy); }
194 
195  if (cutEnergy < tmax) {
196 
197  G4double tau = kineticEnergy/mass;
198  G4double gam = tau + 1.0;
199  G4double bg2 = tau * (tau+2.0);
200  G4double beta2 = bg2/(gam*gam);
201  G4double x = cutEnergy/tmax;
202 
203  dedx += chargeSquare*( G4Log(x) + (1.0 - x)*beta2 ) * twopi_mc2_rcl2
204  * material->GetElectronDensity()/beta2;
205  }
206  if(dedx < 0.0) { dedx = 0.0; }
207  return dedx;
208 }
int twopi_mc2_rcl2
Definition: hepunit.py:294
G4double DEDX(const G4Material *material, G4double kineticEnergy)
void SetParticle(const G4ParticleDefinition *p)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double lowestKinEnergy
G4double GetElectronDensity() const
Definition: G4Material.hh:217
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CorrectionsAlongStep()

void G4ICRU73QOModel::CorrectionsAlongStep ( const G4MaterialCutsCouple ,
const G4DynamicParticle ,
G4double eloss,
G4double niel,
G4double  length 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 364 of file G4ICRU73QOModel.cc.

369 {}

◆ CrossSectionPerVolume()

G4double G4ICRU73QOModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Reimplemented from G4VEmModel.

Reimplemented in G4ICRU73NoDeltaModel.

Definition at line 168 of file G4ICRU73QOModel.cc.

174 {
175  G4double eDensity = material->GetElectronDensity();
177  (p,kineticEnergy,cutEnergy,maxEnergy);
178  return cross;
179 }
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double GetElectronDensity() const
Definition: G4Material.hh:217
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ DEDX()

G4double G4ICRU73QOModel::DEDX ( const G4Material material,
G4double  kineticEnergy 
)
private

Definition at line 212 of file G4ICRU73QOModel.cc.

214 {
215  G4double eloss = 0.0;
216  const G4int numberOfElements = material->GetNumberOfElements();
217  const G4double* theAtomicNumDensityVector =
218  material->GetAtomicNumDensityVector();
219 
220  // Bragg's rule calculation
221  const G4ElementVector* theElementVector =
222  material->GetElementVector() ;
223 
224  // loop for the elements in the material
225  for (G4int i=0; i<numberOfElements; ++i)
226  {
227  const G4Element* element = (*theElementVector)[i] ;
228  eloss += DEDXPerElement(G4lrint(element->GetZ()), kineticEnergy)
229  * theAtomicNumDensityVector[i] * G4int(element->GetZ());
230  }
231  return eloss;
232 }
G4double DEDXPerElement(G4int Z, G4double kineticEnergy)
std::vector< G4Element * > G4ElementVector
int G4int
Definition: G4Types.hh:78
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:216
int G4lrint(double ad)
Definition: templates.hh:163
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
double G4double
Definition: G4Types.hh:76
G4double GetZ() const
Definition: G4Element.hh:131
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DEDXPerElement()

G4double G4ICRU73QOModel::DEDXPerElement ( G4int  Z,
G4double  kineticEnergy 
)
private

Definition at line 236 of file G4ICRU73QOModel.cc.

238 {
239  G4int Z = AtomicNumber;
240  if(Z > 97) { Z = 97; }
241  G4int nbOfShells = GetNumberOfShells(Z);
242  if(nbOfShells < 1) { nbOfShells = 1; }
243 
244  G4double v = CLHEP::c_light * std::sqrt( 2.0*kineticEnergy/proton_mass_c2 );
245 
247 
248  G4double tau = kineticEnergy/proton_mass_c2;
249  G4double gam = tau + 1.0;
250  G4double bg2 = tau * (tau+2.0);
251  G4double beta2 = bg2/(gam*gam);
252 
253  G4double l0Term = 0, l1Term = 0, l2Term = 0;
254 
255  for (G4int nos = 0; nos < nbOfShells; ++nos){
256 
257  G4double NormalizedEnergy = (2.0*CLHEP::electron_mass_c2*beta2) /
258  GetShellEnergy(Z,nos);
259 
260  G4double shStrength = GetShellStrength(Z,nos);
261 
262  G4double l0 = GetL0(NormalizedEnergy);
263  l0Term += shStrength * l0;
264 
265  G4double l1 = GetL1(NormalizedEnergy);
266  l1Term += shStrength * l1;
267 
268  G4double l2 = GetL2(NormalizedEnergy);
269  l2Term += shStrength * l2;
270 
271  }
273  (l0Term + charge*fBetheVelocity*l1Term
274  + chargeSquare*fBetheVelocity*fBetheVelocity*l2Term)/beta2;
275  return dedx;
276 }
G4double GetL0(G4double normEnergy) const
static const G4double factorBethe[99]
G4double GetL2(G4double normEnergy) const
int G4int
Definition: G4Types.hh:78
static const double fine_structure_const
Float_t Z
float proton_mass_c2
Definition: hepunit.py:275
static const double c_light
G4double GetL1(G4double normEnergy) const
static const double twopi_mc2_rcl2
G4double GetShellStrength(G4int Z, G4int nbOfTheShell) const
G4int GetNumberOfShells(G4int Z) const
G4double GetShellEnergy(G4int Z, G4int nbOfTheShell) const
double G4double
Definition: G4Types.hh:76
static const double electron_mass_c2
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetL0()

G4double G4ICRU73QOModel::GetL0 ( G4double  normEnergy) const
private

Definition at line 305 of file G4ICRU73QOModel.cc.

306 {
307  G4int n;
308 
309  for(n = 0; n < sizeL0; n++) {
310  if( normEnergy < L0[n][0] ) break;
311  }
312  if(0 == n) { n = 1; }
313  if(n >= sizeL0) { n = sizeL0 - 1; }
314 
315  G4double l0 = L0[n][1];
316  G4double l0p = L0[n-1][1];
317  G4double bethe = l0p + (l0 - l0p) * ( normEnergy - L0[n-1][0]) /
318  (L0[n][0] - L0[n-1][0]);
319 
320  return bethe ;
321 }
static const G4double L0[67][2]
int G4int
Definition: G4Types.hh:78
Char_t n[5]
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ GetL1()

G4double G4ICRU73QOModel::GetL1 ( G4double  normEnergy) const
private

Definition at line 325 of file G4ICRU73QOModel.cc.

326 {
327  G4int n;
328 
329  for(n = 0; n < sizeL1; n++) {
330  if( normEnergy < L1[n][0] ) break;
331  }
332  if(0 == n) n = 1 ;
333  if(n >= sizeL1) n = sizeL1 - 1 ;
334 
335  G4double l1 = L1[n][1];
336  G4double l1p = L1[n-1][1];
337  G4double barkas= l1p + (l1 - l1p) * ( normEnergy - L1[n-1][0]) /
338  (L1[n][0] - L1[n-1][0]);
339 
340  return barkas;
341 }
static const G4double L1[22][2]
int G4int
Definition: G4Types.hh:78
Char_t n[5]
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ GetL2()

G4double G4ICRU73QOModel::GetL2 ( G4double  normEnergy) const
private

Definition at line 345 of file G4ICRU73QOModel.cc.

346 {
347  G4int n;
348  for(n = 0; n < sizeL2; n++) {
349  if( normEnergy < L2[n][0] ) break;
350  }
351  if(0 == n) n = 1 ;
352  if(n >= sizeL2) n = sizeL2 - 1 ;
353 
354  G4double l2 = L2[n][1];
355  G4double l2p = L2[n-1][1];
356  G4double bloch = l2p + (l2 - l2p) * ( normEnergy - L2[n-1][0]) /
357  (L2[n][0] - L2[n-1][0]);
358 
359  return bloch;
360 }
static const G4double L2[14][2]
int G4int
Definition: G4Types.hh:78
Char_t n[5]
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ GetNumberOfShells()

G4int G4ICRU73QOModel::GetNumberOfShells ( G4int  Z) const
inlineprivate

Definition at line 209 of file G4ICRU73QOModel.hh.

210 {
211  G4int nShell = 0;
212 
213  if(indexZ[Z] >= 0) { nShell = nbofShellsForElement[indexZ[Z]];
214  } else { nShell = G4AtomicShells::GetNumberOfShells(Z); }
215 
216  return nShell;
217 }
static const G4int nbofShellsForElement[NQOELEM]
int G4int
Definition: G4Types.hh:78
Float_t Z
static G4int GetNumberOfShells(G4int Z)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetOscillatorEnergy()

G4double G4ICRU73QOModel::GetOscillatorEnergy ( G4int  Z,
G4int  nbOfTheShell 
) const
private

Definition at line 281 of file G4ICRU73QOModel.cc.

283 {
285  if(idx == -1) { idx = denEffData->GetElementIndex(Z-1, kStateUndefined); }
286  G4double PlasmaEnergy = denEffData->GetPlasmaEnergy(idx);
287 
288  G4double PlasmaEnergy2 = PlasmaEnergy * PlasmaEnergy;
289 
290  G4double plasmonTerm = 0.66667
291  * G4AtomicShells::GetNumberOfElectrons(Z,nbOfTheShell)
292  * PlasmaEnergy2 / (Z*Z) ;
293 
294  G4double ionTerm = G4Exp(0.5) *
295  (G4AtomicShells::GetBindingEnergy(Z,nbOfTheShell)) ;
296  G4double ionTerm2 = ionTerm*ionTerm ;
297 
298  G4double oscShellEnergy = std::sqrt( ionTerm2 + plasmonTerm );
299 
300  return oscShellEnergy;
301 }
G4int GetElementIndex(G4int Z, G4State mState)
int G4int
Definition: G4Types.hh:78
G4DensityEffectData * denEffData
Float_t Z
static G4double GetBindingEnergy(G4int Z, G4int SubshellNb)
G4double GetPlasmaEnergy(G4int idx)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static G4int GetNumberOfElectrons(G4int Z, G4int SubshellNb)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetShellEnergy()

G4double G4ICRU73QOModel::GetShellEnergy ( G4int  Z,
G4int  nbOfTheShell 
) const
inlineprivate

Definition at line 222 of file G4ICRU73QOModel.hh.

223 {
224  G4double shellEnergy = 0.;
225 
226  G4int idx = indexZ[Z];
227 
228  if(idx >= 0) { shellEnergy = ShellEnergy[startElemIndex[idx] + nbOfTheShell]*CLHEP::eV;
229  } else { shellEnergy = GetOscillatorEnergy(Z, nbOfTheShell); }
230 
231  return shellEnergy;
232 }
int G4int
Definition: G4Types.hh:78
Float_t Z
G4double GetOscillatorEnergy(G4int Z, G4int nbOfTheShell) const
double G4double
Definition: G4Types.hh:76
static const G4double ShellEnergy[NQODATA]
static const double eV
static const G4int startElemIndex[NQOELEM]
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetShellStrength()

G4double G4ICRU73QOModel::GetShellStrength ( G4int  Z,
G4int  nbOfTheShell 
) const
inlineprivate

Definition at line 237 of file G4ICRU73QOModel.hh.

238 {
239  G4double shellStrength = 0.;
240 
241  G4int idx = indexZ[Z];
242 
243  if(idx >= 0) { shellStrength = SubShellOccupation[startElemIndex[idx] + nbOfTheShell] / Z;
244  } else { shellStrength = G4double(G4AtomicShells::GetNumberOfElectrons(Z,nbOfTheShell))/Z; }
245 
246  return shellStrength;
247 }
static const G4double SubShellOccupation[NQODATA]
int G4int
Definition: G4Types.hh:78
Float_t Z
static G4int GetNumberOfElectrons(G4int Z, G4int SubshellNb)
double G4double
Definition: G4Types.hh:76
static const G4int startElemIndex[NQOELEM]
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Initialise()

void G4ICRU73QOModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 105 of file G4ICRU73QOModel.cc.

107 {
108  if(p != particle) SetParticle(p);
109 
110  // always false before the run
111  SetDeexcitationFlag(false);
112 
113  if(!isInitialised) {
114  isInitialised = true;
115 
118  }
119 
123  denEffData = (*mtab)[0]->GetIonisation()->GetDensityEffectData();
124  }
125 }
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:120
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:697
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:589
std::vector< G4Material * > G4MaterialTable
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:617
G4ParticleChangeForLoss * fParticleChange
const G4String & GetParticleName() const
G4DensityEffectData * denEffData
void SetParticle(const G4ParticleDefinition *p)
string pname
Definition: eplot.py:33
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:624
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:781
const G4ParticleDefinition * particle
Here is the call graph for this function:

◆ MaxSecondaryEnergy()

G4double G4ICRU73QOModel::MaxSecondaryEnergy ( const G4ParticleDefinition pd,
G4double  kinEnergy 
)
protectedvirtual

Reimplemented from G4VEmModel.

Definition at line 450 of file G4ICRU73QOModel.cc.

452 {
453  if(pd != particle) { SetParticle(pd); }
454  G4double tau = kinEnergy/mass;
455  G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
456  (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
457  return tmax;
458 }
float electron_mass_c2
Definition: hepunit.py:274
void SetParticle(const G4ParticleDefinition *p)
double G4double
Definition: G4Types.hh:76
const G4ParticleDefinition * particle
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

G4ICRU73QOModel& G4ICRU73QOModel::operator= ( const G4ICRU73QOModel right)
private

◆ SampleSecondaries()

void G4ICRU73QOModel::SampleSecondaries ( std::vector< G4DynamicParticle *> *  vdp,
const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 373 of file G4ICRU73QOModel.cc.

378 {
379  G4double tmax = MaxSecondaryKinEnergy(dp);
380  G4double xmax = std::min(tmax, maxEnergy);
381  if(xmin >= xmax) { return; }
382 
383  G4double kineticEnergy = dp->GetKineticEnergy();
384  G4double energy = kineticEnergy + mass;
385  G4double energy2 = energy*energy;
386  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
387  G4double grej = 1.0;
388  G4double deltaKinEnergy, f;
389 
390  G4ThreeVector direction = dp->GetMomentumDirection();
391 
392  // sampling follows ...
393  do {
395  deltaKinEnergy = xmin*xmax/(xmin*(1.0 - x) + xmax*x);
396 
397  f = 1.0 - beta2*deltaKinEnergy/tmax;
398 
399  if(f > grej) {
400  G4cout << "G4ICRU73QOModel::SampleSecondary Warning! "
401  << "Majorant " << grej << " < "
402  << f << " for e= " << deltaKinEnergy
403  << G4endl;
404  }
405 
406  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
407  } while( grej*G4UniformRand() >= f );
408 
409  G4ThreeVector deltaDirection;
410 
412  const G4Material* mat = couple->GetMaterial();
414 
415  deltaDirection =
416  GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
417 
418  } else {
419 
420  G4double deltaMomentum =
421  sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
422  G4double totMomentum = energy*sqrt(beta2);
423  G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
424  (deltaMomentum * totMomentum);
425  if(cost > 1.0) { cost = 1.0; }
426  G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
427 
428  G4double phi = twopi * G4UniformRand() ;
429 
430  deltaDirection.set(sint*cos(phi),sint*sin(phi), cost) ;
431  deltaDirection.rotateUz(direction);
432  }
433  // create G4DynamicParticle object for delta ray
434  G4DynamicParticle* delta =
435  new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
436 
437  // Change kinematics of primary particle
438  kineticEnergy -= deltaKinEnergy;
439  G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
440  finalP = finalP.unit();
441 
442  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
443  fParticleChange->SetProposedMomentumDirection(finalP);
444 
445  vdp->push_back(delta);
446 }
void set(double x, double y, double z)
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:482
const G4Material * GetMaterial() const
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:697
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:617
G4ParticleChangeForLoss * fParticleChange
int G4int
Definition: G4Types.hh:78
Float_t mat
G4ParticleDefinition * theElectron
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
double energy
Definition: plottest35.C:25
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
Float_t Z
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
static const double twopi
Definition: G4SIunits.hh:75
float electron_mass_c2
Definition: hepunit.py:274
const G4ThreeVector & GetMomentumDirection() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.hh:564
Here is the call graph for this function:

◆ SetLowestKinEnergy()

void G4ICRU73QOModel::SetLowestKinEnergy ( const G4double  val)
inlineprivate

Definition at line 251 of file G4ICRU73QOModel.hh.

252 {
253  lowestKinEnergy = val;
254 }
G4double lowestKinEnergy

◆ SetParticle()

void G4ICRU73QOModel::SetParticle ( const G4ParticleDefinition p)
inlineprivate

Definition at line 197 of file G4ICRU73QOModel.hh.

198 {
199  particle = p;
200  mass = particle->GetPDGMass();
205 }
static const double proton_mass_c2
static const double eplus
static const double electron_mass_c2
const G4ParticleDefinition * particle
G4double GetPDGCharge() const
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ charge

G4double G4ICRU73QOModel::charge
private

Definition at line 144 of file G4ICRU73QOModel.hh.

◆ chargeSquare

G4double G4ICRU73QOModel::chargeSquare
private

Definition at line 145 of file G4ICRU73QOModel.hh.

◆ denEffData

G4DensityEffectData* G4ICRU73QOModel::denEffData
private

Definition at line 141 of file G4ICRU73QOModel.hh.

◆ factorBethe

const G4double G4ICRU73QOModel::factorBethe
staticprivate
Initial value:
= { 1.0,
0.9637, 0.9872, 0.9469, 0.9875, 0.91, 0.989, 0.9507, 0.9773, 0.8621, 0.979,
0.8357, 0.868, 0.9417, 0.9466, 0.8911, 0.905, 0.944, 0.9607, 0.928, 0.96,
0.9098, 0.976, 0.8425, 0.8099, 0.7858, 0.947, 0.7248, 0.9106, 0.9246, 0.6821,
0.7223, 0.9784, 0.774, 0.7953, 0.829, 0.9405, 0.8318, 0.8583, 0.8563, 0.8481,
0.8207, 0.9033, 0.8063, 0.7837, 0.7818, 0.744, 0.875, 0.7693, 0.7871, 0.8459,
0.8231, 0.8462, 0.853, 0.8736, 0.856, 0.8762, 0.8629, 0.8323, 0.8064, 0.7828,
0.7533, 0.7273, 0.7093, 0.7157, 0.6823, 0.6612, 0.6418, 0.6395, 0.6323, 0.6221,
0.6497, 0.6746, 0.8568, 0.8541, 0.6958, 0.6962, 0.7051, 0.863, 0.8588, 0.7226,
0.7454, 0.78, 0.7783, 0.7996, 0.8216, 0.8632, 0.8558, 0.8792, 0.8745, 0.8676,
0.8321, 0.8272, 0.7999, 0.7934, 0.7787, 0.7851, 0.7692, 0.7598}

Definition at line 191 of file G4ICRU73QOModel.hh.

◆ fParticleChange

G4ParticleChangeForLoss* G4ICRU73QOModel::fParticleChange
private

Definition at line 140 of file G4ICRU73QOModel.hh.

◆ indexZ

G4int G4ICRU73QOModel::indexZ[100]
private

Definition at line 180 of file G4ICRU73QOModel.hh.

◆ isInitialised

G4bool G4ICRU73QOModel::isInitialised
private

Definition at line 150 of file G4ICRU73QOModel.hh.

◆ L0

const G4double G4ICRU73QOModel::L0
staticprivate

Definition at line 183 of file G4ICRU73QOModel.hh.

◆ L1

const G4double G4ICRU73QOModel::L1
staticprivate
Initial value:
=
{
{0.00, -0.000001},
{0.10, -0.00001},
{0.20, -0.00049},
{0.30, -0.00084},
{0.40, 0.00085},
{0.50, 0.00519},
{0.60, 0.01198},
{0.70, 0.02074},
{0.80, 0.03133},
{0.90, 0.04369},
{1.00, 0.06035},
{2.00, 0.24023},
{3.00, 0.44284},
{4.00, 0.62012},
{5.00, 0.77031},
{6.00, 0.90390},
{7.00, 1.02705},
{8.00, 1.10867},
{9.00, 1.17546},
{10.00, 1.21599},
{15.00, 1.24349},
{20.00, 1.16752}
}

Definition at line 184 of file G4ICRU73QOModel.hh.

◆ L2

const G4double G4ICRU73QOModel::L2
staticprivate
Initial value:
=
{
{0.00, 0.000001},
{0.10, 0.00001},
{0.20, 0.00000},
{0.40, -0.00120},
{0.60, -0.00036},
{0.80, 0.00372},
{1.00, 0.01298},
{2.00, 0.08296},
{4.00, 0.21953},
{6.00, 0.23903},
{8.00, 0.20893},
{10.00, 0.10879},
{20.00, -0.88409},
{40.00, -1.13902}
}

Definition at line 185 of file G4ICRU73QOModel.hh.

◆ lowestKinEnergy

G4double G4ICRU73QOModel::lowestKinEnergy
private

Definition at line 148 of file G4ICRU73QOModel.hh.

◆ mass

G4double G4ICRU73QOModel::mass
private

Definition at line 143 of file G4ICRU73QOModel.hh.

◆ massRate

G4double G4ICRU73QOModel::massRate
private

Definition at line 146 of file G4ICRU73QOModel.hh.

◆ nbofShellsForElement

const G4int G4ICRU73QOModel::nbofShellsForElement
staticprivate
Initial value:
= {1,1,2,3,3,3,3,4,5,4,
5,5,5,5,6,4,6,6,
7,6,6,8,7,7,9,9}

Definition at line 176 of file G4ICRU73QOModel.hh.

◆ NQODATA

const G4int G4ICRU73QOModel::NQODATA = 130
staticprivate

Definition at line 170 of file G4ICRU73QOModel.hh.

◆ NQOELEM

const G4int G4ICRU73QOModel::NQOELEM = 26
staticprivate

Definition at line 169 of file G4ICRU73QOModel.hh.

◆ particle

const G4ParticleDefinition* G4ICRU73QOModel::particle
private

Definition at line 138 of file G4ICRU73QOModel.hh.

◆ ratio

G4double G4ICRU73QOModel::ratio
private

Definition at line 147 of file G4ICRU73QOModel.hh.

◆ ShellEnergy

const G4double G4ICRU73QOModel::ShellEnergy
staticprivate
Initial value:
=
{
19.2,
41.8,
209.11, 21.68,
486.2, 60.95, 23.43,
732.61, 100.646, 23.550,
965.1, 129.85, 31.60,
1525.9, 234.9, 56.18,
2701, 476.5, 150.42, 16.89,
3206.1, 586.4, 186.8, 23.52, 14.91,
5551.6, 472.43, 124.85, 22.332,
8554.6, 850.58, 93.47, 39.19, 19.46,
12254.7, 1279.29, 200.35, 49.19, 17.66,
14346.9, 1532.28, 262.71, 74.37, 23.03,
15438.5, 1667.96, 294.1, 70.69, 16.447,
19022.1, 2150.79, 455.79, 179.87, 57.89, 20.95,
24643, 2906.4, 366.85, 22.24,
34394, 4365.3, 589.36, 129.42, 35.59, 18.42,
43664.3, 5824.91, 909.79, 175.47, 54.89, 19.63,
49948, 6818.2, 1036.1, 172.65, 70.89, 33.87, 14.54,
58987, 8159, 1296.6, 356.75, 101.03, 16.52,
88926, 18012, 3210, 575, 108.7, 30.8,
115025.9, 17827.44, 3214.36, 750.41, 305.21, 105.50, 38.09, 21.25,
128342, 20254, 3601.8, 608.1, 115.0, 42.75, 17.04,
131872, 20903, 3757.4, 682.1, 105.2, 44.89, 17.575,
154449, 25067, 5105.0, 987.44, 247.59, 188.1, 40.61, 19.2, 15.17,
167282, 27868, 6022.7, 1020.4, 244.81, 51.33, 13, 11.06, 14.43
}

Definition at line 177 of file G4ICRU73QOModel.hh.

◆ sizeL0

G4int G4ICRU73QOModel::sizeL0
private

Definition at line 187 of file G4ICRU73QOModel.hh.

◆ sizeL1

G4int G4ICRU73QOModel::sizeL1
private

Definition at line 188 of file G4ICRU73QOModel.hh.

◆ sizeL2

G4int G4ICRU73QOModel::sizeL2
private

Definition at line 189 of file G4ICRU73QOModel.hh.

◆ startElemIndex

const G4int G4ICRU73QOModel::startElemIndex
staticprivate
Initial value:
= {0,1,2,4,7,10,13,16,20,25,
29,34,39,44,49,55,59,65,
71,78,84,90,98,105,112,121}

Definition at line 175 of file G4ICRU73QOModel.hh.

◆ SubShellOccupation

const G4double G4ICRU73QOModel::SubShellOccupation
staticprivate
Initial value:
=
{
1.000,
2.000,
1.930, 2.070,
1.992, 1.841, 2.167,
1.741, 1.680, 3.579,
1.802, 1.849, 4.349,
1.788, 2.028, 6.184,
1.623, 2.147, 6.259, 2.971,
1.631, 2.094, 6.588, 2.041, 1.646,
1.535, 8.655, 1.706, 6.104,
1.581, 8.358, 8.183, 2.000, 1.878,
1.516, 8.325, 8.461, 6.579, 1.119,
1.422, 7.81, 8.385, 8.216, 2.167,
1.458, 8.049, 8.79, 9.695, 1.008,
1.442, 7.791, 7.837, 10.122, 2.463, 2.345,
1.645, 7.765, 19.192, 7.398,
1.313, 6.409, 19.229, 8.633, 5.036, 1.380,
1.295, 6.219, 18.751, 8.748, 10.184, 1.803,
1.277, 6.099, 20.386, 8.011, 10.007, 2.272, 1.948,
1.563, 6.312, 21.868, 5.762, 11.245, 7.250,
0.9198, 6.5408, 18.9727, 24.9149, 15.0161, 6.6284,
1.202, 5.582, 19.527, 18.741, 8.411, 14.387, 4.042, 2.108,
1.159, 5.467, 18.802, 33.905, 8.300, 9.342, 1.025,
1.124, 5.331, 18.078, 34.604, 8.127, 10.414, 1.322,
2.000, 8.000, 18.000, 18.000, 14.000, 8.000, 10.000, 2.000, 2.000,
2.000, 8.000, 18.000, 32.000, 18.000, 8.000, 2.000, 1.000, 3.000
}

Definition at line 178 of file G4ICRU73QOModel.hh.

◆ theElectron

G4ParticleDefinition* G4ICRU73QOModel::theElectron
private

Definition at line 139 of file G4ICRU73QOModel.hh.

◆ ZElementAvailable

const G4int G4ICRU73QOModel::ZElementAvailable
staticprivate
Initial value:
= {1,2,4,6,7,8,10,13,14,-18,
22,26,28,29,32,36,42,47,
50,54,73,74,78,79,82,92}

Definition at line 171 of file G4ICRU73QOModel.hh.


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