Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4eeToTwoGammaModel.cc
이 파일의 문서화 페이지로 가기
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 109177 2018-04-03 06:55:14Z 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),
91 {
93  fParticleChange = nullptr;
94 }
95 
96 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
97 
99 {}
100 
101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
102 
104  const G4DataVector&)
105 {
106  if(fParticleChange) { return; }
108 }
109 
110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
111 
112 G4double
114 {
115  // Calculates the cross section per electron of annihilation into two photons
116  // from the Heilter formula.
117 
118  G4double ekin = std::max(eV,kineticEnergy);
119 
120  G4double tau = ekin/electron_mass_c2;
121  G4double gam = tau + 1.0;
122  G4double gamma2= gam*gam;
123  G4double bg2 = tau * (tau+2.0);
124  G4double bg = sqrt(bg2);
125 
126  G4double cross = pi_rcl2*((gamma2+4*gam+1.)*G4Log(gam+bg) - (gam+3.)*bg)
127  / (bg2*(gam+1.));
128  return cross;
129 }
130 
131 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
132 
134  const G4ParticleDefinition*,
137 {
138  // Calculates the cross section per atom of annihilation into two photons
139 
140  G4double cross = Z*ComputeCrossSectionPerElectron(kineticEnergy);
141  return cross;
142 }
143 
144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145 
147  const G4Material* material,
148  const G4ParticleDefinition*,
151 {
152  // Calculates the cross section per volume of annihilation into two photons
153 
154  G4double eDensity = material->GetElectronDensity();
155  G4double cross = eDensity*ComputeCrossSectionPerElectron(kineticEnergy);
156  return cross;
157 }
158 
159 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
160 
161 // Polarisation of gamma according to M.H.L.Pryce and J.C.Ward,
162 // Nature 4065 (1947) 435.
163 
164 void G4eeToTwoGammaModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
165  const G4MaterialCutsCouple*,
166  const G4DynamicParticle* dp,
167  G4double,
168  G4double)
169 {
170  G4double posiKinEnergy = dp->GetKineticEnergy();
171  G4DynamicParticle *aGamma1, *aGamma2;
172 
173  CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
174 
175  // Case at rest
176  if(posiKinEnergy == 0.0) {
177  G4double cost = 2.*rndmEngine->flat()-1.;
178  G4double sint = sqrt((1. - cost)*(1. + cost));
179  G4double phi = twopi * rndmEngine->flat();
180  G4ThreeVector dir(sint*cos(phi), sint*sin(phi), cost);
181  phi = twopi * rndmEngine->flat();
182  G4double cosphi = cos(phi);
183  G4double sinphi = sin(phi);
184  G4ThreeVector pol(cosphi, sinphi, 0.0);
185  pol.rotateUz(dir);
186  aGamma1 = new G4DynamicParticle(theGamma, dir, electron_mass_c2);
187  aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
188  aGamma2 = new G4DynamicParticle(theGamma,-dir, electron_mass_c2);
189  pol.set(-sinphi, cosphi, 0.0);
190  pol.rotateUz(dir);
191  aGamma2->SetPolarization(pol.x(),pol.y(),pol.z());
192  /*
193  G4cout << "Annihilation at rest fly: e0= " << " dir= " << dir
194  << G4endl;
195  */
196  } else {
197 
198  G4ThreeVector posiDirection = dp->GetMomentumDirection();
199 
200  G4double tau = posiKinEnergy/electron_mass_c2;
201  G4double gam = tau + 1.0;
202  G4double tau2 = tau + 2.0;
203  G4double sqgrate = sqrt(tau/tau2)*0.5;
204  G4double sqg2m1 = sqrt(tau*tau2);
205 
206  // limits of the energy sampling
207  G4double epsilmin = 0.5 - sqgrate;
208  G4double epsilmax = 0.5 + sqgrate;
209  G4double epsilqot = epsilmax/epsilmin;
210 
211  //
212  // sample the energy rate of the created gammas
213  //
214  G4double epsil, greject;
215 
216  do {
217  epsil = epsilmin*G4Exp(G4Log(epsilqot)*rndmEngine->flat());
218  greject = 1. - epsil + (2.*gam*epsil-1.)/(epsil*tau2*tau2);
219  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
220  } while( greject < rndmEngine->flat());
221 
222  //
223  // scattered Gamma angles. ( Z - axis along the parent positron)
224  //
225 
226  G4double cost = (epsil*tau2-1.)/(epsil*sqg2m1);
227  if(std::abs(cost) > 1.0) {
228  G4cout << "### G4eeToTwoGammaModel WARNING cost= " << cost
229  << " positron Ekin(MeV)= " << posiKinEnergy
230  << " gamma epsil= " << epsil
231  << G4endl;
232  if(cost > 1.0) cost = 1.0;
233  else cost = -1.0;
234  }
235  G4double sint = sqrt((1.+cost)*(1.-cost));
236  G4double phi = twopi * rndmEngine->flat();
237 
238  //
239  // kinematic of the created pair
240  //
241 
242  G4double totalEnergy = posiKinEnergy + 2.0*electron_mass_c2;
243  G4double phot1Energy = epsil*totalEnergy;
244 
245  G4ThreeVector phot1Direction(sint*cos(phi), sint*sin(phi), cost);
246  phot1Direction.rotateUz(posiDirection);
247  aGamma1 = new G4DynamicParticle (theGamma,phot1Direction, phot1Energy);
248  phi = twopi * rndmEngine->flat();
249  G4double cosphi = cos(phi);
250  G4double sinphi = sin(phi);
251  G4ThreeVector pol(cosphi, sinphi, 0.0);
252  pol.rotateUz(phot1Direction);
253  aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
254 
255  G4double phot2Energy =(1.-epsil)*totalEnergy;
256  G4double posiP= sqrt(posiKinEnergy*(posiKinEnergy+2.*electron_mass_c2));
257  G4ThreeVector dir = posiDirection*posiP - phot1Direction*phot1Energy;
258  G4ThreeVector phot2Direction = dir.unit();
259 
260  // create G4DynamicParticle object for the particle2
261  aGamma2 = new G4DynamicParticle (theGamma, phot2Direction, phot2Energy);
262 
264  pol.set(-sinphi, cosphi, 0.0);
265  pol.rotateUz(phot1Direction);
266  cost = pol*phot2Direction;
267  pol -= cost*phot2Direction;
268  pol = pol.unit();
269  aGamma2->SetPolarization(pol.x(),pol.y(),pol.z());
270  /*
271  G4cout << "Annihilation on fly: e0= " << posiKinEnergy
272  << " m= " << electron_mass_c2
273  << " e1= " << phot1Energy
274  << " e2= " << phot2Energy << " dir= " << dir
275  << " -> " << phot1Direction << " "
276  << phot2Direction << G4endl;
277  */
278  }
279 
280  vdp->push_back(aGamma1);
281  vdp->push_back(aGamma2);
282 
283  // kill primary positron
286 }
287 
288 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
void set(double x, double y, double z)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void SetPolarization(G4double polX, G4double polY, G4double polZ)
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) override
G4ParticleChangeForGamma * fParticleChange
#define G4endl
Definition: G4ios.hh:61
const G4ThreeVector & GetMomentumDirection() const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void SetProposedKineticEnergy(G4double proposedKinEnergy)
double z() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
virtual G4double ComputeCrossSectionPerElectron(G4double kinEnergy)
static constexpr double classic_electr_radius
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
Float_t Z
double G4double
Definition: G4Types.hh:76
virtual double flat()=0
static constexpr double electron_mass_c2
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double eV
Definition: G4SIunits.hh:215
Hep3Vector unit() const
TDirectory * dir
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0., G4double maxEnergy=DBL_MAX) override
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
double x() const
static constexpr double pi
Definition: G4SIunits.hh:75
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:131
double y() const
G4ParticleDefinition * theGamma
void ProposeTrackStatus(G4TrackStatus status)
double flat()
Definition: G4AblaRandom.cc:48
G4double GetElectronDensity() const
Definition: G4Material.hh:218
G4eeToTwoGammaModel(const G4ParticleDefinition *p=0, const G4String &nam="eplus2gg")