Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4MuPairProductionModel Class Reference

#include <G4MuPairProductionModel.hh>

Inheritance diagram for G4MuPairProductionModel:
Collaboration 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 *)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
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 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 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)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, 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 *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
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)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () 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)
 

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
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Static Protected Attributes

static const G4double xgi [8]
 
static const G4double wgi [8]
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 73 of file G4MuPairProductionModel.hh.

Constructor & Destructor Documentation

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

Definition at line 106 of file G4MuPairProductionModel.cc.

108  : G4VEmModel(nam),
109  particle(nullptr),
112  sqrte(sqrt(G4Exp(1.))),
113  currentZ(0),
114  fParticleChange(nullptr),
115  minPairEnergy(4.*electron_mass_c2),
116  lowestKinEnergy(1.0*GeV),
117  nzdat(5),
118  nYBinPerDecade(4),
119  nbiny(1000),
120  nbine(0),
121  ymin(-5.),
122  dy(0.005)
123 {
125 
126  theElectron = G4Electron::Electron();
127  thePositron = G4Positron::Positron();
128 
129  particleMass = lnZ = z13 = z23 = 0.;
130 
131  // setup lowest limit dependent on particle mass
132  if(p) {
133  SetParticle(p);
134  lowestKinEnergy = std::max(lowestKinEnergy,p->GetPDGMass()*8.0);
135  }
136  emin = lowestKinEnergy;
137  emax = 10.*TeV;
138 }
const G4ParticleDefinition * particle
void SetParticle(const G4ParticleDefinition *)
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
static constexpr double electron_mass_c2
static constexpr double TeV
Definition: G4SIunits.hh:218
static constexpr double classic_electr_radius
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4double GetPDGMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static constexpr double GeV
Definition: G4SIunits.hh:217
static G4Electron * Electron()
Definition: G4Electron.cc:94
static constexpr double pi
Definition: G4SIunits.hh:75
static constexpr double fine_structure_const

Here is the call graph for this function:

G4MuPairProductionModel::~G4MuPairProductionModel ( )
virtual

Definition at line 142 of file G4MuPairProductionModel.cc.

143 {}

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 416 of file G4MuPairProductionModel.cc.

422 {
423  G4double cross = 0.0;
424  if (kineticEnergy <= lowestKinEnergy) { return cross; }
425 
426  G4double maxPairEnergy = MaxSecondaryEnergyForElement(kineticEnergy, Z);
427  G4double tmax = std::min(maxEnergy, maxPairEnergy);
428  G4double cut = std::max(cutEnergy, minPairEnergy);
429  if (cut >= tmax) { return cross; }
430 
431  cross = ComputeMicroscopicCrossSection(kineticEnergy, Z, cut);
432  if(tmax < kineticEnergy) {
433  cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
434  }
435  return cross;
436 }
G4double MaxSecondaryEnergyForElement(G4double kineticEnergy, G4double Z)
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 199 of file G4MuPairProductionModel.cc.

204 {
205  G4double dedx = 0.0;
206  if (cutEnergy <= minPairEnergy || kineticEnergy <= lowestKinEnergy)
207  { return dedx; }
208 
209  const G4ElementVector* theElementVector = material->GetElementVector();
210  const G4double* theAtomicNumDensityVector =
211  material->GetAtomicNumDensityVector();
212 
213  // loop for elements in the material
214  for (size_t i=0; i<material->GetNumberOfElements(); ++i) {
215  G4double Z = (*theElementVector)[i]->GetZ();
216  G4double tmax = MaxSecondaryEnergyForElement(kineticEnergy, Z);
217  G4double loss = ComputMuPairLoss(Z, kineticEnergy, cutEnergy, tmax);
218  dedx += loss*theAtomicNumDensityVector[i];
219  }
220  dedx = std::max(dedx, 0.0);
221  return dedx;
222 }
G4double MaxSecondaryEnergyForElement(G4double kineticEnergy, G4double Z)
std::vector< G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
G4double ComputMuPairLoss(G4double Z, G4double tkin, G4double cut, G4double tmax)
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:216
T max(const T t1, const T t2)
brief Return the largest of the two arguments
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Reimplemented in G4hPairProductionModel.

Definition at line 303 of file G4MuPairProductionModel.cc.

310 {
311  static const G4double bbbtf= 183. ;
312  static const G4double bbbh = 202.4 ;
313  static const G4double g1tf = 1.95e-5 ;
314  static const G4double g2tf = 5.3e-5 ;
315  static const G4double g1h = 4.4e-5 ;
316  static const G4double g2h = 4.8e-5 ;
317 
318  G4double totalEnergy = tkin + particleMass;
319  G4double residEnergy = totalEnergy - pairEnergy;
320  G4double massratio = particleMass/electron_mass_c2;
321  G4double massratio2 = massratio*massratio;
322  G4double cross = 0.;
323 
324  G4double c3 = 0.75*sqrte*particleMass;
325  if (residEnergy <= c3*z13) { return cross; }
326 
328  G4double c8 = 6.*particleMass*particleMass;
329  G4double alf = c7/pairEnergy;
330  G4double a3 = 1. - alf;
331  if (a3 <= 0.) { return cross; }
332 
333  // zeta calculation
334  G4double bbb,g1,g2;
335  if( Z < 1.5 ) { bbb = bbbh ; g1 = g1h ; g2 = g2h ; }
336  else { bbb = bbbtf; g1 = g1tf; g2 = g2tf; }
337 
338  G4double zeta = 0;
339  G4double zeta1 =
340  0.073*G4Log(totalEnergy/(particleMass+g1*z23*totalEnergy))-0.26;
341  if ( zeta1 > 0.)
342  {
343  G4double zeta2 =
344  0.058*G4Log(totalEnergy/(particleMass+g2*z13*totalEnergy))-0.14;
345  zeta = zeta1/zeta2 ;
346  }
347 
348  G4double z2 = Z*(Z+zeta);
349  G4double screen0 = 2.*electron_mass_c2*sqrte*bbb/(z13*pairEnergy);
350  G4double a0 = totalEnergy*residEnergy;
351  G4double a1 = pairEnergy*pairEnergy/a0;
352  G4double bet = 0.5*a1;
353  G4double xi0 = 0.25*massratio2*a1;
354  G4double del = c8/a0;
355 
356  G4double rta3 = sqrt(a3);
357  G4double tmnexp = alf/(1. + rta3) + del*rta3;
358  if(tmnexp >= 1.0) { return cross; }
359 
360  G4double tmn = G4Log(tmnexp);
361  G4double sum = 0.;
362 
363  // Gaussian integration in ln(1-ro) ( with 8 points)
364  for (G4int i=0; i<8; ++i)
365  {
366  G4double a4 = G4Exp(tmn*xgi[i]); // a4 = (1.-asymmetry)
367  G4double a5 = a4*(2.-a4) ;
368  G4double a6 = 1.-a5 ;
369  G4double a7 = 1.+a6 ;
370  G4double a9 = 3.+a6 ;
371  G4double xi = xi0*a5 ;
372  G4double xii = 1./xi ;
373  G4double xi1 = 1.+xi ;
374  G4double screen = screen0*xi1/a5 ;
375  G4double yeu = 5.-a6+4.*bet*a7 ;
376  G4double yed = 2.*(1.+3.*bet)*G4Log(3.+xii)-a6-a1*(2.-a6) ;
377  G4double ye1 = 1.+yeu/yed ;
378  G4double ale = G4Log(bbb/z13*sqrt(xi1*ye1)/(1.+screen*ye1)) ;
379  G4double cre = 0.5*G4Log(1.+2.25*z23*xi1*ye1/massratio2) ;
380  G4double be;
381 
382  if (xi <= 1.e3) {
383  be = ((2.+a6)*(1.+bet)+xi*a9)*G4Log(1.+xii)+(a5-bet)/xi1-a9;
384  } else {
385  be = (3.-a6+a1*a7)/(2.*xi);
386  }
387  G4double fe = (ale-cre)*be;
388  if ( fe < 0.) fe = 0. ;
389 
390  G4double ymu = 4.+a6 +3.*bet*a7 ;
391  G4double ymd = a7*(1.5+a1)*G4Log(3.+xi)+1.-1.5*a6 ;
392  G4double ym1 = 1.+ymu/ymd ;
393  G4double alm_crm = G4Log(bbb*massratio/(1.5*z23*(1.+screen*ym1)));
394  G4double a10,bm;
395  if ( xi >= 1.e-3)
396  {
397  a10 = (1.+a1)*a5 ;
398  bm = (a7*(1.+1.5*bet)-a10*xii)*G4Log(xi1)+xi*(a5-bet)/xi1+a10;
399  } else {
400  bm = (5.-a6+bet*a9)*(xi/2.);
401  }
402 
403  G4double fm = alm_crm*bm;
404  if ( fm < 0.) { fm = 0.; }
405 
406  sum += wgi[i]*a4*(fe+fm/massratio2);
407  }
408 
409  cross = -tmn*sum*factorForCross*z2*residEnergy/(totalEnergy*pairEnergy);
410  cross = std::max(cross, 0.0);
411  return cross;
412 }
const G4double a0
static const G4double xgi[8]
int G4int
Definition: G4Types.hh:78
static constexpr double electron_mass_c2
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static const G4double wgi[8]
double G4double
Definition: G4Types.hh:76
G4fissionEvent * fe

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 265 of file G4MuPairProductionModel.cc.

269 {
270  G4double cross = 0.;
272  G4double cut = std::max(cutEnergy, minPairEnergy);
273  if (tmax <= cut) { return cross; }
274 
275  // G4double ak1=6.9 ;
276  // G4double ak2=1.0 ;
277  G4double aaa = G4Log(cut);
278  G4double bbb = G4Log(tmax);
279  G4int kkk = (G4int)((bbb-aaa)/ak1 + ak2);
280  if(kkk > 8) { kkk = 8; }
281  else if (kkk < 1) { kkk = 1; }
282 
283  G4double hhh = (bbb-aaa)/G4double(kkk);
284  G4double x = aaa;
285 
286  for(G4int l=0; l<kkk; ++l)
287  {
288  for(G4int i=0; i<8; ++i)
289  {
290  G4double ep = G4Exp(x + xgi[i]*hhh);
291  cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
292  }
293  x += hhh;
294  }
295 
296  cross *= hhh;
297  cross = std::max(cross, 0.0);
298  return cross;
299 }
G4double MaxSecondaryEnergyForElement(G4double kineticEnergy, G4double Z)
static const G4double xgi[8]
static const G4double ak2
int G4int
Definition: G4Types.hh:78
static const G4double ak1
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
T max(const T t1, const T t2)
brief Return the largest of the two arguments
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double pairEnergy)
static const G4double wgi[8]
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 226 of file G4MuPairProductionModel.cc.

230 {
231  G4double loss = 0.0;
232 
233  G4double cut = std::min(cutEnergy,tmax);
234  if(cut <= minPairEnergy) { return loss; }
235 
236  // calculate the rectricted loss
237  // numerical integration in log(PairEnergy)
238  G4double aaa = G4Log(minPairEnergy);
239  G4double bbb = G4Log(cut);
240 
241  G4int kkk = (G4int)((bbb-aaa)/ak1+ak2);
242  if (kkk > 8) { kkk = 8; }
243  else if (kkk < 1) { kkk = 1; }
244 
245  G4double hhh = (bbb-aaa)/(G4double)kkk;
246  G4double x = aaa;
247 
248  for (G4int l=0 ; l<kkk; l++)
249  {
250 
251  for (G4int ll=0; ll<8; ll++)
252  {
253  G4double ep = G4Exp(x+xgi[ll]*hhh);
254  loss += wgi[ll]*ep*ep*ComputeDMicroscopicCrossSection(tkin, Z, ep);
255  }
256  x += hhh;
257  }
258  loss *= hhh;
259  loss = std::max(loss, 0.0);
260  return loss;
261 }
static const G4double xgi[8]
static const G4double ak2
int G4int
Definition: G4Types.hh:78
static const G4double ak1
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
T max(const T t1, const T t2)
brief Return the largest of the two arguments
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double pairEnergy)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static const G4double wgi[8]
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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

Implements G4VEmModel.

Definition at line 156 of file G4MuPairProductionModel.cc.

158 {
159  SetParticle(p);
160  if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
161 
162  // for low-energy application this process should not work
163  if(lowestKinEnergy >= HighEnergyLimit()) { return; }
164 
165  // define scale of internal table for each thread only once
166  if(0 == nbine) {
167  emin = std::max(lowestKinEnergy, LowEnergyLimit());
168  emax = std::max(HighEnergyLimit(), emin*2);
169  nbine = size_t(nYBinPerDecade*std::log10(emax/emin));
170  if(nbine < 3) { nbine = 3; }
171 
172  ymin = G4Log(minPairEnergy/emin);
173  dy = -ymin/G4double(nbiny);
174  }
175 
176  if(IsMaster() && p == particle) {
177 
178  if(!fElementData) {
179  fElementData = new G4ElementData();
180  MakeSamplingTables();
181  }
182  InitialiseElementSelectors(p, cuts);
183  }
184 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
const G4ParticleDefinition * particle
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
void SetParticle(const G4ParticleDefinition *)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
T max(const T t1, const T t2)
brief Return the largest of the two arguments
double G4double
Definition: G4Types.hh:76
G4ElementData * fElementData
Definition: G4VEmModel.hh:421

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 188 of file G4MuPairProductionModel.cc.

190 {
191  if(p == particle && lowestKinEnergy < HighEnergyLimit()) {
193  fElementData = masterModel->GetElementData();
194  }
195 }
const G4ParticleDefinition * particle
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
G4ElementData * GetElementData()
Definition: G4VEmModel.hh:818
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:801
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:809
G4ElementData * fElementData
Definition: G4VEmModel.hh:421

Here is the call graph for this function:

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

Definition at line 200 of file G4MuPairProductionModel.hh.

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

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 147 of file G4MuPairProductionModel.cc.

150 {
151  return std::max(lowestKinEnergy,cut);
152 }
T max(const T t1, const T t2)
brief Return the largest of the two arguments

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 509 of file G4MuPairProductionModel.cc.

515 {
516  G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
517  //G4cout << "------- G4MuPairProductionModel::SampleSecondaries E(MeV)= "
518  // << kineticEnergy << " "
519  // << aDynamicParticle->GetDefinition()->GetParticleName() << G4endl;
520  G4double totalEnergy = kineticEnergy + particleMass;
521  G4double totalMomentum =
522  sqrt(kineticEnergy*(kineticEnergy + 2.0*particleMass));
523 
524  G4ThreeVector partDirection = aDynamicParticle->GetMomentumDirection();
525 
526  // select randomly one element constituing the material
527  const G4Element* anElement = SelectRandomAtom(couple,particle,kineticEnergy);
528 
529  // define interval of energy transfer
530  G4double maxPairEnergy = MaxSecondaryEnergyForElement(kineticEnergy,
531  anElement->GetZ());
532  G4double maxEnergy = std::min(tmax, maxPairEnergy);
533  G4double minEnergy = std::max(tmin, minPairEnergy);
534 
535  if(minEnergy >= maxEnergy) { return; }
536  //G4cout << "emin= " << minEnergy << " emax= " << maxEnergy
537  // << " minPair= " << minPairEnergy << " maxpair= " << maxPairEnergy
538  // << " ymin= " << ymin << " dy= " << dy << G4endl;
539 
540  G4double coeff = G4Log(minPairEnergy/kineticEnergy)/ymin;
541 
542  // compute limits
543  G4double yymin = G4Log(minEnergy/kineticEnergy)/coeff;
544  G4double yymax = G4Log(maxEnergy/kineticEnergy)/coeff;
545 
546  //G4cout << "yymin= " << yymin << " yymax= " << yymax << G4endl;
547 
548  // units should not be used, bacause table was built without
549  G4double logTkin = G4Log(kineticEnergy/MeV);
550 
551  // sample e-e+ energy, pair energy first
552 
553  // select sample table via Z
554  G4int iz1(0), iz2(0);
555  for(G4int iz=0; iz<nzdat; ++iz) {
556  if(currentZ == zdat[iz]) {
557  iz1 = iz2 = currentZ;
558  break;
559  } else if(currentZ < zdat[iz]) {
560  iz2 = zdat[iz];
561  if(iz > 0) { iz1 = zdat[iz-1]; }
562  else { iz1 = iz2; }
563  break;
564  }
565  }
566  if(0 == iz1) { iz1 = iz2 = zdat[nzdat-1]; }
567 
568  G4double PairEnergy = 0.0;
569  G4int count = 0;
570  //G4cout << "start loop Z1= " << iz1 << " Z2= " << iz2 << G4endl;
571  do {
572  ++count;
573  // sampling using only one random number
574  G4double rand = G4UniformRand();
575 
576  G4double x = FindScaledEnergy(iz1, rand, logTkin, yymin, yymax);
577  if(iz1 != iz2) {
578  G4double x2 = FindScaledEnergy(iz2, rand, logTkin, yymin, yymax);
579  G4double lz1= nist->GetLOGZ(iz1);
580  G4double lz2= nist->GetLOGZ(iz2);
581  //G4cout << count << ". x= " << x << " x2= " << x2
582  // << " Z1= " << iz1 << " Z2= " << iz2 << G4endl;
583  x += (x2 - x)*(lnZ - lz1)/(lz2 - lz1);
584  }
585  //G4cout << "x= " << x << " coeff= " << coeff << G4endl;
586  PairEnergy = kineticEnergy*G4Exp(x*coeff);
587 
588  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
589  } while((PairEnergy < minEnergy || PairEnergy > maxEnergy) && 10 > count);
590 
591  //G4cout << "## PairEnergy(GeV)= " << PairEnergy/GeV
592  // << " Etot(GeV)= " << totalEnergy/GeV << G4endl;
593 
594  // sample r=(E+-E-)/PairEnergy ( uniformly .....)
595  G4double rmax =
596  (1.-6.*particleMass*particleMass/(totalEnergy*(totalEnergy-PairEnergy)))
597  *sqrt(1.-minPairEnergy/PairEnergy);
598  G4double r = rmax * (-1.+2.*G4UniformRand()) ;
599 
600  // compute energies from PairEnergy,r
601  G4double ElectronEnergy = (1.-r)*PairEnergy*0.5;
602  G4double PositronEnergy = PairEnergy - ElectronEnergy;
603 
604  // The angle of the emitted virtual photon is sampled
605  // according to the muon bremsstrahlung model
606 
607  G4double gam = totalEnergy/particleMass;
608  G4double gmax = gam*std::min(1.0, totalEnergy/PairEnergy - 1.0);
609  G4double gmax2= gmax*gmax;
610  G4double x = G4UniformRand()*gmax2/(1.0 + gmax2);
611 
612  G4double theta = sqrt(x/(1.0 - x))/gam;
613  G4double sint = sin(theta);
614  G4double phi = twopi * G4UniformRand() ;
615  G4double dirx = sint*cos(phi), diry = sint*sin(phi), dirz = cos(theta) ;
616 
617  G4ThreeVector gDirection(dirx, diry, dirz);
618  gDirection.rotateUz(partDirection);
619 
620  // the angles of e- and e+ assumed to be the same as virtual gamma
621 
622  // create G4DynamicParticle object for the particle1
623  G4DynamicParticle* aParticle1 =
624  new G4DynamicParticle(theElectron, gDirection,
625  ElectronEnergy - electron_mass_c2);
626 
627  // create G4DynamicParticle object for the particle2
628  G4DynamicParticle* aParticle2 =
629  new G4DynamicParticle(thePositron, gDirection,
630  PositronEnergy - electron_mass_c2);
631 
632  // primary change
633  kineticEnergy -= (ElectronEnergy + PositronEnergy);
634  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
635 
636  partDirection *= totalMomentum;
637  partDirection -= (aParticle1->GetMomentum() + aParticle2->GetMomentum());
638  partDirection = partDirection.unit();
639  fParticleChange->SetProposedMomentumDirection(partDirection);
640 
641  // add secondary
642  vdp->push_back(aParticle1);
643  vdp->push_back(aParticle2);
644  //G4cout << "-- G4MuPairProductionModel::SampleSecondaries done" << G4endl;
645 }
G4double MaxSecondaryEnergyForElement(G4double kineticEnergy, G4double Z)
G4double GetKineticEnergy() const
const G4ParticleDefinition * particle
G4double GetZ() const
Definition: G4Element.hh:131
int G4int
Definition: G4Types.hh:78
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double electron_mass_c2
static const G4int zdat[5]
#define G4UniformRand()
Definition: Randomize.hh:97
const G4ThreeVector & GetMomentumDirection() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4double GetLOGZ(G4int Z) const
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
T max(const T t1, const T t2)
brief Return the largest of the two arguments
Hep3Vector unit() const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:541
G4ThreeVector GetMomentum() const

Here is the call graph for this function:

void G4MuPairProductionModel::SetLowestKineticEnergy ( G4double  e)
inline

Definition at line 181 of file G4MuPairProductionModel.hh.

182 {
183  lowestKinEnergy = e;
184 }

Here is the caller graph for this function:

void G4MuPairProductionModel::SetParticle ( const G4ParticleDefinition p)
inline

Definition at line 189 of file G4MuPairProductionModel.hh.

190 {
191  if(!particle) {
192  particle = p;
194  }
195 }
const G4ParticleDefinition * particle
const char * p
Definition: xmltok.h:285
G4double GetPDGMass() const

Here is the call graph for this function:

Here is the caller graph for this function:

Member Data Documentation

G4int G4MuPairProductionModel::currentZ
protected

Definition at line 154 of file G4MuPairProductionModel.hh.

G4double G4MuPairProductionModel::factorForCross
protected

Definition at line 148 of file G4MuPairProductionModel.hh.

G4double G4MuPairProductionModel::lnZ
protected

Definition at line 153 of file G4MuPairProductionModel.hh.

G4NistManager* G4MuPairProductionModel::nist
protected

Definition at line 146 of file G4MuPairProductionModel.hh.

const G4ParticleDefinition* G4MuPairProductionModel::particle
protected

Definition at line 145 of file G4MuPairProductionModel.hh.

G4double G4MuPairProductionModel::particleMass
protected

Definition at line 150 of file G4MuPairProductionModel.hh.

G4double G4MuPairProductionModel::sqrte
protected

Definition at line 149 of file G4MuPairProductionModel.hh.

const G4double G4MuPairProductionModel::wgi
staticprotected
Initial value:
=
{ 0.0506, 0.1112, 0.1569, 0.1813, 0.1813, 0.1569, 0.1112, 0.0506 }

Definition at line 156 of file G4MuPairProductionModel.hh.

const G4double G4MuPairProductionModel::xgi
staticprotected
Initial value:
=
{ 0.0199, 0.1017, 0.2372, 0.4083, 0.5917, 0.7628, 0.8983, 0.9801 }

Definition at line 156 of file G4MuPairProductionModel.hh.

G4double G4MuPairProductionModel::z13
protected

Definition at line 151 of file G4MuPairProductionModel.hh.

G4double G4MuPairProductionModel::z23
protected

Definition at line 152 of file G4MuPairProductionModel.hh.


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