Geant4  10.02.p03
G4LowEnergyGammaConversion Class Reference

#include <G4LowEnergyGammaConversion.hh>

Inheritance diagram for G4LowEnergyGammaConversion:
Collaboration diagram for G4LowEnergyGammaConversion:

Public Member Functions

 G4LowEnergyGammaConversion (const G4String &processName="LowEnConversion")
 
 ~G4LowEnergyGammaConversion ()
 
G4bool IsApplicable (const G4ParticleDefinition &photon)
 
void BuildPhysicsTable (const G4ParticleDefinition &photon)
 
G4VParticleChange * PostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
G4double DumpMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChange * AtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChange * AlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Private Member Functions

G4LowEnergyGammaConversionoperator= (const G4LowEnergyGammaConversion &right)
 
 G4LowEnergyGammaConversion (const G4LowEnergyGammaConversion &)
 
G4double ScreenFunction1 (G4double screenVariable)
 
G4double ScreenFunction2 (G4double screenVariable)
 

Private Attributes

G4double lowEnergyLimit
 
G4double highEnergyLimit
 
G4RDVEMDataSetmeanFreePathTable
 
G4RDVCrossSectionHandlercrossSectionHandler
 
G4RDVRangeTestrangeTest
 
const G4double intrinsicLowEnergyLimit
 
const G4double intrinsicHighEnergyLimit
 
const G4double smallEnergy
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChange * pParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 62 of file G4LowEnergyGammaConversion.hh.

Constructor & Destructor Documentation

◆ G4LowEnergyGammaConversion() [1/2]

G4LowEnergyGammaConversion::G4LowEnergyGammaConversion ( const G4String processName = "LowEnConversion")

Definition at line 78 of file G4LowEnergyGammaConversion.cc.

79  : G4VDiscreteProcess(processName),
80  lowEnergyLimit(1.022000*MeV),
81  highEnergyLimit(100*GeV),
82  intrinsicLowEnergyLimit(1.022000*MeV),
84  smallEnergy(2.*MeV)
85 
86 {
89  {
90  G4Exception("G4LowEnergyGammaConversion::G4LowEnergyGammaConversion()",
91  "OutOfRange", FatalException,
92  "Energy limit outside intrinsic process validity range!");
93  }
94 
95  // The following pointer is owned by G4DataHandler
96 
98  crossSectionHandler->Initialise(0,1.0220*MeV,100.*GeV,400);
100  rangeTest = new G4RDRangeTest;
101 
102  if (verboseLevel > 0)
103  {
104  G4cout << GetProcessName() << " is created " << G4endl
105  << "Energy range: "
106  << lowEnergyLimit / MeV << " MeV - "
107  << highEnergyLimit / GeV << " GeV"
108  << G4endl;
109  }
110 }
static const double MeV
Definition: G4SIunits.hh:211
G4int verboseLevel
Definition: G4VProcess.hh:368
void Initialise(G4RDVDataSetAlgorithm *interpolation=0, G4double minE=250 *CLHEP::eV, G4double maxE=100 *CLHEP::GeV, G4int numberOfBins=200, G4double unitE=CLHEP::MeV, G4double unitData=CLHEP::barn, G4int minZ=1, G4int maxZ=99)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
G4RDVCrossSectionHandler * crossSectionHandler
static const double GeV
Definition: G4SIunits.hh:214
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ~G4LowEnergyGammaConversion()

G4LowEnergyGammaConversion::~G4LowEnergyGammaConversion ( )

Definition at line 112 of file G4LowEnergyGammaConversion.cc.

113 {
114  delete meanFreePathTable;
115  delete crossSectionHandler;
116  delete rangeTest;
117 }
G4RDVCrossSectionHandler * crossSectionHandler

◆ G4LowEnergyGammaConversion() [2/2]

G4LowEnergyGammaConversion::G4LowEnergyGammaConversion ( const G4LowEnergyGammaConversion )
private

Member Function Documentation

◆ BuildPhysicsTable()

void G4LowEnergyGammaConversion::BuildPhysicsTable ( const G4ParticleDefinition photon)
virtual

Reimplemented from G4VProcess.

Definition at line 119 of file G4LowEnergyGammaConversion.cc.

120 {
121 
123  G4String crossSectionFile = "pair/pp-cs-";
124  crossSectionHandler->LoadData(crossSectionFile);
125 
126  delete meanFreePathTable;
128 }
G4RDVEMDataSet * BuildMeanFreePathForMaterials(const G4DataVector *energyCuts=0)
G4RDVCrossSectionHandler * crossSectionHandler
void LoadData(const G4String &dataFile)
Here is the call graph for this function:

◆ DumpMeanFreePath()

G4double G4LowEnergyGammaConversion::DumpMeanFreePath ( const G4Track &  aTrack,
G4double  previousStepSize,
G4ForceCondition *  condition 
)
inline

Definition at line 78 of file G4LowEnergyGammaConversion.hh.

81  { return GetMeanFreePath(aTrack, previousStepSize, condition); }
G4double condition(const G4ErrorSymMatrix &m)
G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
Here is the call graph for this function:

◆ GetMeanFreePath()

G4double G4LowEnergyGammaConversion::GetMeanFreePath ( const G4Track &  aTrack,
G4double  previousStepSize,
G4ForceCondition *  condition 
)
protectedvirtual

Implements G4VDiscreteProcess.

Definition at line 319 of file G4LowEnergyGammaConversion.cc.

322 {
323  const G4DynamicParticle* photon = track.GetDynamicParticle();
324  G4double energy = photon->GetKineticEnergy();
325  const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
326  size_t materialIndex = couple->GetIndex();
327 
328  G4double meanFreePath;
329  if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex);
330  else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX;
331  else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex);
332  return meanFreePath;
333 }
G4double GetKineticEnergy() const
double energy
Definition: plottest35.C:25
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
Here is the call graph for this function:
Here is the caller graph for this function:

◆ IsApplicable()

G4bool G4LowEnergyGammaConversion::IsApplicable ( const G4ParticleDefinition photon)
virtual

Reimplemented from G4VProcess.

Definition at line 314 of file G4LowEnergyGammaConversion.cc.

315 {
316  return ( &particle == G4Gamma::Gamma() );
317 }
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
Here is the call graph for this function:

◆ operator=()

G4LowEnergyGammaConversion& G4LowEnergyGammaConversion::operator= ( const G4LowEnergyGammaConversion right)
private
Here is the caller graph for this function:

◆ PostStepDoIt()

G4VParticleChange * G4LowEnergyGammaConversion::PostStepDoIt ( const G4Track &  aTrack,
const G4Step &  aStep 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 130 of file G4LowEnergyGammaConversion.cc.

132 {
133 // The energies of the e+ e- secondaries are sampled using the Bethe - Heitler
134 // cross sections with Coulomb correction. A modified version of the random
135 // number techniques of Butcher & Messel is used (Nuc Phys 20(1960),15).
136 
137 // Note 1 : Effects due to the breakdown of the Born approximation at low
138 // energy are ignored.
139 // Note 2 : The differential cross section implicitly takes account of
140 // pair creation in both nuclear and atomic electron fields. However triplet
141 // prodution is not generated.
142 
143  aParticleChange.Initialize(aTrack);
144 
145  const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
146 
147  const G4DynamicParticle* incidentPhoton = aTrack.GetDynamicParticle();
148  G4double photonEnergy = incidentPhoton->GetKineticEnergy();
149  G4ParticleMomentum photonDirection = incidentPhoton->GetMomentumDirection();
150 
151  G4double epsilon ;
152  G4double epsilon0 = electron_mass_c2 / photonEnergy ;
153 
154  // Do it fast if photon energy < 2. MeV
155  if (photonEnergy < smallEnergy )
156  {
157  epsilon = epsilon0 + (0.5 - epsilon0) * G4UniformRand();
158  }
159  else
160  {
161  // Select randomly one element in the current material
162  const G4Element* element = crossSectionHandler->SelectRandomElement(couple,photonEnergy);
163 
164  if (element == 0)
165  {
166  G4cout << "G4LowEnergyGammaConversion::PostStepDoIt - element = 0" << G4endl;
167  }
168  G4IonisParamElm* ionisation = element->GetIonisation();
169  if (ionisation == 0)
170  {
171  G4cout << "G4LowEnergyGammaConversion::PostStepDoIt - ionisation = 0" << G4endl;
172  }
173 
174  // Extract Coulomb factor for this Element
175  G4double fZ = 8. * (ionisation->GetlogZ3());
176  if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
177 
178  // Limits of the screening variable
179  G4double screenFactor = 136. * epsilon0 / (element->GetIonisation()->GetZ3()) ;
180  G4double screenMax = std::exp ((42.24 - fZ)/8.368) - 0.952 ;
181  G4double screenMin = std::min(4.*screenFactor,screenMax) ;
182 
183  // Limits of the energy sampling
184  G4double epsilon1 = 0.5 - 0.5 * std::sqrt(1. - screenMin / screenMax) ;
185  G4double epsilonMin = std::max(epsilon0,epsilon1);
186  G4double epsilonRange = 0.5 - epsilonMin ;
187 
188  // Sample the energy rate of the created electron (or positron)
189  G4double screen;
190  G4double gReject ;
191 
192  G4double f10 = ScreenFunction1(screenMin) - fZ;
193  G4double f20 = ScreenFunction2(screenMin) - fZ;
194  G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
195  G4double normF2 = std::max(1.5 * f20,0.);
196 
197  do {
198  if (normF1 / (normF1 + normF2) > G4UniformRand() )
199  {
200  epsilon = 0.5 - epsilonRange * std::pow(G4UniformRand(), 0.3333) ;
201  screen = screenFactor / (epsilon * (1. - epsilon));
202  gReject = (ScreenFunction1(screen) - fZ) / f10 ;
203  }
204  else
205  {
206  epsilon = epsilonMin + epsilonRange * G4UniformRand();
207  screen = screenFactor / (epsilon * (1 - epsilon));
208  gReject = (ScreenFunction2(screen) - fZ) / f20 ;
209  }
210  } while ( gReject < G4UniformRand() );
211 
212  } // End of epsilon sampling
213 
214  // Fix charges randomly
215 
216  G4double electronTotEnergy;
217  G4double positronTotEnergy;
218 
220  {
221  electronTotEnergy = (1. - epsilon) * photonEnergy;
222  positronTotEnergy = epsilon * photonEnergy;
223  }
224  else
225  {
226  positronTotEnergy = (1. - epsilon) * photonEnergy;
227  electronTotEnergy = epsilon * photonEnergy;
228  }
229 
230  // Scattered electron (positron) angles. ( Z - axis along the parent photon)
231  // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
232  // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
233 
234  G4double u;
235  const G4double a1 = 0.625;
236  G4double a2 = 3. * a1;
237  // G4double d = 27. ;
238 
239  // if (9. / (9. + d) > G4UniformRand())
240  if (0.25 > G4UniformRand())
241  {
242  u = - std::log(G4UniformRand() * G4UniformRand()) / a1 ;
243  }
244  else
245  {
246  u = - std::log(G4UniformRand() * G4UniformRand()) / a2 ;
247  }
248 
249  G4double thetaEle = u*electron_mass_c2/electronTotEnergy;
250  G4double thetaPos = u*electron_mass_c2/positronTotEnergy;
251  G4double phi = twopi * G4UniformRand();
252 
253  G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
254  G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
255 
256 
257  // Kinematics of the created pair:
258  // the electron and positron are assumed to have a symetric angular
259  // distribution with respect to the Z axis along the parent photon
260 
261  G4double localEnergyDeposit = 0. ;
262 
263  aParticleChange.SetNumberOfSecondaries(2) ;
264  G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
265 
266  // Generate the electron only if with large enough range w.r.t. cuts and safety
267 
268  G4double safety = aStep.GetPostStepPoint()->GetSafety();
269 
270  if (rangeTest->Escape(G4Electron::Electron(),couple,electronKineEnergy,safety))
271  {
272  G4ThreeVector electronDirection (dxEle, dyEle, dzEle);
273  electronDirection.rotateUz(photonDirection);
274 
276  electronDirection,
277  electronKineEnergy);
278  aParticleChange.AddSecondary(particle1) ;
279  }
280  else
281  {
282  localEnergyDeposit += electronKineEnergy ;
283  }
284 
285  // The e+ is always created (even with kinetic energy = 0) for further annihilation
286  G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
287 
288  // Is the local energy deposit correct, if the positron is always created?
289  if (! (rangeTest->Escape(G4Positron::Positron(),couple,positronKineEnergy,safety)))
290  {
291  localEnergyDeposit += positronKineEnergy ;
292  positronKineEnergy = 0. ;
293  }
294 
295  G4ThreeVector positronDirection (dxPos, dyPos, dzPos);
296  positronDirection.rotateUz(photonDirection);
297 
298  // Create G4DynamicParticle object for the particle2
300  positronDirection, positronKineEnergy);
301  aParticleChange.AddSecondary(particle2) ;
302 
303  aParticleChange.ProposeLocalEnergyDeposit(localEnergyDeposit) ;
304 
305  // Kill the incident photon
306  aParticleChange.ProposeMomentumDirection(0.,0.,0.) ;
307  aParticleChange.ProposeEnergy(0.) ;
308  aParticleChange.ProposeTrackStatus(fStopAndKill) ;
309 
310  // Reset NbOfInteractionLengthLeft and return aParticleChange
311  return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
312 }
static const double MeV
Definition: G4SIunits.hh:211
G4double ScreenFunction1(G4double screenVariable)
G4double GetlogZ3() const
static const G4double a1
G4double GetZ3() const
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4RDVCrossSectionHandler * crossSectionHandler
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
virtual G4bool Escape(const G4ParticleDefinition *particle, const G4MaterialCutsCouple *couple, G4double energy, G4double safety) const =0
static G4Positron * Positron()
Definition: G4Positron.cc:94
const G4ThreeVector & GetMomentumDirection() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double ScreenFunction2(G4double screenVariable)
static int shootBit()
const G4Element * SelectRandomElement(const G4MaterialCutsCouple *material, G4double e) const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
double epsilon(double density, double temperature)
static const G4double a2
Here is the call graph for this function:

◆ ScreenFunction1()

G4double G4LowEnergyGammaConversion::ScreenFunction1 ( G4double  screenVariable)
private

Definition at line 335 of file G4LowEnergyGammaConversion.cc.

336 {
337  // Compute the value of the screening function 3*phi1 - phi2
338 
339  G4double value;
340 
341  if (screenVariable > 1.)
342  value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
343  else
344  value = 42.392 - screenVariable * (7.796 - 1.961 * screenVariable);
345 
346  return value;
347 }
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ ScreenFunction2()

G4double G4LowEnergyGammaConversion::ScreenFunction2 ( G4double  screenVariable)
private

Definition at line 349 of file G4LowEnergyGammaConversion.cc.

350 {
351  // Compute the value of the screening function 1.5*phi1 - 0.5*phi2
352 
353  G4double value;
354 
355  if (screenVariable > 1.)
356  value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
357  else
358  value = 41.405 - screenVariable * (5.828 - 0.8945 * screenVariable);
359 
360  return value;
361 }
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

Member Data Documentation

◆ crossSectionHandler

G4RDVCrossSectionHandler* G4LowEnergyGammaConversion::crossSectionHandler
private

Definition at line 102 of file G4LowEnergyGammaConversion.hh.

◆ highEnergyLimit

G4double G4LowEnergyGammaConversion::highEnergyLimit
private

Definition at line 99 of file G4LowEnergyGammaConversion.hh.

◆ intrinsicHighEnergyLimit

const G4double G4LowEnergyGammaConversion::intrinsicHighEnergyLimit
private

Definition at line 107 of file G4LowEnergyGammaConversion.hh.

◆ intrinsicLowEnergyLimit

const G4double G4LowEnergyGammaConversion::intrinsicLowEnergyLimit
private

Definition at line 106 of file G4LowEnergyGammaConversion.hh.

◆ lowEnergyLimit

G4double G4LowEnergyGammaConversion::lowEnergyLimit
private

Definition at line 98 of file G4LowEnergyGammaConversion.hh.

◆ meanFreePathTable

G4RDVEMDataSet* G4LowEnergyGammaConversion::meanFreePathTable
private

Definition at line 101 of file G4LowEnergyGammaConversion.hh.

◆ rangeTest

G4RDVRangeTest* G4LowEnergyGammaConversion::rangeTest
private

Definition at line 104 of file G4LowEnergyGammaConversion.hh.

◆ smallEnergy

const G4double G4LowEnergyGammaConversion::smallEnergy
private

Definition at line 109 of file G4LowEnergyGammaConversion.hh.


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