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

#include <G4MuPairProductionModel.hh>

+ Inheritance diagram for G4MuPairProductionModel:

Public Member Functions

 G4MuPairProductionModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="muPairProd")
 
virtual ~G4MuPairProductionModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double) override
 
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 *, 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

G4double ComputMuPairLoss (G4double Z, G4double tkin, G4double cut, G4double tmax)
 
G4double ComputeMicroscopicCrossSection (G4double tkin, G4double Z, G4double cut)
 
virtual G4double ComputeDMicroscopicCrossSection (G4double tkin, G4double Z, G4double pairEnergy)
 
G4double MaxSecondaryEnergyForElement (G4double kineticEnergy, G4double Z)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

const G4ParticleDefinitionparticle
 
G4NistManagernist
 
G4double factorForCross
 
G4double sqrte
 
G4double particleMass
 
G4double z13
 
G4double z23
 
G4double lnZ
 
G4int currentZ
 
G4ParticleDefinitiontheElectron
 
G4ParticleDefinitionthePositron
 
G4ParticleChangeForLossfParticleChange
 
G4double minPairEnergy
 
G4double lowestKinEnergy
 
G4int nYBinPerDecade
 
size_t nbiny
 
size_t nbine
 
G4double ymin
 
G4double dy
 
G4double emin
 
G4double emax
 
G4bool fTableToFile
 
- 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 71 of file G4MuPairProductionModel.hh.

Constructor & Destructor Documentation

◆ G4MuPairProductionModel() [1/2]

G4MuPairProductionModel::G4MuPairProductionModel ( const G4ParticleDefinition p = nullptr,
const G4String nam = "muPairProd" 
)
explicit

Definition at line 112 of file G4MuPairProductionModel.cc.

114 : G4VEmModel(nam),
115 particle(nullptr),
116 factorForCross(4.*fine_structure_const*fine_structure_const
117 *classic_electr_radius*classic_electr_radius/(3.*pi)),
118 sqrte(sqrt(G4Exp(1.))),
119 currentZ(0),
120 fParticleChange(nullptr),
121 minPairEnergy(4.*electron_mass_c2),
122 lowestKinEnergy(1.0*GeV),
124 nbiny(1000),
125 nbine(0),
126 ymin(-5.),
127 dy(0.005),
128 fTableToFile(false)
129{
131
134
135 particleMass = lnZ = z13 = z23 = 0.;
136
137 // setup lowest limit dependent on particle mass
138 if(p) {
139 SetParticle(p);
140 lowestKinEnergy = std::max(lowestKinEnergy,p->GetPDGMass()*8.0);
141 }
143 emax = 10.*TeV;
145}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
static G4Electron * Electron()
Definition: G4Electron.cc:93
void SetParticle(const G4ParticleDefinition *)
G4ParticleDefinition * thePositron
const G4ParticleDefinition * particle
G4ParticleChangeForLoss * fParticleChange
G4ParticleDefinition * theElectron
static G4NistManager * Instance()
static G4Positron * Positron()
Definition: G4Positron.cc:93
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618

◆ ~G4MuPairProductionModel()

G4MuPairProductionModel::~G4MuPairProductionModel ( )
virtual

Definition at line 149 of file G4MuPairProductionModel.cc.

150{}

◆ G4MuPairProductionModel() [2/2]

G4MuPairProductionModel::G4MuPairProductionModel ( const G4MuPairProductionModel )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 447 of file G4MuPairProductionModel.cc.

453{
454 G4double cross = 0.0;
455 if (kineticEnergy <= lowestKinEnergy) { return cross; }
456
457 G4double maxPairEnergy = MaxSecondaryEnergyForElement(kineticEnergy, Z);
458 G4double tmax = std::min(maxEnergy, maxPairEnergy);
459 G4double cut = std::max(cutEnergy, minPairEnergy);
460 if (cut >= tmax) { return cross; }
461
462 cross = ComputeMicroscopicCrossSection(kineticEnergy, Z, cut);
463 if(tmax < kineticEnergy) {
464 cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
465 }
466 return cross;
467}
double G4double
Definition: G4Types.hh:83
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)
G4double MaxSecondaryEnergyForElement(G4double kineticEnergy, G4double Z)

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 208 of file G4MuPairProductionModel.cc.

213{
214 G4double dedx = 0.0;
215 if (cutEnergy <= minPairEnergy || kineticEnergy <= lowestKinEnergy)
216 { return dedx; }
217
218 const G4ElementVector* theElementVector = material->GetElementVector();
219 const G4double* theAtomicNumDensityVector =
220 material->GetAtomicNumDensityVector();
221
222 // loop for elements in the material
223 for (size_t i=0; i<material->GetNumberOfElements(); ++i) {
224 G4double Z = (*theElementVector)[i]->GetZ();
225 G4double tmax = MaxSecondaryEnergyForElement(kineticEnergy, Z);
226 G4double loss = ComputMuPairLoss(Z, kineticEnergy, cutEnergy, tmax);
227 dedx += loss*theAtomicNumDensityVector[i];
228 }
229 dedx = std::max(dedx, 0.0);
230 return dedx;
231}
std::vector< G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:214
G4double ComputMuPairLoss(G4double Z, G4double tkin, G4double cut, G4double tmax)

◆ ComputeDMicroscopicCrossSection()

G4double G4MuPairProductionModel::ComputeDMicroscopicCrossSection ( G4double  tkin,
G4double  Z,
G4double  pairEnergy 
)
protectedvirtual

Definition at line 310 of file G4MuPairProductionModel.cc.

317{
318 static const G4double bbbtf= 183. ;
319 static const G4double bbbh = 202.4 ;
320 static const G4double g1tf = 1.95e-5 ;
321 static const G4double g2tf = 5.3e-5 ;
322 static const G4double g1h = 4.4e-5 ;
323 static const G4double g2h = 4.8e-5 ;
324
325 if (pairEnergy <= 4.0 * electron_mass_c2)
326 return 0.0;
327
328 G4double totalEnergy = tkin + particleMass;
329 G4double residEnergy = totalEnergy - pairEnergy;
330
331 if (residEnergy <= 0.75*sqrte*z13*particleMass)
332 return 0.0;
333
334 G4double a0 = 1.0 / (totalEnergy * residEnergy);
335 G4double alf = 4.0 * electron_mass_c2 / pairEnergy;
336 G4double rt = sqrt(1.0 - alf);
337 G4double delta = 6.0 * particleMass * particleMass * a0;
338 G4double tmnexp = alf/(1.0 + rt) + delta*rt;
339
340 if(tmnexp >= 1.0) { return 0.0; }
341
342 G4double tmn = G4Log(tmnexp);
343
344 G4double massratio = particleMass/electron_mass_c2;
345 G4double massratio2 = massratio*massratio;
346 G4double inv_massratio2 = 1.0 / massratio2;
347
348 // zeta calculation
349 G4double bbb,g1,g2;
350 if( Z < 1.5 ) { bbb = bbbh ; g1 = g1h ; g2 = g2h ; }
351 else { bbb = bbbtf; g1 = g1tf; g2 = g2tf; }
352
353 G4double zeta = 0.0;
354 G4double z1exp = totalEnergy / (particleMass + g1*z23*totalEnergy);
355
356 // 35.221047195922 is the root of zeta1(x) = 0.073 * log(x) - 0.26, so the
357 // condition below is the same as zeta1 > 0.0, but without calling log(x)
358 if (z1exp > 35.221047195922)
359 {
360 G4double z2exp = totalEnergy / (particleMass + g2*z13*totalEnergy);
361 zeta = (0.073 * G4Log(z1exp) - 0.26) / (0.058 * G4Log(z2exp) - 0.14);
362 }
363
364 G4double z2 = Z*(Z+zeta);
365 G4double screen0 = 2.*electron_mass_c2*sqrte*bbb/(z13*pairEnergy);
366 G4double beta = 0.5*pairEnergy*pairEnergy*a0;
367 G4double xi0 = 0.5*massratio2*beta;
368
369 // Gaussian integration in ln(1-ro) ( with 8 points)
370 G4double rho[8];
371 G4double rho2[8];
372 G4double xi[8];
373 G4double xi1[8];
374 G4double xii[8];
375
376 for (G4int i = 0; i < 8; ++i)
377 {
378 rho[i] = G4Exp(tmn*xgi[i]) - 1.0; // rho = -asymmetry
379 rho2[i] = rho[i] * rho[i];
380 xi[i] = xi0*(1.0-rho2[i]);
381 xi1[i] = 1.0 + xi[i];
382 xii[i] = 1.0 / xi[i];
383 }
384
385 G4double ye1[8];
386 G4double ym1[8];
387
388 G4double b40 = 4.0 * beta;
389 G4double b62 = 6.0 * beta + 2.0;
390
391 for (G4int i = 0; i < 8; ++i)
392 {
393 G4double yeu = (b40 + 5.0) + (b40 - 1.0) * rho2[i];
394 G4double yed = b62 * G4Log(3.0 + xii[i]) + (2.0 * beta - 1.0) * rho2[i] - b40;
395
396 G4double ymu = b62 * (1.0 + rho2[i]) + 6.0;
397 G4double ymd = (b40 + 3.0) * (1.0 + rho2[i]) * G4Log(3.0 + xi[i]) + 2.0 - 3.0 * rho2[i];
398
399 ye1[i] = 1.0 + yeu / yed;
400 ym1[i] = 1.0 + ymu / ymd;
401 }
402
403 G4double be[8];
404 G4double bm[8];
405
406 for(G4int i = 0; i < 8; ++i)
407 {
408 if(xi[i] <= 1000.0) {
409 be[i] = ((2.0 + rho2[i]) * (1.0 + beta) + xi[i] * (3.0 + rho2[i])) *
410 G4Log(1.0 + xii[i]) + (1.0 - rho2[i] - beta) / xi1[i] - (3.0 + rho2[i]);
411 } else {
412 be[i] = 0.5 * (3.0 - rho2[i] + 2.0 * beta * (1.0 + rho2[i])) * xii[i];
413 }
414
415 if(xi[i] >= 0.001) {
416 G4double a10 = (1.0 + 2.0 * beta) * (1.0 - rho2[i]);
417 bm[i] = ((1.0 + rho2[i]) * (1.0 + 1.5 * beta) - a10 * xii[i]) * G4Log(xi1[i]) +
418 xi[i] * (1.0 - rho2[i] - beta) / xi1[i] + a10;
419 } else {
420 bm[i] = 0.5 * (5.0 - rho2[i] + beta * (3.0 + rho2[i])) * xi[i];
421 }
422 }
423
424 G4double sum = 0.0;
425
426 for (G4int i = 0; i < 8; ++i)
427 {
428 G4double screen = screen0*xi1[i]/(1.0-rho2[i]) ;
429 G4double ale = G4Log(bbb/z13*sqrt(xi1[i]*ye1[i])/(1.+screen*ye1[i])) ;
430 G4double cre = 0.5*G4Log(1.+2.25*z23*xi1[i]*ye1[i]*inv_massratio2) ;
431
432 G4double fe = (ale-cre)*be[i];
433 fe *= (fe > 0.0);
434
435 G4double alm_crm = G4Log(bbb*massratio/(1.5*z23*(1.+screen*ym1[i])));
436 G4double fm = alm_crm*bm[i];
437 fm *= (fm > 0.0) * inv_massratio2;
438
439 sum += wgi[i]*(1.0 + rho[i])*(fe+fm);
440 }
441
442 return -tmn*sum*factorForCross*z2*residEnergy/(totalEnergy*pairEnergy);
443}
const G4double a0
G4fissionEvent * fe
G4double G4Log(G4double x)
Definition: G4Log.hh:226
int G4int
Definition: G4Types.hh:85

Referenced by ComputeMicroscopicCrossSection(), and ComputMuPairLoss().

◆ ComputeMicroscopicCrossSection()

G4double G4MuPairProductionModel::ComputeMicroscopicCrossSection ( G4double  tkin,
G4double  Z,
G4double  cut 
)
protected

Definition at line 272 of file G4MuPairProductionModel.cc.

276{
277 G4double cross = 0.;
279 G4double cut = std::max(cutEnergy, minPairEnergy);
280 if (tmax <= cut) { return cross; }
281
282 // G4double ak1=6.9 ;
283 // G4double ak2=1.0 ;
284 G4double aaa = G4Log(cut);
285 G4double bbb = G4Log(tmax);
286 G4int kkk = (G4int)((bbb-aaa)/ak1 + ak2);
287 if(kkk > 8) { kkk = 8; }
288 else if (kkk < 1) { kkk = 1; }
289
290 G4double hhh = (bbb-aaa)/G4double(kkk);
291 G4double x = aaa;
292
293 for(G4int l=0; l<kkk; ++l)
294 {
295 for(G4int i=0; i<8; ++i)
296 {
297 G4double ep = G4Exp(x + xgi[i]*hhh);
298 cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
299 }
300 x += hhh;
301 }
302
303 cross *= hhh;
304 cross = std::max(cross, 0.0);
305 return cross;
306}
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double pairEnergy)

Referenced by ComputeCrossSectionPerAtom().

◆ ComputMuPairLoss()

G4double G4MuPairProductionModel::ComputMuPairLoss ( G4double  Z,
G4double  tkin,
G4double  cut,
G4double  tmax 
)
protected

Definition at line 235 of file G4MuPairProductionModel.cc.

239{
240 G4double loss = 0.0;
241
242 G4double cut = std::min(cutEnergy,tmax);
243 if(cut <= minPairEnergy) { return loss; }
244
245 // calculate the rectricted loss
246 // numerical integration in log(PairEnergy)
248 G4double bbb = G4Log(cut);
249
250 G4int kkk = (G4int)((bbb-aaa)/ak1+ak2);
251 if (kkk > 8) { kkk = 8; }
252 else if (kkk < 1) { kkk = 1; }
253
254 G4double hhh = (bbb-aaa)/(G4double)kkk;
255 G4double x = aaa;
256
257 for (G4int l=0 ; l<kkk; ++l) {
258 for (G4int ll=0; ll<8; ++ll)
259 {
260 G4double ep = G4Exp(x+xgi[ll]*hhh);
261 loss += wgi[ll]*ep*ep*ComputeDMicroscopicCrossSection(tkin, Z, ep);
262 }
263 x += hhh;
264 }
265 loss *= hhh;
266 loss = std::max(loss, 0.0);
267 return loss;
268}

Referenced by ComputeDEDXPerVolume().

◆ Initialise()

void G4MuPairProductionModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
overridevirtual

Implements G4VEmModel.

Definition at line 163 of file G4MuPairProductionModel.cc.

165{
166 SetParticle(p);
168
169 // for low-energy application this process should not work
170 if(lowestKinEnergy >= HighEnergyLimit()) { return; }
171
172 // define scale of internal table for each thread only once
173 if(0 == nbine) {
174 emin = std::max(lowestKinEnergy, LowEnergyLimit());
175 emax = std::max(HighEnergyLimit(), emin*2);
176 nbine = size_t(nYBinPerDecade*std::log10(emax/emin));
177 if(nbine < 3) { nbine = 3; }
178
180 dy = -ymin/G4double(nbiny);
181 }
182
183 if(IsMaster() && p == particle) {
184 if(!fElementData) {
187 if(dataFile) { dataFile = RetrieveTables(); }
188 if(!dataFile) { MakeSamplingTables(); }
189 if(fTableToFile) { StoreTables(); }
190 }
192 }
193}
bool G4bool
Definition: G4Types.hh:86
static G4EmParameters * Instance()
G4bool RetrieveMuDataFromFile() const
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
G4ElementData * fElementData
Definition: G4VEmModel.hh:438
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118

◆ InitialiseLocal()

void G4MuPairProductionModel::InitialiseLocal ( const G4ParticleDefinition p,
G4VEmModel masterModel 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 197 of file G4MuPairProductionModel.cc.

199{
202 fElementData = masterModel->GetElementData();
203 }
204}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:842
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:834
G4ElementData * GetElementData()
Definition: G4VEmModel.hh:853

◆ MaxSecondaryEnergyForElement()

G4double G4MuPairProductionModel::MaxSecondaryEnergyForElement ( G4double  kineticEnergy,
G4double  Z 
)
inlineprotected

Definition at line 199 of file G4MuPairProductionModel.hh.

201{
202 G4int Z = G4lrint(ZZ);
203 if(Z != currentZ) {
204 currentZ = Z;
205 z13 = nist->GetZ13(Z);
206 z23 = z13*z13;
207 lnZ = nist->GetLOGZ(Z);
208 }
209 return kineticEnergy + particleMass*(1.0 - 0.75*sqrte*z13);
210}
G4double GetZ13(G4double Z) const
G4double GetLOGZ(G4int Z) const
int G4lrint(double ad)
Definition: templates.hh:134

Referenced by ComputeCrossSectionPerAtom(), ComputeDEDXPerVolume(), ComputeMicroscopicCrossSection(), and SampleSecondaries().

◆ MinPrimaryEnergy()

G4double G4MuPairProductionModel::MinPrimaryEnergy ( const G4Material ,
const G4ParticleDefinition ,
G4double  cut 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 154 of file G4MuPairProductionModel.cc.

157{
158 return std::max(lowestKinEnergy,cut);
159}

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 540 of file G4MuPairProductionModel.cc.

546{
547 G4double kinEnergy = aDynamicParticle->GetKineticEnergy();
548 //G4cout << "------- G4MuPairProductionModel::SampleSecondaries E(MeV)= "
549 // << kinEnergy << " "
550 // << aDynamicParticle->GetDefinition()->GetParticleName() << G4endl;
551 G4double totalEnergy = kinEnergy + particleMass;
552 G4double totalMomentum =
553 sqrt(kinEnergy*(kinEnergy + 2.0*particleMass));
554
555 G4ThreeVector partDirection = aDynamicParticle->GetMomentumDirection();
556
557 // select randomly one element constituing the material
558 const G4Element* anElement = SelectRandomAtom(couple,particle,kinEnergy);
559
560 // define interval of energy transfer
561 G4double maxPairEnergy = MaxSecondaryEnergyForElement(kinEnergy,
562 anElement->GetZ());
563 G4double maxEnergy = std::min(tmax, maxPairEnergy);
564 G4double minEnergy = std::max(tmin, minPairEnergy);
565
566 if(minEnergy >= maxEnergy) { return; }
567 //G4cout << "emin= " << minEnergy << " emax= " << maxEnergy
568 // << " minPair= " << minPairEnergy << " maxpair= " << maxPairEnergy
569 // << " ymin= " << ymin << " dy= " << dy << G4endl;
570
571 G4double coeff = G4Log(minPairEnergy/kinEnergy)/ymin;
572
573 // compute limits
574 G4double yymin = G4Log(minEnergy/kinEnergy)/coeff;
575 G4double yymax = G4Log(maxEnergy/kinEnergy)/coeff;
576
577 //G4cout << "yymin= " << yymin << " yymax= " << yymax << G4endl;
578
579 // units should not be used, bacause table was built without
580 G4double logTkin = G4Log(kinEnergy/MeV);
581
582 // sample e-e+ energy, pair energy first
583
584 // select sample table via Z
585 G4int iz1(0), iz2(0);
586 for(G4int iz=0; iz<nzdat; ++iz) {
587 if(currentZ == zdat[iz]) {
588 iz1 = iz2 = currentZ;
589 break;
590 } else if(currentZ < zdat[iz]) {
591 iz2 = zdat[iz];
592 if(iz > 0) { iz1 = zdat[iz-1]; }
593 else { iz1 = iz2; }
594 break;
595 }
596 }
597 if(0 == iz1) { iz1 = iz2 = zdat[nzdat-1]; }
598
599 G4double pairEnergy = 0.0;
600 G4int count = 0;
601 //G4cout << "start loop Z1= " << iz1 << " Z2= " << iz2 << G4endl;
602 do {
603 ++count;
604 // sampling using only one random number
605 G4double rand = G4UniformRand();
606
607 G4double x = FindScaledEnergy(iz1, rand, logTkin, yymin, yymax);
608 if(iz1 != iz2) {
609 G4double x2 = FindScaledEnergy(iz2, rand, logTkin, yymin, yymax);
610 G4double lz1= nist->GetLOGZ(iz1);
611 G4double lz2= nist->GetLOGZ(iz2);
612 //G4cout << count << ". x= " << x << " x2= " << x2
613 // << " Z1= " << iz1 << " Z2= " << iz2 << G4endl;
614 x += (x2 - x)*(lnZ - lz1)/(lz2 - lz1);
615 }
616 //G4cout << "x= " << x << " coeff= " << coeff << G4endl;
617 pairEnergy = kinEnergy*G4Exp(x*coeff);
618
619 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
620 } while((pairEnergy < minEnergy || pairEnergy > maxEnergy) && 10 > count);
621
622 //G4cout << "## pairEnergy(GeV)= " << pairEnergy/GeV
623 // << " Etot(GeV)= " << totalEnergy/GeV << G4endl;
624
625 // sample r=(E+-E-)/pairEnergy ( uniformly .....)
626 G4double rmax =
627 (1.-6.*particleMass*particleMass/(totalEnergy*(totalEnergy-pairEnergy)))
628 *sqrt(1.-minPairEnergy/pairEnergy);
629 G4double r = rmax * (-1.+2.*G4UniformRand()) ;
630
631 // compute energies from pairEnergy,r
632 G4double eEnergy = (1.-r)*pairEnergy*0.5;
633 G4double pEnergy = pairEnergy - eEnergy;
634
635 // Sample angles
636 G4ThreeVector eDirection, pDirection;
637 //
638 GetAngularDistribution()->SamplePairDirections(aDynamicParticle,
639 eEnergy, pEnergy,
640 eDirection, pDirection);
641 // create G4DynamicParticle object for e+e-
642 eEnergy = std::max(eEnergy - CLHEP::electron_mass_c2, 0.0);
643 pEnergy = std::max(pEnergy - CLHEP::electron_mass_c2, 0.0);
644 G4DynamicParticle* aParticle1 =
645 new G4DynamicParticle(theElectron,eDirection,eEnergy);
646 G4DynamicParticle* aParticle2 =
647 new G4DynamicParticle(thePositron,pDirection,pEnergy);
648 // Fill output vector
649 vdp->push_back(aParticle1);
650 vdp->push_back(aParticle2);
651
652 // primary change
653 kinEnergy -= pairEnergy;
654 partDirection *= totalMomentum;
655 partDirection -= (aParticle1->GetMomentum() + aParticle2->GetMomentum());
656 partDirection = partDirection.unit();
657
658 // if energy transfer is higher than threshold (very high by default)
659 // then stop tracking the primary particle and create a new secondary
660 if (pairEnergy > SecondaryThreshold()) {
663 G4DynamicParticle* newdp =
664 new G4DynamicParticle(particle, partDirection, kinEnergy);
665 vdp->push_back(newdp);
666 } else { // continue tracking the primary e-/e+ otherwise
669 }
670 //G4cout << "-- G4MuPairProductionModel::SampleSecondaries done" << G4endl;
671}
@ 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:130
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:611
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:570
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:680
void ProposeTrackStatus(G4TrackStatus status)

◆ SetLowestKineticEnergy()

void G4MuPairProductionModel::SetLowestKineticEnergy ( G4double  e)
inline

◆ SetParticle()

void G4MuPairProductionModel::SetParticle ( const G4ParticleDefinition p)
inline

Definition at line 188 of file G4MuPairProductionModel.hh.

189{
190 if(!particle) {
191 particle = p;
193 }
194}

Referenced by G4MuPairProductionModel(), and Initialise().

Member Data Documentation

◆ currentZ

G4int G4MuPairProductionModel::currentZ
protected

Definition at line 157 of file G4MuPairProductionModel.hh.

Referenced by MaxSecondaryEnergyForElement(), and SampleSecondaries().

◆ dy

G4double G4MuPairProductionModel::dy
protected

Definition at line 171 of file G4MuPairProductionModel.hh.

Referenced by Initialise().

◆ emax

G4double G4MuPairProductionModel::emax
protected

Definition at line 173 of file G4MuPairProductionModel.hh.

Referenced by G4MuPairProductionModel(), and Initialise().

◆ emin

G4double G4MuPairProductionModel::emin
protected

Definition at line 172 of file G4MuPairProductionModel.hh.

Referenced by G4MuPairProductionModel(), and Initialise().

◆ factorForCross

G4double G4MuPairProductionModel::factorForCross
protected

Definition at line 151 of file G4MuPairProductionModel.hh.

Referenced by ComputeDMicroscopicCrossSection().

◆ fParticleChange

G4ParticleChangeForLoss* G4MuPairProductionModel::fParticleChange
protected

Definition at line 161 of file G4MuPairProductionModel.hh.

Referenced by Initialise(), and SampleSecondaries().

◆ fTableToFile

G4bool G4MuPairProductionModel::fTableToFile
protected

Definition at line 175 of file G4MuPairProductionModel.hh.

Referenced by Initialise().

◆ lnZ

G4double G4MuPairProductionModel::lnZ
protected

◆ lowestKinEnergy

◆ minPairEnergy

◆ nbine

size_t G4MuPairProductionModel::nbine
protected

Definition at line 169 of file G4MuPairProductionModel.hh.

Referenced by Initialise().

◆ nbiny

size_t G4MuPairProductionModel::nbiny
protected

Definition at line 168 of file G4MuPairProductionModel.hh.

Referenced by Initialise().

◆ nist

G4NistManager* G4MuPairProductionModel::nist
protected

◆ nYBinPerDecade

G4int G4MuPairProductionModel::nYBinPerDecade
protected

Definition at line 167 of file G4MuPairProductionModel.hh.

Referenced by Initialise().

◆ particle

const G4ParticleDefinition* G4MuPairProductionModel::particle
protected

◆ particleMass

G4double G4MuPairProductionModel::particleMass
protected

◆ sqrte

G4double G4MuPairProductionModel::sqrte
protected

◆ theElectron

G4ParticleDefinition* G4MuPairProductionModel::theElectron
protected

Definition at line 159 of file G4MuPairProductionModel.hh.

Referenced by G4MuPairProductionModel(), and SampleSecondaries().

◆ thePositron

G4ParticleDefinition* G4MuPairProductionModel::thePositron
protected

Definition at line 160 of file G4MuPairProductionModel.hh.

Referenced by G4MuPairProductionModel(), and SampleSecondaries().

◆ ymin

G4double G4MuPairProductionModel::ymin
protected

Definition at line 170 of file G4MuPairProductionModel.hh.

Referenced by Initialise(), and SampleSecondaries().

◆ z13

G4double G4MuPairProductionModel::z13
protected

◆ z23

G4double G4MuPairProductionModel::z23
protected

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