Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4LowEnergyPhotoElectric Class Reference

#include <G4LowEnergyPhotoElectric.hh>

Inheritance diagram for G4LowEnergyPhotoElectric:
Collaboration diagram for G4LowEnergyPhotoElectric:

Public Member Functions

 G4LowEnergyPhotoElectric (const G4String &processName="LowEnPhotoElec")
 
 ~G4LowEnergyPhotoElectric ()
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &photon)
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
void SetCutForLowEnSecPhotons (G4double)
 
void SetCutForLowEnSecElectrons (G4double)
 
void ActivateAuger (G4bool)
 
void SetAngularGenerator (G4RDVPhotoElectricAngularDistribution *distribution)
 
void SetAngularGenerator (const G4String &name)
 
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 G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (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 ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
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 65 of file G4LowEnergyPhotoElectric.hh.

Constructor & Destructor Documentation

G4LowEnergyPhotoElectric::G4LowEnergyPhotoElectric ( const G4String processName = "LowEnPhotoElec")

Definition at line 97 of file G4LowEnergyPhotoElectric.cc.

98  : G4VDiscreteProcess(processName), lowEnergyLimit(250*eV), highEnergyLimit(100*GeV),
99  intrinsicLowEnergyLimit(10*eV),
100  intrinsicHighEnergyLimit(100*GeV),
101  cutForLowEnergySecondaryPhotons(250.*eV),
102  cutForLowEnergySecondaryElectrons(250.*eV)
103 {
104  if (lowEnergyLimit < intrinsicLowEnergyLimit ||
105  highEnergyLimit > intrinsicHighEnergyLimit)
106  {
107  G4Exception("G4LowEnergyPhotoElectric::G4LowEnergyPhotoElectric()",
108  "OutOfRange", FatalException,
109  "Energy limit outside intrinsic process validity range!");
110  }
111 
112  crossSectionHandler = new G4RDCrossSectionHandler();
113  shellCrossSectionHandler = new G4RDCrossSectionHandler();
114  meanFreePathTable = 0;
115  rangeTest = new G4RDRangeNoTest;
116  generatorName = "geant4.6.2";
117  ElectronAngularGenerator = new G4RDPhotoElectricAngularGeneratorSimple("GEANTSimpleGenerator"); // default generator
118 
119 
120  if (verboseLevel > 0)
121  {
122  G4cout << GetProcessName() << " is created " << G4endl
123  << "Energy range: "
124  << lowEnergyLimit / keV << " keV - "
125  << highEnergyLimit / GeV << " GeV"
126  << G4endl;
127  }
128 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4GLOB_DLL std::ostream G4cout
static constexpr double eV
Definition: G4SIunits.hh:215
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static constexpr double GeV
Definition: G4SIunits.hh:217
#define G4endl
Definition: G4ios.hh:61
static constexpr double keV
Definition: G4SIunits.hh:216

Here is the call graph for this function:

G4LowEnergyPhotoElectric::~G4LowEnergyPhotoElectric ( )

Definition at line 130 of file G4LowEnergyPhotoElectric.cc.

131 {
132  delete crossSectionHandler;
133  delete shellCrossSectionHandler;
134  delete meanFreePathTable;
135  delete rangeTest;
136  delete ElectronAngularGenerator;
137 }

Member Function Documentation

void G4LowEnergyPhotoElectric::ActivateAuger ( G4bool  val)

Definition at line 359 of file G4LowEnergyPhotoElectric.cc.

360 {
361  deexcitationManager.ActivateAugerElectronProduction(val);
362 }
void ActivateAugerElectronProduction(G4bool val)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4LowEnergyPhotoElectric::BuildPhysicsTable ( const G4ParticleDefinition photon)
virtual

Reimplemented from G4VProcess.

Definition at line 139 of file G4LowEnergyPhotoElectric.cc.

140 {
141 
142  crossSectionHandler->Clear();
143  G4String crossSectionFile = "phot/pe-cs-";
144  crossSectionHandler->LoadData(crossSectionFile);
145 
146  shellCrossSectionHandler->Clear();
147  G4String shellCrossSectionFile = "phot/pe-ss-cs-";
148  shellCrossSectionHandler->LoadShellData(shellCrossSectionFile);
149 
150  delete meanFreePathTable;
151  meanFreePathTable = crossSectionHandler->BuildMeanFreePathForMaterials();
152 }
void LoadShellData(const G4String &dataFile)
G4RDVEMDataSet * BuildMeanFreePathForMaterials(const G4DataVector *energyCuts=0)
void LoadData(const G4String &dataFile)

Here is the call graph for this function:

G4double G4LowEnergyPhotoElectric::DumpMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
inline

Definition at line 91 of file G4LowEnergyPhotoElectric.hh.

94  { 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:

G4double G4LowEnergyPhotoElectric::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
protectedvirtual

Implements G4VDiscreteProcess.

Definition at line 325 of file G4LowEnergyPhotoElectric.cc.

328 {
329  const G4DynamicParticle* photon = track.GetDynamicParticle();
330  G4double energy = photon->GetKineticEnergy();
331  G4Material* material = track.GetMaterial();
332  // size_t materialIndex = material->GetIndex();
333 
334  G4double meanFreePath = DBL_MAX;
335 
336  // if (energy > highEnergyLimit)
337  // meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex);
338  // else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX;
339  // else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex);
340 
341  G4double cross = shellCrossSectionHandler->ValueForMaterial(material,energy);
342  if(cross > 0.0) meanFreePath = 1.0/cross;
343 
344  return meanFreePath;
345 }
G4double ValueForMaterial(const G4Material *material, G4double e) const
G4double GetKineticEnergy() const
string material
Definition: eplot.py:19
G4double energy(const ThreeVector &p, const G4double m)
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:

G4bool G4LowEnergyPhotoElectric::IsApplicable ( const G4ParticleDefinition particle)
virtual

Reimplemented from G4VProcess.

Definition at line 320 of file G4LowEnergyPhotoElectric.cc.

321 {
322  return ( &particle == G4Gamma::Gamma() );
323 }
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86

Here is the call graph for this function:

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

Reimplemented from G4VDiscreteProcess.

Definition at line 154 of file G4LowEnergyPhotoElectric.cc.

156 {
157  // Fluorescence generated according to:
158  // J. Stepanek ,"A program to determine the radiation spectra due to a single atomic
159  // subshell ionisation by a particle or due to deexcitation or decay of radionuclides",
160  // Comp. Phys. Comm. 1206 pp 1-1-9 (1997)
161 
162  aParticleChange.Initialize(aTrack);
163 
164  const G4DynamicParticle* incidentPhoton = aTrack.GetDynamicParticle();
165  G4double photonEnergy = incidentPhoton->GetKineticEnergy();
166  if (photonEnergy <= lowEnergyLimit)
167  {
171  return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
172  }
173 
174  G4ThreeVector photonDirection = incidentPhoton->GetMomentumDirection(); // Returns the normalized direction of the momentum
175 
176  // Select randomly one element in the current material
177  const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
178  G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy);
179 
180  // Select the ionised shell in the current atom according to shell cross sections
181  size_t shellIndex = shellCrossSectionHandler->SelectRandomShell(Z,photonEnergy);
182 
183  // Retrieve the corresponding identifier and binding energy of the selected shell
185  const G4RDAtomicShell* shell = transitionManager->Shell(Z,shellIndex);
187  G4int shellId = shell->ShellId();
188 
189  // Create lists of pointers to DynamicParticles (photons and electrons)
190  // (Is the electron vector necessary? To be checked)
191  std::vector<G4DynamicParticle*>* photonVector = 0;
192  std::vector<G4DynamicParticle*> electronVector;
193 
194  G4double energyDeposit = 0.0;
195 
196  // Primary outcoming electron
197  G4double eKineticEnergy = photonEnergy - bindingEnergy;
198 
199  // There may be cases where the binding energy of the selected shell is > photon energy
200  // In such cases do not generate secondaries
201  if (eKineticEnergy > 0.)
202  {
203  // Generate the electron only if with large enough range w.r.t. cuts and safety
204  G4double safety = aStep.GetPostStepPoint()->GetSafety();
205 
206  if (rangeTest->Escape(G4Electron::Electron(),couple,eKineticEnergy,safety))
207  {
208 
209  // Calculate direction of the photoelectron
210  G4ThreeVector gammaPolarization = incidentPhoton->GetPolarization();
211  G4ThreeVector electronDirection = ElectronAngularGenerator->GetPhotoElectronDirection(photonDirection,eKineticEnergy,gammaPolarization,shellId);
212 
213  // The electron is created ...
215  electronDirection,
216  eKineticEnergy);
217  electronVector.push_back(electron);
218  }
219  else
220  {
221  energyDeposit += eKineticEnergy;
222  }
223  }
224  else
225  {
226  bindingEnergy = photonEnergy;
227  }
228 
229  G4int nElectrons = electronVector.size();
230  size_t nTotPhotons = 0;
231  G4int nPhotons=0;
232  const G4ProductionCutsTable* theCoupleTable=
234 
235  size_t index = couple->GetIndex();
236  G4double cutg = (*(theCoupleTable->GetEnergyCutsVector(0)))[index];
237  cutg = std::min(cutForLowEnergySecondaryPhotons,cutg);
238 
239  G4double cute = (*(theCoupleTable->GetEnergyCutsVector(1)))[index];
240  cute = std::min(cutForLowEnergySecondaryPhotons,cute);
241 
242  G4DynamicParticle* aPhoton;
243 
244  // Generation of fluorescence
245  // Data in EADL are available only for Z > 5
246  // Protection to avoid generating photons in the unphysical case of
247  // shell binding energy > photon energy
248  if (Z > 5 && (bindingEnergy > cutg || bindingEnergy > cute))
249  {
250  photonVector = deexcitationManager.GenerateParticles(Z,shellId);
251  nTotPhotons = photonVector->size();
252  for (size_t k=0; k<nTotPhotons; k++)
253  {
254  aPhoton = (*photonVector)[k];
255  if (aPhoton)
256  {
257  G4double itsCut = cutg;
258  if(aPhoton->GetDefinition() == G4Electron::Electron()) itsCut = cute;
259  G4double itsEnergy = aPhoton->GetKineticEnergy();
260 
261  if (itsEnergy > itsCut && itsEnergy <= bindingEnergy)
262  {
263  nPhotons++;
264  // Local energy deposit is given as the sum of the
265  // energies of incident photons minus the energies
266  // of the outcoming fluorescence photons
267  bindingEnergy -= itsEnergy;
268 
269  }
270  else
271  {
272  delete aPhoton;
273  (*photonVector)[k] = 0;
274  }
275  }
276  }
277  }
278 
279  energyDeposit += bindingEnergy;
280 
281  G4int nSecondaries = nElectrons + nPhotons;
283 
284  for (G4int l = 0; l<nElectrons; l++ )
285  {
286  aPhoton = electronVector[l];
287  if(aPhoton) {
288  aParticleChange.AddSecondary(aPhoton);
289  }
290  }
291  for ( size_t ll = 0; ll < nTotPhotons; ll++)
292  {
293  aPhoton = (*photonVector)[ll];
294  if(aPhoton) {
295  aParticleChange.AddSecondary(aPhoton);
296  }
297  }
298 
299  delete photonVector;
300 
301  if (energyDeposit < 0)
302  {
303  G4cout << "WARNING - "
304  << "G4LowEnergyPhotoElectric::PostStepDoIt - Negative energy deposit"
305  << G4endl;
306  energyDeposit = 0;
307  }
308 
309  // Kill the incident photon
312 
315 
316  // Reset NbOfInteractionLengthLeft and return aParticleChange
317  return G4VDiscreteProcess::PostStepDoIt( aTrack, aStep );
318 }
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
std::vector< G4DynamicParticle * > * GenerateParticles(G4int Z, G4int shellId)
G4double GetKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static G4RDAtomicTransitionManager * Instance()
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
G4double BindingEnergy() const
G4int SelectRandomShell(G4int Z, G4double e) const
virtual G4bool Escape(const G4ParticleDefinition *particle, const G4MaterialCutsCouple *couple, G4double energy, G4double safety) const =0
virtual void Initialize(const G4Track &)
static G4ProductionCutsTable * GetProductionCutsTable()
void SetNumberOfSecondaries(G4int totSecondaries)
G4StepPoint * GetPostStepPoint() const
G4int SelectRandomAtom(const G4MaterialCutsCouple *couple, G4double e) const
const G4ThreeVector & GetPolarization() const
void ProposeEnergy(G4double finalEnergy)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4RDAtomicShell * Shell(G4int Z, size_t shellIndex) const
G4double GetSafety() const
void AddSecondary(G4Track *aSecondary)
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
virtual G4ThreeVector GetPhotoElectronDirection(const G4ThreeVector &direction, const G4double kineticEnergy, const G4ThreeVector &polarization, const G4int shellID) const =0
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4double bindingEnergy(G4int A, G4int Z)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4int ShellId() const

Here is the call graph for this function:

void G4LowEnergyPhotoElectric::SetAngularGenerator ( G4RDVPhotoElectricAngularDistribution distribution)

Definition at line 364 of file G4LowEnergyPhotoElectric.cc.

365 {
366  ElectronAngularGenerator = distribution;
367  ElectronAngularGenerator->PrintGeneratorInformation();
368 }

Here is the call graph for this function:

void G4LowEnergyPhotoElectric::SetAngularGenerator ( const G4String name)

Definition at line 370 of file G4LowEnergyPhotoElectric.cc.

371 {
372  if (name == "default")
373  {
374  delete ElectronAngularGenerator;
375  ElectronAngularGenerator = new G4RDPhotoElectricAngularGeneratorSimple("GEANT4LowEnergySimpleGenerator");
376  generatorName = name;
377  }
378  else if (name == "standard")
379  {
380  delete ElectronAngularGenerator;
381  ElectronAngularGenerator = new G4RDPhotoElectricAngularGeneratorSauterGavrila("GEANT4SauterGavrilaGenerator");
382  generatorName = name;
383  }
384  else if (name == "polarized")
385  {
386  delete ElectronAngularGenerator;
387  ElectronAngularGenerator = new G4RDPhotoElectricAngularGeneratorPolarized("GEANT4LowEnergyPolarizedGenerator");
388  generatorName = name;
389  }
390  else
391  {
392  G4Exception("G4LowEnergyPhotoElectric::SetAngularGenerator()",
393  "InvalidSetup", FatalException,
394  "Generator does not exist!");
395  }
396 
397  ElectronAngularGenerator->PrintGeneratorInformation();
398 }
const XML_Char * name
Definition: expat.h:151
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

void G4LowEnergyPhotoElectric::SetCutForLowEnSecElectrons ( G4double  cut)

Definition at line 353 of file G4LowEnergyPhotoElectric.cc.

354 {
355  cutForLowEnergySecondaryElectrons = cut;
356  deexcitationManager.SetCutForAugerElectrons(cut);
357 }
void SetCutForAugerElectrons(G4double cut)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4LowEnergyPhotoElectric::SetCutForLowEnSecPhotons ( G4double  cut)

Definition at line 347 of file G4LowEnergyPhotoElectric.cc.

348 {
349  cutForLowEnergySecondaryPhotons = cut;
350  deexcitationManager.SetCutForSecondaryPhotons(cut);
351 }
void SetCutForSecondaryPhotons(G4double cut)

Here is the call graph for this function:

Here is the caller graph for this function:


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