Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ICRU73QOModel Class Reference

#include <G4ICRU73QOModel.hh>

+ Inheritance 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 Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
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
 
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 ForceBuildTable (G4bool val)
 
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 ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Additional Inherited Members

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

Detailed Description

Definition at line 71 of file G4ICRU73QOModel.hh.

Constructor & Destructor Documentation

◆ G4ICRU73QOModel()

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

Definition at line 64 of file G4ICRU73QOModel.cc.

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

◆ ~G4ICRU73QOModel()

G4ICRU73QOModel::~G4ICRU73QOModel ( )
virtual

Definition at line 97 of file G4ICRU73QOModel.cc.

98{}

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 146 of file G4ICRU73QOModel.cc.

152{
154 (p,kineticEnergy,cutEnergy,maxEnergy);
155 return cross;
156}
double G4double
Definition: G4Types.hh:64
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)

◆ ComputeCrossSectionPerElectron()

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

Definition at line 122 of file G4ICRU73QOModel.cc.

127{
128 G4double cross = 0.0;
129 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
130 G4double maxEnergy = std::min(tmax,maxKinEnergy);
131 if(cutEnergy < maxEnergy) {
132
133 G4double energy = kineticEnergy + mass;
134 G4double energy2 = energy*energy;
135 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
136 cross = 1.0/cutEnergy - 1.0/maxEnergy - beta2*log(maxEnergy/cutEnergy)/tmax;
137
138 cross *= CLHEP::twopi_mc2_rcl2*chargeSquare/beta2;
139 }
140
141 return cross;
142}
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Reimplemented in G4ICRU73NoDeltaModel.

Definition at line 175 of file G4ICRU73QOModel.cc.

179{
180 SetParticle(p);
181 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
182 G4double tkin = kineticEnergy/massRate;
183 G4double dedx = 0.0;
184 if(tkin > lowestKinEnergy) { dedx = DEDX(material, tkin); }
185 else { dedx = DEDX(material, lowestKinEnergy)*sqrt(tkin/lowestKinEnergy); }
186
187 if (cutEnergy < tmax) {
188
189 G4double tau = kineticEnergy/mass;
190 G4double gam = tau + 1.0;
191 G4double bg2 = tau * (tau+2.0);
192 G4double beta2 = bg2/(gam*gam);
193 G4double x = cutEnergy/tmax;
194
195 dedx += chargeSquare*( log(x) + (1.0 - x)*beta2 ) * twopi_mc2_rcl2
196 * material->GetElectronDensity()/beta2;
197 }
198 if(dedx < 0.0) { dedx = 0.0; }
199 return dedx;
200}
G4double GetElectronDensity() const
Definition: G4Material.hh:216

Referenced by G4ICRU73NoDeltaModel::ComputeDEDXPerVolume().

◆ CorrectionsAlongStep()

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

Reimplemented from G4VEmModel.

Definition at line 354 of file G4ICRU73QOModel.cc.

359{}

◆ 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 160 of file G4ICRU73QOModel.cc.

166{
167 G4double eDensity = material->GetElectronDensity();
169 (p,kineticEnergy,cutEnergy,maxEnergy);
170 return cross;
171}

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 102 of file G4ICRU73QOModel.cc.

104{
105 if(p != particle) SetParticle(p);
106
107 // always false before the run
108 SetDeexcitationFlag(false);
109
110 if(!isInitialised) {
111 isInitialised = true;
112
113 G4String pname = particle->GetParticleName();
114 fParticleChange = GetParticleChangeForLoss();
116 denEffData = (*mtab)[0]->GetIonisation()->GetDensityEffectData();
117 }
118}
std::vector< G4Material * > G4MaterialTable
static const G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:562
const G4String & GetParticleName() const
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:641
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:95

◆ MaxSecondaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 428 of file G4ICRU73QOModel.cc.

430{
431 if(pd != particle) { SetParticle(pd); }
432 G4double tau = kinEnergy/mass;
433 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
434 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
435 return tmax;
436}

Referenced by ComputeCrossSectionPerElectron(), and ComputeDEDXPerVolume().

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 363 of file G4ICRU73QOModel.cc.

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

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