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

#include <G4PEEffectFluoModel.hh>

Inheritance diagram for G4PEEffectFluoModel:
Collaboration diagram for G4PEEffectFluoModel:

Public Member Functions

 G4PEEffectFluoModel (const G4String &nam="PhotoElectric")
 
virtual ~G4PEEffectFluoModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double, G4double) override
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
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 MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
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 * > *)
 
virtual 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=nullptr)
 
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 SetFluctuationFlag (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)
 

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
 
G4bool lossFlucFlag
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 58 of file G4PEEffectFluoModel.hh.

Constructor & Destructor Documentation

G4PEEffectFluoModel::G4PEEffectFluoModel ( const G4String nam = "PhotoElectric")
explicit

Definition at line 66 of file G4PEEffectFluoModel.cc.

67  : G4VEmModel(nam)
68 {
69  theGamma = G4Gamma::Gamma();
70  theElectron = G4Electron::Electron();
71  fminimalEnergy = 1.0*eV;
72  SetDeexcitationFlag(true);
73  fParticleChange = nullptr;
74  fAtomDeexcitation = nullptr;
75 
76  fSandiaCof.resize(4,0.0);
77 
78  // default generator
80 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
static constexpr double eV
Definition: G4SIunits.hh:215
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:626
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:788

Here is the call graph for this function:

G4PEEffectFluoModel::~G4PEEffectFluoModel ( )
virtual

Definition at line 84 of file G4PEEffectFluoModel.cc.

85 {}

Member Function Documentation

G4double G4PEEffectFluoModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A,
G4double  ,
G4double   
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 107 of file G4PEEffectFluoModel.cc.

111 {
112  // This method may be used only if G4MaterialCutsCouple pointer
113  // has been set properly
115  ->GetSandiaTable()->GetSandiaCofPerAtom((G4int)Z, energy, fSandiaCof);
116 
117  G4double energy2 = energy*energy;
118  G4double energy3 = energy*energy2;
119  G4double energy4 = energy2*energy2;
120 
121  return fSandiaCof[0]/energy + fSandiaCof[1]/energy2 +
122  fSandiaCof[2]/energy3 + fSandiaCof[3]/energy4;
123 }
int G4int
Definition: G4Types.hh:78
G4SandiaTable * GetSandiaTable() const
Definition: G4Material.hh:229
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:454
G4double energy(const ThreeVector &p, const G4double m)
double G4double
Definition: G4Types.hh:76
void GetSandiaCofPerAtom(G4int Z, G4double energy, std::vector< G4double > &coeff) const
const G4Material * GetMaterial() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4PEEffectFluoModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition ,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 128 of file G4PEEffectFluoModel.cc.

132 {
133  // This method may be used only if G4MaterialCutsCouple pointer
134  // has been set properly
135  energy = std::max(energy, fMatEnergyTh[material->GetIndex()]);
136  const G4double* SandiaCof =
138 
139  G4double energy2 = energy*energy;
140  G4double energy3 = energy*energy2;
141  G4double energy4 = energy2*energy2;
142 
143  return SandiaCof[0]/energy + SandiaCof[1]/energy2 +
144  SandiaCof[2]/energy3 + SandiaCof[3]/energy4;
145 }
size_t GetIndex() const
Definition: G4Material.hh:262
G4double GetSandiaCofForMaterial(G4int, G4int) const
G4SandiaTable * GetSandiaTable() const
Definition: G4Material.hh:229
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double energy(const ThreeVector &p, const G4double m)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4PEEffectFluoModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
overridevirtual

Implements G4VEmModel.

Reimplemented in G4PolarizedPEEffectModel.

Definition at line 89 of file G4PEEffectFluoModel.cc.

91 {
92  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
93  if(nullptr == fParticleChange) { fParticleChange = GetParticleChangeForGamma(); }
94  size_t nmat = G4Material::GetNumberOfMaterials();
95  fMatEnergyTh.resize(nmat, 0.0);
96  for(size_t i=0; i<nmat; ++i) {
97  fMatEnergyTh[i] = (*(G4Material::GetMaterialTable()))[i]
98  ->GetSandiaTable()->GetSandiaCofForMaterial(0, 0);
99  //G4cout << "G4PEEffectFluoModel::Initialise Eth(eV)= "
100  // << fMatEnergyTh[i]/eV << G4endl;
101  }
102 }
static G4LossTableManager * Instance()
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:594
G4VAtomDeexcitation * AtomDeexcitation()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132

Here is the call graph for this function:

Here is the caller graph for this function:

void G4PEEffectFluoModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicPhoton,
G4double  tmin,
G4double  maxEnergy 
)
overridevirtual

Implements G4VEmModel.

Reimplemented in G4PolarizedPEEffectModel.

Definition at line 150 of file G4PEEffectFluoModel.cc.

155 {
156  SetCurrentCouple(couple);
157  const G4Material* aMaterial = couple->GetMaterial();
158 
159  G4double energy = aDynamicPhoton->GetKineticEnergy();
160 
161  // select randomly one element constituing the material.
162  const G4Element* anElement = SelectRandomAtom(aMaterial,theGamma,energy);
163 
164  //
165  // Photo electron
166  //
167 
168  // Select atomic shell
169  G4int nShells = anElement->GetNbOfAtomicShells();
170  G4int i = 0;
171  for(; i<nShells; ++i) {
172  /*
173  G4cout << "i= " << i << " E(eV)= " << energy/eV
174  << " Eb(eV)= " << anElement->GetAtomicShell(i)/eV
175  << " " << anElement->GetName()
176  << G4endl;
177  */
178  if(energy >= anElement->GetAtomicShell(i)) { break; }
179  }
180 
181  G4double edep = energy;
182 
183  // Normally one shell is available
184  if (i < nShells) {
185 
186  G4double bindingEnergy = anElement->GetAtomicShell(i);
187  edep = bindingEnergy;
188  G4double esec = 0.0;
189 
190  // sample deexcitation
191  //
192  if(fAtomDeexcitation) {
193  G4int index = couple->GetIndex();
194  if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
195  G4int Z = G4lrint(anElement->GetZ());
197  const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
198  G4double eshell = shell->BindingEnergy();
199  if(eshell > bindingEnergy && eshell <= energy) {
200  bindingEnergy = eshell;
201  edep = eshell;
202  }
203  G4int nbefore = fvect->size();
204  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
205  G4int nafter = fvect->size();
206  for (G4int j=nbefore; j<nafter; ++j) {
207  G4double e = ((*fvect)[j])->GetKineticEnergy();
208  if(esec + e > edep) {
209  // correct energy in order to have energy balance
210  e = edep - esec;
211  ((*fvect)[j])->SetKineticEnergy(e);
212  esec += e;
213  /*
214  G4cout << "### G4PEffectFluoModel Edep(eV)= " << edep/eV
215  << " Esec(eV)= " << esec/eV
216  << " E["<< j << "](eV)= " << e/eV
217  << " N= " << nafter
218  << " Z= " << Z << " shell= " << i
219  << " Ebind(keV)= " << bindingEnergy/keV
220  << " Eshell(keV)= " << shell->BindingEnergy()/keV
221  << G4endl;
222  */
223  // delete the rest of secondaries (should not happens)
224  for (G4int jj=nafter-1; jj>j; --jj) {
225  delete (*fvect)[jj];
226  fvect->pop_back();
227  }
228  break;
229  }
230  esec += e;
231  }
232  edep -= esec;
233  }
234  }
235  // create photo electron
236  //
237  G4double elecKineEnergy = energy - bindingEnergy;
238  if (elecKineEnergy > fminimalEnergy) {
239  G4DynamicParticle* aParticle = new G4DynamicParticle(theElectron,
240  GetAngularDistribution()->SampleDirection(aDynamicPhoton,
241  elecKineEnergy,
242  i, couple->GetMaterial()),
243  elecKineEnergy);
244  fvect->push_back(aParticle);
245  } else {
246  edep += elecKineEnergy;
247  elecKineEnergy = 0.0;
248  }
249  if(fabs(energy - elecKineEnergy - esec - edep) > eV) {
250  G4cout << "### G4PEffectFluoModel dE(eV)= "
251  << (energy - elecKineEnergy - esec - edep)/eV
252  << " shell= " << i
253  << " E(keV)= " << energy/keV
254  << " Ebind(keV)= " << bindingEnergy/keV
255  << " Ee(keV)= " << elecKineEnergy/keV
256  << " Esec(keV)= " << esec/keV
257  << " Edep(keV)= " << edep/keV
258  << G4endl;
259  }
260  }
261 
262  // kill primary photon
263  fParticleChange->SetProposedKineticEnergy(0.);
264  fParticleChange->ProposeTrackStatus(fStopAndKill);
265  if(edep > 0.0) {
266  fParticleChange->ProposeLocalEnergyDeposit(edep);
267  }
268 }
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
G4int GetNbOfAtomicShells() const
Definition: G4Element.hh:147
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:619
G4double BindingEnergy() const
int G4int
Definition: G4Types.hh:78
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4GLOB_DLL std::ostream G4cout
static constexpr double eV
Definition: G4SIunits.hh:215
int G4lrint(double ad)
Definition: templates.hh:163
G4double energy(const ThreeVector &p, const G4double m)
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:447
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
G4double GetAtomicShell(G4int index) const
Definition: G4Element.cc:367
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
G4double bindingEnergy(G4int A, G4int Z)
static constexpr double keV
Definition: G4SIunits.hh:216
G4AtomicShellEnumerator
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:544
const G4Material * GetMaterial() const

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: