Geant4 11.1.1
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")
 
 ~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

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

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 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 factorForCross(CLHEP::fine_structure_const*CLHEP::fine_structure_const*
116 CLHEP::classic_electr_radius*CLHEP::classic_electr_radius*
117 4./(3.*CLHEP::pi)),
118 sqrte(sqrt(G4Exp(1.))),
119 minPairEnergy(4.*CLHEP::electron_mass_c2),
120 lowestKinEnergy(0.85*CLHEP::GeV)
121{
123
124 theElectron = G4Electron::Electron();
125 thePositron = G4Positron::Positron();
126
127 if(nullptr != p) {
128 SetParticle(p);
129 lowestKinEnergy = std::max(lowestKinEnergy, p->GetPDGMass()*8.0);
130 }
132 emax = emin*10000.;
134}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
static G4Electron * Electron()
Definition: G4Electron.cc:93
void SetParticle(const G4ParticleDefinition *)
static G4NistManager * Instance()
static G4Positron * Positron()
Definition: G4Positron.cc:93
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607

◆ ~G4MuPairProductionModel()

G4MuPairProductionModel::~G4MuPairProductionModel ( )
default

◆ 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 427 of file G4MuPairProductionModel.cc.

433{
434 G4double cross = 0.0;
435 if (kineticEnergy <= lowestKinEnergy) { return cross; }
436
437 G4double maxPairEnergy = MaxSecondaryEnergyForElement(kineticEnergy, Z);
438 G4double tmax = std::min(maxEnergy, maxPairEnergy);
439 G4double cut = std::max(cutEnergy, minPairEnergy);
440 if (cut >= tmax) { return cross; }
441
442 cross = ComputeMicroscopicCrossSection(kineticEnergy, Z, cut);
443 if(tmax < kineticEnergy) {
444 cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
445 }
446 return cross;
447}
double G4double
Definition: G4Types.hh:83
const G4int Z[17]
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 195 of file G4MuPairProductionModel.cc.

200{
201 G4double dedx = 0.0;
202 if (cutEnergy <= minPairEnergy || kineticEnergy <= lowestKinEnergy)
203 { return dedx; }
204
205 const G4ElementVector* theElementVector = material->GetElementVector();
206 const G4double* theAtomicNumDensityVector =
207 material->GetAtomicNumDensityVector();
208
209 // loop for elements in the material
210 for (size_t i=0; i<material->GetNumberOfElements(); ++i) {
211 G4double Z = (*theElementVector)[i]->GetZ();
212 G4double tmax = MaxSecondaryEnergyForElement(kineticEnergy, Z);
213 G4double loss = ComputMuPairLoss(Z, kineticEnergy, cutEnergy, tmax);
214 dedx += loss*theAtomicNumDensityVector[i];
215 }
216 dedx = std::max(dedx, 0.0);
217 return dedx;
218}
std::vector< const G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:211
G4double ComputMuPairLoss(G4double Z, G4double tkin, G4double cut, G4double tmax)

◆ ComputeDMicroscopicCrossSection()

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

Reimplemented in G4MuonToMuonPairProductionModel.

Definition at line 291 of file G4MuPairProductionModel.cc.

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

Referenced by ComputeMicroscopicCrossSection(), ComputMuPairLoss(), and MakeSamplingTables().

◆ ComputeMicroscopicCrossSection()

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

Definition at line 257 of file G4MuPairProductionModel.cc.

261{
262 G4double cross = 0.;
264 G4double cut = std::max(cutEnergy, minPairEnergy);
265 if (tmax <= cut) { return cross; }
266
267 G4double aaa = G4Log(cut);
268 G4double bbb = G4Log(tmax);
269 G4int kkk = G4lrint((bbb-aaa)/ak1 + ak2);
270 if(kkk > 8) { kkk = 8; }
271 else if (kkk < 1) { kkk = 1; }
272
273 G4double hhh = (bbb-aaa)/(kkk);
274 G4double x = aaa;
275
276 for(G4int l=0; l<kkk; ++l) {
277 for(G4int i=0; i<8; ++i) {
278 G4double ep = G4Exp(x + xgi[i]*hhh);
279 cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
280 }
281 x += hhh;
282 }
283
284 cross *= hhh;
285 cross = std::max(cross, 0.0);
286 return cross;
287}
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double pairEnergy)
int G4lrint(double ad)
Definition: templates.hh:134

Referenced by ComputeCrossSectionPerAtom().

◆ ComputMuPairLoss()

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

Definition at line 222 of file G4MuPairProductionModel.cc.

226{
227 G4double loss = 0.0;
228
229 G4double cut = std::min(cutEnergy, tmax);
230 if(cut <= minPairEnergy) { return loss; }
231
232 // calculate the rectricted loss
233 // numerical integration in log(PairEnergy)
235 G4double bbb = G4Log(cut);
236
237 G4int kkk = G4lrint((bbb-aaa)/ak1+ak2);
238 if(kkk > 8) { kkk = 8; }
239 else if (kkk < 1) { kkk = 1; }
240 G4double hhh = (bbb-aaa)/kkk;
241 G4double x = aaa;
242
243 for (G4int l=0 ; l<kkk; ++l) {
244 for (G4int ll=0; ll<8; ++ll) {
245 G4double ep = G4Exp(x+xgi[ll]*hhh);
246 loss += wgi[ll]*ep*ep*ComputeDMicroscopicCrossSection(tkin, Z, ep);
247 }
248 x += hhh;
249 }
250 loss *= hhh;
251 loss = std::max(loss, 0.0);
252 return loss;
253}

Referenced by ComputeDEDXPerVolume().

◆ DataCorrupted()

void G4MuPairProductionModel::DataCorrupted ( G4int  Z,
G4double  logTkin 
) const
protectedvirtual

Reimplemented in G4MuonToMuonPairProductionModel.

Definition at line 676 of file G4MuPairProductionModel.cc.

677{
679 ed << "G4ElementData is not properly initialized Z= " << Z
680 << " Ekin(MeV)= " << G4Exp(logTkin)
681 << " IsMasterThread= " << IsMaster()
682 << " Model " << GetName();
683 G4Exception("G4MuPairProductionModel::()", "em0033", FatalException, ed, "");
684}
@ 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

Referenced by FindScaledEnergy(), and StoreTables().

◆ FindScaledEnergy()

G4double G4MuPairProductionModel::FindScaledEnergy ( G4int  Z,
G4double  rand,
G4double  logTkin,
G4double  yymin,
G4double  yymax 
)
protected

Definition at line 656 of file G4MuPairProductionModel.cc.

659{
660 G4double res = yymin;
662 if(nullptr != pv) {
663 G4double pmin = pv->Value(yymin, logTkin);
664 G4double pmax = pv->Value(yymax, logTkin);
665 G4double p0 = pv->Value(0.0, logTkin);
666 if(p0 <= 0.0) { DataCorrupted(Z, logTkin); }
667 else { res = pv->FindLinearX((pmin + rand*(pmax - pmin))/p0, logTkin); }
668 } else {
669 DataCorrupted(Z, logTkin);
670 }
671 return res;
672}
G4Physics2DVector * GetElement2DData(G4int Z)
virtual void DataCorrupted(G4int Z, G4double logTkin) const
G4double Value(G4double x, G4double y, std::size_t &lastidx, std::size_t &lastidy) const
G4double FindLinearX(G4double rand, G4double y, std::size_t &lastidy) const
G4ElementData * fElementData
Definition: G4VEmModel.hh:420

Referenced by G4MuonToMuonPairProductionModel::SampleSecondaries(), and SampleSecondaries().

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 147 of file G4MuPairProductionModel.cc.

149{
150 SetParticle(p);
151
152 if(nullptr == fParticleChange) {
154 }
155
156 // for low-energy application this process should not work
157 if(lowestKinEnergy >= HighEnergyLimit()) { return; }
158
159 // define scale of internal table for each thread only once
160 if(0 == nbine) {
161 emin = std::max(lowestKinEnergy, LowEnergyLimit());
162 emax = std::max(HighEnergyLimit(), emin*2);
163 nbine = size_t(nYBinPerDecade*std::log10(emax/emin));
164 if(nbine < 3) { nbine = 3; }
165
167 dy = -ymin/G4double(nbiny);
168 }
169
170 if(IsMaster() && p == particle) {
171 if(nullptr == fElementData) {
174 if(dataFile) { dataFile = RetrieveTables(); }
175 if(!dataFile) { MakeSamplingTables(); }
176 if(fTableToFile) { StoreTables(); }
177 }
179 }
180}
bool G4bool
Definition: G4Types.hh:86
static G4EmParameters * Instance()
G4bool RetrieveMuDataFromFile() const
const G4ParticleDefinition * particle
G4ParticleChangeForLoss * fParticleChange
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:139
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:109

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 184 of file G4MuPairProductionModel.cc.

186{
189 fElementData = masterModel->GetElementData();
190 }
191}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:831
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:823
G4ElementData * GetElementData()
Definition: G4VEmModel.hh:842

◆ MakeSamplingTables()

void G4MuPairProductionModel::MakeSamplingTables ( )
protected

Definition at line 451 of file G4MuPairProductionModel.cc.

452{
454
455 for (G4int iz=0; iz<nzdat; ++iz) {
456
457 G4double Z = zdat[iz];
459 G4double kinEnergy = emin;
460
461 for (size_t it=0; it<=nbine; ++it) {
462
463 pv->PutY(it, G4Log(kinEnergy/CLHEP::MeV));
464 G4double maxPairEnergy = MaxSecondaryEnergyForElement(kinEnergy, Z);
465 /*
466 G4cout << "it= " << it << " E= " << kinEnergy
467 << " " << particle->GetParticleName()
468 << " maxE= " << maxPairEnergy << " minE= " << minPairEnergy
469 << " ymin= " << ymin << G4endl;
470 */
471 G4double coef = G4Log(minPairEnergy/kinEnergy)/ymin;
472 G4double ymax = G4Log(maxPairEnergy/kinEnergy)/coef;
473 G4double fac = (ymax - ymin)/dy;
474 size_t imax = (size_t)fac;
475 fac -= (G4double)imax;
476
477 G4double xSec = 0.0;
478 G4double x = ymin;
479 /*
480 G4cout << "Z= " << currentZ << " z13= " << z13
481 << " mE= " << maxPairEnergy << " ymin= " << ymin
482 << " dy= " << dy << " c= " << coef << G4endl;
483 */
484 // start from zero
485 pv->PutValue(0, it, 0.0);
486 if(0 == it) { pv->PutX(nbiny, 0.0); }
487
488 for (size_t i=0; i<nbiny; ++i) {
489
490 if(0 == it) { pv->PutX(i, x); }
491
492 if(i < imax) {
493 G4double ep = kinEnergy*G4Exp(coef*(x + dy*0.5));
494
495 // not multiplied by interval, because table
496 // will be used only for sampling
497 //G4cout << "i= " << i << " x= " << x << "E= " << kinEnergy
498 // << " Egamma= " << ep << G4endl;
499 xSec += ep*ComputeDMicroscopicCrossSection(kinEnergy, Z, ep);
500
501 // last bin before the kinematic limit
502 } else if(i == imax) {
503 G4double ep = kinEnergy*G4Exp(coef*(x + fac*dy*0.5));
504 xSec += ep*fac*ComputeDMicroscopicCrossSection(kinEnergy, Z, ep);
505 }
506 pv->PutValue(i + 1, it, xSec);
507 x += dy;
508 }
509 kinEnergy *= factore;
510
511 // to avoid precision lost
512 if(it+1 == nbine) { kinEnergy = emax; }
513 }
514 fElementData->InitialiseForElement(zdat[iz], pv);
515 }
516}
void InitialiseForElement(G4int Z, G4PhysicsVector *v)
void PutY(std::size_t idy, G4double value)
void PutValue(std::size_t idx, std::size_t idy, G4double value)
void PutX(std::size_t idx, G4double value)

Referenced by Initialise().

◆ MaxSecondaryEnergyForElement()

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

Definition at line 195 of file G4MuPairProductionModel.hh.

197{
198 G4int Z = G4lrint(ZZ);
199 if(Z != currentZ) {
200 currentZ = Z;
201 z13 = nist->GetZ13(Z);
202 z23 = z13*z13;
203 lnZ = nist->GetLOGZ(Z);
204 }
205 return kineticEnergy + particleMass*(1.0 - 0.75*sqrte*z13);
206}
G4double GetZ13(G4double Z) const
G4double GetLOGZ(G4int Z) const

Referenced by ComputeCrossSectionPerAtom(), ComputeDEDXPerVolume(), ComputeMicroscopicCrossSection(), MakeSamplingTables(), G4MuonToMuonPairProductionModel::SampleSecondaries(), and SampleSecondaries().

◆ MinPrimaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 138 of file G4MuPairProductionModel.cc.

141{
142 return std::max(lowestKinEnergy, cut);
143}

◆ operator=()

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

◆ RetrieveTables()

G4bool G4MuPairProductionModel::RetrieveTables ( )
protected

Definition at line 706 of file G4MuPairProductionModel.cc.

707{
708 const char* path = G4FindDataDir("G4LEDATA");
709 G4String dir("");
710 if (path) {
711 std::ostringstream ost;
712 ost << path << "/mupair/";
713 dir = ost.str();
714 } else {
715 dir = "./mupair/";
716 }
717
718 for (G4int iz=0; iz<nzdat; ++iz) {
719 G4double Z = zdat[iz];
721 std::ostringstream ss;
722 ss << dir << particle->GetParticleName() << Z << ".dat";
723 std::ifstream infile(ss.str(), std::ios::in);
724 if(!pv->Retrieve(infile)) {
725 delete pv;
726 return false;
727 }
729 }
730 return true;
731}
const char * G4FindDataDir(const char *)
const G4String & GetParticleName() const
G4bool Retrieve(std::ifstream &fIn)

Referenced by Initialise().

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 520 of file G4MuPairProductionModel.cc.

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

◆ SetLowestKineticEnergy()

void G4MuPairProductionModel::SetLowestKineticEnergy ( G4double  e)
inline

Definition at line 176 of file G4MuPairProductionModel.hh.

177{
178 lowestKinEnergy = e;
179}

Referenced by G4ePairProduction::InitialiseEnergyLossProcess().

◆ SetParticle()

void G4MuPairProductionModel::SetParticle ( const G4ParticleDefinition p)
inline

Definition at line 184 of file G4MuPairProductionModel.hh.

185{
186 if(nullptr == particle) {
187 particle = p;
189 }
190}

Referenced by G4MuPairProductionModel(), and Initialise().

◆ StoreTables()

void G4MuPairProductionModel::StoreTables ( ) const
protected

Definition at line 688 of file G4MuPairProductionModel.cc.

689{
690 for (G4int iz=0; iz<nzdat; ++iz) {
691 G4int Z = zdat[iz];
693 if(nullptr == pv) {
694 DataCorrupted(Z, 1.0);
695 return;
696 }
697 std::ostringstream ss;
698 ss << "mupair/" << particle->GetParticleName() << Z << ".dat";
699 std::ofstream outfile(ss.str());
700 pv->Store(outfile);
701 }
702}
void Store(std::ofstream &fOut) const

Referenced by Initialise().

Member Data Documentation

◆ currentZ

G4int G4MuPairProductionModel::currentZ = 0
protected

◆ dy

G4double G4MuPairProductionModel::dy = 0.005
protected

Definition at line 159 of file G4MuPairProductionModel.hh.

Referenced by Initialise(), and MakeSamplingTables().

◆ emax

G4double G4MuPairProductionModel::emax
protected

◆ emin

G4double G4MuPairProductionModel::emin
protected

◆ factorForCross

G4double G4MuPairProductionModel::factorForCross
protected

Definition at line 146 of file G4MuPairProductionModel.hh.

Referenced by ComputeDMicroscopicCrossSection().

◆ fParticleChange

G4ParticleChangeForLoss* G4MuPairProductionModel::fParticleChange = nullptr
protected

◆ fTableToFile

G4bool G4MuPairProductionModel::fTableToFile = false
protected

Definition at line 166 of file G4MuPairProductionModel.hh.

Referenced by Initialise().

◆ lnZ

G4double G4MuPairProductionModel::lnZ = 0.0
protected

◆ lowestKinEnergy

◆ minPairEnergy

◆ nbine

size_t G4MuPairProductionModel::nbine = 0
protected

Definition at line 164 of file G4MuPairProductionModel.hh.

Referenced by Initialise(), MakeSamplingTables(), and RetrieveTables().

◆ nbiny

size_t G4MuPairProductionModel::nbiny = 1000
protected

Definition at line 163 of file G4MuPairProductionModel.hh.

Referenced by Initialise(), MakeSamplingTables(), and RetrieveTables().

◆ nist

◆ nYBinPerDecade

G4int G4MuPairProductionModel::nYBinPerDecade = 4
protected

Definition at line 162 of file G4MuPairProductionModel.hh.

Referenced by Initialise().

◆ particle

◆ particleMass

◆ sqrte

G4double G4MuPairProductionModel::sqrte
protected

◆ ymin

G4double G4MuPairProductionModel::ymin = -5.0
protected

◆ z13

G4double G4MuPairProductionModel::z13 = 0.0
protected

◆ z23

G4double G4MuPairProductionModel::z23 = 0.0
protected

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