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

Last change on this file since 846 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 6.0 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.6 2007/05/22 17:34:36 vnivanch Exp $
27// GEANT4 tag $Name:  $
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//
44// Class Description:
45//
46// -------------------------------------------------------------------
47//
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
50
51#include "G4PEEffectModel.hh"
52#include "G4Electron.hh"
53#include "G4Gamma.hh"
54#include "Randomize.hh"
55#include "G4DataVector.hh"
56#include "G4ParticleChangeForGamma.hh"
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
59
60using namespace std;
61
62G4PEEffectModel::G4PEEffectModel(const G4ParticleDefinition*,
63                                         const G4String& nam)
64  : G4VEmModel(nam),isInitialized(false)
65{
66  theGamma    = G4Gamma::Gamma();
67  theElectron = G4Electron::Electron();
68}
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
71
72G4PEEffectModel::~G4PEEffectModel()
73{
74}
75
76//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77
78void G4PEEffectModel::Initialise(const G4ParticleDefinition*,
79                                 const G4DataVector&)
80{
81 if (isInitialized) return;
82 if (pParticleChange)
83   fParticleChange =
84                  reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
85  else
86   fParticleChange = new G4ParticleChangeForGamma();
87
88 fminimalEnergy = 1.0*eV;
89}
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92
93void G4PEEffectModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
94                                        const G4MaterialCutsCouple* couple,
95                                        const G4DynamicParticle* aDynamicPhoton,
96                                        G4double,
97                                        G4double)
98{
99  const G4Material* aMaterial = couple->GetMaterial();
100
101  G4double energy = aDynamicPhoton->GetKineticEnergy();
102  G4ParticleMomentum PhotonDirection = aDynamicPhoton->GetMomentumDirection();
103
104  // select randomly one element constituing the material.
105  const G4Element* anElement = SelectRandomAtom(aMaterial,theGamma,energy);
106 
107  //
108  // Photo electron
109  //
110
111  // Select atomic shell
112  G4int nShells = anElement->GetNbOfAtomicShells();
113  G4int i  = 0; 
114  while ((i<nShells) && (energy<anElement->GetAtomicShell(i))) i++;
115
116  // no shell available
117  if (i == nShells) return;
118 
119  G4double bindingEnergy  = anElement->GetAtomicShell(i);
120  G4double ElecKineEnergy = energy - bindingEnergy;
121
122  if (ElecKineEnergy > fminimalEnergy)
123    {
124     // direction of the photo electron
125     //
126     G4double cosTeta = ElecCosThetaDistribution(ElecKineEnergy);
127     G4double sinTeta = sqrt(1.-cosTeta*cosTeta);
128     G4double Phi     = twopi * G4UniformRand();
129     G4double dirx = sinTeta*cos(Phi),diry = sinTeta*sin(Phi),dirz = cosTeta;
130     G4ThreeVector ElecDirection(dirx,diry,dirz);
131     ElecDirection.rotateUz(PhotonDirection);
132     //
133     G4DynamicParticle* aParticle = new G4DynamicParticle (
134                       theElectron,ElecDirection, ElecKineEnergy);
135     fvect->push_back(aParticle);
136    }
137   
138  fParticleChange->SetProposedKineticEnergy(0.);
139  fParticleChange->ProposeTrackStatus(fStopAndKill);
140  fParticleChange->ProposeLocalEnergyDeposit(bindingEnergy);
141}
142
143//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
144
145G4double G4PEEffectModel::ElecCosThetaDistribution(G4double kineEnergy)
146{
147  // Compute Theta distribution of the emitted electron, with respect to the
148  // incident Gamma.
149  // The Sauter-Gavrila distribution for the K-shell is used.
150  //
151  G4double costeta = 1.;
152  G4double gamma   = 1. + kineEnergy/electron_mass_c2;
153  if (gamma > 5.) return costeta;
154  G4double beta  = sqrt(gamma*gamma-1.)/gamma;
155  G4double b     = 0.5*gamma*(gamma-1.)*(gamma-2);
156
157  G4double rndm,term,greject,grejsup;
158  if (gamma < 2.) grejsup = gamma*gamma*(1.+b-beta*b);
159  else            grejsup = gamma*gamma*(1.+b+beta*b);
160
161  do { rndm = 1.-2*G4UniformRand();
162       costeta = (rndm+beta)/(rndm*beta+1.);
163       term = 1.-beta*costeta;
164       greject = (1.-costeta*costeta)*(1.+b*term)/(term*term);
165  } while(greject < G4UniformRand()*grejsup);
166
167  return costeta;
168}
169
170//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.