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

#include <G4LivermorePolarizedPhotoElectricModel.hh>

Inheritance diagram for G4LivermorePolarizedPhotoElectricModel:
G4VEmModel

Public Member Functions

 G4LivermorePolarizedPhotoElectricModel (const G4String &nam="LivermorePolarizedPhotoElectric")
 
virtual ~G4LivermorePolarizedPhotoElectricModel ()
 
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 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 CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, 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 *)
 
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
 

Protected Attributes

G4ParticleChangeForGammafParticleChange
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

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 *)
 

Detailed Description

Definition at line 41 of file G4LivermorePolarizedPhotoElectricModel.hh.

Constructor & Destructor Documentation

G4LivermorePolarizedPhotoElectricModel::G4LivermorePolarizedPhotoElectricModel ( const G4String nam = "LivermorePolarizedPhotoElectric")

Definition at line 48 of file G4LivermorePolarizedPhotoElectricModel.cc.

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

50  :G4VEmModel(nam),fParticleChange(0),
51  crossSectionHandler(0), shellCrossSectionHandler(0)
52 {
53  lowEnergyLimit = 10 * eV; // SI - Could be 10 eV ?
54  highEnergyLimit = 100 * GeV;
55 
56  verboseLevel= 0;
57  // Verbosity scale:
58  // 0 = nothing
59  // 1 = warning for energy non-conservation
60  // 2 = details of energy budget
61  // 3 = calculation of cross sections, file openings, sampling of atoms
62  // 4 = entering in methods
63 
64  theGamma = G4Gamma::Gamma();
65  theElectron = G4Electron::Electron();
66 
67  // use default generator
69 
70  // polarized generator needs fix
71  //SetAngularDistribution(new G4PhotoElectricAngularGeneratorPolarized());
72  SetDeexcitationFlag(true);
73  fAtomDeexcitation = 0;
74  fDeexcitationActive = false;
75 
76  if (verboseLevel > 1) {
77  G4cout << "Livermore Polarized PhotoElectric is constructed " << G4endl
78  << "Energy range: "
79  << lowEnergyLimit / keV << " keV - "
80  << highEnergyLimit / GeV << " GeV"
81  << G4endl;
82  }
83 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
G4GLOB_DLL std::ostream G4cout
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:585
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:739
G4LivermorePolarizedPhotoElectricModel::~G4LivermorePolarizedPhotoElectricModel ( )
virtual

Definition at line 87 of file G4LivermorePolarizedPhotoElectricModel.cc.

88 {
89  delete crossSectionHandler;
90  delete shellCrossSectionHandler;
91 }

Member Function Documentation

G4double G4LivermorePolarizedPhotoElectricModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 145 of file G4LivermorePolarizedPhotoElectricModel.cc.

References python.hepunit::barn, G4VCrossSectionHandler::FindValue(), G4cout, G4endl, G4lrint(), and python.hepunit::keV.

150 {
151  G4double cs = crossSectionHandler->FindValue(G4lrint(Z), GammaEnergy);
152 
153  if (verboseLevel > 3) {
154  G4cout << "G4LivermorePolarizedPhotoElectricModel::ComputeCrossSectionPerAtom()"
155  << G4endl;
156  G4cout << " E(keV)= " << GammaEnergy/keV << " Z= " << Z
157  << " CrossSection(barn)= "
158  << cs/barn << G4endl;
159  }
160  return cs;
161 }
G4double FindValue(G4int Z, G4double e) const
G4GLOB_DLL std::ostream G4cout
int G4lrint(double ad)
Definition: templates.hh:163
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void G4LivermorePolarizedPhotoElectricModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 96 of file G4LivermorePolarizedPhotoElectricModel.cc.

References G4LossTableManager::AtomDeexcitation(), G4VCrossSectionHandler::Clear(), fParticleChange, G4cout, G4endl, G4VEmModel::GetParticleChangeForGamma(), python.hepunit::GeV, G4VEmModel::HighEnergyLimit(), G4LossTableManager::Instance(), G4VAtomDeexcitation::IsFluoActive(), python.hepunit::keV, G4VCrossSectionHandler::LoadData(), G4VCrossSectionHandler::LoadShellData(), and G4VEmModel::LowEnergyLimit().

99 {
100  if (verboseLevel > 3) {
101  G4cout << "Calling G4LivermorePolarizedPhotoElectricModel::Initialise()" << G4endl;
102  }
103  if (crossSectionHandler)
104  {
105  crossSectionHandler->Clear();
106  delete crossSectionHandler;
107  }
108 
109  if (shellCrossSectionHandler)
110  {
111  shellCrossSectionHandler->Clear();
112  delete shellCrossSectionHandler;
113  }
114 
115 
116  // Reading of data files - all materials are read
117  crossSectionHandler = new G4CrossSectionHandler;
118  crossSectionHandler->Clear();
119  G4String crossSectionFile = "phot/pe-cs-";
120  crossSectionHandler->LoadData(crossSectionFile);
121 
122  shellCrossSectionHandler = new G4CrossSectionHandler();
123  shellCrossSectionHandler->Clear();
124  G4String shellCrossSectionFile = "phot/pe-ss-cs-";
125  shellCrossSectionHandler->LoadShellData(shellCrossSectionFile);
126 
128 
129  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
130  if(fAtomDeexcitation) {
131  fDeexcitationActive = fAtomDeexcitation->IsFluoActive();
132  }
133  if (verboseLevel > 1) {
134  G4cout << "Livermore Polarized PhotoElectric model is initialized "
135  << G4endl
136  << "Energy range: "
137  << LowEnergyLimit() / keV << " keV - "
138  << HighEnergyLimit() / GeV << " GeV"
139  << G4endl;
140  }
141 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
static G4LossTableManager * Instance()
G4bool IsFluoActive() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
G4GLOB_DLL std::ostream G4cout
void LoadShellData(const G4String &dataFile)
void LoadData(const G4String &dataFile)
#define G4endl
Definition: G4ios.hh:61
G4VAtomDeexcitation * AtomDeexcitation()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
void G4LivermorePolarizedPhotoElectricModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 165 of file G4LivermorePolarizedPhotoElectricModel.cc.

References G4InuclSpecialFunctions::bindingEnergy(), G4AtomicShell::BindingEnergy(), G4VAtomDeexcitation::CheckDeexcitationActiveRegion(), G4InuclParticleNames::electron, fParticleChange, fStopAndKill, G4cout, G4endl, G4lrint(), G4VAtomDeexcitation::GenerateParticles(), G4VEmModel::GetAngularDistribution(), G4VAtomDeexcitation::GetAtomicShell(), G4Element::GetAtomicShell(), G4MaterialCutsCouple::GetIndex(), G4DynamicParticle::GetKineticEnergy(), G4MaterialCutsCouple::GetMaterial(), G4Element::GetNbOfAtomicShells(), G4Element::GetZ(), G4INCL::Math::max(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4VParticleChange::ProposeTrackStatus(), G4VEmAngularDistribution::SampleDirection(), G4VEmModel::SelectRandomAtom(), G4VCrossSectionHandler::SelectRandomShell(), and G4ParticleChangeForGamma::SetProposedKineticEnergy().

171 {
172  if (verboseLevel > 3) {
173  G4cout << "Calling SampleSecondaries() of G4LivermorePolarizedPhotoElectricModel"
174  << G4endl;
175  }
176 
177  G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
178 
179  // kill incident photon
182 
183  // low-energy gamma is absorpted by this process
184 
185  if (photonEnergy <= lowEnergyLimit)
186  {
188  return;
189  }
190 
191  // Select randomly one element in the current material
192  //G4cout << "Select random atom Egamma(keV)= " << photonEnergy/keV << G4endl;
193  const G4Element* elm =
194  SelectRandomAtom(couple->GetMaterial(),theGamma,photonEnergy);
195  G4int Z = G4lrint(elm->GetZ());
196 
197  // Select the ionised shell in the current atom according to shell cross sections
198  //G4cout << "Select random shell Z= " << Z << G4endl;
199  size_t shellIndex = shellCrossSectionHandler->SelectRandomShell(Z,photonEnergy);
200  //G4cout << "Shell index= " << shellIndex
201  // << " nShells= " << elm->GetNbOfAtomicShells() << G4endl;
203  const G4AtomicShell* shell = 0;
204  if(fDeexcitationActive) {
206  shell = fAtomDeexcitation->GetAtomicShell(Z, as);
207  bindingEnergy = shell->BindingEnergy();
208  } else {
209  G4int nshells = elm->GetNbOfAtomicShells() - 1;
210  if(G4int(shellIndex) > nshells) { shellIndex = std::max(0, nshells); }
211  bindingEnergy = elm->GetAtomicShell(shellIndex);
212  }
213 
214  // There may be cases where the binding energy of the selected
215  // shell is > photon energy
216  // In such cases do not generate secondaries
217  if(photonEnergy < bindingEnergy) {
219  return;
220  }
221  //G4cout << "Z= " << Z << " shellIndex= " << shellIndex
222  // << " Ebind(keV)= " << bindingEnergy/keV << G4endl;
223 
224 
225  // Primary outcoming electron
226  G4double eKineticEnergy = photonEnergy - bindingEnergy;
227  G4double edep = bindingEnergy;
228 
229  // Calculate direction of the photoelectron
230  G4ThreeVector electronDirection =
231  GetAngularDistribution()->SampleDirection(aDynamicGamma,
232  eKineticEnergy,
233  shellIndex,
234  couple->GetMaterial());
235 
236  // The electron is created ...
237  G4DynamicParticle* electron = new G4DynamicParticle (theElectron,
238  electronDirection,
239  eKineticEnergy);
240  fvect->push_back(electron);
241 
242  // Sample deexcitation
243  if(fDeexcitationActive) {
244  G4int index = couple->GetIndex();
245  if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
246  size_t nbefore = fvect->size();
247 
248  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
249  size_t nafter = fvect->size();
250  if(nafter > nbefore) {
251  for (size_t j=nbefore; j<nafter; ++j) {
252  edep -= ((*fvect)[j])->GetKineticEnergy();
253  }
254  }
255  }
256  }
257 
258  // energy balance - excitation energy left
259  if(edep > 0.0) {
261  }
262 }
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
G4int GetNbOfAtomicShells() const
Definition: G4Element.hh:146
G4double GetKineticEnergy() const
G4int SelectRandomShell(G4int Z, G4double e) 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
G4GLOB_DLL std::ostream G4cout
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
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

Field Documentation

G4ParticleChangeForGamma* G4LivermorePolarizedPhotoElectricModel::fParticleChange
protected

Definition at line 69 of file G4LivermorePolarizedPhotoElectricModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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