Geant4  10.02.p03
G4LivermoreNuclearGammaConversionModel Class Reference

#include <G4LivermoreNuclearGammaConversionModel.hh>

Inheritance diagram for G4LivermoreNuclearGammaConversionModel:
Collaboration diagram for G4LivermoreNuclearGammaConversionModel:

Public Member Functions

 G4LivermoreNuclearGammaConversionModel (const G4ParticleDefinition *p=0, const G4String &nam="LivermoreNuclearConversion")
 
virtual ~G4LivermoreNuclearGammaConversionModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
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)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
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 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

void ReadData (size_t Z, const char *path=0)
 
G4double ScreenFunction1 (G4double screenVariable)
 
G4double ScreenFunction2 (G4double screenVariable)
 
G4LivermoreNuclearGammaConversionModeloperator= (const G4LivermoreNuclearGammaConversionModel &right)
 
 G4LivermoreNuclearGammaConversionModel (const G4LivermoreNuclearGammaConversionModel &)
 

Private Attributes

G4bool isInitialised
 
G4int verboseLevel
 
G4double lowEnergyLimit
 
G4double smallEnergy
 
G4ParticleChangeForGamma * fParticleChange
 

Static Private Attributes

static G4int maxZ = 100
 
static G4LPhysicsFreeVectordata [100] = {0}
 

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 38 of file G4LivermoreNuclearGammaConversionModel.hh.

Constructor & Destructor Documentation

◆ G4LivermoreNuclearGammaConversionModel() [1/2]

G4LivermoreNuclearGammaConversionModel::G4LivermoreNuclearGammaConversionModel ( const G4ParticleDefinition p = 0,
const G4String nam = "LivermoreNuclearConversion" 
)

Definition at line 47 of file G4LivermoreNuclearGammaConversionModel.cc.

48 :G4VEmModel(nam),isInitialised(false),smallEnergy(2.*MeV)
49 {
50  fParticleChange = 0;
51 
53 
54  verboseLevel= 0;
55  // Verbosity scale for debugging purposes:
56  // 0 = nothing
57  // 1 = calculation of cross sections, file openings...
58  // 2 = entering in methods
59 
60  if(verboseLevel > 0)
61  {
62  G4cout << "G4LivermoreNuclearGammaConversionModel is constructed " << G4endl;
63  }
64 }
static const double MeV
Definition: G4SIunits.hh:211
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
G4GLOB_DLL std::ostream G4cout
float electron_mass_c2
Definition: hepunit.py:274
#define G4endl
Definition: G4ios.hh:61

◆ ~G4LivermoreNuclearGammaConversionModel()

G4LivermoreNuclearGammaConversionModel::~G4LivermoreNuclearGammaConversionModel ( )
virtual

Definition at line 68 of file G4LivermoreNuclearGammaConversionModel.cc.

69 {
70  if(IsMaster()) {
71  for(G4int i=0; i<maxZ; ++i) {
72  if(data[i]) {
73  delete data[i];
74  data[i] = 0;
75  }
76  }
77  }
78 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
int G4int
Definition: G4Types.hh:78

◆ G4LivermoreNuclearGammaConversionModel() [2/2]

G4LivermoreNuclearGammaConversionModel::G4LivermoreNuclearGammaConversionModel ( const G4LivermoreNuclearGammaConversionModel )
private

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 217 of file G4LivermoreNuclearGammaConversionModel.cc.

221 {
222  if (verboseLevel > 1)
223  {
224  G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreNuclearGammaConversionModel"
225  << G4endl;
226  }
227 
228  if (GammaEnergy < lowEnergyLimit) { return 0.0; }
229 
230  G4double xs = 0.0;
231 
232  G4int intZ=G4int(Z);
233 
234  if(intZ < 1 || intZ > maxZ) { return xs; }
235 
236  G4LPhysicsFreeVector* pv = data[intZ];
237 
238  // if element was not initialised
239  // do initialisation safely for MT mode
240  if(!pv)
241  {
242  InitialiseForElement(0, intZ);
243  pv = data[intZ];
244  if(!pv) { return xs; }
245  }
246  // x-section is taken from the table
247  xs = pv->Value(GammaEnergy);
248 
249  if(verboseLevel > 0)
250  {
251  G4int n = pv->GetVectorLength() - 1;
252  G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)="
253  << GammaEnergy/MeV << G4endl;
254  G4cout << " cs (Geant4 internal unit)=" << xs << G4endl;
255  G4cout << " -> first cs value in EADL data file (iu) =" << (*pv)[0] << G4endl;
256  G4cout << " -> last cs value in EADL data file (iu) =" << (*pv)[n] << G4endl;
257  G4cout << "*********************************************************" << G4endl;
258  }
259 
260  return xs;
261 
262 }
static const double MeV
Definition: G4SIunits.hh:211
int G4int
Definition: G4Types.hh:78
Char_t n[5]
G4GLOB_DLL std::ostream G4cout
Float_t Z
size_t GetVectorLength() const
G4double Value(G4double theEnergy, size_t &lastidx) const
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ Initialise()

void G4LivermoreNuclearGammaConversionModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 82 of file G4LivermoreNuclearGammaConversionModel.cc.

85 {
86 
87  if (verboseLevel > 1)
88  {
89  G4cout << "Calling Initialise() of G4LivermoreNuclearGammaConversionModel."
90  << G4endl
91  << "Energy range: "
92  << LowEnergyLimit() / MeV << " MeV - "
93  << HighEnergyLimit() / GeV << " GeV"
94  << G4endl;
95  }
96 
97  if(IsMaster())
98  {
99 
100  // Initialise element selector
101 
102  InitialiseElementSelectors(particle, cuts);
103 
104  // Access to elements
105 
106  char* path = getenv("G4LEDATA");
107 
108  G4ProductionCutsTable* theCoupleTable =
110 
111  G4int numOfCouples = theCoupleTable->GetTableSize();
112 
113  for(G4int i=0; i<numOfCouples; ++i)
114  {
115  const G4Material* material =
116  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
117  const G4ElementVector* theElementVector = material->GetElementVector();
118  G4int nelm = material->GetNumberOfElements();
119 
120  for (G4int j=0; j<nelm; ++j)
121  {
122  G4int Z = (G4int)(*theElementVector)[j]->GetZ();
123  if(Z < 1) { Z = 1; }
124  else if(Z > maxZ) { Z = maxZ; }
125  if(!data[Z]) { ReadData(Z, path); }
126  }
127  }
128  }
129  if(isInitialised) { return; }
131  isInitialised = true;
132 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
static const double MeV
Definition: G4SIunits.hh:211
std::vector< G4Element * > G4ElementVector
const G4Material * GetMaterial() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
int G4int
Definition: G4Types.hh:78
string material
Definition: eplot.py:19
G4GLOB_DLL std::ostream G4cout
Float_t Z
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
static const double GeV
Definition: G4SIunits.hh:214
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
#define G4endl
Definition: G4ios.hh:61
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
Here is the call graph for this function:

◆ InitialiseForElement()

void G4LivermoreNuclearGammaConversionModel::InitialiseForElement ( const G4ParticleDefinition ,
G4int  Z 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 476 of file G4LivermoreNuclearGammaConversionModel.cc.

479 {
480  G4AutoLock l(&LivermoreNuclearGammaConversionModelMutex);
481  // G4cout << "G4LivermoreNuclearGammaConversionModel::InitialiseForElement Z= "
482  // << Z << G4endl;
483  if(!data[Z]) { ReadData(Z); }
484  l.unlock();
485 }
Float_t Z
Here is the call graph for this function:

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 136 of file G4LivermoreNuclearGammaConversionModel.cc.

138 {
140 }
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:

◆ MinPrimaryEnergy()

G4double G4LivermoreNuclearGammaConversionModel::MinPrimaryEnergy ( const G4Material ,
const G4ParticleDefinition ,
G4double   
)
virtual

Reimplemented from G4VEmModel.

Definition at line 145 of file G4LivermoreNuclearGammaConversionModel.cc.

◆ operator=()

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

◆ ReadData()

void G4LivermoreNuclearGammaConversionModel::ReadData ( size_t  Z,
const char *  path = 0 
)
private

Definition at line 154 of file G4LivermoreNuclearGammaConversionModel.cc.

155 {
156  if (verboseLevel > 1)
157  {
158  G4cout << "Calling ReadData() of G4LivermoreNuclearGammaConversionModel"
159  << G4endl;
160  }
161 
162 
163  if(data[Z]) { return; }
164 
165  const char* datadir = path;
166 
167  if(!datadir)
168  {
169  datadir = getenv("G4LEDATA");
170  if(!datadir)
171  {
172  G4Exception("G4LivermoreNuclearGammaConversionModel::ReadData()",
173  "em0006",FatalException,
174  "Environment variable G4LEDATA not defined");
175  return;
176  }
177  }
178 
179  //
180 
181  data[Z] = new G4LPhysicsFreeVector();
182 
183  //
184 
185  std::ostringstream ost;
186  ost << datadir << "livermore/pairdata/pp-pair-cs-" << Z <<".dat";
187  std::ifstream fin(ost.str().c_str());
188 
189  if( !fin.is_open())
190  {
192  ed << "G4LivermoreNuclearGammaConversionModel data file <" << ost.str().c_str()
193  << "> is not opened!" << G4endl;
194  G4Exception("G4LivermoreNuclearGammaConversionModel::ReadData()",
195  "em0003",FatalException,
196  ed,"G4LEDATA version should be G4EMLOW6.27 or later.");
197  return;
198  }
199 
200  else
201  {
202 
203  if(verboseLevel > 3) { G4cout << "File " << ost.str()
204  << " is opened by G4LivermoreNuclearGammaConversionModel" << G4endl;}
205 
206  data[Z]->Retrieve(fin, true);
207  }
208 
209  // Activation of spline interpolation
210  data[Z] ->SetSpline(true);
211 
212 }
TString fin
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void SetSpline(G4bool)
G4GLOB_DLL std::ostream G4cout
Float_t Z
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 266 of file G4LivermoreNuclearGammaConversionModel.cc.

271 {
272 
273 // The energies of the e+ e- secondaries are sampled using the Bethe - Heitler
274 // cross sections with Coulomb correction. A modified version of the random
275 // number techniques of Butcher & Messel is used (Nuc Phys 20(1960),15).
276 
277 // Note 1 : Effects due to the breakdown of the Born approximation at low
278 // energy are ignored.
279 // Note 2 : The differential cross section implicitly takes account of
280 // pair creation in both nuclear and atomic electron fields. However triplet
281 // prodution is not generated.
282 
283  if (verboseLevel > 1) {
284  G4cout << "Calling SampleSecondaries() of G4LivermoreNuclearGammaConversionModel"
285  << G4endl;
286  }
287 
288  G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
289  G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
290 
291  G4double epsilon ;
292  G4double epsilon0Local = electron_mass_c2 / photonEnergy ;
293 
294  // Do it fast if photon energy < 2. MeV
295  if (photonEnergy < smallEnergy )
296  {
297  epsilon = epsilon0Local + (0.5 - epsilon0Local) * G4UniformRand();
298  }
299  else
300  {
301  // Select randomly one element in the current material
302 
303  const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
304  const G4Element* element = SelectRandomAtom(couple,particle,photonEnergy);
305 
306  if (element == 0)
307  {
308  G4cout << "G4LivermoreNuclearGammaConversionModel::SampleSecondaries - element = 0"
309  << G4endl;
310  return;
311  }
312  G4IonisParamElm* ionisation = element->GetIonisation();
313  if (ionisation == 0)
314  {
315  G4cout << "G4LivermoreNuclearGammaConversionModel::SampleSecondaries - ionisation = 0"
316  << G4endl;
317  return;
318  }
319 
320  // Extract Coulomb factor for this Elements
321  G4double fZ = 8. * (ionisation->GetlogZ3());
322  if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
323 
324  // Limits of the screening variable
325  G4double screenFactor = 136. * epsilon0Local / (element->GetIonisation()->GetZ3()) ;
326  G4double screenMax = G4Exp ((42.24 - fZ)/8.368) - 0.952 ;
327  G4double screenMin = std::min(4.*screenFactor,screenMax) ;
328 
329  // Limits of the energy sampling
330  G4double epsilon1 = 0.5 - 0.5 * std::sqrt(1. - screenMin / screenMax) ;
331  G4double epsilonMin = std::max(epsilon0Local,epsilon1);
332  G4double epsilonRange = 0.5 - epsilonMin ;
333 
334  // Sample the energy rate of the created electron (or positron)
335  G4double screen;
336  G4double gReject ;
337 
338  G4double f10 = ScreenFunction1(screenMin) - fZ;
339  G4double f20 = ScreenFunction2(screenMin) - fZ;
340  G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
341  G4double normF2 = std::max(1.5 * f20,0.);
342 
343  do
344  {
345  if (normF1 / (normF1 + normF2) > G4UniformRand() )
346  {
347  epsilon = 0.5 - epsilonRange * std::pow(G4UniformRand(), 0.333333) ;
348  screen = screenFactor / (epsilon * (1. - epsilon));
349  gReject = (ScreenFunction1(screen) - fZ) / f10 ;
350  }
351  else
352  {
353  epsilon = epsilonMin + epsilonRange * G4UniformRand();
354  screen = screenFactor / (epsilon * (1 - epsilon));
355  gReject = (ScreenFunction2(screen) - fZ) / f20 ;
356  }
357  } while ( gReject < G4UniformRand() );
358 
359  } // End of epsilon sampling
360 
361  // Fix charges randomly
362 
363  G4double electronTotEnergy;
364  G4double positronTotEnergy;
365 
366  if (G4UniformRand() > 0.5)
367  {
368  electronTotEnergy = (1. - epsilon) * photonEnergy;
369  positronTotEnergy = epsilon * photonEnergy;
370  }
371  else
372  {
373  positronTotEnergy = (1. - epsilon) * photonEnergy;
374  electronTotEnergy = epsilon * photonEnergy;
375  }
376 
377  // Scattered electron (positron) angles. ( Z - axis along the parent photon)
378  // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
379  // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
380 
381  G4double u;
382  const G4double a1 = 0.625;
383  G4double a2 = 3. * a1;
384  // G4double d = 27. ;
385 
386  // if (9. / (9. + d) > G4UniformRand())
387  if (0.25 > G4UniformRand())
388  {
389  u = - G4Log(G4UniformRand() * G4UniformRand()) / a1 ;
390  }
391  else
392  {
393  u = - G4Log(G4UniformRand() * G4UniformRand()) / a2 ;
394  }
395 
396  G4double thetaEle = u*electron_mass_c2/electronTotEnergy;
397  G4double thetaPos = u*electron_mass_c2/positronTotEnergy;
398  G4double phi = twopi * G4UniformRand();
399 
400  G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
401  G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
402 
403 
404  // Kinematics of the created pair:
405  // the electron and positron are assumed to have a symetric angular
406  // distribution with respect to the Z axis along the parent photon
407 
408  G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
409 
410  G4ThreeVector electronDirection (dxEle, dyEle, dzEle);
411  electronDirection.rotateUz(photonDirection);
412 
414  electronDirection,
415  electronKineEnergy);
416 
417  // The e+ is always created
418  G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
419 
420  G4ThreeVector positronDirection (dxPos, dyPos, dzPos);
421  positronDirection.rotateUz(photonDirection);
422 
423  // Create G4DynamicParticle object for the particle2
425  positronDirection,
426  positronKineEnergy);
427  // Fill output vector
428  fvect->push_back(particle1);
429  fvect->push_back(particle2);
430 
431  // kill incident photon
432  fParticleChange->SetProposedKineticEnergy(0.);
433  fParticleChange->ProposeTrackStatus(fStopAndKill);
434 
435 }
static const double MeV
Definition: G4SIunits.hh:211
G4double GetlogZ3() const
static const G4double a1
G4double GetZ3() const
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
static const double twopi
Definition: G4SIunits.hh:75
G4IonisParamElm * GetIonisation() const
Definition: G4Element.hh:198
float electron_mass_c2
Definition: hepunit.py:274
G4double GetfCoulomb() const
Definition: G4Element.hh:190
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static G4Positron * Positron()
Definition: G4Positron.cc:94
const G4ThreeVector & GetMomentumDirection() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
double epsilon(double density, double temperature)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:544
static const G4double a2
Here is the call graph for this function:

◆ ScreenFunction1()

G4double G4LivermoreNuclearGammaConversionModel::ScreenFunction1 ( G4double  screenVariable)
private

Definition at line 440 of file G4LivermoreNuclearGammaConversionModel.cc.

441 {
442  // Compute the value of the screening function 3*phi1 - phi2
443 
444  G4double value;
445 
446  if (screenVariable > 1.)
447  value = 42.24 - 8.368 * G4Log(screenVariable + 0.952);
448  else
449  value = 42.392 - screenVariable * (7.796 - 1.961 * screenVariable);
450 
451  return value;
452 }
G4double G4Log(G4double x)
Definition: G4Log.hh:230
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ ScreenFunction2()

G4double G4LivermoreNuclearGammaConversionModel::ScreenFunction2 ( G4double  screenVariable)
private

Definition at line 457 of file G4LivermoreNuclearGammaConversionModel.cc.

458 {
459  // Compute the value of the screening function 1.5*phi1 - 0.5*phi2
460 
461  G4double value;
462 
463  if (screenVariable > 1.)
464  value = 42.24 - 8.368 * G4Log(screenVariable + 0.952);
465  else
466  value = 41.405 - screenVariable * (5.828 - 0.8945 * screenVariable);
467 
468  return value;
469 }
G4double G4Log(G4double x)
Definition: G4Log.hh:230
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

Member Data Documentation

◆ data

G4LPhysicsFreeVector * G4LivermoreNuclearGammaConversionModel::data = {0}
staticprivate

Definition at line 94 of file G4LivermoreNuclearGammaConversionModel.hh.

◆ fParticleChange

G4ParticleChangeForGamma* G4LivermoreNuclearGammaConversionModel::fParticleChange
private

Definition at line 99 of file G4LivermoreNuclearGammaConversionModel.hh.

◆ isInitialised

G4bool G4LivermoreNuclearGammaConversionModel::isInitialised
private

Definition at line 86 of file G4LivermoreNuclearGammaConversionModel.hh.

◆ lowEnergyLimit

G4double G4LivermoreNuclearGammaConversionModel::lowEnergyLimit
private

Definition at line 89 of file G4LivermoreNuclearGammaConversionModel.hh.

◆ maxZ

G4int G4LivermoreNuclearGammaConversionModel::maxZ = 100
staticprivate

Definition at line 93 of file G4LivermoreNuclearGammaConversionModel.hh.

◆ smallEnergy

G4double G4LivermoreNuclearGammaConversionModel::smallEnergy
private

Definition at line 90 of file G4LivermoreNuclearGammaConversionModel.hh.

◆ verboseLevel

G4int G4LivermoreNuclearGammaConversionModel::verboseLevel
private

Definition at line 87 of file G4LivermoreNuclearGammaConversionModel.hh.


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