source: trunk/source/processes/electromagnetic/standard/src/G4PEEffectFluoModel.cc @ 1350

Last change on this file since 1350 was 1350, checked in by garnier, 13 years ago

update to last version 4.9.4

File size: 7.4 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer                                           *
4// *                                                                  *
5// * The  Geant4 software  is  copyright of the Copyright Holders  of *
6// * the Geant4 Collaboration.  It is provided  under  the terms  and *
7// * conditions of the Geant4 Software License,  included in the file *
8// * LICENSE and available at  http://cern.ch/geant4/license .  These *
9// * include a list of copyright holders.                             *
10// *                                                                  *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work  make  any representation or  warranty, express or implied, *
14// * regarding  this  software system or assume any liability for its *
15// * use.  Please see the license in the file  LICENSE  and URL above *
16// * for the full disclaimer and the limitation of liability.         *
17// *                                                                  *
18// * This  code  implementation is the result of  the  scientific and *
19// * technical work of the GEANT4 collaboration.                      *
20// * By using,  copying,  modifying or  distributing the software (or *
21// * any work based  on the software)  you  agree  to acknowledge its *
22// * use  in  resulting  scientific  publications,  and indicate your *
23// * acceptance of all terms of the Geant4 Software license.          *
24// ********************************************************************
25//
26// $Id: G4PEEffectFluoModel.cc,v 1.1 2010/09/03 14:11:16 vnivanch Exp $
27// GEANT4 tag $Name: emstand-V09-03-24 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name:     G4PEEffectFluoModel
35//
36// Author:        Vladimir Ivanchenko on base of G4PEEffectModel
37//
38// Creation date: 13.06.2010
39//
40// Modifications:
41//
42// Class Description:
43// Implementation of the photo-electric effect with deexcitation
44//
45// -------------------------------------------------------------------
46//
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49
50#include "G4PEEffectFluoModel.hh"
51#include "G4Electron.hh"
52#include "G4Gamma.hh"
53#include "Randomize.hh"
54#include "G4DataVector.hh"
55#include "G4ParticleChangeForGamma.hh"
56#include "G4VAtomDeexcitation.hh"
57#include "G4LossTableManager.hh"
58
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
60
61using namespace std;
62
63G4PEEffectFluoModel::G4PEEffectFluoModel(const G4String& nam)
64  : G4VEmModel(nam),isInitialized(false)
65{
66  theGamma    = G4Gamma::Gamma();
67  theElectron = G4Electron::Electron();
68  fminimalEnergy = 1.0*eV;
69}
70
71//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
72
73G4PEEffectFluoModel::~G4PEEffectFluoModel()
74{}
75
76//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77
78void G4PEEffectFluoModel::Initialise(const G4ParticleDefinition*,
79                                     const G4DataVector&)
80{
81  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
82
83  if (isInitialized) return;
84  fParticleChange = GetParticleChangeForGamma();
85  isInitialized = true;
86}
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
89
90G4double
91G4PEEffectFluoModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
92                                                G4double energy,
93                                                G4double Z, G4double,
94                                                G4double, G4double)
95{
96  G4double* SandiaCof = G4SandiaTable::GetSandiaCofPerAtom((G4int)Z, energy);
97
98  G4double energy2 = energy*energy;
99  G4double energy3 = energy*energy2;
100  G4double energy4 = energy2*energy2;
101
102  return SandiaCof[0]/energy  + SandiaCof[1]/energy2 +
103    SandiaCof[2]/energy3 + SandiaCof[3]/energy4;
104}
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107
108G4double
109G4PEEffectFluoModel::CrossSectionPerVolume(const G4Material* material,
110                                           const G4ParticleDefinition*,
111                                           G4double energy,
112                                           G4double, G4double)
113{
114  G4double* SandiaCof = 
115    material->GetSandiaTable()->GetSandiaCofForMaterial(energy);
116                               
117  G4double energy2 = energy*energy;
118  G4double energy3 = energy*energy2;
119  G4double energy4 = energy2*energy2;
120         
121  return SandiaCof[0]/energy  + SandiaCof[1]/energy2 +
122    SandiaCof[2]/energy3 + SandiaCof[3]/energy4; 
123}
124
125//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
126
127void 
128G4PEEffectFluoModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
129                                       const G4MaterialCutsCouple* couple,
130                                       const G4DynamicParticle* aDynamicPhoton,
131                                       G4double,
132                                       G4double)
133{
134  const G4Material* aMaterial = couple->GetMaterial();
135
136  G4double energy = aDynamicPhoton->GetKineticEnergy();
137  G4ParticleMomentum PhotonDirection = aDynamicPhoton->GetMomentumDirection();
138
139  // select randomly one element constituing the material.
140  const G4Element* anElement = SelectRandomAtom(aMaterial,theGamma,energy);
141 
142  //
143  // Photo electron
144  //
145
146  // Select atomic shell
147  G4int nShells = anElement->GetNbOfAtomicShells();
148  G4int i  = 0; 
149  while ((i<nShells) && (energy<anElement->GetAtomicShell(i))) { ++i; }
150
151  // no shell available
152  if (i == nShells) { return; }
153 
154  G4double bindingEnergy  = anElement->GetAtomicShell(i);
155  G4double ElecKineEnergy = energy - bindingEnergy;
156
157  // create photo electron
158  //
159  if (ElecKineEnergy > fminimalEnergy) {
160    G4double cosTeta = ElecCosThetaDistribution(ElecKineEnergy);
161    G4double sinTeta = sqrt(1.-cosTeta*cosTeta);
162    G4double Phi     = twopi * G4UniformRand();
163    G4double dirx = sinTeta*cos(Phi),diry = sinTeta*sin(Phi),dirz = cosTeta;
164    G4ThreeVector ElecDirection(dirx,diry,dirz);
165    ElecDirection.rotateUz(PhotonDirection);
166    //
167    G4DynamicParticle* aParticle = new G4DynamicParticle (
168                       theElectron,ElecDirection, ElecKineEnergy);
169    fvect->push_back(aParticle);
170  } 
171
172  // sample deexcitation
173  //
174  if(fAtomDeexcitation) {
175    G4int Z = (G4int)anElement->GetZ();
176    G4AtomicShellEnumerator as = G4AtomicShellEnumerator(i);
177    const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);   
178    fAtomDeexcitation->GenerateParticles(fvect, shell, Z, couple->GetIndex());
179  }
180
181  // kill primary photon
182  fParticleChange->SetProposedKineticEnergy(0.);
183  fParticleChange->ProposeTrackStatus(fStopAndKill);
184  fParticleChange->ProposeLocalEnergyDeposit(bindingEnergy);
185}
186
187//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
188
189G4double G4PEEffectFluoModel::ElecCosThetaDistribution(G4double kineEnergy)
190{
191  // Compute Theta distribution of the emitted electron, with respect to the
192  // incident Gamma.
193  // The Sauter-Gavrila distribution for the K-shell is used.
194  //
195  G4double costeta = 1.;
196  G4double gamma   = 1. + kineEnergy/electron_mass_c2;
197  if (gamma > 5.) return costeta;
198  G4double beta  = sqrt(gamma*gamma-1.)/gamma;
199  G4double b     = 0.5*gamma*(gamma-1.)*(gamma-2);
200
201  G4double rndm,term,greject,grejsup;
202  if (gamma < 2.) grejsup = gamma*gamma*(1.+b-beta*b);
203  else            grejsup = gamma*gamma*(1.+b+beta*b);
204
205  do { rndm = 1.-2*G4UniformRand();
206       costeta = (rndm+beta)/(rndm*beta+1.);
207       term = 1.-beta*costeta;
208       greject = (1.-costeta*costeta)*(1.+b*term)/(term*term);
209  } while(greject < G4UniformRand()*grejsup);
210
211  return costeta;
212}
213
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.