Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4eplusTo3GammaOKVIModel.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: G4eplusTo3GammaOKVIModel.cc 101193 2016-11-08 18:02:50Z vnivanch $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4eplusTo3GammaOKVIModel
34 //
35 // Author: Vladimir Ivanchenko and Omrame Kadri
36 //
37 // Creation date: 29.03.2018
38 //
39 // -------------------------------------------------------------------
40 //
41 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
43 
45 #include "G4PhysicalConstants.hh"
46 #include "G4SystemOfUnits.hh"
47 #include "G4EmParameters.hh"
48 #include "G4TrackStatus.hh"
49 #include "G4Electron.hh"
50 #include "G4Positron.hh"
51 #include "G4Gamma.hh"
52 #include "Randomize.hh"
54 #include "G4Log.hh"
55 #include "G4Exp.hh"
56 
57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
58 
59 using namespace std;
60 
62  const G4String& nam)
63  : G4VEmModel(nam),
65  energyTh(10*MeV)
66 {
68  fParticleChange = nullptr;
69 }
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
72 
74 {}
75 
76 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77 
79  const G4DataVector&)
80 {
82 }
83 
84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
85 
86 G4double
88 {
89  // Calculates the cross section per electron of annihilation into two photons
90  // from the Heilter formula.
91 
92  G4double ekin = std::max(eV,kineticEnergy);
93 
94  G4double tau = ekin/electron_mass_c2;
95  G4double gam = tau + 1.0;
96  G4double gamma2= gam*gam;
97  G4double bg2 = tau * (tau+2.0);
98  G4double bg = sqrt(bg2);
99 
100  G4double cross = pi_rcl2*((gamma2+4*gam+1.)*G4Log(gam+bg) - (gam+3.)*bg)
101  / (bg2*(gam+1.));
102  return cross;
103 }
104 
105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
106 
108  const G4ParticleDefinition*,
111 {
112  // Calculates the cross section per atom of annihilation into two photons
113 
114  G4double cross = Z*ComputeCrossSectionPerElectron(kineticEnergy);
115  return cross;
116 }
117 
118 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
119 
121  const G4Material* material,
122  const G4ParticleDefinition*,
125 {
126  // Calculates the cross section per volume of annihilation into two photons
127 
128  G4double eDensity = material->GetElectronDensity();
129  G4double cross = eDensity*ComputeCrossSectionPerElectron(kineticEnergy);
130  return cross;
131 }
132 
133 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134 
135 // Polarisation of gamma according to M.H.L.Pryce and J.C.Ward,
136 // Nature 4065 (1947) 435.
137 
138 void
139 G4eplusTo3GammaOKVIModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
140  const G4MaterialCutsCouple*,
141  const G4DynamicParticle* dp,
143 {
144  G4double posiKinEnergy = dp->GetKineticEnergy();
145  G4DynamicParticle *aGamma1, *aGamma2, *aGamma3;
146 
147  CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
148 
149  // Case at rest
150  if(posiKinEnergy == 0.0) {
151  G4double cost = 2.*rndmEngine->flat()-1.;
152  G4double sint = sqrt((1. - cost)*(1. + cost));
153  G4double phi = twopi * rndmEngine->flat();
154  G4ThreeVector dir(sint*cos(phi), sint*sin(phi), cost);
155  phi = twopi * rndmEngine->flat();
156  G4double cosphi = cos(phi);
157  G4double sinphi = sin(phi);
158  G4ThreeVector pol(cosphi, sinphi, 0.0);
159  pol.rotateUz(dir);
160  aGamma1 = new G4DynamicParticle(theGamma, dir, electron_mass_c2);
161  aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
162  aGamma2 = new G4DynamicParticle(theGamma,-dir, electron_mass_c2);
163  pol.set(-sinphi, cosphi, 0.0);
164  pol.rotateUz(dir);
165  aGamma2->SetPolarization(pol.x(),pol.y(),pol.z());
166 
167  } else {
168 
169  G4ThreeVector posiDirection = dp->GetMomentumDirection();
170 
171  G4double tau = posiKinEnergy/electron_mass_c2;
172  G4double gam = tau + 1.0;
173  G4double tau2 = tau + 2.0;
174  G4double sqgrate = sqrt(tau/tau2)*0.5;
175  G4double sqg2m1 = sqrt(tau*tau2);
176 
177  // limits of the energy sampling
178  G4double epsilmin = 0.5 - sqgrate;
179  G4double epsilmax = 0.5 + sqgrate;
180  G4double epsilqot = epsilmax/epsilmin;
181 
182  //
183  // sample the energy rate of the created gammas
184  //
185  G4double epsil, greject;
186 
187  do {
188  epsil = epsilmin*G4Exp(G4Log(epsilqot)*rndmEngine->flat());
189  greject = 1. - epsil + (2.*gam*epsil-1.)/(epsil*tau2*tau2);
190  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
191  } while( greject < rndmEngine->flat());
192 
193  //
194  // scattered Gamma angles. ( Z - axis along the parent positron)
195  //
196 
197  G4double cost = (epsil*tau2-1.)/(epsil*sqg2m1);
198  if(std::abs(cost) > 1.0) {
199  G4cout << "### G4eplusTo3GammaOKVIModel WARNING cost= " << cost
200  << " positron Ekin(MeV)= " << posiKinEnergy
201  << " gamma epsil= " << epsil
202  << G4endl;
203  if(cost > 1.0) cost = 1.0;
204  else cost = -1.0;
205  }
206  G4double sint = sqrt((1.+cost)*(1.-cost));
207  G4double phi = twopi * rndmEngine->flat();
208 
209  //
210  // kinematic of the created pair
211  //
212 
213  G4double TotalAvailableEnergy = posiKinEnergy + 2.0*electron_mass_c2;
214  G4double phot1Energy = epsil*TotalAvailableEnergy;
215 
216  G4ThreeVector phot1Direction(sint*cos(phi), sint*sin(phi), cost);
217  phot1Direction.rotateUz(posiDirection);
218  aGamma1 = new G4DynamicParticle (theGamma,phot1Direction, phot1Energy);
219  phi = twopi * rndmEngine->flat();
220  G4double cosphi = cos(phi);
221  G4double sinphi = sin(phi);
222  G4ThreeVector pol(cosphi, sinphi, 0.0);
223  pol.rotateUz(phot1Direction);
224  aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
225 
226  G4double phot2Energy =(1.-epsil)*TotalAvailableEnergy;
227  G4double posiP= sqrt(posiKinEnergy*(posiKinEnergy+2.*electron_mass_c2));
228  G4ThreeVector dir = posiDirection*posiP - phot1Direction*phot1Energy;
229  G4ThreeVector phot2Direction = dir.unit();
230 
231  // create G4DynamicParticle object for the particle2
232  aGamma2 = new G4DynamicParticle (theGamma,phot2Direction, phot2Energy);
233 
235  pol.set(-sinphi, cosphi, 0.0);
236  pol.rotateUz(phot1Direction);
237  cost = pol*phot2Direction;
238  pol -= cost*phot2Direction;
239  pol = pol.unit();
240  aGamma2->SetPolarization(pol.x(),pol.y(),pol.z());
241  }
242  /*
243  G4cout << "Annihilation in fly: e0= " << posiKinEnergy
244  << " m= " << electron_mass_c2
245  << " e1= " << phot1Energy
246  << " e2= " << phot2Energy << " dir= " << dir
247  << " -> " << phot1Direction << " "
248  << phot2Direction << G4endl;
249  */
250 
251  vdp->push_back(aGamma1);
252  vdp->push_back(aGamma2);
253  vdp->push_back(aGamma3);
254 
255  // kill primary positron
258 }
259 
260 //....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)
static constexpr double MeV
Definition: G4SIunits.hh:214
#define G4endl
Definition: G4ios.hh:61
const G4ThreeVector & GetMomentumDirection() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
double z() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
G4ParticleChangeForGamma * fParticleChange
static constexpr double classic_electr_radius
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
Float_t Z
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) final
double G4double
Definition: G4Types.hh:76
virtual double flat()=0
static constexpr double electron_mass_c2
G4double ComputeCrossSectionPerElectron(G4double kinEnergy)
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double eV
Definition: G4SIunits.hh:215
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double maxEnergy=DBL_MAX) final
G4ParticleDefinition * theGamma
Hep3Vector unit() const
G4double LowestTripletEnergy() const
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0., G4double maxEnergy=DBL_MAX) final
TDirectory * dir
G4eplusTo3GammaOKVIModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="eplus3ggOKVI")
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
double x() const
static constexpr double pi
Definition: G4SIunits.hh:75
double y() const
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) final
void ProposeTrackStatus(G4TrackStatus status)
double flat()
Definition: G4AblaRandom.cc:48
G4double GetElectronDensity() const
Definition: G4Material.hh:218
static G4EmParameters * Instance()