source: trunk/source/processes/electromagnetic/standard/src/G4eeToTwoGammaModel.cc @ 1006

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

fichiers oublies

File size: 9.2 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: G4eeToTwoGammaModel.cc,v 1.14 2007/05/23 08:47:35 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name:   G4eeToTwoGammaModel
35//
36// Author:        Vladimir Ivanchenko on base of Michel Maire code
37//
38// Creation date: 02.08.2004
39//
40// Modifications:
41// 08-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
42// 18-04-05 Compute CrossSectionPerVolume (V.Ivanchenko)
43// 06-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
44// 29-06-06 Fix problem for zero energy incident positron (V.Ivanchenko)
45// 20-10-06 Add theGamma as a member (V.Ivanchenko)
46//
47//
48// Class Description:
49//
50// Implementation of e+ annihilation into 2 gamma
51//
52// The secondaries Gamma energies are sampled using the Heitler cross section.
53//
54// A modified version of the random number techniques of Butcher & Messel
55// is used (Nuc Phys 20(1960),15).
56//
57// GEANT4 internal units.
58//
59// Note 1: The initial electron is assumed free and at rest.
60//
61// Note 2: The annihilation processes producing one or more than two photons are
62//         ignored, as negligible compared to the two photons process.
63
64
65
66//
67// -------------------------------------------------------------------
68//
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
71
72#include "G4eeToTwoGammaModel.hh"
73#include "G4TrackStatus.hh"
74#include "G4Electron.hh"
75#include "G4Positron.hh"
76#include "G4Gamma.hh"
77#include "Randomize.hh"
78#include "G4ParticleChangeForGamma.hh"
79
80//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
81
82using namespace std;
83
84G4eeToTwoGammaModel::G4eeToTwoGammaModel(const G4ParticleDefinition*,
85                                         const G4String& nam)
86  : G4VEmModel(nam),
87    pi_rcl2(pi*classic_electr_radius*classic_electr_radius),
88    isInitialised(false)
89{
90  theGamma = G4Gamma::Gamma();
91}
92
93//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
94
95G4eeToTwoGammaModel::~G4eeToTwoGammaModel()
96{}
97
98//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
99
100void G4eeToTwoGammaModel::Initialise(const G4ParticleDefinition*,
101                                     const G4DataVector&)
102{
103  if(isInitialised) return;
104
105  if(pParticleChange)
106    fParticleChange =
107      reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
108  else
109    fParticleChange = new G4ParticleChangeForGamma();
110
111  isInitialised = true;
112}
113
114//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
115
116G4double G4eeToTwoGammaModel::ComputeCrossSectionPerElectron(
117                                       const G4ParticleDefinition*,
118                                       G4double kineticEnergy,
119                                       G4double, G4double)
120{
121  // Calculates the cross section per electron of annihilation into two photons
122  // from the Heilter formula.
123 
124  G4double tau   = kineticEnergy/electron_mass_c2;
125  G4double gam   = tau + 1.0;
126  G4double gamma2= gam*gam;
127  G4double bg2   = tau * (tau+2.0);
128  G4double bg    = sqrt(bg2);
129
130  G4double cross = pi_rcl2*((gamma2+4*gam+1.)*log(gam+bg) - (gam+3.)*bg)
131                 / (bg2*(gam+1.));
132  return cross; 
133}
134
135//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
136
137G4double G4eeToTwoGammaModel::ComputeCrossSectionPerAtom(
138                                    const G4ParticleDefinition* p,
139                                    G4double kineticEnergy, G4double Z,
140                                    G4double, G4double, G4double)
141{
142  // Calculates the cross section per atom of annihilation into two photons
143 
144  G4double cross = Z*ComputeCrossSectionPerElectron(p,kineticEnergy);
145  return cross; 
146}
147
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
149
150G4double G4eeToTwoGammaModel::CrossSectionPerVolume(
151                                        const G4Material* material,
152                                        const G4ParticleDefinition* p,
153                                              G4double kineticEnergy,
154                                              G4double, G4double)
155{
156  // Calculates the cross section per volume of annihilation into two photons
157 
158  G4double eDensity = material->GetElectronDensity();
159  G4double cross = eDensity*ComputeCrossSectionPerElectron(p,kineticEnergy);
160  return cross;
161}
162
163//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
164
165void G4eeToTwoGammaModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
166                                            const G4MaterialCutsCouple*,
167                                            const G4DynamicParticle* dp,
168                                            G4double,
169                                            G4double)
170{
171  G4double PositKinEnergy = dp->GetKineticEnergy();
172
173  // Case at rest
174  if(PositKinEnergy == 0.0) {
175    G4double cost = 2.*G4UniformRand()-1.;
176    G4double sint = sqrt((1. - cost)*(1. + cost));
177    G4double phi  = twopi * G4UniformRand();
178    G4ThreeVector dir (sint*cos(phi), sint*sin(phi), cost);
179    G4DynamicParticle* aGamma1 = new G4DynamicParticle(theGamma,
180                                                       dir, electron_mass_c2);
181    G4DynamicParticle* aGamma2 = new G4DynamicParticle(theGamma,
182                                                       -dir, electron_mass_c2);
183    vdp->push_back(aGamma1);
184    vdp->push_back(aGamma2);
185
186  } else {
187
188    G4ThreeVector PositDirection = dp->GetMomentumDirection();
189
190    G4double tau     = PositKinEnergy/electron_mass_c2;
191    G4double gam     = tau + 1.0;
192    G4double tau2    = tau + 2.0;
193    G4double sqgrate = sqrt(tau/tau2)*0.5;
194    G4double sqg2m1  = sqrt(tau*tau2);
195
196    // limits of the energy sampling
197    G4double epsilmin = 0.5 - sqgrate;
198    G4double epsilmax = 0.5 + sqgrate;
199    G4double epsilqot = epsilmax/epsilmin;
200
201    //
202    // sample the energy rate of the created gammas
203    //
204    G4double epsil, greject;
205
206    do {
207      epsil = epsilmin*pow(epsilqot,G4UniformRand());
208      greject = 1. - epsil + (2.*gam*epsil-1.)/(epsil*tau2*tau2);
209    } while( greject < G4UniformRand() );
210
211    //
212    // scattered Gamma angles. ( Z - axis along the parent positron)
213    //
214
215    G4double cost = (epsil*tau2-1.)/(epsil*sqg2m1);
216    if(std::abs(cost) > 1.0) {
217      G4cout << "### G4eeToTwoGammaModel WARNING cost= " << cost
218             << " positron Ekin(MeV)= " << PositKinEnergy
219             << " gamma epsil= " << epsil
220             << G4endl;
221      if(cost > 1.0) cost = 1.0;
222      else cost = -1.0; 
223    }
224    G4double sint = sqrt((1.+cost)*(1.-cost));
225    G4double phi  = twopi * G4UniformRand();
226
227    G4double dirx = sint*cos(phi) , diry = sint*sin(phi) , dirz = cost;
228
229    //
230    // kinematic of the created pair
231    //
232
233    G4double TotalAvailableEnergy = PositKinEnergy + 2.0*electron_mass_c2;
234    G4double Phot1Energy = epsil*TotalAvailableEnergy;
235
236    G4ThreeVector Phot1Direction (dirx, diry, dirz);
237    Phot1Direction.rotateUz(PositDirection);
238    G4DynamicParticle* aGamma1 = 
239      new G4DynamicParticle (theGamma,Phot1Direction, Phot1Energy);
240    vdp->push_back(aGamma1);
241
242    G4double Phot2Energy =(1.-epsil)*TotalAvailableEnergy;
243    G4double PositP= sqrt(PositKinEnergy*(PositKinEnergy+2.*electron_mass_c2));
244    G4ThreeVector dir = PositDirection*PositP - Phot1Direction*Phot1Energy;
245    G4ThreeVector Phot2Direction = dir.unit();
246
247    // create G4DynamicParticle object for the particle2
248    G4DynamicParticle* aGamma2= 
249      new G4DynamicParticle (theGamma,Phot2Direction, Phot2Energy);
250    vdp->push_back(aGamma2);
251    /*
252      G4cout << "Annihilation in fly: e0= " << PositKinEnergy
253      << " m= " << electron_mass_c2
254      << " e1= " << Phot1Energy
255      << " e2= " << Phot2Energy << " dir= " <<  dir
256      << " -> " << Phot1Direction << " "
257      << Phot2Direction << G4endl;
258    */
259  }
260  fParticleChange->SetProposedKineticEnergy(0.);
261  fParticleChange->ProposeTrackStatus(fStopAndKill);
262}
263
264//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
Note: See TracBrowser for help on using the repository browser.