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

#include <G4eplusTo2or3GammaModel.hh>

+ Inheritance diagram for G4eplusTo2or3GammaModel:

Public Member Functions

 G4eplusTo2or3GammaModel ()
 
 ~G4eplusTo2or3GammaModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double ComputeCrossSectionPerElectron (G4double kinEnergy)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0., G4double maxEnergy=DBL_MAX) override
 
G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double maxEnergy=DBL_MAX) override
 
void SetDelta (G4double val)
 
G4eplusTo2or3GammaModeloperator= (const G4eplusTo2or3GammaModel &right)=delete
 
 G4eplusTo2or3GammaModel (const G4eplusTo2or3GammaModel &)=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 ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
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 G4MaterialCutsCouple *, const G4DynamicParticle *, 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
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- 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 66 of file G4eplusTo2or3GammaModel.hh.

Constructor & Destructor Documentation

◆ G4eplusTo2or3GammaModel() [1/2]

G4eplusTo2or3GammaModel::G4eplusTo2or3GammaModel ( )

Definition at line 68 of file G4eplusTo2or3GammaModel.cc.

69 : G4VEmModel("eplusTo2or3gamma"),
70 fDeltaMin(0.001),
71 fDelta(fDeltaMin),
72 fGammaTh(CLHEP::MeV)
73{
74 theGamma = G4Gamma::Gamma();
75 fParticleChange = nullptr;
76 f3GModel = new G4eplusTo3GammaOKVIModel();
77 SetTripletModel(f3GModel);
78
79 // instantiate vectors once
80 if (nullptr == fCrossSection) {
81 G4double emin = 10*CLHEP::eV;
82 G4double emax = 100*CLHEP::TeV;
83 G4int nbins = 20*G4lrint(std::log10(emax/emin));
84 fCrossSection = new G4PhysicsLogVector(emin, emax, nbins, true);
85 f3GProbability = new G4PhysicsLogVector(emin, emax, nbins, true);
86 }
87}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
void SetTripletModel(G4VEmModel *)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
int G4lrint(double ad)
Definition templates.hh:134

Referenced by G4eplusTo2or3GammaModel(), and operator=().

◆ ~G4eplusTo2or3GammaModel()

G4eplusTo2or3GammaModel::~G4eplusTo2or3GammaModel ( )
override

Definition at line 91 of file G4eplusTo2or3GammaModel.cc.

92{
93 if (IsMaster()) {
94 delete fCrossSection;
95 delete f3GProbability;
96 fCrossSection = nullptr;
97 f3GProbability = nullptr;
98 }
99}
G4bool IsMaster() const

◆ G4eplusTo2or3GammaModel() [2/2]

G4eplusTo2or3GammaModel::G4eplusTo2or3GammaModel ( const G4eplusTo2or3GammaModel & )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4eplusTo2or3GammaModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition * ,
G4double kinEnergy,
G4double Z,
G4double A = 0.,
G4double cutEnergy = 0.,
G4double maxEnergy = DBL_MAX )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 161 of file G4eplusTo2or3GammaModel.cc.

165{
166 // Calculates the cross section per atom of annihilation into two photons
167 G4double cross = Z*fCrossSection->Value(kineticEnergy);
168 return cross;
169}

◆ ComputeCrossSectionPerElectron()

G4double G4eplusTo2or3GammaModel::ComputeCrossSectionPerElectron ( G4double kinEnergy)

Definition at line 134 of file G4eplusTo2or3GammaModel.cc.

135{
136 // Calculates the cross section per electron of annihilation into two
137 // photons from the Heilter formula with the radiation correction to 3 gamma
138 // annihilation channel. (A.A.) rho is changed
139
140 G4double ekin = std::max(CLHEP::eV, kinEnergy);
141 G4double tau = ekin/CLHEP::electron_mass_c2;
142 G4double gam = tau + 1.0;
143 G4double gamma2 = gam*gam;
144 G4double bg2 = tau * (tau+2.0);
145 G4double bg = std::sqrt(bg2);
146 G4double rho = (gamma2+4.*gam+1.)*G4Log(gam+bg)/(gamma2-1.)
147 - (gam+3.)/(std::sqrt(gam*gam - 1.));
148 G4double eGammaCMS = CLHEP::electron_mass_c2 * std::sqrt(0.5*(tau + 2.0));
149 fDelta = std::max(fDeltaMin, fGammaTh/eGammaCMS);
150 f3GModel->SetDelta(fDelta);
151
152 static const G4double pir2 =
153 CLHEP::pi*CLHEP::classic_electr_radius*CLHEP::classic_electr_radius;
154 G4double cross = (pir2*rho + alpha_rcl2*2.*G4Log(fDelta)*rho*rho)/(gam+1.);
155
156 return cross;
157}
G4double G4Log(G4double x)
Definition G4Log.hh:227

Referenced by Initialise().

◆ CrossSectionPerVolume()

G4double G4eplusTo2or3GammaModel::CrossSectionPerVolume ( const G4Material * material,
const G4ParticleDefinition * ,
G4double kineticEnergy,
G4double cutEnergy = 0.0,
G4double maxEnergy = DBL_MAX )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 173 of file G4eplusTo2or3GammaModel.cc.

178{
179 // Calculates the cross section per volume of annihilation into two photons
180 G4double eDensity = material->GetElectronDensity();
181 G4double cross = eDensity*fCrossSection->Value(kineticEnergy);
182 return cross;
183}
G4double GetElectronDensity() const

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 103 of file G4eplusTo2or3GammaModel.cc.

105{
106 // here particle change is set for the triplet model
107 if (nullptr == fParticleChange) {
108 fParticleChange = GetParticleChangeForGamma();
109 }
110 // initialialise 3-gamma model before new run
111 f3GModel->Initialise(p, cuts);
113
114 // initialise vectors
115 if (IsMaster()) {
116 std::size_t num = fCrossSection->GetVectorLength();
117 for (std::size_t i=0; i<num; ++i) {
118 G4double e = fCrossSection->Energy(i);
120 G4double cs3 = f3GModel->ComputeCrossSectionPerElectron(e);
121 cs2 += cs3;
122 fCrossSection->PutValue(i, cs2);
123 G4double y = (cs2 > 0.0) ? cs3/cs2 : 0.0;
124 f3GProbability->PutValue(i, y);
125 }
126 fCrossSection->FillSecondDerivatives();
127 f3GProbability->FillSecondDerivatives();
128 }
129}
static G4EmParameters * Instance()
G4double LowestTripletEnergy() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double ComputeCrossSectionPerElectron(G4double kinEnergy)

◆ operator=()

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

◆ SampleSecondaries()

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

!! boost of polarisation vector is not yet implemented

Implements G4VEmModel.

Definition at line 190 of file G4eplusTo2or3GammaModel.cc.

195{
196 // kill primary positron
197 fParticleChange->SetProposedKineticEnergy(0.0);
198 fParticleChange->ProposeTrackStatus(fStopAndKill);
199
200 // Case at rest not considered anymore
201 G4double posiKinEnergy = dp->GetKineticEnergy();
203 posiKinEnergy + 2*CLHEP::electron_mass_c2);
204 G4double eGammaCMS = 0.5 * lv.mag();
205
206 if (G4UniformRand() < f3GProbability->Value(posiKinEnergy)) {
207 fDelta = std::max(fDeltaMin, fGammaTh/eGammaCMS);
208 f3GModel->SetDelta(fDelta);
209 f3GModel->SampleSecondaries(vdp, couple, dp);
210 return;
211 }
212
214 G4double phi = CLHEP::twopi * G4UniformRand();
215 G4double cosphi = std::cos(phi);
216 G4double sinphi = std::sin(phi);
217 G4ThreeVector pol1(cosphi, sinphi, 0.0);
218 pol1.rotateUz(dir1);
219 G4LorentzVector lv1(eGammaCMS*dir1, eGammaCMS);
220
221 G4ThreeVector pol2(-sinphi, cosphi, 0.0);
222 pol2.rotateUz(dir1);
223
224 // transformation to lab system
225 lv1.boost(lv.boostVector());
226 lv -= lv1;
227
228 //!!! boost of polarisation vector is not yet implemented
229
230 // use constructors optimal for massless particle
231 auto aGamma1 = new G4DynamicParticle(G4Gamma::Gamma(), lv1.vect());
232 aGamma1->SetPolarization(pol1);
233 auto aGamma2 = new G4DynamicParticle(G4Gamma::Gamma(), lv.vect());
234 aGamma2->SetPolarization(pol2);
235
236 vdp->push_back(aGamma1);
237 vdp->push_back(aGamma2);
238}
CLHEP::HepLorentzVector G4LorentzVector
G4ThreeVector G4RandomDirection()
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
#define G4UniformRand()
Definition Randomize.hh:52
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const

◆ SetDelta()

void G4eplusTo2or3GammaModel::SetDelta ( G4double val)
inline

Definition at line 99 of file G4eplusTo2or3GammaModel.hh.

99{ if(val > 0.0) { fDeltaMin = val; } };

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