Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4eplusTo2GammaOKVIModel.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: G4eplusTo2GammaOKVIModel.cc 101193 2016-11-08 18:02:50Z vnivanch $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4eplusTo2GammaOKVIModel
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 
46 #include "G4PhysicalConstants.hh"
47 #include "G4SystemOfUnits.hh"
48 #include "G4EmParameters.hh"
49 #include "G4TrackStatus.hh"
50 #include "G4Electron.hh"
51 #include "G4Positron.hh"
52 #include "G4Gamma.hh"
53 #include "G4PhysicsVector.hh"
54 #include "G4PhysicsLogVector.hh"
55 #include "Randomize.hh"
57 #include "G4Log.hh"
58 #include "G4Exp.hh"
59 
60 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
61 
62 using namespace std;
63 
66 
68  const G4String& nam)
69  : G4VEmModel(nam),
71  energyTh(10*MeV)
72 {
74  fParticleChange = nullptr;
77 }
78 
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80 
82 {}
83 
84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85 
87  const G4DataVector& cuts)
88 {
90  f3GModel->Initialise(p, cuts);
91 
92  if(IsMaster()) {
93  if(!fCrossSection) {
94  G4double emin = 10*eV;
95  G4double emax = 100*TeV;
96  G4int nbins = 20*G4lrint(std::log10(emax/emin));
97  fCrossSection = new G4PhysicsLogVector(emin, emax, nbins);
98  f3GProbability= new G4PhysicsLogVector(emin, emax, nbins);
99  fCrossSection->SetSpline(true);
100  f3GProbability->SetSpline(true);
101  for(G4int i=0; i<= nbins; ++i) {
105  cs2 += cs3;
106  fCrossSection->PutValue(i, cs2);
107  f3GProbability->PutValue(i, cs3/cs2);
108  }
109  }
110  }
111  // here particle change is set for the triplet model
112  if(fParticleChange) { return; }
114 }
115 
116 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
117 
118 G4double
120 {
121  // Calculates the cross section per electron of annihilation into two photons
122  // from the Heilter formula.
123 
124  G4double ekin = std::max(eV,kineticEnergy);
125 
126  G4double tau = ekin/electron_mass_c2;
127  G4double gam = tau + 1.0;
128  G4double gamma2= gam*gam;
129  G4double bg2 = tau * (tau+2.0);
130  G4double bg = sqrt(bg2);
131 
132  G4double cross = pi_rcl2*((gamma2+4*gam+1.)*G4Log(gam+bg) - (gam+3.)*bg)
133  / (bg2*(gam+1.));
134 
135  return cross;
136 }
137 
138 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139 
141  const G4ParticleDefinition*,
144 {
145  // Calculates the cross section per atom of annihilation into two photons
146  G4double cross = Z*fCrossSection->Value(kineticEnergy);
147  return cross;
148 }
149 
150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
151 
153  const G4Material* material,
154  const G4ParticleDefinition*,
157 {
158  // Calculates the cross section per volume of annihilation into two photons
159  G4double eDensity = material->GetElectronDensity();
160  G4double cross = eDensity*fCrossSection->Value(kineticEnergy);
161  return cross;
162 }
163 
164 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
165 
166 // Polarisation of gamma according to M.H.L.Pryce and J.C.Ward,
167 // Nature 4065 (1947) 435.
168 
169 void
170 G4eplusTo2GammaOKVIModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
171  const G4MaterialCutsCouple* mcc,
172  const G4DynamicParticle* dp,
174 {
175  G4double posiKinEnergy = dp->GetKineticEnergy();
176  CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
177 
178  if(rndmEngine->flat() < f3GProbability->Value(posiKinEnergy)) {
179  f3GModel->SampleSecondaries(vdp, mcc, dp);
180  return;
181  }
182 
183  G4DynamicParticle *aGamma1, *aGamma2;
184 
185  // Case at rest
186  if(posiKinEnergy == 0.0) {
187  G4double cost = 2.*rndmEngine->flat()-1.;
188  G4double sint = sqrt((1. - cost)*(1. + cost));
189  G4double phi = twopi * rndmEngine->flat();
190  G4ThreeVector dir(sint*cos(phi), sint*sin(phi), cost);
191  phi = twopi * rndmEngine->flat();
192  G4double cosphi = cos(phi);
193  G4double sinphi = sin(phi);
194  G4ThreeVector pol(cosphi, sinphi, 0.0);
195  pol.rotateUz(dir);
196  aGamma1 = new G4DynamicParticle(theGamma, dir, electron_mass_c2);
197  aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
198  aGamma2 = new G4DynamicParticle(theGamma,-dir, electron_mass_c2);
199  pol.set(-sinphi, cosphi, 0.0);
200  pol.rotateUz(dir);
201  aGamma2->SetPolarization(pol.x(),pol.y(),pol.z());
202 
203  } else {
204 
205  G4ThreeVector posiDirection = dp->GetMomentumDirection();
206 
207  G4double tau = posiKinEnergy/electron_mass_c2;
208  G4double gam = tau + 1.0;
209  G4double tau2 = tau + 2.0;
210  G4double sqgrate = sqrt(tau/tau2)*0.5;
211  G4double sqg2m1 = sqrt(tau*tau2);
212 
213  // limits of the energy sampling
214  G4double epsilmin = 0.5 - sqgrate;
215  G4double epsilmax = 0.5 + sqgrate;
216  G4double epsilqot = epsilmax/epsilmin;
217 
218  //
219  // sample the energy rate of the created gammas
220  //
221  G4double epsil, greject;
222 
223  do {
224  epsil = epsilmin*G4Exp(G4Log(epsilqot)*rndmEngine->flat());
225  greject = 1. - epsil + (2.*gam*epsil-1.)/(epsil*tau2*tau2);
226  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
227  } while( greject < rndmEngine->flat());
228 
229  //
230  // scattered Gamma angles. ( Z - axis along the parent positron)
231  //
232 
233  G4double cost = (epsil*tau2-1.)/(epsil*sqg2m1);
234  if(std::abs(cost) > 1.0) {
235  G4cout << "### G4eplusTo2GammaOKVIModel WARNING cost= " << cost
236  << " positron Ekin(MeV)= " << posiKinEnergy
237  << " gamma epsil= " << epsil
238  << G4endl;
239  if(cost > 1.0) cost = 1.0;
240  else cost = -1.0;
241  }
242  G4double sint = sqrt((1.+cost)*(1.-cost));
243  G4double phi = twopi * rndmEngine->flat();
244 
245  //
246  // kinematic of the created pair
247  //
248 
249  G4double TotalAvailableEnergy = posiKinEnergy + 2.0*electron_mass_c2;
250  G4double phot1Energy = epsil*TotalAvailableEnergy;
251 
252  G4ThreeVector phot1Direction(sint*cos(phi), sint*sin(phi), cost);
253  phot1Direction.rotateUz(posiDirection);
254  aGamma1 = new G4DynamicParticle (theGamma,phot1Direction, phot1Energy);
255  phi = twopi * rndmEngine->flat();
256  G4double cosphi = cos(phi);
257  G4double sinphi = sin(phi);
258  G4ThreeVector pol(cosphi, sinphi, 0.0);
259  pol.rotateUz(phot1Direction);
260  aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
261 
262  G4double phot2Energy =(1.-epsil)*TotalAvailableEnergy;
263  G4double posiP= sqrt(posiKinEnergy*(posiKinEnergy+2.*electron_mass_c2));
264  G4ThreeVector dir = posiDirection*posiP - phot1Direction*phot1Energy;
265  G4ThreeVector phot2Direction = dir.unit();
266 
267  // create G4DynamicParticle object for the particle2
268  aGamma2 = new G4DynamicParticle (theGamma,phot2Direction, phot2Energy);
269 
271  pol.set(-sinphi, cosphi, 0.0);
272  pol.rotateUz(phot1Direction);
273  cost = pol*phot2Direction;
274  pol -= cost*phot2Direction;
275  pol = pol.unit();
276  aGamma2->SetPolarization(pol.x(),pol.y(),pol.z());
277  }
278  /*
279  G4cout << "Annihilation in fly: e0= " << posiKinEnergy
280  << " m= " << electron_mass_c2
281  << " e1= " << phot1Energy
282  << " e2= " << phot2Energy << " dir= " << dir
283  << " -> " << phot1Direction << " "
284  << phot2Direction << G4endl;
285  */
286 
287  vdp->push_back(aGamma1);
288  vdp->push_back(aGamma2);
289 
290  // kill primary positron
293 }
294 
295 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double Energy(size_t index) const
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
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) final
#define G4endl
Definition: G4ios.hh:61
const G4ThreeVector & GetMomentumDirection() const
const char * p
Definition: xmltok.h:285
G4double ComputeCrossSectionPerElectron(G4double kinEnergy)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
double z() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
static const G4double emax
G4double Value(G4double theEnergy, size_t &lastidx) const
G4ParticleChangeForGamma * fParticleChange
static constexpr double classic_electr_radius
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4eplusTo2GammaOKVIModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="eplus2ggOKVI")
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
Float_t Z
G4eplusTo3GammaOKVIModel * f3GModel
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) final
static constexpr double TeV
Definition: G4SIunits.hh:218
double G4double
Definition: G4Types.hh:76
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0., G4double maxEnergy=DBL_MAX) final
virtual double flat()=0
static constexpr double electron_mass_c2
G4double ComputeCrossSectionPerElectron(G4double kinEnergy)
static constexpr double twopi
Definition: G4SIunits.hh:76
void SetSpline(G4bool)
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
Hep3Vector unit() const
static G4PhysicsVector * f3GProbability
G4double LowestTripletEnergy() const
int G4lrint(double ad)
Definition: templates.hh:151
int G4int
Definition: G4Types.hh:78
TDirectory * dir
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
double x() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4bool IsMaster() const
Definition: G4VEmModel.hh:700
static G4PhysicsVector * fCrossSection
static constexpr double pi
Definition: G4SIunits.hh:75
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:131
double y() const
void ProposeTrackStatus(G4TrackStatus status)
double flat()
Definition: G4AblaRandom.cc:48
void SetTripletModel(G4VEmModel *)
Definition: G4VEmModel.hh:599
G4double GetElectronDensity() const
Definition: G4Material.hh:218
static G4EmParameters * Instance()
void PutValue(size_t index, G4double theValue)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double maxEnergy=DBL_MAX) final