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

#include <G4MuonToMuonPairProductionModel.hh>

+ Inheritance diagram for G4MuonToMuonPairProductionModel:

Public Member Functions

 G4MuonToMuonPairProductionModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="muToMuonPairProd")
 
 ~G4MuonToMuonPairProductionModel ()=default
 
G4MuonToMuonPairProductionModeloperator= (const G4MuonToMuonPairProductionModel &right)=delete
 
 G4MuonToMuonPairProductionModel (const G4MuonToMuonPairProductionModel &)=delete
 
G4double ComputeDMicroscopicCrossSection (G4double tkin, G4double Z, G4double pairEnergy) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4double U_func (G4double Z, G4double rho2, G4double xi, G4double Y, G4double pairEnergy, const G4double B=183)
 
- Public Member Functions inherited from G4MuPairProductionModel
 G4MuPairProductionModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="muPairProd")
 
 ~G4MuPairProductionModel ()=default
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double) override
 
virtual G4double ComputeDMicroscopicCrossSection (G4double tkin, G4double Z, G4double pairEnergy)
 
void SetLowestKineticEnergy (G4double e)
 
void SetParticle (const G4ParticleDefinition *)
 
G4MuPairProductionModeloperator= (const G4MuPairProductionModel &right)=delete
 
 G4MuPairProductionModel (const G4MuPairProductionModel &)=delete
 
- 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 *, 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 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 *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Member Functions

void DataCorrupted (G4int Z, G4double logTkin) const override
 
- Protected Member Functions inherited from G4MuPairProductionModel
G4double ComputMuPairLoss (G4double Z, G4double tkin, G4double cut, G4double tmax)
 
G4double ComputeMicroscopicCrossSection (G4double tkin, G4double Z, G4double cut)
 
G4double FindScaledEnergy (G4int Z, G4double rand, G4double logTkin, G4double yymin, G4double yymax)
 
G4double MaxSecondaryEnergyForElement (G4double kineticEnergy, G4double Z)
 
void MakeSamplingTables ()
 
void StoreTables () const
 
G4bool RetrieveTables ()
 
virtual void DataCorrupted (G4int Z, G4double logTkin) const
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

G4ParticleDefinitiontheMuonMinus = nullptr
 
G4ParticleDefinitiontheMuonPlus = nullptr
 
G4double factorForCross
 
G4double minPairEnergy
 
G4double muonMass
 
G4double mueRatio
 
- Protected Attributes inherited from G4MuPairProductionModel
G4ParticleChangeForLossfParticleChange = nullptr
 
const G4ParticleDefinitionparticle = nullptr
 
G4NistManagernist = nullptr
 
G4double factorForCross
 
G4double sqrte
 
G4double particleMass = 0.0
 
G4double z13 = 0.0
 
G4double z23 = 0.0
 
G4double lnZ = 0.0
 
G4double minPairEnergy
 
G4double lowestKinEnergy
 
G4double emin
 
G4double emax
 
G4double ymin = -5.0
 
G4double dy = 0.005
 
G4int currentZ = 0
 
G4int nYBinPerDecade = 4
 
size_t nbiny = 1000
 
size_t nbine = 0
 
G4bool fTableToFile = false
 
- 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
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 61 of file G4MuonToMuonPairProductionModel.hh.

Constructor & Destructor Documentation

◆ G4MuonToMuonPairProductionModel() [1/2]

G4MuonToMuonPairProductionModel::G4MuonToMuonPairProductionModel ( const G4ParticleDefinition p = nullptr,
const G4String nam = "muToMuonPairProd" 
)
explicit

◆ ~G4MuonToMuonPairProductionModel()

G4MuonToMuonPairProductionModel::~G4MuonToMuonPairProductionModel ( )
default

◆ G4MuonToMuonPairProductionModel() [2/2]

G4MuonToMuonPairProductionModel::G4MuonToMuonPairProductionModel ( const G4MuonToMuonPairProductionModel )
delete

Member Function Documentation

◆ ComputeDMicroscopicCrossSection()

G4double G4MuonToMuonPairProductionModel::ComputeDMicroscopicCrossSection ( G4double  tkin,
G4double  Z,
G4double  pairEnergy 
)
overridevirtual

Reimplemented from G4MuPairProductionModel.

Definition at line 94 of file G4MuonToMuonPairProductionModel.cc.

101{
102 if (pairEnergy <= minPairEnergy)
103 return 0.0;
104
105 G4double totalEnergy = tkin + particleMass;
106 G4double residEnergy = totalEnergy - pairEnergy;
107
108 if (residEnergy <= muonMass) { return 0.0; }
109
110 G4double a0 = 1.0 / (totalEnergy * residEnergy);
111 G4double rhomax = 1.0 - 2*muonMass/pairEnergy;
112 G4double tmnexp = 1. - rhomax;
113
114 if(tmnexp >= 1.0) { return 0.0; }
115
116 G4double tmn = G4Log(tmnexp);
117
118 G4double z2 = Z*Z;
119 G4double beta = 0.5*pairEnergy*pairEnergy*a0;
120 G4double xi0 = 0.5*beta;
121
122 // Gaussian integration in ln(1-ro) ( with 8 points)
123 G4double rho[8];
124 G4double rho2[8];
125 G4double xi[8];
126 G4double xi1[8];
127 G4double xii[8];
128
129 for (G4int i = 0; i < 8; ++i)
130 {
131 rho[i] = G4Exp(tmn*xgi[i]) - 1.0; // rho = -asymmetry
132 rho2[i] = rho[i] * rho[i];
133 xi[i] = xi0*(1.0-rho2[i]);
134 xi1[i] = 1.0 + xi[i];
135 xii[i] = 1.0 / xi[i];
136 }
137
138 G4double ximax = xi0*(1. - rhomax*rhomax);
139
140 G4double Y = 10 * sqrt(particleMass/totalEnergy);
141 G4double U[8];
142
143 for (G4int i = 0; i < 8; ++i)
144 {
145 U[i] = U_func(Z, rho2[i], xi[i], Y, pairEnergy);
146 }
147
148 G4double UMax = U_func(Z, rhomax*rhomax, ximax, Y, pairEnergy);
149
150 G4double sum = 0.0;
151 for (G4int i = 0; i < 8; ++i)
152 {
153 G4double X = 1 + U[i] - UMax;
154 G4double lnX = G4Log(X);
155 G4double phi = ((2 + rho2[i])*(1 + beta) + xi[i]*(3 + rho2[i]))*
156 G4Log(1 + xii[i]) - 1 - 3*rho2[i] + beta*(1 - 2*rho2[i])
157 + ((1 + rho2[i])*(1 + 1.5*beta) - xii[i]*(1 + 2*beta)
158 *(1 - rho2[i]))*G4Log(xi1[i]);
159 sum += wgi[i]*(1.0 + rho[i])*phi*lnX;
160 }
161
162 return -tmn*sum*factorForCross*z2*residEnergy/(totalEnergy*pairEnergy);
163
164}
G4double Y(G4double density)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
const G4double a0
G4double G4Log(G4double x)
Definition: G4Log.hh:227
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
G4double U_func(G4double Z, G4double rho2, G4double xi, G4double Y, G4double pairEnergy, const G4double B=183)

◆ DataCorrupted()

void G4MuonToMuonPairProductionModel::DataCorrupted ( G4int  Z,
G4double  logTkin 
) const
overrideprotectedvirtual

Reimplemented from G4MuPairProductionModel.

Definition at line 314 of file G4MuonToMuonPairProductionModel.cc.

315{
317 ed << "G4ElementData is not properly initialized Z= " << Z
318 << " Ekin(MeV)= " << G4Exp(logTkin)
319 << " IsMasterThread= " << IsMaster()
320 << " Model " << GetName();
321 G4Exception("G4MuonToMuonPairProductionModel","em0033",FatalException,ed,"");
322}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
const G4String & GetName() const
Definition: G4VEmModel.hh:816

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 183 of file G4MuonToMuonPairProductionModel.cc.

189{
190 G4double kinEnergy = aDynamicParticle->GetKineticEnergy();
191 //G4cout << "--- G4MuonToMuonPairProductionModel::SampleSecondaries E(MeV)= "
192 // << kinEnergy << " "
193 // << aDynamicParticle->GetDefinition()->GetParticleName() << G4endl;
194 G4double totalMomentum = std::sqrt(kinEnergy*(kinEnergy + 2.0*muonMass));
195
196 G4ThreeVector partDirection = aDynamicParticle->GetMomentumDirection();
197
198 // select randomly one element constituing the material
199 const G4Element* anElement = SelectRandomAtom(couple,particle,kinEnergy);
200
201 // define interval of energy transfer
202 G4double maxPairEnergy = MaxSecondaryEnergyForElement(kinEnergy,
203 anElement->GetZ());
204 G4double maxEnergy = std::min(tmax, maxPairEnergy);
205 G4double minEnergy = std::max(tmin, minPairEnergy);
206
207 if(minEnergy >= maxEnergy) { return; }
208 //G4cout << "emin= " << minEnergy << " emax= " << maxEnergy
209 // << " minPair= " << minPairEnergy << " maxpair= " << maxPairEnergy
210 // << " ymin= " << ymin << " dy= " << dy << G4endl;
211
212 G4double coeff = G4Log(minPairEnergy/kinEnergy)/ymin;
213
214 // compute limits
215 G4double yymin = G4Log(minEnergy/kinEnergy)/coeff;
216 G4double yymax = G4Log(maxEnergy/kinEnergy)/coeff;
217
218 //G4cout << "yymin= " << yymin << " yymax= " << yymax << G4endl;
219
220 // units should not be used, bacause table was built without
221 G4double logTkin = G4Log(kinEnergy/CLHEP::MeV);
222
223 // sample mu-mu+ pair energy first
224
225 // select sample table via Z
226 G4int iz1(0), iz2(0);
227 for(G4int iz=0; iz<nzdat; ++iz) {
228 if(currentZ == zdat[iz]) {
229 iz1 = iz2 = currentZ;
230 break;
231 } else if(currentZ < zdat[iz]) {
232 iz2 = zdat[iz];
233 if(iz > 0) { iz1 = zdat[iz-1]; }
234 else { iz1 = iz2; }
235 break;
236 }
237 }
238 if(0 == iz1) { iz1 = iz2 = zdat[nzdat-1]; }
239
240 G4double pairEnergy = 0.0;
241 G4int count = 0;
242 //G4cout << "start loop Z1= " << iz1 << " Z2= " << iz2 << G4endl;
243 do {
244 ++count;
245 // sampling using only one random number
246 G4double rand = G4UniformRand();
247
248 G4double x = FindScaledEnergy(iz1, rand, logTkin, yymin, yymax);
249 if(iz1 != iz2) {
250 G4double x2 = FindScaledEnergy(iz2, rand, logTkin, yymin, yymax);
251 G4double lz1= nist->GetLOGZ(iz1);
252 G4double lz2= nist->GetLOGZ(iz2);
253 //G4cout << count << ". x= " << x << " x2= " << x2
254 // << " Z1= " << iz1 << " Z2= " << iz2 << G4endl;
255 x += (x2 - x)*(lnZ - lz1)/(lz2 - lz1);
256 }
257 //G4cout << "x= " << x << " coeff= " << coeff << G4endl;
258 pairEnergy = kinEnergy*G4Exp(x*coeff);
259
260 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
261 } while((pairEnergy < minEnergy || pairEnergy > maxEnergy) && 10 > count);
262
263 //G4cout << "## pairEnergy(GeV)= " << pairEnergy/GeV
264 // << " Etot(GeV)= " << totalEnergy/GeV << G4endl;
265
266 // sample r=(mu+mu-)/pairEnergy ( uniformly .....)
267 G4double rmax = 1 - 2*muonMass/(pairEnergy);
268 G4double r = rmax * (-1.+2.*G4UniformRand()) ;
269
270 // compute energies from pairEnergy,r
271 G4double muMinusEnergy = (1.-r)*pairEnergy*0.5;
272 G4double muPlusEnergy = pairEnergy - muMinusEnergy;
273
274 // Sample angles
275 G4ThreeVector muMinusDirection, muPlusDirection;
276 //
277 GetAngularDistribution()->SamplePairDirections(aDynamicParticle,
278 muMinusEnergy, muPlusEnergy,
279 muMinusDirection, muPlusDirection);
280 // create G4DynamicParticle object for mu+mu-
281 muMinusEnergy = std::max(muMinusEnergy - muonMass, 0.0);
282 muPlusEnergy = std::max(muPlusEnergy - muonMass, 0.0);
283 G4DynamicParticle* aParticle1 =
284 new G4DynamicParticle(theMuonMinus,muMinusDirection,muMinusEnergy);
285 G4DynamicParticle* aParticle2 =
286 new G4DynamicParticle(theMuonPlus,muPlusDirection,muPlusEnergy);
287 // Fill output vector
288 vdp->push_back(aParticle1);
289 vdp->push_back(aParticle2);
290
291 // primary change
292 kinEnergy -= pairEnergy;
293 partDirection *= totalMomentum;
294 partDirection -= (aParticle1->GetMomentum() + aParticle2->GetMomentum());
295 partDirection = partDirection.unit();
296
297 // if energy transfer is higher than threshold (very high by default)
298 // then stop tracking the primary particle and create a new secondary
299 if (pairEnergy > SecondaryThreshold()) {
302 G4DynamicParticle* newdp =
303 new G4DynamicParticle(particle, partDirection, kinEnergy);
304 vdp->push_back(newdp);
305 } else { // continue tracking the primary e-/e+ otherwise
308 }
309 //G4cout << "--- G4MuonToMuonPairProductionModel::SampleSecondaries done" << G4endl;
310}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetZ() const
Definition: G4Element.hh:131
G4double FindScaledEnergy(G4int Z, G4double rand, G4double logTkin, G4double yymin, G4double yymax)
const G4ParticleDefinition * particle
G4ParticleChangeForLoss * fParticleChange
G4double MaxSecondaryEnergyForElement(G4double kineticEnergy, G4double Z)
G4double GetLOGZ(G4int Z) const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
virtual void SamplePairDirections(const G4DynamicParticle *dp, G4double elecKinEnergy, G4double posiKinEnergy, G4ThreeVector &dirElectron, G4ThreeVector &dirPositron, G4int Z=0, const G4Material *mat=nullptr)
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:600
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:561
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:669
void ProposeTrackStatus(G4TrackStatus status)

◆ U_func()

G4double G4MuonToMuonPairProductionModel::U_func ( G4double  Z,
G4double  rho2,
G4double  xi,
G4double  Y,
G4double  pairEnergy,
const G4double  B = 183 
)

Definition at line 166 of file G4MuonToMuonPairProductionModel.cc.

170{
171 G4int Z = G4lrint(ZZ);
172 G4double A27 = nist->GetA27(Z);
173 G4double Z13 = nist->GetZ13(Z);
174 static const G4double sqe = std::sqrt(G4Exp(1.0));
175 G4double res = (0.65 * B / (A27*Z13) * mueRatio)/
176 (1 + (2*sqe*muonMass*muonMass*(B/Z13)*(1 + xi)*(1 + Y))
177 /(CLHEP::electron_mass_c2*pairEnergy*(1 - rho2)));
178 return res;
179}
G4double B(G4double temperature)
G4double GetA27(G4int Z) const
G4double GetZ13(G4double Z) const
int G4lrint(double ad)
Definition: templates.hh:134

Referenced by ComputeDMicroscopicCrossSection().

Member Data Documentation

◆ factorForCross

G4double G4MuonToMuonPairProductionModel::factorForCross
protected

◆ minPairEnergy

G4double G4MuonToMuonPairProductionModel::minPairEnergy
protected

◆ mueRatio

G4double G4MuonToMuonPairProductionModel::mueRatio
protected

Definition at line 98 of file G4MuonToMuonPairProductionModel.hh.

Referenced by G4MuonToMuonPairProductionModel(), and U_func().

◆ muonMass

G4double G4MuonToMuonPairProductionModel::muonMass
protected

◆ theMuonMinus

G4ParticleDefinition* G4MuonToMuonPairProductionModel::theMuonMinus = nullptr
protected

◆ theMuonPlus

G4ParticleDefinition* G4MuonToMuonPairProductionModel::theMuonPlus = nullptr
protected

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