Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4eplusTo3GammaOKVIModel.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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4eplusTo3GammaOKVIModel
33//
34// Authors: Andrei Alkin, Vladimir Ivanchenko, Omrame Kadri
35//
36// Creation date: 29.03.2018
37//
38// -------------------------------------------------------------------
39//
40//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
42
45#include "G4SystemOfUnits.hh"
46#include "G4EmParameters.hh"
47#include "G4TrackStatus.hh"
48#include "G4Electron.hh"
49#include "G4Positron.hh"
50#include "G4Gamma.hh"
51#include "Randomize.hh"
53#include "G4Log.hh"
54#include "G4Exp.hh"
55
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
57
58using namespace std;
59
61 const G4String& nam)
62 : G4VEmModel(nam), fDelta(0.001)
63{
64 theGamma = G4Gamma::Gamma();
65 fParticleChange = nullptr;
66}
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
69
71{}
72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74
76 const G4DataVector&)
77{
78 // here particle change is set for the triplet model
79 if(fParticleChange) { return; }
80 fParticleChange = GetParticleChangeForGamma();
81}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
84
85// (A.A.) F_{ijk} calculation method
87 G4double fr3, G4double kinEnergy)
88{
89 G4double ekin = std::max(eV,kinEnergy);
90 G4double tau = ekin/electron_mass_c2;
91 G4double gam = tau + 1.0;
92 G4double gamma2 = gam*gam;
93 G4double bg2 = tau * (tau+2.0);
94 G4double bg = sqrt(bg2);
95
96 G4double rho = (gamma2+4.*gam+1.)*G4Log(gam+bg)/(gamma2-1.)
97 - (gam+3.)/(sqrt(gam*gam - 1.)) + 1.;
98 G4double border;
99
100 if(ekin < 500*MeV) {
101 border = 1. - (electron_mass_c2)/(2*(ekin + electron_mass_c2));
102 } else {
103 border = 1. - (100*electron_mass_c2)/(2*(ekin + electron_mass_c2));
104 }
105
106 border = std::min(border, 0.9999);
107
108 if (fr1>border) { fr1 = border; }
109 if (fr2>border) { fr2 = border; }
110 if (fr3>border) { fr3 = border; }
111
112 G4double fr1s = fr1*fr1; // "s" for "squared"
113 G4double fr2s = fr2*fr2;
114 G4double fr3s = fr3*fr3;
115
116 G4double aa = (1.-fr1)*(1.-fr2);
117 G4double ab = fr3s + (fr1-fr2)*(fr1-fr2);
118 G4double add= ((1.-fr1)*(1.-fr1) + (1.-fr2)*(1.-fr2))/(fr3s*aa);
119
120 G4double fres = -rho*(1./fr1s + 1./fr2s)
121 + (ab/(2.*(fr1*fr2*aa)))*(G4Log(2.*gam*aa/(fr1*fr2)))
122 + (ab/(2.*fr1*fr2*(1-fr3)))*G4Log(2.*gam*(1.-fr3)/(fr1*fr2)) - add;
123
124 return fres;
125}
126
127//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
128
129// (A.A.) F_{ijk} calculation method
131 G4double fr3)
132{
133 G4double tau = 0.0;
134 G4double gam = tau + 1.0;
135 G4double gamma2 = gam*gam;
136 G4double bg2 = tau * (tau+2.0);
137 G4double bg = sqrt(bg2);
138
139 G4double rho = (gamma2+4.*gam+1.)*G4Log(gam+bg)/(gamma2-1.)
140 - (gam+3.)/(sqrt(gam*gam - 1.)) + 1.;
141 G4double border = 0.5;
142
143 if (fr1>border) { fr1 = border; }
144 if (fr2>border) { fr2 = border; }
145 if (fr3>border) { fr3 = border; }
146
147 G4double fr1s = fr1*fr1; // "s" for "squared"
148 G4double fr2s = fr2*fr2;
149 G4double fr3s = fr3*fr3;
150
151 G4double aa = (1.-fr1)*(1.-fr2);
152 G4double ab = fr3s + (fr1-fr2)*(fr1-fr2);
153 G4double add= ((1.-fr1)*(1.-fr1) + (1.-fr2)*(1.-fr2))/(fr3s*aa);
154
155 G4double fres = -rho*(1./fr1s + 1./fr2s)
156 + (ab/(2.*(fr1*fr2*aa)))*(G4Log(2.*gam*aa/(fr1*fr2)))
157 + (ab/(2.*fr1*fr2*(1-fr3)))*G4Log(2.*gam*(1.-fr3)/(fr1*fr2)) - add;
158
159 return fres;
160}
161
162//(A.A.) diff x-sections for maximum search and rejection
164 G4double fr2, G4double fr3, G4double kinEnergy)
165{
166 G4double ekin = std::max(eV,kinEnergy);
167 G4double tau = ekin/electron_mass_c2;
168 G4double gam = tau + 1.0;
169
170 G4double fsum = fr1*fr1*(ComputeF(fr1,fr2,fr3,ekin) +
171 ComputeF(fr3,fr1,fr2,ekin) +
172 ComputeF(fr2,fr3,fr1,ekin));
173
174 G4double dcross = fsum/((3*fr1*fr1*(gam+1.)));
175
176 return dcross;
177}
178
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
180
183{
184 // Calculates the cross section per electron of annihilation into 3 photons
185 // from the Heilter formula.
186
187 G4double ekin = std::max(eV,kinEnergy);
188 G4double tau = ekin/electron_mass_c2;
189 G4double gam = tau + 1.0;
190 G4double gamma2 = gam*gam;
191 G4double bg2 = tau * (tau+2.0);
192 G4double bg = sqrt(bg2);
193
194 G4double rho = (gamma2+4*gam+1.)*G4Log(gam+bg)/(gamma2-1.)
195 - (gam+3.)/(sqrt(gam*gam - 1.));
196
197 G4double cross = alpha_rcl2*(4.2 - (2.*G4Log(fDelta)+1.)*rho*rho)/(gam+1.);
198 return cross;
199}
200
201//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
202
205 G4double kineticEnergy, G4double Z,
207{
208 // Calculates the cross section per atom of annihilation into two photons
209
210
211
212 G4double cross = Z*ComputeCrossSectionPerElectron(kineticEnergy);
213 return cross;
214}
215
216//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
217
219 const G4Material* material,
221 G4double kineticEnergy,
223{
224 // Calculates the cross section per volume of annihilation into two photons
225
226 G4double eDensity = material->GetElectronDensity();
227 G4double cross = eDensity*ComputeCrossSectionPerElectron(kineticEnergy);
228 return cross;
229}
230
231//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
232
233// Polarisation of gamma according to M.H.L.Pryce and J.C.Ward,
234// Nature 4065 (1947) 435.
235
236void
237G4eplusTo3GammaOKVIModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
239 const G4DynamicParticle* dp,
241{
242
243 G4double posiKinEnergy = dp->GetKineticEnergy();
244 G4DynamicParticle *aGamma1, *aGamma2;
245 G4DynamicParticle* aGamma3 = nullptr;
246 G4double border;
247
248 if(posiKinEnergy < 500*MeV) {
249 border = 1. - (electron_mass_c2)/(2*(posiKinEnergy + electron_mass_c2));
250 } else {
251 border = 1. - (100*electron_mass_c2)/(2*(posiKinEnergy + electron_mass_c2));
252 }
253 border = std::min(border, 0.9999);
254
255 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
256
257 // Case at rest
258 if(posiKinEnergy == 0.0) {
259 G4double cost = 2.*rndmEngine->flat()-1.;
260 G4double sint = sqrt((1. - cost)*(1. + cost));
261 G4double phi = twopi * rndmEngine->flat();
262 G4ThreeVector dir(sint*cos(phi), sint*sin(phi), cost);
263 phi = twopi * rndmEngine->flat();
264 G4double cosphi = cos(phi);
265 G4double sinphi = sin(phi);
266 G4ThreeVector pol(cosphi, sinphi, 0.0);
267 pol.rotateUz(dir);
268 aGamma1 = new G4DynamicParticle(theGamma, dir, electron_mass_c2);
269 aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
270 aGamma2 = new G4DynamicParticle(theGamma,-dir, electron_mass_c2);
271 pol.set(-sinphi, cosphi, 0.0);
272 pol.rotateUz(dir);
273 aGamma2->SetPolarization(pol.x(),pol.y(),pol.z());
274
275 } else {
276
277 G4ThreeVector posiDirection = dp->GetMomentumDirection();
278
279 // (A.A.) LIMITS FOR 1st GAMMA
280 G4double xmin = 0.01;
281 G4double xmax = 0.667; // CHANGE to 3/2
282
283 G4double d1, d0, x1, x2, dmax, x2min;
284
285 // (A.A.) sampling of x1 x2 x3 (whole cycle of rejection)
286 do {
287 x1 = 1/((1/xmin) - ((1/xmin)-(1/xmax))*rndmEngine->flat());
288 dmax = ComputeFS(posiKinEnergy, x1,1.-x1,border);
289 x2min = 1.-x1;
290 x2 = 1 - rndmEngine->flat()*(1-x2min);
291 d1 = dmax*rndmEngine->flat();
292 d0 = ComputeFS(posiKinEnergy,x1,x2,2-x1-x2);
293 }
294 while(d0 < d1);
295
296 G4double x3 = 2 - x1 - x2;
297
298 //
299 // angles between Gammas
300 //
301
302 G4double psi13 = 2*asin(sqrt(std::abs((x1+x3-1)/(x1*x3))));
303 G4double psi12 = 2*asin(sqrt(std::abs((x1+x2-1)/(x1*x2))));
304
305 // sin^t
306
307 //G4double phi = twopi * rndmEngine->flat();
308 //G4double psi = acos(x3); // Angle of the plane
309
310 //
311 // kinematic of the created pair
312 //
313
314 G4double TotalAvailableEnergy = posiKinEnergy + 2.0*electron_mass_c2;
315
316 G4double phot1Energy = 0.5*x1*TotalAvailableEnergy;
317 G4double phot2Energy = 0.5*x2*TotalAvailableEnergy;
318 G4double phot3Energy = 0.5*x3*TotalAvailableEnergy;
319
320
321 // DIRECTIONS
322
323 // The azimuthal angles of ql and q3 with respect to some plane
324 // through the beam axis are generated at random.
325
326 G4ThreeVector phot1Direction(0, 0, 1);
327 G4ThreeVector phot2Direction(0, sin(psi12), cos(psi12));
328 G4ThreeVector phot3Direction(0, sin(psi13), cos(psi13));
329
330 phot1Direction.rotateUz(posiDirection);
331 phot2Direction.rotateUz(posiDirection);
332 phot3Direction.rotateUz(posiDirection);
333
334 aGamma1 = new G4DynamicParticle (theGamma,phot1Direction, phot1Energy);
335 aGamma2 = new G4DynamicParticle (theGamma,phot2Direction, phot2Energy);
336 aGamma3 = new G4DynamicParticle (theGamma,phot3Direction, phot3Energy);
337
338
339 //POLARIZATION - ???
340 /*
341
342
343 phi = twopi * rndmEngine->flat();
344 G4double cosphi = cos(phi);
345 G4double sinphi = sin(phi);
346 G4ThreeVector pol(cosphi, sinphi, 0.0);
347 pol.rotateUz(phot1Direction);
348 aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
349
350 G4double phot2Energy =(1.-epsil)*TotalAvailableEnergy;
351 G4double posiP= sqrt(posiKinEnergy*(posiKinEnergy+2.*electron_mass_c2));
352 G4ThreeVector dir = posiDirection*posiP - phot1Direction*phot1Energy;
353 G4ThreeVector phot2Direction = dir.unit();
354
355 // create G4DynamicParticle object for the particle2
356 aGamma2 = new G4DynamicParticle (theGamma,phot2Direction, phot2Energy);
357
358 //!!! likely problematic direction to be checked
359 pol.set(-sinphi, cosphi, 0.0);
360 pol.rotateUz(phot1Direction);
361 cost = pol*phot2Direction;
362 pol -= cost*phot2Direction;
363 pol = pol.unit();
364 aGamma2->SetPolarization(pol.x(),pol.y(),pol.z());
365
366 */
367
368 }
369 /*
370 G4cout << "Annihilation in fly: e0= " << posiKinEnergy
371 << " m= " << electron_mass_c2
372 << " e1= " << phot1Energy
373 << " e2= " << phot2Energy << " dir= " << dir
374 << " -> " << phot1Direction << " "
375 << phot2Direction << G4endl;
376 */
377
378 vdp->push_back(aGamma1);
379 vdp->push_back(aGamma2);
380 if(aGamma3 != nullptr) { vdp->push_back(aGamma3); }
381
382 // kill primary positron
383 fParticleChange->SetProposedKineticEnergy(0.0);
384 fParticleChange->ProposeTrackStatus(fStopAndKill);
385}
386
387//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double G4Log(G4double x)
Definition: G4Log.hh:226
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
double z() const
double x() const
double y() const
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
virtual double flat()=0
void SetPolarization(const G4ThreeVector &)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
G4double GetElectronDensity() const
Definition: G4Material.hh:215
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
void ProposeTrackStatus(G4TrackStatus status)
G4double ComputeF0(G4double fr1, G4double fr2, G4double fr3)
G4double ComputeCrossSectionPerElectron(G4double kinEnergy)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0., G4double maxEnergy=DBL_MAX) final
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double maxEnergy=DBL_MAX) final
G4double ComputeF(G4double fr1, G4double fr2, G4double fr3, G4double kinEnergy)
G4eplusTo3GammaOKVIModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="eplus3ggOKVI")
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) final
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) final
G4double ComputeFS(G4double fr1, G4double fr2, G4double fr3, G4double kinEnergy)