source: trunk/source/processes/electromagnetic/standard/src/G4PEEffectModel.cc @ 961

Last change on this file since 961 was 961, checked in by garnier, 15 years ago

update processes

File size: 7.3 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: G4PEEffectModel.cc,v 1.7 2009/02/20 12:06:37 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name:     G4PEEffectModel
35//
36// Author:        Vladimir Ivanchenko on base of Michel Maire code
37//
38// Creation date: 21.03.2005
39//
40// Modifications:
41//
42// 04.12.05 : SetProposedKineticEnergy(0.) for the killed photon (mma)
43// 20.02.09 : Added initialisation of deexcitation flag and method
44//            CrossSectionPerVolume instead of mfp (V.Ivanchenko)
45//
46// Class Description:
47//
48// -------------------------------------------------------------------
49//
50//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
52
53#include "G4PEEffectModel.hh"
54#include "G4Electron.hh"
55#include "G4Gamma.hh"
56#include "Randomize.hh"
57#include "G4DataVector.hh"
58#include "G4ParticleChangeForGamma.hh"
59
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
61
62using namespace std;
63
64G4PEEffectModel::G4PEEffectModel(const G4ParticleDefinition*,
65                                         const G4String& nam)
66  : G4VEmModel(nam),isInitialized(false)
67{
68  theGamma    = G4Gamma::Gamma();
69  theElectron = G4Electron::Electron();
70  fminimalEnergy = 1.0*eV;
71}
72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74
75G4PEEffectModel::~G4PEEffectModel()
76{}
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79
80void G4PEEffectModel::Initialise(const G4ParticleDefinition*,
81                                 const G4DataVector&)
82{
83  // always false before the run
84  SetDeexcitationFlag(false);
85
86  if (isInitialized) return;
87  if (pParticleChange) {
88    fParticleChange = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
89  } else {
90    fParticleChange = new G4ParticleChangeForGamma();
91  }
92  isInitialized = true;
93}
94
95//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
96
97G4double G4PEEffectModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
98                                                     G4double energy,
99                                                     G4double Z, G4double,
100                                                     G4double, G4double)
101{
102  G4double* SandiaCof = G4SandiaTable::GetSandiaCofPerAtom((G4int)Z, energy);
103
104  G4double energy2 = energy*energy;
105  G4double energy3 = energy*energy2;
106  G4double energy4 = energy2*energy2;
107
108  return SandiaCof[0]/energy  + SandiaCof[1]/energy2 +
109    SandiaCof[2]/energy3 + SandiaCof[3]/energy4;
110}
111
112//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
113
114G4double G4PEEffectModel::CrossSectionPerVolume(const G4Material* material,
115                                                const G4ParticleDefinition*,
116                                                G4double energy,
117                                                G4double, G4double)
118{
119  G4double* SandiaCof = 
120    material->GetSandiaTable()->GetSandiaCofForMaterial(energy);
121                               
122  G4double energy2 = energy*energy;
123  G4double energy3 = energy*energy2;
124  G4double energy4 = energy2*energy2;
125         
126  return SandiaCof[0]/energy  + SandiaCof[1]/energy2 +
127    SandiaCof[2]/energy3 + SandiaCof[3]/energy4; 
128}
129
130//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
131
132void G4PEEffectModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
133                                        const G4MaterialCutsCouple* couple,
134                                        const G4DynamicParticle* aDynamicPhoton,
135                                        G4double,
136                                        G4double)
137{
138  const G4Material* aMaterial = couple->GetMaterial();
139
140  G4double energy = aDynamicPhoton->GetKineticEnergy();
141  G4ParticleMomentum PhotonDirection = aDynamicPhoton->GetMomentumDirection();
142
143  // select randomly one element constituing the material.
144  const G4Element* anElement = SelectRandomAtom(aMaterial,theGamma,energy);
145 
146  //
147  // Photo electron
148  //
149
150  // Select atomic shell
151  G4int nShells = anElement->GetNbOfAtomicShells();
152  G4int i  = 0; 
153  while ((i<nShells) && (energy<anElement->GetAtomicShell(i))) i++;
154
155  // no shell available
156  if (i == nShells) return;
157 
158  G4double bindingEnergy  = anElement->GetAtomicShell(i);
159  G4double ElecKineEnergy = energy - bindingEnergy;
160
161  if (ElecKineEnergy > fminimalEnergy)
162    {
163     // direction of the photo electron
164     //
165     G4double cosTeta = ElecCosThetaDistribution(ElecKineEnergy);
166     G4double sinTeta = sqrt(1.-cosTeta*cosTeta);
167     G4double Phi     = twopi * G4UniformRand();
168     G4double dirx = sinTeta*cos(Phi),diry = sinTeta*sin(Phi),dirz = cosTeta;
169     G4ThreeVector ElecDirection(dirx,diry,dirz);
170     ElecDirection.rotateUz(PhotonDirection);
171     //
172     G4DynamicParticle* aParticle = new G4DynamicParticle (
173                       theElectron,ElecDirection, ElecKineEnergy);
174     fvect->push_back(aParticle);
175    }
176   
177  fParticleChange->SetProposedKineticEnergy(0.);
178  fParticleChange->ProposeTrackStatus(fStopAndKill);
179  fParticleChange->ProposeLocalEnergyDeposit(bindingEnergy);
180}
181
182//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
183
184G4double G4PEEffectModel::ElecCosThetaDistribution(G4double kineEnergy)
185{
186  // Compute Theta distribution of the emitted electron, with respect to the
187  // incident Gamma.
188  // The Sauter-Gavrila distribution for the K-shell is used.
189  //
190  G4double costeta = 1.;
191  G4double gamma   = 1. + kineEnergy/electron_mass_c2;
192  if (gamma > 5.) return costeta;
193  G4double beta  = sqrt(gamma*gamma-1.)/gamma;
194  G4double b     = 0.5*gamma*(gamma-1.)*(gamma-2);
195
196  G4double rndm,term,greject,grejsup;
197  if (gamma < 2.) grejsup = gamma*gamma*(1.+b-beta*b);
198  else            grejsup = gamma*gamma*(1.+b+beta*b);
199
200  do { rndm = 1.-2*G4UniformRand();
201       costeta = (rndm+beta)/(rndm*beta+1.);
202       term = 1.-beta*costeta;
203       greject = (1.-costeta*costeta)*(1.+b*term)/(term*term);
204  } while(greject < G4UniformRand()*grejsup);
205
206  return costeta;
207}
208
209//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.