Geant4 11.4.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=nullptr, const G4String &nam="ICRU73QO")
 ~G4ICRU73QOModel ()=default
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double) override
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4ICRU73QOModeloperator= (const G4ICRU73QOModel &right)=delete
 G4ICRU73QOModel (const G4ICRU73QOModel &)=delete
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 void CorrectionsAlongStep (const G4Material *, const G4ParticleDefinition *, const G4double kinEnergy, const G4double cutEnergy, const G4double &length, G4double &eloss)
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 FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
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)
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)
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
G4int SelectRandomAtomNumber (const G4Material *) const
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
G4int SelectIsotopeNumber (const G4Element *) const
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 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 SetDeexcitationFlag (G4bool val)
void SetForceBuildTable (G4bool val)
void SetFluctuationFlag (G4bool val)
G4bool IsMaster () const
void SetUseBaseMaterials (G4bool val)
G4bool UseBaseMaterials () const
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
const G4StringGetName () const
void SetCurrentCouple (const G4MaterialCutsCouple *)
G4bool IsLocked () const
void SetLocked (G4bool)
void SetLPMFlag (G4bool)
void SetMasterThread (G4bool)
G4VEmModeloperator= (const G4VEmModel &right)=delete
 G4VEmModel (const G4VEmModel &)=delete

Protected Member Functions

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 = nullptr
G4VParticleChangepParticleChange = nullptr
G4PhysicsTablexSectionTable = nullptr
const G4MaterialpBaseMaterial = nullptr
const std::vector< G4double > * theDensityFactor = nullptr
const std::vector< G4int > * theDensityIdx = nullptr
G4double inveplus
G4double pFactor = 1.0
std::size_t currentCoupleIndex = 0
std::size_t basedCoupleIndex = 0
G4bool lossFlucFlag = true

Detailed Description

Definition at line 70 of file G4ICRU73QOModel.hh.

Constructor & Destructor Documentation

◆ G4ICRU73QOModel() [1/2]

G4ICRU73QOModel::G4ICRU73QOModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "ICRU73QO" )
explicit

Definition at line 64 of file G4ICRU73QOModel.cc.

65 : G4VEmModel(nam),
66 particle(nullptr),
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 = nullptr;
92 denEffData = nullptr;
93}
int G4int
Definition G4Types.hh:85
static G4Electron * Electron()
Definition G4Electron.cc:91
void SetHighEnergyLimit(G4double)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

Referenced by G4ICRU73NoDeltaModel::G4ICRU73NoDeltaModel(), G4ICRU73QOModel(), and operator=().

◆ ~G4ICRU73QOModel()

G4ICRU73QOModel::~G4ICRU73QOModel ( )
default

◆ G4ICRU73QOModel() [2/2]

G4ICRU73QOModel::G4ICRU73QOModel ( const G4ICRU73QOModel & )
delete

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

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

◆ ComputeCrossSectionPerElectron()

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

Definition at line 121 of file G4ICRU73QOModel.cc.

126{
127 G4double cross = 0.0;
128 const G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
129 const G4double maxEnergy = std::min(tmax, maxKinEnergy);
130 const G4double cutEnergy = std::max(cut, lowestKinEnergy*massRate);
131 if(cutEnergy < maxEnergy) {
132
133 const G4double energy = kineticEnergy + mass;
134 const G4double energy2 = energy*energy;
135 const G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
136 cross = (maxEnergy - cutEnergy)/(cutEnergy*maxEnergy)
137 - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
138
139 cross *= CLHEP::twopi_mc2_rcl2*chargeSquare/beta2;
140 }
141
142 return cross;
143}
G4double G4Log(G4double x)
Definition G4Log.hh:169
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 cut )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 176 of file G4ICRU73QOModel.cc.

180{
181 SetParticle(p);
182 const G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
183 const G4double tkin = kineticEnergy/massRate;
184 const G4double cutEnergy = std::max(cut, lowestKinEnergy*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 const G4double tau = kineticEnergy/mass;
192 const G4double x = cutEnergy/tmax;
193 dedx += (G4Log(x)*(tau + 1.)*(tau + 1.)/(tau * (tau + 2.0)) + 1.0 - x) *
194 CLHEP::twopi_mc2_rcl2 *chargeSquare * material->GetElectronDensity();
195 }
196 dedx = std::max(dedx, 0.0);
197 return dedx;
198}
G4double GetElectronDensity() const

Referenced by G4ICRU73NoDeltaModel::ComputeDEDXPerVolume().

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 161 of file G4ICRU73QOModel.cc.

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

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 97 of file G4ICRU73QOModel.cc.

99{
100 if(p != particle) SetParticle(p);
101
102 // always false before the run
103 SetDeexcitationFlag(false);
104
105 if(!isInitialised) {
106 isInitialised = true;
107
109 SetAngularDistribution(new G4DeltaAngle());
110 }
111
112 G4String pname = particle->GetParticleName();
113 fParticleChange = GetParticleChangeForLoss();
115 denEffData = (*mtab)[0]->GetIonisation()->GetDensityEffectData();
116 }
117}
std::vector< G4Material * > G4MaterialTable
static G4MaterialTable * GetMaterialTable()
G4VEmAngularDistribution * GetAngularDistribution()
void SetDeexcitationFlag(G4bool val)
void SetAngularDistribution(G4VEmAngularDistribution *)
G4bool UseAngularGeneratorFlag() const
G4ParticleChangeForLoss * GetParticleChangeForLoss()

◆ MaxSecondaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 478 of file G4ICRU73QOModel.cc.

480{
481 if(pd != particle) { SetParticle(pd); }
482 G4double tau = kinEnergy/mass;
483 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
484 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
485 return tmax;
486}

Referenced by ComputeCrossSectionPerElectron(), and ComputeDEDXPerVolume().

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 401 of file G4ICRU73QOModel.cc.

406{
407 const G4double tmax = MaxSecondaryKinEnergy(dp);
408 const G4double xmax = std::min(tmax, maxEnergy);
409 const G4double xmin = std::max(minEnergy, lowestKinEnergy*massRate);
410 if(xmin >= xmax) { return; }
411
412 G4double kineticEnergy = dp->GetKineticEnergy();
413 const G4double energy = kineticEnergy + mass;
414 const G4double energy2 = energy*energy;
415 const G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
416 G4double grej = 1.0;
417 G4double deltaKinEnergy, f;
418
419 G4ThreeVector direction = dp->GetMomentumDirection();
420
421 // sampling follows ...
422 do {
424 deltaKinEnergy = xmin*xmax/(xmin*(1.0 - x) + xmax*x);
425
426 f = 1.0 - beta2*deltaKinEnergy/tmax;
427
428 if(f > grej) {
429 G4cout << "G4ICRU73QOModel::SampleSecondary Warning! "
430 << "Majorant " << grej << " < "
431 << f << " for e= " << deltaKinEnergy
432 << G4endl;
433 }
434
435 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
436 } while( grej*G4UniformRand() >= f );
437
438 G4ThreeVector deltaDirection;
439
441 const G4Material* mat = couple->GetMaterial();
443
444 deltaDirection =
445 GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
446
447 } else {
448
449 G4double deltaMomentum =
450 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
451 G4double totMomentum = energy*sqrt(beta2);
452 G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
453 (deltaMomentum * totMomentum);
454 if(cost > 1.0) { cost = 1.0; }
455 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
456
457 G4double phi = twopi * G4UniformRand() ;
458
459 deltaDirection.set(sint*cos(phi),sint*sin(phi), cost) ;
460 deltaDirection.rotateUz(direction);
461 }
462 // create G4DynamicParticle object for delta ray
463 auto delta = new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
464
465 // Change kinematics of primary particle
466 kineticEnergy -= deltaKinEnergy;
467 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
468 finalP = finalP.unit();
469
470 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
471 fParticleChange->SetProposedMomentumDirection(finalP);
472
473 vdp->push_back(delta);
474}
CLHEP::Hep3Vector G4ThreeVector
#define G4endl
Definition G4ios.hh:67
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 &)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
const G4Material * GetMaterial() const
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4int SelectRandomAtomNumber(const G4Material *) const
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)

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