Geant4  10.02.p03
G4BetheHeitlerModel Class Reference

#include <G4BetheHeitlerModel.hh>

Inheritance diagram for G4BetheHeitlerModel:
Collaboration diagram for G4BetheHeitlerModel:

Public Member Functions

 G4BetheHeitlerModel (const G4ParticleDefinition *p=0, const G4String &nam="BetheHeitler")
 
virtual ~G4BetheHeitlerModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle *> *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
- 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 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, const G4ParticleDefinition *, G4double)
 
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 *> *)
 
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=0)
 
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)
 

Private Member Functions

G4double ScreenFunction1 (G4double ScreenVariable)
 
G4double ScreenFunction2 (G4double ScreenVariable)
 
G4BetheHeitlerModeloperator= (const G4BetheHeitlerModel &right)
 
 G4BetheHeitlerModel (const G4BetheHeitlerModel &)
 

Private Attributes

G4Powg4pow
 
G4ParticleDefinitiontheGamma
 
G4ParticleDefinitiontheElectron
 
G4ParticleDefinitionthePositron
 
G4ParticleChangeForGamma * fParticleChange
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLoss * GetParticleChangeForLoss ()
 
G4ParticleChangeForGamma * GetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChange * pParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 60 of file G4BetheHeitlerModel.hh.

Constructor & Destructor Documentation

◆ G4BetheHeitlerModel() [1/2]

G4BetheHeitlerModel::G4BetheHeitlerModel ( const G4ParticleDefinition p = 0,
const G4String nam = "BetheHeitler" 
)

Definition at line 85 of file G4BetheHeitlerModel.cc.

87  : G4VEmModel(nam)
88 {
89  fParticleChange = 0;
94 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4ParticleChangeForGamma * fParticleChange
G4ParticleDefinition * thePositron
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
G4ParticleDefinition * theGamma
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4ParticleDefinition * theElectron
static G4Electron * Electron()
Definition: G4Electron.cc:94
Here is the call graph for this function:

◆ ~G4BetheHeitlerModel()

G4BetheHeitlerModel::~G4BetheHeitlerModel ( )
virtual

Definition at line 98 of file G4BetheHeitlerModel.cc.

99 {}

◆ G4BetheHeitlerModel() [2/2]

G4BetheHeitlerModel::G4BetheHeitlerModel ( const G4BetheHeitlerModel )
private

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4BetheHeitlerModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0.,
G4double  cut = 0.,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 121 of file G4BetheHeitlerModel.cc.

130 {
131  G4double xSection = 0.0 ;
132  if ( Z < 0.9 || GammaEnergy <= 2.0*electron_mass_c2 ) { return xSection; }
133 
134 
135  G4double GammaEnergySave = GammaEnergy;
136  if (GammaEnergy < GammaEnergyLimit) { GammaEnergy = GammaEnergyLimit; }
137 
138  G4double X=G4Log(GammaEnergy/electron_mass_c2), X2=X*X, X3=X2*X, X4=X3*X, X5=X4*X;
139 
140  G4double F1 = a0 + a1*X + a2*X2 + a3*X3 + a4*X4 + a5*X5,
141  F2 = b0 + b1*X + b2*X2 + b3*X3 + b4*X4 + b5*X5,
142  F3 = c0 + c1*X + c2*X2 + c3*X3 + c4*X4 + c5*X5;
143 
144  xSection = (Z + 1.)*(F1*Z + F2*Z*Z + F3);
145 
146  if (GammaEnergySave < GammaEnergyLimit) {
147 
148  X = (GammaEnergySave - 2.*electron_mass_c2)
150  xSection *= X*X;
151  }
152 
153  if (xSection < 0.) { xSection = 0.; }
154  return xSection;
155 }
static const G4double a0
static const G4double c2
static const G4double a1
static const G4double a4
Double_t X2
Float_t X
static const G4double b3
Float_t Z
static const G4double c3
static const G4double c4
static const G4double b2
float electron_mass_c2
Definition: hepunit.py:274
static const G4double c1
static const G4double a3
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static const G4double c0
static const G4double b1
static const G4double c5
static const G4double a5
static const G4double GammaEnergyLimit
double G4double
Definition: G4Types.hh:76
static const G4double b0
static const G4double b5
static const G4double b4
static const G4double a2
Here is the call graph for this function:

◆ Initialise()

void G4BetheHeitlerModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Reimplemented in G4PolarizedGammaConversionModel.

Definition at line 103 of file G4BetheHeitlerModel.cc.

105 {
107  if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
108 }
G4ParticleChangeForGamma * fParticleChange
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
Here is the call graph for this function:
Here is the caller graph for this function:

◆ InitialiseLocal()

void G4BetheHeitlerModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 112 of file G4BetheHeitlerModel.cc.

114 {
116 }
void SetElementSelectors(std::vector< G4EmElementSelector *> *)
Definition: G4VEmModel.hh:810
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:802
Here is the call graph for this function:

◆ operator=()

G4BetheHeitlerModel& G4BetheHeitlerModel::operator= ( const G4BetheHeitlerModel right)
private

◆ SampleSecondaries()

void G4BetheHeitlerModel::SampleSecondaries ( std::vector< G4DynamicParticle *> *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Reimplemented in G4PolarizedGammaConversionModel.

Definition at line 159 of file G4BetheHeitlerModel.cc.

176 {
177  const G4Material* aMaterial = couple->GetMaterial();
178 
179  G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
180  G4ParticleMomentum GammaDirection = aDynamicGamma->GetMomentumDirection();
181 
182  G4double epsil ;
183  G4double epsil0 = electron_mass_c2/GammaEnergy ;
184  if(epsil0 > 1.0) { return; }
185 
186  // do it fast if GammaEnergy < Egsmall
187  // select randomly one element constituing the material
188  const G4Element* anElement = SelectRandomAtom(aMaterial, theGamma, GammaEnergy);
189 
190  CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
191 
192  if (GammaEnergy < Egsmall) {
193 
194  epsil = epsil0 + (0.5-epsil0)*rndmEngine->flat();
195 
196  } else {
197  // now comes the case with GammaEnergy >= 2. MeV
198 
199  // Extract Coulomb factor for this Element
200  G4double FZ = 8.*(anElement->GetIonisation()->GetlogZ3());
201  if (GammaEnergy > 50.*MeV) { FZ += 8.*(anElement->GetfCoulomb()); }
202 
203  // limits of the screening variable
204  G4double screenfac = 136.*epsil0/(anElement->GetIonisation()->GetZ3());
205  G4double screenmax = exp ((42.24 - FZ)/8.368) - 0.952 ;
206  G4double screenmin = min(4.*screenfac,screenmax);
207 
208  // limits of the energy sampling
209  G4double epsil1 = 0.5 - 0.5*sqrt(1. - screenmin/screenmax) ;
210  G4double epsilmin = max(epsil0,epsil1) , epsilrange = 0.5 - epsilmin;
211 
212  //
213  // sample the energy rate of the created electron (or positron)
214  //
215  //G4double epsil, screenvar, greject ;
216  G4double screenvar, greject ;
217 
218  G4double F10 = ScreenFunction1(screenmin) - FZ;
219  G4double F20 = ScreenFunction2(screenmin) - FZ;
220  G4double NormF1 = max(F10*epsilrange*epsilrange,0.);
221  G4double NormF2 = max(1.5*F20,0.);
222 
223  do {
224  if ( NormF1/(NormF1+NormF2) > rndmEngine->flat()) {
225  epsil = 0.5 - epsilrange*g4pow->A13(rndmEngine->flat());
226  screenvar = screenfac/(epsil*(1-epsil));
227  greject = (ScreenFunction1(screenvar) - FZ)/F10;
228 
229  } else {
230  epsil = epsilmin + epsilrange*rndmEngine->flat();
231  screenvar = screenfac/(epsil*(1-epsil));
232  greject = (ScreenFunction2(screenvar) - FZ)/F20;
233  }
234 
235  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
236  } while( greject < rndmEngine->flat());
237 
238  } // end of epsil sampling
239 
240  //
241  // fixe charges randomly
242  //
243 
244  G4double ElectTotEnergy, PositTotEnergy;
245  if (rndmEngine->flat() > 0.5) {
246 
247  ElectTotEnergy = (1.-epsil)*GammaEnergy;
248  PositTotEnergy = epsil*GammaEnergy;
249 
250  } else {
251 
252  PositTotEnergy = (1.-epsil)*GammaEnergy;
253  ElectTotEnergy = epsil*GammaEnergy;
254  }
255 
256  //
257  // scattered electron (positron) angles. ( Z - axis along the parent photon)
258  //
259  // universal distribution suggested by L. Urban
260  // (Geant3 manual (1993) Phys211),
261  // derived from Tsai distribution (Rev Mod Phys 49,421(1977))
262 
263  G4double u= - G4Log(rndmEngine->flat()*rndmEngine->flat());
264 
265  if (9. > 36.*rndmEngine->flat()) { u *= 1.6; }
266  else { u *= 0.53333; }
267 
268  G4double TetEl = u*electron_mass_c2/ElectTotEnergy;
269  G4double TetPo = u*electron_mass_c2/PositTotEnergy;
270  G4double Phi = twopi * rndmEngine->flat();
271  G4double dxEl= sin(TetEl)*cos(Phi),dyEl= sin(TetEl)*sin(Phi),dzEl=cos(TetEl);
272  G4double dxPo=-sin(TetPo)*cos(Phi),dyPo=-sin(TetPo)*sin(Phi),dzPo=cos(TetPo);
273 
274  //
275  // kinematic of the created pair
276  //
277  // the electron and positron are assumed to have a symetric
278  // angular distribution with respect to the Z axis along the parent photon.
279 
280  G4double ElectKineEnergy = max(0.,ElectTotEnergy - electron_mass_c2);
281 
282  G4ThreeVector ElectDirection (dxEl, dyEl, dzEl);
283  ElectDirection.rotateUz(GammaDirection);
284 
285  // create G4DynamicParticle object for the particle1
286  G4DynamicParticle* aParticle1= new G4DynamicParticle(
287  theElectron,ElectDirection,ElectKineEnergy);
288 
289  // the e+ is always created (even with Ekine=0) for further annihilation.
290 
291  G4double PositKineEnergy = max(0.,PositTotEnergy - electron_mass_c2);
292 
293  G4ThreeVector PositDirection (dxPo, dyPo, dzPo);
294  PositDirection.rotateUz(GammaDirection);
295 
296  // create G4DynamicParticle object for the particle2
297  G4DynamicParticle* aParticle2= new G4DynamicParticle(
298  thePositron,PositDirection,PositKineEnergy);
299 
300  // Fill output vector
301  fvect->push_back(aParticle1);
302  fvect->push_back(aParticle2);
303 
304  // kill incident photon
305  fParticleChange->SetProposedKineticEnergy(0.);
306  fParticleChange->ProposeTrackStatus(fStopAndKill);
307 }
static const double MeV
Definition: G4SIunits.hh:211
G4ParticleChangeForGamma * fParticleChange
const G4Material * GetMaterial() const
G4double GetlogZ3() const
G4ParticleDefinition * thePositron
virtual double flat()=0
G4double ScreenFunction2(G4double ScreenVariable)
#define F20
G4double GetZ3() const
static const G4double Egsmall
G4double GetKineticEnergy() const
G4double ScreenFunction1(G4double ScreenVariable)
static const double twopi
Definition: G4SIunits.hh:75
G4IonisParamElm * GetIonisation() const
Definition: G4Element.hh:198
float electron_mass_c2
Definition: hepunit.py:274
G4ParticleDefinition * theGamma
#define F10
double flat()
Definition: G4AblaRandom.cc:47
G4double GetfCoulomb() const
Definition: G4Element.hh:190
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4ParticleDefinition * theElectron
G4double A13(G4double A) const
Definition: G4Pow.hh:132
const G4ThreeVector & GetMomentumDirection() const
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:544
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ScreenFunction1()

G4double G4BetheHeitlerModel::ScreenFunction1 ( G4double  ScreenVariable)
inlineprivate

Definition at line 108 of file G4BetheHeitlerModel.hh.

112 {
113  G4double screenVal;
114 
115  if (ScreenVariable > 1.)
116  screenVal = 42.24 - 8.368*G4Log(ScreenVariable+0.952);
117  else
118  screenVal = 42.392 - ScreenVariable*(7.796 - 1.961*ScreenVariable);
119 
120  return screenVal;
121 }
G4double G4Log(G4double x)
Definition: G4Log.hh:230
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ScreenFunction2()

G4double G4BetheHeitlerModel::ScreenFunction2 ( G4double  ScreenVariable)
inlineprivate

Definition at line 125 of file G4BetheHeitlerModel.hh.

129 {
130  G4double screenVal;
131 
132  if (ScreenVariable > 1.)
133  screenVal = 42.24 - 8.368*G4Log(ScreenVariable+0.952);
134  else
135  screenVal = 41.405 - ScreenVariable*(5.828 - 0.8945*ScreenVariable);
136 
137  return screenVal;
138 }
G4double G4Log(G4double x)
Definition: G4Log.hh:230
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForGamma* G4BetheHeitlerModel::fParticleChange
private

Definition at line 103 of file G4BetheHeitlerModel.hh.

◆ g4pow

G4Pow* G4BetheHeitlerModel::g4pow
private

Definition at line 99 of file G4BetheHeitlerModel.hh.

◆ theElectron

G4ParticleDefinition* G4BetheHeitlerModel::theElectron
private

Definition at line 101 of file G4BetheHeitlerModel.hh.

◆ theGamma

G4ParticleDefinition* G4BetheHeitlerModel::theGamma
private

Definition at line 100 of file G4BetheHeitlerModel.hh.

◆ thePositron

G4ParticleDefinition* G4BetheHeitlerModel::thePositron
private

Definition at line 102 of file G4BetheHeitlerModel.hh.


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