Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4LivermoreGammaConversionModel Class Reference

#include <G4LivermoreGammaConversionModel.hh>

Inheritance diagram for G4LivermoreGammaConversionModel:
G4VEmModel

Public Member Functions

 G4LivermoreGammaConversionModel (const G4ParticleDefinition *p=0, const G4String &nam="LivermoreConversion")
 
virtual ~G4LivermoreGammaConversionModel ()
 
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 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 *)
 
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
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Detailed Description

Definition at line 41 of file G4LivermoreGammaConversionModel.hh.

Constructor & Destructor Documentation

G4LivermoreGammaConversionModel::G4LivermoreGammaConversionModel ( const G4ParticleDefinition p = 0,
const G4String nam = "LivermoreConversion" 
)

Definition at line 47 of file G4LivermoreGammaConversionModel.cc.

References python.hepunit::electron_mass_c2, G4cout, and G4endl.

48 :G4VEmModel(nam),isInitialised(false),smallEnergy(2.*MeV)
49 {
50  fParticleChange = 0;
51 
52  lowEnergyLimit = 2.0*electron_mass_c2;
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 << "G4LivermoreGammaConversionModel is constructed " << G4endl;
63  }
64 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
G4GLOB_DLL std::ostream G4cout
float electron_mass_c2
Definition: hepunit.py:274
#define G4endl
Definition: G4ios.hh:61
G4LivermoreGammaConversionModel::~G4LivermoreGammaConversionModel ( )
virtual

Definition at line 68 of file G4LivermoreGammaConversionModel.cc.

69 {}

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 209 of file G4LivermoreGammaConversionModel.cc.

References G4cout, G4endl, G4PhysicsVector::GetVectorLength(), python.hepunit::MeV, n, and G4PhysicsVector::Value().

213 {
214  if (verboseLevel > 1)
215  {
216  G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreGammaConversionModel"
217  << G4endl;
218  }
219 
220  if (GammaEnergy < lowEnergyLimit) { return 0.0; }
221 
222  G4double xs = 0.0;
223 
224  G4int intZ=G4int(Z);
225 
226  if(intZ < 1 || intZ > maxZ) { return xs; }
227 
228  G4LPhysicsFreeVector* pv = data[intZ];
229 
230  // if element was not initialised
231  // do initialisation safely for MT mode
232  if(!pv)
233  {
234  InitialiseForElement(0, intZ);
235  pv = data[intZ];
236  if(!pv) { return xs; }
237  }
238  // x-section is taken from the table
239  xs = pv->Value(GammaEnergy);
240 
241  if(verboseLevel > 0)
242  {
243  G4int n = pv->GetVectorLength() - 1;
244  G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)="
245  << GammaEnergy/MeV << G4endl;
246  G4cout << " cs (Geant4 internal unit)=" << xs << G4endl;
247  G4cout << " -> first cs value in EADL data file (iu) =" << (*pv)[0] << G4endl;
248  G4cout << " -> last cs value in EADL data file (iu) =" << (*pv)[n] << G4endl;
249  G4cout << "*********************************************************" << G4endl;
250  }
251 
252  return xs;
253 
254 }
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
size_t GetVectorLength() const
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
const G4int n
G4double Value(G4double theEnergy, size_t &lastidx) const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
const XML_Char const XML_Char * data
void G4LivermoreGammaConversionModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 73 of file G4LivermoreGammaConversionModel.cc.

References G4cout, G4endl, G4Material::GetElementVector(), G4MaterialCutsCouple::GetMaterial(), G4ProductionCutsTable::GetMaterialCutsCouple(), G4Material::GetNumberOfElements(), G4ProductionCutsTable::GetProductionCutsTable(), G4ProductionCutsTable::GetTableSize(), python.hepunit::GeV, eplot::material, and python.hepunit::MeV.

76 {
77  if (verboseLevel > 1)
78  {
79  G4cout << "Calling Initialise() of G4LivermoreGammaConversionModel."
80  << G4endl
81  << "Energy range: "
82  << LowEnergyLimit() / MeV << " MeV - "
83  << HighEnergyLimit() / GeV << " GeV"
84  << G4endl;
85  }
86 
87  if(IsMaster())
88  {
89 
90  // Initialise element selector
91 
92  InitialiseElementSelectors(particle, cuts);
93 
94  // Access to elements
95 
96  char* path = getenv("G4LEDATA");
97 
98  G4ProductionCutsTable* theCoupleTable =
100 
101  G4int numOfCouples = theCoupleTable->GetTableSize();
102 
103  for(G4int i=0; i<numOfCouples; ++i)
104  {
105  const G4Material* material =
106  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
107  const G4ElementVector* theElementVector = material->GetElementVector();
108  G4int nelm = material->GetNumberOfElements();
109 
110  for (G4int j=0; j<nelm; ++j)
111  {
112  G4int Z = (G4int)(*theElementVector)[j]->GetZ();
113  if(Z < 1) { Z = 1; }
114  else if(Z > maxZ) { Z = maxZ; }
115  if(!data[Z]) { ReadData(Z, path); }
116  }
117  }
118  }
119 
120  //
121 
122  if(isInitialised) { return; }
123  fParticleChange = GetParticleChangeForGamma();
124  isInitialised = true;
125 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
std::vector< G4Element * > G4ElementVector
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:135
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
string material
Definition: eplot.py:19
G4GLOB_DLL std::ostream G4cout
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const XML_Char const XML_Char * data
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
const G4Material * GetMaterial() const
void G4LivermoreGammaConversionModel::InitialiseForElement ( const G4ParticleDefinition ,
G4int  Z 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 468 of file G4LivermoreGammaConversionModel.cc.

References G4TemplateAutoLock< M, L, U >::unlock().

471 {
472  G4AutoLock l(&LivermoreGammaConversionModelMutex);
473  // G4cout << "G4LivermoreGammaConversionModel::InitialiseForElement Z= "
474  // << Z << G4endl;
475  if(!data[Z]) { ReadData(Z); }
476  l.unlock();
477 }
const XML_Char const XML_Char * data
void G4LivermoreGammaConversionModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 129 of file G4LivermoreGammaConversionModel.cc.

References G4VEmModel::GetElementSelectors().

131 {
133 }
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:760
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:768
G4double G4LivermoreGammaConversionModel::MinPrimaryEnergy ( const G4Material ,
const G4ParticleDefinition ,
G4double   
)
virtual

Reimplemented from G4VEmModel.

Definition at line 138 of file G4LivermoreGammaConversionModel.cc.

141 {
142  return lowEnergyLimit;
143 }
void G4LivermoreGammaConversionModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 258 of file G4LivermoreGammaConversionModel.cc.

References G4Electron::Electron(), python.hepunit::electron_mass_c2, fStopAndKill, G4cout, G4endl, G4Exp(), G4Log(), G4UniformRand, G4DynamicParticle::GetDefinition(), G4Element::GetfCoulomb(), G4Element::GetIonisation(), G4DynamicParticle::GetKineticEnergy(), G4IonisParamElm::GetlogZ3(), G4DynamicParticle::GetMomentumDirection(), G4IonisParamElm::GetZ3(), G4INCL::Math::max(), python.hepunit::MeV, G4INCL::Math::min(), G4Positron::Positron(), CLHEP::Hep3Vector::rotateUz(), and python.hepunit::twopi.

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

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