Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions
G4ICRU73QOModel Class Reference

#include <G4ICRU73QOModel.hh>

Inheritance diagram for G4ICRU73QOModel:
G4VEmModel G4ICRU73NoDeltaModel

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 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 *)
 
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
 

Protected Member Functions

virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy)
 
- 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
 

Detailed Description

Definition at line 71 of file G4ICRU73QOModel.hh.

Constructor & Destructor Documentation

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

Definition at line 66 of file G4ICRU73QOModel.cc.

References G4Electron::Electron(), python.hepunit::keV, python.hepunit::MeV, and G4VEmModel::SetHighEnergyLimit().

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

Definition at line 99 of file G4ICRU73QOModel.cc.

100 {}

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 149 of file G4ICRU73QOModel.cc.

References ComputeCrossSectionPerElectron().

155 {
157  (p,kineticEnergy,cutEnergy,maxEnergy);
158  return cross;
159 }
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
double G4double
Definition: G4Types.hh:76
G4double G4ICRU73QOModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Definition at line 124 of file G4ICRU73QOModel.cc.

References energy(), G4Log(), MaxSecondaryEnergy(), and G4INCL::Math::min().

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

129 {
130  G4double cross = 0.0;
131  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
132  G4double maxEnergy = std::min(tmax,maxKinEnergy);
133  if(cutEnergy < maxEnergy) {
134 
135  G4double energy = kineticEnergy + mass;
136  G4double energy2 = energy*energy;
137  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
138  cross = 1.0/cutEnergy - 1.0/maxEnergy
139  - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
140 
141  cross *= CLHEP::twopi_mc2_rcl2*chargeSquare/beta2;
142  }
143 
144  return cross;
145 }
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4double G4Log(G4double x)
Definition: G4Log.hh:227
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
double G4double
Definition: G4Types.hh:76
G4double G4ICRU73QOModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy 
)
virtual

Reimplemented from G4VEmModel.

Reimplemented in G4ICRU73NoDeltaModel.

Definition at line 178 of file G4ICRU73QOModel.cc.

References G4Log(), G4InuclParticleNames::gam, G4Material::GetElectronDensity(), MaxSecondaryEnergy(), python.hepunit::twopi_mc2_rcl2, and test::x.

Referenced by G4ICRU73NoDeltaModel::ComputeDEDXPerVolume().

182 {
183  SetParticle(p);
184  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
185  G4double tkin = kineticEnergy/massRate;
186  G4double dedx = 0.0;
187  if(tkin > lowestKinEnergy) { dedx = DEDX(material, tkin); }
188  else { dedx = DEDX(material, lowestKinEnergy)*sqrt(tkin/lowestKinEnergy); }
189 
190  if (cutEnergy < tmax) {
191 
192  G4double tau = kineticEnergy/mass;
193  G4double gam = tau + 1.0;
194  G4double bg2 = tau * (tau+2.0);
195  G4double beta2 = bg2/(gam*gam);
196  G4double x = cutEnergy/tmax;
197 
198  dedx += chargeSquare*( G4Log(x) + (1.0 - x)*beta2 ) * twopi_mc2_rcl2
199  * material->GetElectronDensity()/beta2;
200  }
201  if(dedx < 0.0) { dedx = 0.0; }
202  return dedx;
203 }
G4double GetElectronDensity() const
Definition: G4Material.hh:215
G4double G4Log(G4double x)
Definition: G4Log.hh:227
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
double G4double
Definition: G4Types.hh:76
void G4ICRU73QOModel::CorrectionsAlongStep ( const G4MaterialCutsCouple ,
const G4DynamicParticle ,
G4double eloss,
G4double niel,
G4double  length 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 359 of file G4ICRU73QOModel.cc.

364 {}
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 163 of file G4ICRU73QOModel.cc.

References ComputeCrossSectionPerElectron(), and G4Material::GetElectronDensity().

169 {
170  G4double eDensity = material->GetElectronDensity();
172  (p,kineticEnergy,cutEnergy,maxEnergy);
173  return cross;
174 }
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double GetElectronDensity() const
Definition: G4Material.hh:215
double G4double
Definition: G4Types.hh:76
void G4ICRU73QOModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 104 of file G4ICRU73QOModel.cc.

References G4Material::GetMaterialTable(), G4VEmModel::GetParticleChangeForLoss(), G4ParticleDefinition::GetParticleName(), eplot::pname, and G4VEmModel::SetDeexcitationFlag().

106 {
107  if(p != particle) SetParticle(p);
108 
109  // always false before the run
110  SetDeexcitationFlag(false);
111 
112  if(!isInitialised) {
113  isInitialised = true;
114 
115  G4String pname = particle->GetParticleName();
116  fParticleChange = GetParticleChangeForLoss();
118  denEffData = (*mtab)[0]->GetIonisation()->GetDensityEffectData();
119  }
120 }
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:564
std::vector< G4Material * > G4MaterialTable
const G4String & GetParticleName() const
string pname
Definition: eplot.py:33
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:739
G4double G4ICRU73QOModel::MaxSecondaryEnergy ( const G4ParticleDefinition pd,
G4double  kinEnergy 
)
protectedvirtual

Reimplemented from G4VEmModel.

Definition at line 433 of file G4ICRU73QOModel.cc.

References python.hepunit::electron_mass_c2.

Referenced by ComputeCrossSectionPerElectron(), and ComputeDEDXPerVolume().

435 {
436  if(pd != particle) { SetParticle(pd); }
437  G4double tau = kinEnergy/mass;
438  G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
439  (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
440  return tmax;
441 }
float electron_mass_c2
Definition: hepunit.py:274
double G4double
Definition: G4Types.hh:76
void G4ICRU73QOModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  vdp,
const G4MaterialCutsCouple ,
const G4DynamicParticle dp,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 368 of file G4ICRU73QOModel.cc.

References python.hepunit::electron_mass_c2, energy(), G4cout, G4endl, G4UniformRand, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4VEmModel::MaxSecondaryKinEnergy(), G4INCL::Math::min(), CLHEP::Hep3Vector::rotateUz(), G4ParticleChangeForLoss::SetProposedKineticEnergy(), G4ParticleChangeForLoss::SetProposedMomentumDirection(), python.hepunit::twopi, CLHEP::Hep3Vector::unit(), and test::x.

373 {
374  G4double tmax = MaxSecondaryKinEnergy(dp);
375  G4double xmax = std::min(tmax, maxEnergy);
376  if(xmin >= xmax) { return; }
377 
378  G4double kineticEnergy = dp->GetKineticEnergy();
379  G4double energy = kineticEnergy + mass;
380  G4double energy2 = energy*energy;
381  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
382  G4double grej = 1.0;
383  G4double deltaKinEnergy, f;
384 
385  G4ThreeVector direction = dp->GetMomentumDirection();
386 
387  // sampling follows ...
388  do {
390  deltaKinEnergy = xmin*xmax/(xmin*(1.0 - x) + xmax*x);
391 
392  f = 1.0 - beta2*deltaKinEnergy/tmax;
393 
394  if(f > grej) {
395  G4cout << "G4ICRU73QOModel::SampleSecondary Warning! "
396  << "Majorant " << grej << " < "
397  << f << " for e= " << deltaKinEnergy
398  << G4endl;
399  }
400 
401  } while( grej*G4UniformRand() >= f );
402 
403  G4double deltaMomentum =
404  sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
405  G4double totMomentum = energy*sqrt(beta2);
406  G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
407  (deltaMomentum * totMomentum);
408  if(cost > 1.0) { cost = 1.0; }
409  G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
410 
411  G4double phi = twopi * G4UniformRand() ;
412 
413  G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost) ;
414  deltaDirection.rotateUz(direction);
415 
416  // Change kinematics of primary particle
417  kineticEnergy -= deltaKinEnergy;
418  G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
419  finalP = finalP.unit();
420 
421  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
422  fParticleChange->SetProposedMomentumDirection(finalP);
423 
424  // create G4DynamicParticle object for delta ray
425  G4DynamicParticle* delta = new G4DynamicParticle(theElectron,deltaDirection,
426  deltaKinEnergy);
427 
428  vdp->push_back(delta);
429 }
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:448
G4double GetKineticEnergy() const
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
float electron_mass_c2
Definition: hepunit.py:274
void SetProposedMomentumDirection(const G4ThreeVector &dir)
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

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