Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 &) override
 
virtual G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) final
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length) override
 
- 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 level, const G4ParticleDefinition *, G4double kineticEnergy)
 
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 * > *)
 
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 Member Functions

virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy) final
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
- 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::G4ICRU73QOModel ( const G4ParticleDefinition p = 0,
const G4String nam = "ICRU73QO" 
)
explicit

Definition at line 67 of file G4ICRU73QOModel.cc.

68  : G4VEmModel(nam),
69  particle(nullptr),
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 
82  theElectron = G4Electron::Electron();
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 = nullptr;
95  denEffData = nullptr;
96 }
int G4int
Definition: G4Types.hh:78
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
static G4Electron * Electron()
Definition: G4Electron.cc:94
static constexpr double MeV
Definition: G4SIunits.hh:214
static constexpr double keV
Definition: G4SIunits.hh:216

Here is the call graph for this function:

G4ICRU73QOModel::~G4ICRU73QOModel ( )
virtual

Definition at line 100 of file G4ICRU73QOModel.cc.

101 {}

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 154 of file G4ICRU73QOModel.cc.

160 {
162  (p,kineticEnergy,cutEnergy,maxEnergy);
163  return cross;
164 }
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) final
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

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 }
static constexpr double twopi_mc2_rcl2
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double energy(const ThreeVector &p, const G4double m)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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

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 }
tuple x
Definition: test.py:50
G4double GetElectronDensity() const
Definition: G4Material.hh:217
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
G4double G4Log(G4double x)
Definition: G4Log.hh:230
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 365 of file G4ICRU73QOModel.cc.

370 {}
G4double G4ICRU73QOModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
overridevirtual

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) final
G4double GetElectronDensity() const
Definition: G4Material.hh:217
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

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 
120  G4String pname = particle->GetParticleName();
121  fParticleChange = GetParticleChangeForLoss();
123  denEffData = (*mtab)[0]->GetIonisation()->GetDensityEffectData();
124  }
125 }
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
std::vector< G4Material * > G4MaterialTable
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:619
const G4String & GetParticleName() const
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:699
string pname
Definition: eplot.py:33
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:626
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:788

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 451 of file G4ICRU73QOModel.cc.

453 {
454  if(pd != particle) { SetParticle(pd); }
455  G4double tau = kinEnergy/mass;
456  G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
457  (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
458  return tmax;
459 }
float electron_mass_c2
Definition: hepunit.py:274
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

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

Implements G4VEmModel.

Definition at line 374 of file G4ICRU73QOModel.cc.

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

Here is the call graph for this function:


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