Geant4 10.7.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")
 
 ~G4ICRU73QOModel ()=default
 
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 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 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 CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=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)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, 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 *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
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) const
 
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
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () 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)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Member Functions

virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy) final
 
- 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
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Detailed Description

Definition at line 70 of file G4ICRU73QOModel.hh.

Constructor & Destructor Documentation

◆ G4ICRU73QOModel()

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

Definition at line 66 of file G4ICRU73QOModel.cc.

67 : G4VEmModel(nam),
68 particle(nullptr),
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 = nullptr;
94 denEffData = nullptr;
95}
int G4int
Definition: G4Types.hh:85
static G4Electron * Electron()
Definition: G4Electron.cc:93
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757

◆ ~G4ICRU73QOModel()

G4ICRU73QOModel::~G4ICRU73QOModel ( )
default

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 148 of file G4ICRU73QOModel.cc.

154{
156 (p,kineticEnergy,cutEnergy,maxEnergy);
157 return cross;
158}
double G4double
Definition: G4Types.hh:83
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) final

◆ ComputeCrossSectionPerElectron()

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

Definition at line 123 of file G4ICRU73QOModel.cc.

128{
129 G4double cross = 0.0;
130 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
131 G4double maxEnergy = std::min(tmax,maxKinEnergy);
132 if(cutEnergy < maxEnergy) {
133
134 G4double energy = kineticEnergy + mass;
135 G4double energy2 = energy*energy;
136 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
137 cross = (maxEnergy - cutEnergy)/(cutEnergy*maxEnergy)
138 - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
139
140 cross *= CLHEP::twopi_mc2_rcl2*chargeSquare/beta2;
141 }
142
143 return cross;
144}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
G4double energy(const ThreeVector &p, const G4double m)

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Reimplemented in G4ICRU73NoDeltaModel.

Definition at line 177 of file G4ICRU73QOModel.cc.

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

Referenced by G4ICRU73NoDeltaModel::ComputeDEDXPerVolume().

◆ CorrectionsAlongStep()

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

Reimplemented from G4VEmModel.

Definition at line 405 of file G4ICRU73QOModel.cc.

410{}

◆ CrossSectionPerVolume()

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

168{
169 G4double eDensity = material->GetElectronDensity();
171 (p,kineticEnergy,cutEnergy,maxEnergy);
172 return cross;
173}

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 99 of file G4ICRU73QOModel.cc.

101{
102 if(p != particle) SetParticle(p);
103
104 // always false before the run
105 SetDeexcitationFlag(false);
106
107 if(!isInitialised) {
108 isInitialised = true;
109
112 }
113
114 G4String pname = particle->GetParticleName();
115 fParticleChange = GetParticleChangeForLoss();
117 denEffData = (*mtab)[0]->GetIonisation()->GetDensityEffectData();
118 }
119}
std::vector< G4Material * > G4MaterialTable
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:637
const G4String & GetParticleName() const
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:611
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:813
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:708
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118

◆ MaxSecondaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 491 of file G4ICRU73QOModel.cc.

493{
494 if(pd != particle) { SetParticle(pd); }
495 G4double tau = kinEnergy/mass;
496 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
497 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
498 return tmax;
499}

Referenced by ComputeCrossSectionPerElectron(), and ComputeDEDXPerVolume().

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 414 of file G4ICRU73QOModel.cc.

419{
421 G4double xmax = std::min(tmax, maxEnergy);
422 if(xmin >= xmax) { return; }
423
424 G4double kineticEnergy = dp->GetKineticEnergy();
425 G4double energy = kineticEnergy + mass;
426 G4double energy2 = energy*energy;
427 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
428 G4double grej = 1.0;
429 G4double deltaKinEnergy, f;
430
431 G4ThreeVector direction = dp->GetMomentumDirection();
432
433 // sampling follows ...
434 do {
436 deltaKinEnergy = xmin*xmax/(xmin*(1.0 - x) + xmax*x);
437
438 f = 1.0 - beta2*deltaKinEnergy/tmax;
439
440 if(f > grej) {
441 G4cout << "G4ICRU73QOModel::SampleSecondary Warning! "
442 << "Majorant " << grej << " < "
443 << f << " for e= " << deltaKinEnergy
444 << G4endl;
445 }
446
447 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
448 } while( grej*G4UniformRand() >= f );
449
450 G4ThreeVector deltaDirection;
451
453 const G4Material* mat = couple->GetMaterial();
455
456 deltaDirection =
457 GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
458
459 } else {
460
461 G4double deltaMomentum =
462 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
463 G4double totMomentum = energy*sqrt(beta2);
464 G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
465 (deltaMomentum * totMomentum);
466 if(cost > 1.0) { cost = 1.0; }
467 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
468
469 G4double phi = twopi * G4UniformRand() ;
470
471 deltaDirection.set(sint*cos(phi),sint*sin(phi), cost) ;
472 deltaDirection.rotateUz(direction);
473 }
474 // create G4DynamicParticle object for delta ray
475 G4DynamicParticle* delta =
476 new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
477
478 // Change kinematics of primary particle
479 kineticEnergy -= deltaKinEnergy;
480 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
481 finalP = finalP.unit();
482
483 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
484 fParticleChange->SetProposedMomentumDirection(finalP);
485
486 vdp->push_back(delta);
487}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.cc:315
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:510

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