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

#include <G4hPairProductionModel.hh>

Inheritance diagram for G4hPairProductionModel:
Collaboration diagram for G4hPairProductionModel:

Public Member Functions

 G4hPairProductionModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="hPairProd")
 
virtual ~G4hPairProductionModel ()
 
- Public Member Functions inherited from G4MuPairProductionModel
 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

virtual G4double ComputeDMicroscopicCrossSection (G4double tkin, G4double Z, G4double pairEnergy) override
 
- Protected Member Functions inherited from G4MuPairProductionModel
G4double ComputMuPairLoss (G4double Z, G4double tkin, G4double cut, G4double tmax)
 
G4double ComputeMicroscopicCrossSection (G4double tkin, G4double Z, G4double cut)
 
G4double 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 *)
 

Additional Inherited Members

- Protected Attributes inherited from G4MuPairProductionModel
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 inherited from G4MuPairProductionModel
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 56 of file G4hPairProductionModel.hh.

Constructor & Destructor Documentation

G4hPairProductionModel::G4hPairProductionModel ( const G4ParticleDefinition p = nullptr,
const G4String nam = "hPairProd" 
)
explicit

Definition at line 58 of file G4hPairProductionModel.cc.

60  : G4MuPairProductionModel(p, nam)
61 {}
G4MuPairProductionModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="muPairProd")
G4hPairProductionModel::~G4hPairProductionModel ( )
virtual

Definition at line 65 of file G4hPairProductionModel.cc.

66 {}

Member Function Documentation

G4double G4hPairProductionModel::ComputeDMicroscopicCrossSection ( G4double  tkin,
G4double  Z,
G4double  pairEnergy 
)
overrideprotectedvirtual

Reimplemented from G4MuPairProductionModel.

Definition at line 70 of file G4hPairProductionModel.cc.

75 {
76  static const G4double bbbtf= 183. ;
77  static const G4double bbbh = 202.4 ;
78  static const G4double g1tf = 1.95e-5 ;
79  static const G4double g2tf = 5.3e-5 ;
80  static const G4double g1h = 4.4e-5 ;
81  static const G4double g2h = 4.8e-5 ;
82 
83  G4double totalEnergy = tkin + particleMass;
84  G4double residEnergy = totalEnergy - pairEnergy;
85  G4double massratio = particleMass/electron_mass_c2 ;
86  G4double massratio2 = massratio*massratio ;
87  G4double cross = 0.;
88 
89  G4double c3 = 0.75*sqrte*particleMass;
90  if (residEnergy <= c3*z13) { return cross; }
91 
93  G4double c8 = 6.*particleMass*particleMass;
94  G4double alf = c7/pairEnergy;
95  G4double a3 = 1. - alf;
96  if (a3 <= 0.) { return cross; }
97 
98  // zeta calculation
99  G4double bbb,g1,g2;
100  if( Z < 1.5 ) { bbb = bbbh ; g1 = g1h ; g2 = g2h ; }
101  else { bbb = bbbtf; g1 = g1tf; g2 = g2tf; }
102 
103  G4double zeta = 0.;
104  G4double zeta1 =
105  0.073*G4Log(totalEnergy/(particleMass+g1*z23*totalEnergy))-0.26;
106  if ( zeta1 > 0.)
107  {
108  G4double zeta2 =
109  0.058*G4Log(totalEnergy/(particleMass+g2*z13*totalEnergy))-0.14;
110  zeta = zeta1/zeta2 ;
111  }
112 
113  G4double z2 = Z*(Z+zeta);
114  G4double screen0 = 2.*electron_mass_c2*sqrte*bbb/(z13*pairEnergy);
115  G4double a0 = totalEnergy*residEnergy;
116  G4double a1 = pairEnergy*pairEnergy/a0;
117  G4double bet = 0.5*a1;
118  G4double xi0 = 0.25*massratio2*a1;
119  G4double del = c8/a0;
120 
121  G4double rta3 = sqrt(a3);
122  G4double tmnexp = alf/(1. + rta3) + del*rta3;
123  if(tmnexp >= 1.0) { return cross; }
124 
125  G4double tmn = G4Log(tmnexp);
126  G4double sum = 0.;
127 
128  // Gaussian integration in ln(1-ro) ( with 8 points)
129  for (G4int i=0; i<8; i++)
130  {
131  G4double a4 = G4Exp(tmn*xgi[i]); // a4 = (1.-asymmetry)
132  G4double a5 = a4*(2.-a4) ;
133  G4double a6 = 1.-a5 ;
134  G4double a7 = 1.+a6 ;
135  G4double a9 = 3.+a6 ;
136  G4double xi = xi0*a5 ;
137  G4double xii = 1./xi ;
138  G4double xi1 = 1.+xi ;
139  G4double screen = screen0*xi1/a5 ;
140  G4double yeu = 5.-a6+4.*bet*a7 ;
141  G4double yed = 2.*(1.+3.*bet)*G4Log(3.+xii)-a6-a1*(2.-a6) ;
142  G4double ye1 = 1.+yeu/yed ;
143  G4double ale=G4Log(bbb/z13*sqrt(xi1*ye1)/(1.+screen*ye1)) ;
144  G4double cre = 0.5*G4Log(1.+2.25*z23*xi1*ye1/massratio2) ;
145  G4double be;
146 
147  if (xi <= 1.e3) {
148  be = ((2.+a6)*(1.+bet)+xi*a9)*G4Log(1.+xii)+(a5-bet)/xi1-a9;
149  } else {
150  be = (3.-a6+a1*a7)/(2.*xi);
151  }
152  G4double fe = (ale-cre)*be;
153  if ( fe < 0.) { fe = 0.; }
154 
155  G4double ymu = 4.+a6 +3.*bet*a7 ;
156  G4double ymd = a7*(1.5+a1)*G4Log(3.+xi)+1.-1.5*a6 ;
157  G4double ym1 = 1.+ymu/ymd ;
158  G4double alm_crm = G4Log(bbb*massratio/(1.5*z23*(1.+screen*ym1)));
159  G4double a10,bm;
160  if ( xi >= 1.e-3)
161  {
162  a10 = (1.+a1)*a5 ;
163  bm = (a7*(1.+1.5*bet)-a10*xii)*G4Log(xi1)+xi*(a5-bet)/xi1+a10;
164  } else {
165  bm = (5.-a6+bet*a9)*(xi/2.);
166  }
167 
168  G4double fm = alm_crm*bm;
169  if ( fm < 0.) { fm = 0.; }
170 
171  sum += wgi[i]*a4*(fe+fm/massratio2);
172  }
173 
174  cross = -tmn*sum*factorForCross*z2*residEnergy/(totalEnergy*pairEnergy);
175 
176  return cross;
177 }
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
static const G4double wgi[8]
double G4double
Definition: G4Types.hh:76
G4fissionEvent * fe

Here is the call graph for this function:


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