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

#include <G4PEEffectFluoModel.hh>

Inheritance diagram for G4PEEffectFluoModel:
G4VEmModel G4PolarizedPEEffectModel

Public Member Functions

 G4PEEffectFluoModel (const G4String &nam="PhotoElectric")
 
virtual ~G4PEEffectFluoModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double, G4double)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
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 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 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 *)
 
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 58 of file G4PEEffectFluoModel.hh.

Constructor & Destructor Documentation

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

Definition at line 66 of file G4PEEffectFluoModel.cc.

References G4Electron::Electron(), python.hepunit::eV, G4Gamma::Gamma(), G4VEmModel::SetAngularDistribution(), and G4VEmModel::SetDeexcitationFlag().

67  : G4VEmModel(nam)
68 {
69  theGamma = G4Gamma::Gamma();
70  theElectron = G4Electron::Electron();
71  fminimalEnergy = 1.0*eV;
72  SetDeexcitationFlag(true);
73  fParticleChange = 0;
74  fAtomDeexcitation = 0;
75 
76  fSandiaCof.resize(4,0.0);
77 
78  // default generator
80 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:585
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:739
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   
)
virtual

Reimplemented from G4VEmModel.

Definition at line 99 of file G4PEEffectFluoModel.cc.

References G4VEmModel::CurrentCouple(), energy(), G4MaterialCutsCouple::GetMaterial(), G4SandiaTable::GetSandiaCofPerAtom(), and G4Material::GetSandiaTable().

Referenced by G4AdjointPhotoElectricModel::AdjointCrossSectionPerAtom().

103 {
104  // This method may be used only if G4MaterialCutsCouple pointer
105  // has been set properly
106 
108  ->GetSandiaTable()->GetSandiaCofPerAtom((G4int)Z, energy, fSandiaCof);
109 
110  G4double energy2 = energy*energy;
111  G4double energy3 = energy*energy2;
112  G4double energy4 = energy2*energy2;
113 
114  return fSandiaCof[0]/energy + fSandiaCof[1]/energy2 +
115  fSandiaCof[2]/energy3 + fSandiaCof[3]/energy4;
116 }
int G4int
Definition: G4Types.hh:78
G4SandiaTable * GetSandiaTable() const
Definition: G4Material.hh:227
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:426
void GetSandiaCofPerAtom(G4int Z, G4double energy, std::vector< G4double > &coeff)
double G4double
Definition: G4Types.hh:76
const G4Material * GetMaterial() const
G4double G4PEEffectFluoModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition ,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 121 of file G4PEEffectFluoModel.cc.

References energy(), G4SandiaTable::GetSandiaCofForMaterial(), and G4Material::GetSandiaTable().

125 {
126  const G4double* SandiaCof =
128 
129  G4double energy2 = energy*energy;
130  G4double energy3 = energy*energy2;
131  G4double energy4 = energy2*energy2;
132 
133  return SandiaCof[0]/energy + SandiaCof[1]/energy2 +
134  SandiaCof[2]/energy3 + SandiaCof[3]/energy4;
135 }
G4SandiaTable * GetSandiaTable() const
Definition: G4Material.hh:227
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4double GetSandiaCofForMaterial(G4int, G4int)
double G4double
Definition: G4Types.hh:76
void G4PEEffectFluoModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Reimplemented in G4PolarizedPEEffectModel.

Definition at line 89 of file G4PEEffectFluoModel.cc.

References G4LossTableManager::AtomDeexcitation(), G4VEmModel::GetParticleChangeForGamma(), and G4LossTableManager::Instance().

Referenced by G4PolarizedPEEffectModel::Initialise().

91 {
92  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
93  if(!fParticleChange) { fParticleChange = GetParticleChangeForGamma(); }
94 }
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
void G4PEEffectFluoModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicPhoton,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Reimplemented in G4PolarizedPEEffectModel.

Definition at line 140 of file G4PEEffectFluoModel.cc.

References G4InuclSpecialFunctions::bindingEnergy(), G4AtomicShell::BindingEnergy(), G4VAtomDeexcitation::CheckDeexcitationActiveRegion(), energy(), python.hepunit::eV, fStopAndKill, G4cout, G4endl, G4lrint(), G4VAtomDeexcitation::GenerateParticles(), G4VEmModel::GetAngularDistribution(), G4VAtomDeexcitation::GetAtomicShell(), G4Element::GetAtomicShell(), G4MaterialCutsCouple::GetIndex(), G4DynamicParticle::GetKineticEnergy(), G4MaterialCutsCouple::GetMaterial(), G4Element::GetNbOfAtomicShells(), G4Element::GetZ(), python.hepunit::keV, G4VParticleChange::ProposeLocalEnergyDeposit(), G4VParticleChange::ProposeTrackStatus(), G4VEmModel::SelectRandomAtom(), G4VEmModel::SetCurrentCouple(), and G4ParticleChangeForGamma::SetProposedKineticEnergy().

Referenced by G4PolarizedPEEffectModel::SampleSecondaries().

145 {
146  SetCurrentCouple(couple);
147  const G4Material* aMaterial = couple->GetMaterial();
148 
149  G4double energy = aDynamicPhoton->GetKineticEnergy();
150 
151  // select randomly one element constituing the material.
152  const G4Element* anElement = SelectRandomAtom(aMaterial,theGamma,energy);
153 
154  //
155  // Photo electron
156  //
157 
158  // Select atomic shell
159  G4int nShells = anElement->GetNbOfAtomicShells();
160  G4int i = 0;
161  for(; i<nShells; ++i) {
162  /*
163  G4cout << "i= " << i << " E(eV)= " << energy/eV
164  << " Eb(eV)= " << anElement->GetAtomicShell(i)/eV
165  << " " << anElement->GetName()
166  << G4endl;
167  */
168  if(energy >= anElement->GetAtomicShell(i)) { break; }
169  }
170 
171  G4double edep = energy;
172 
173  // Normally one shell is available
174  if (i < nShells) {
175 
176  G4double bindingEnergy = anElement->GetAtomicShell(i);
177  edep = bindingEnergy;
178  G4double esec = 0.0;
179 
180  // sample deexcitation
181  //
182  if(fAtomDeexcitation) {
183  G4int index = couple->GetIndex();
184  if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
185  G4int Z = G4lrint(anElement->GetZ());
187  const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
188  G4double eshell = shell->BindingEnergy();
189  if(eshell > bindingEnergy && eshell <= energy) {
190  bindingEnergy = eshell;
191  edep = eshell;
192  }
193  size_t nbefore = fvect->size();
194  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
195  size_t nafter = fvect->size();
196  if(nafter > nbefore) {
197  for (size_t j=nbefore; j<nafter; ++j) {
198  G4double e = ((*fvect)[j])->GetKineticEnergy();
199  if(esec + e > edep) {
200  /*
201  G4cout << "### G4PEffectFluoModel Edep(eV)= " << edep/eV
202  << " Esec(eV)= " << esec/eV
203  << " E["<< j << "](eV)= " << e/eV
204  << " N= " << nafter
205  << " Z= " << Z << " shell= " << i
206  << " Ebind(keV)= " << bindingEnergy/keV
207  << " Eshell(keV)= " << shell->BindingEnergy()/keV
208  << G4endl;
209  */
210  for (size_t jj=j; jj<nafter; ++jj) { delete (*fvect)[jj]; }
211  for (size_t jj=j; jj<nafter; ++jj) { fvect->pop_back(); }
212  break;
213  }
214  esec += e;
215  }
216  }
217  edep -= esec;
218  }
219  }
220  // create photo electron
221  //
222  G4double elecKineEnergy = energy - bindingEnergy;
223  if (elecKineEnergy > fminimalEnergy) {
224  G4DynamicParticle* aParticle = new G4DynamicParticle(theElectron,
225  GetAngularDistribution()->SampleDirection(aDynamicPhoton,
226  elecKineEnergy,
227  i, couple->GetMaterial()),
228  elecKineEnergy);
229  fvect->push_back(aParticle);
230  } else {
231  edep += elecKineEnergy;
232  elecKineEnergy = 0.0;
233  }
234  if(fabs(energy - elecKineEnergy - esec - edep) > eV) {
235  G4cout << "### G4PEffectFluoModel dE(eV)= "
236  << (energy - elecKineEnergy - esec - edep)/eV
237  << " shell= " << i
238  << " E(keV)= " << energy/keV
239  << " Ebind(keV)= " << bindingEnergy/keV
240  << " Ee(keV)= " << elecKineEnergy/keV
241  << " Esec(keV)= " << esec/keV
242  << " Edep(keV)= " << edep/keV
243  << G4endl;
244  }
245  }
246 
247  // kill primary photon
248  fParticleChange->SetProposedKineticEnergy(0.);
249  fParticleChange->ProposeTrackStatus(fStopAndKill);
250  if(edep > 0.0) {
251  fParticleChange->ProposeLocalEnergyDeposit(edep);
252  }
253 }
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
G4int GetNbOfAtomicShells() const
Definition: G4Element.hh:146
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:578
G4double BindingEnergy() const
int G4int
Definition: G4Types.hh:78
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4GLOB_DLL std::ostream G4cout
int G4lrint(double ad)
Definition: templates.hh:163
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:419
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
G4double GetAtomicShell(G4int index) const
Definition: G4Element.cc:365
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)
G4AtomicShellEnumerator
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:510
const G4Material * GetMaterial() const

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