Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LivermoreGammaConversionModel Class Reference

#include <G4LivermoreGammaConversionModel.hh>

+ Inheritance diagram for G4LivermoreGammaConversionModel:

Public Member Functions

 G4LivermoreGammaConversionModel (const G4ParticleDefinition *p=0, const G4String &nam="LivermoreConversion")
 
virtual ~G4LivermoreGammaConversionModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
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 Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
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
 
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 ForceBuildTable (G4bool val)
 
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
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Detailed Description

Definition at line 40 of file G4LivermoreGammaConversionModel.hh.

Constructor & Destructor Documentation

◆ G4LivermoreGammaConversionModel()

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

Definition at line 40 of file G4LivermoreGammaConversionModel.cc.

42:G4VEmModel(nam),smallEnergy(2.*MeV),isInitialised(false),maxZ(99)
43{
44 fParticleChange = 0;
45
46 lowEnergyLimit = 2.0*electron_mass_c2;
47 data.resize(maxZ+1,0);
48
49 verboseLevel= 0;
50 // Verbosity scale for debugging purposes:
51 // 0 = nothing
52 // 1 = calculation of cross sections, file openings...
53 // 2 = entering in methods
54
55 if(verboseLevel > 0)
56 {
57 G4cout << "G4LivermoreGammaConversionModel is constructed " << G4endl;
58 }
59}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout

◆ ~G4LivermoreGammaConversionModel()

G4LivermoreGammaConversionModel::~G4LivermoreGammaConversionModel ( )
virtual

Definition at line 63 of file G4LivermoreGammaConversionModel.cc.

64{
65 for(G4int i=0; i<=maxZ; ++i) { delete data[i]; }
66}
int G4int
Definition: G4Types.hh:66

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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 182 of file G4LivermoreGammaConversionModel.cc.

186{
187 if (verboseLevel > 1)
188 {
189 G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreGammaConversionModel"
190 << G4endl;
191 }
192
193 if (GammaEnergy < lowEnergyLimit) { return 0.0; }
194
195 G4double xs = 0.0;
196
197 G4int intZ=G4int(Z);
198
199 if(intZ < 1 || intZ > maxZ) { return xs; }
200
201 G4LPhysicsFreeVector* pv = data[intZ];
202
203 // element was not initialised
204 if(!pv)
205 {
206 char* path = getenv("G4LEDATA");
207 ReadData(intZ, path);
208 pv = data[intZ];
209 if(!pv) { return xs; }
210 }
211 // x-section is taken from the table
212 xs = pv->Value(GammaEnergy);
213
214 if(verboseLevel > 0)
215 {
216 G4int n = pv->GetVectorLength() - 1;
217 G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)=" << GammaEnergy/MeV << G4endl;
218 G4cout << " cs (Geant4 internal unit)=" << xs << G4endl;
219 G4cout << " -> first cs value in EADL data file (iu) =" << (*pv)[0] << G4endl;
220 G4cout << " -> last cs value in EADL data file (iu) =" << (*pv)[n] << G4endl;
221 G4cout << "*********************************************************" << G4endl;
222 }
223
224 return xs;
225
226}
double G4double
Definition: G4Types.hh:64
G4double Value(G4double theEnergy)
size_t GetVectorLength() const

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 71 of file G4LivermoreGammaConversionModel.cc.

73{
74 if (verboseLevel > 1)
75 {
76 G4cout << "Calling Initialise() of G4LivermoreGammaConversionModel." << G4endl
77 << "Energy range: "
78 << LowEnergyLimit() / MeV << " MeV - "
79 << HighEnergyLimit() / GeV << " GeV"
80 << G4endl;
81 }
82
83 // Initialise element selector
84
85 InitialiseElementSelectors(particle, cuts);
86
87 // Access to elements
88
89 char* path = getenv("G4LEDATA");
90
91 G4ProductionCutsTable* theCoupleTable =
93 G4int numOfCouples = theCoupleTable->GetTableSize();
94
95 for(G4int i=0; i<numOfCouples; ++i)
96 {
97 const G4Material* material =
98 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
99 const G4ElementVector* theElementVector = material->GetElementVector();
100 G4int nelm = material->GetNumberOfElements();
101
102 for (G4int j=0; j<nelm; ++j)
103 {
104
105 G4int Z = (G4int)(*theElementVector)[j]->GetZ();
106 if(Z < 1) { Z = 1; }
107 else if(Z > maxZ) { Z = maxZ; }
108 if(!data[Z]) { ReadData(Z, path); }
109 }
110 }
111 //
112
113 if(isInitialised) { return; }
114 fParticleChange = GetParticleChangeForGamma();
115 isInitialised = true;
116}
std::vector< G4Element * > G4ElementVector
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:123

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 230 of file G4LivermoreGammaConversionModel.cc.

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

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