Geant4  10.01.p02
G4eeToTwoGammaModel.cc
Go to the documentation of this file.
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 74581 2013-10-15 12:03:25Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4eeToTwoGammaModel
34 //
35 // Author: Vladimir Ivanchenko on base of Michel Maire code
36 //
37 // Creation date: 02.08.2004
38 //
39 // Modifications:
40 // 08-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
41 // 18-04-05 Compute CrossSectionPerVolume (V.Ivanchenko)
42 // 06-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
43 // 29-06-06 Fix problem for zero energy incident positron (V.Ivanchenko)
44 // 20-10-06 Add theGamma as a member (V.Ivanchenko)
45 //
46 //
47 // Class Description:
48 //
49 // Implementation of e+ annihilation into 2 gamma
50 //
51 // The secondaries Gamma energies are sampled using the Heitler cross section.
52 //
53 // A modified version of the random number techniques of Butcher & Messel
54 // is used (Nuc Phys 20(1960),15).
55 //
56 // GEANT4 internal units.
57 //
58 // Note 1: The initial electron is assumed free and at rest.
59 //
60 // Note 2: The annihilation processes producing one or more than two photons are
61 // ignored, as negligible compared to the two photons process.
62 
63 
64 
65 //
66 // -------------------------------------------------------------------
67 //
68 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70 
71 #include "G4eeToTwoGammaModel.hh"
72 #include "G4PhysicalConstants.hh"
73 #include "G4SystemOfUnits.hh"
74 #include "G4TrackStatus.hh"
75 #include "G4Electron.hh"
76 #include "G4Positron.hh"
77 #include "G4Gamma.hh"
78 #include "Randomize.hh"
80 #include "G4Log.hh"
81 #include "G4Exp.hh"
82 
83 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
84 
85 using namespace std;
86 
88  const G4String& nam)
89  : G4VEmModel(nam),
90  pi_rcl2(pi*classic_electr_radius*classic_electr_radius),
91  isInitialised(false)
92 {
94  fParticleChange = 0;
95 }
96 
97 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
98 
100 {}
101 
102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
103 
105  const G4DataVector&)
106 {
107  if(isInitialised) { return; }
109  isInitialised = true;
110 }
111 
112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
113 
115  const G4ParticleDefinition*,
116  G4double kineticEnergy,
118 {
119  // Calculates the cross section per electron of annihilation into two photons
120  // from the Heilter formula.
121 
122  G4double ekin = std::max(eV,kineticEnergy);
123 
124  G4double tau = ekin/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.)*G4Log(gam+bg) - (gam+3.)*bg)
131  / (bg2*(gam+1.));
132  return cross;
133 }
134 
135 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
136 
138  const G4ParticleDefinition* p,
139  G4double kineticEnergy, G4double Z,
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 
151  const G4Material* material,
152  const G4ParticleDefinition* p,
153  G4double kineticEnergy,
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 
165 void G4eeToTwoGammaModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
166  const G4MaterialCutsCouple*,
167  const G4DynamicParticle* dp,
168  G4double,
169  G4double)
170 {
171  G4double PositKinEnergy = dp->GetKineticEnergy();
172  G4DynamicParticle *aGamma1, *aGamma2;
173 
174  // Case at rest
175  if(PositKinEnergy == 0.0) {
176  G4double cost = 2.*G4UniformRand()-1.;
177  G4double sint = sqrt((1. - cost)*(1. + cost));
178  G4double phi = twopi * G4UniformRand();
179  G4ThreeVector dir(sint*cos(phi), sint*sin(phi), cost);
180  phi = twopi * G4UniformRand();
181  G4ThreeVector pol(cos(phi), sin(phi), 0.0);
182  pol.rotateUz(dir);
183  aGamma1 = new G4DynamicParticle(theGamma, dir, electron_mass_c2);
184  aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
185  aGamma2 = new G4DynamicParticle(theGamma,-dir, electron_mass_c2);
186  aGamma1->SetPolarization(-pol.x(),-pol.y(),-pol.z());
187 
188  } else {
189 
190  G4ThreeVector PositDirection = dp->GetMomentumDirection();
191 
192  G4double tau = PositKinEnergy/electron_mass_c2;
193  G4double gam = tau + 1.0;
194  G4double tau2 = tau + 2.0;
195  G4double sqgrate = sqrt(tau/tau2)*0.5;
196  G4double sqg2m1 = sqrt(tau*tau2);
197 
198  // limits of the energy sampling
199  G4double epsilmin = 0.5 - sqgrate;
200  G4double epsilmax = 0.5 + sqgrate;
201  G4double epsilqot = epsilmax/epsilmin;
202 
203  //
204  // sample the energy rate of the created gammas
205  //
206  G4double epsil, greject;
207 
208  do {
209  epsil = epsilmin*G4Exp(G4Log(epsilqot)*G4UniformRand());
210  greject = 1. - epsil + (2.*gam*epsil-1.)/(epsil*tau2*tau2);
211  } while( greject < G4UniformRand() );
212 
213  //
214  // scattered Gamma angles. ( Z - axis along the parent positron)
215  //
216 
217  G4double cost = (epsil*tau2-1.)/(epsil*sqg2m1);
218  if(std::abs(cost) > 1.0) {
219  G4cout << "### G4eeToTwoGammaModel WARNING cost= " << cost
220  << " positron Ekin(MeV)= " << PositKinEnergy
221  << " gamma epsil= " << epsil
222  << G4endl;
223  if(cost > 1.0) cost = 1.0;
224  else cost = -1.0;
225  }
226  G4double sint = sqrt((1.+cost)*(1.-cost));
227  G4double phi = twopi * G4UniformRand();
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(sint*cos(phi), sint*sin(phi), cost);
237  Phot1Direction.rotateUz(PositDirection);
238  aGamma1 = new G4DynamicParticle (theGamma,Phot1Direction, Phot1Energy);
239  phi = twopi * G4UniformRand();
240  G4ThreeVector pol1(cos(phi), sin(phi), 0.0);
241  pol1.rotateUz(Phot1Direction);
242  aGamma1->SetPolarization(pol1.x(),pol1.y(),pol1.z());
243 
244  G4double Phot2Energy =(1.-epsil)*TotalAvailableEnergy;
245  G4double PositP= sqrt(PositKinEnergy*(PositKinEnergy+2.*electron_mass_c2));
246  G4ThreeVector dir = PositDirection*PositP - Phot1Direction*Phot1Energy;
247  G4ThreeVector Phot2Direction = dir.unit();
248 
249  // create G4DynamicParticle object for the particle2
250  aGamma2 = new G4DynamicParticle (theGamma,Phot2Direction, Phot2Energy);
251 
253  aGamma2->SetPolarization(-pol1.x(),-pol1.y(),-pol1.z());
254  }
255  /*
256  G4cout << "Annihilation in fly: e0= " << PositKinEnergy
257  << " m= " << electron_mass_c2
258  << " e1= " << Phot1Energy
259  << " e2= " << Phot2Energy << " dir= " << dir
260  << " -> " << Phot1Direction << " "
261  << Phot2Direction << G4endl;
262  */
263 
264  vdp->push_back(aGamma1);
265  vdp->push_back(aGamma2);
268 }
269 
270 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4ParticleChangeForGamma * fParticleChange
const G4double pi
#define G4UniformRand()
Definition: Randomize.hh:93
G4GLOB_DLL std::ostream G4cout
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kinEnergy, G4double cutEnergy=0., G4double maxEnergy=DBL_MAX)
G4double GetElectronDensity() const
Definition: G4Material.hh:217
const G4ThreeVector & GetMomentumDirection() const
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0., G4double maxEnergy=DBL_MAX)
void SetPolarization(G4double polX, G4double polY, G4double polZ)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static const double eV
Definition: G4SIunits.hh:194
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4eeToTwoGammaModel(const G4ParticleDefinition *p=0, const G4String &nam="eplus2gg")
G4ParticleDefinition * theGamma
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134