Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RToEConvForPositron.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// G4RToEConvForPositron class implementation
27//
28// Author: H.Kurashige, 05 October 2002 - First implementation
29// --------------------------------------------------------------------
30
33#include "G4ParticleTable.hh"
34
36#include "G4SystemOfUnits.hh"
37#include "G4Pow.hh"
38#include "G4Log.hh"
39#include "G4Exp.hh"
40
41// --------------------------------------------------------------------
44{
46 if (theParticle == nullptr)
47 {
48#ifdef G4VERBOSE
49 if (GetVerboseLevel()>0)
50 {
51 G4cout << "G4RToEConvForPositron::G4RToEConvForPositron() - ";
52 G4cout << "Positron is not defined !!" << G4endl;
53 }
54#endif
55 }
56 else
57 {
59 }
60}
61
62// --------------------------------------------------------------------
64{}
65
66// --------------------------------------------------------------------
68 const G4double kinEnergy)
69{
70 const G4double cbr1=0.02, cbr2=-5.7e-5, cbr3=1., cbr4=0.072;
71 const G4double Tlow=10.*CLHEP::keV, Thigh=1.*CLHEP::GeV;
72 const G4double taul = Tlow/CLHEP::electron_mass_c2;
73 const G4double logtaul = G4Log(taul);
74 const G4double taul12 = std::sqrt(taul);
75 const G4double bremfactor = 0.1;
76
78 G4double ionpot =
79 1.6e-5*CLHEP::MeV*G4Exp(0.9*Zlog)/CLHEP::electron_mass_c2;
80 G4double ionpotlog = G4Log(ionpot);
81
82 G4double tau = kinEnergy/CLHEP::electron_mass_c2;
83 G4double dEdx = 0.0;
84
85
86 if(tau<taul)
87 {
88 G4double t1 = taul+1.;
89 G4double t2 = taul+2.;
90 G4double tsq = taul*taul;
91 G4double beta2 = taul*t2/(t1*t1);
92 G4double f = 2.*logtaul -
93 (6.*taul+1.5*tsq-taul*(1.-tsq/3.)/t2
94 -tsq*(0.5-tsq/12.)/(t2*t2))/(t1*t1);
95 dEdx = (G4Log(2.*taul+4.)-2.*ionpotlog+f)/beta2;
96 dEdx *= Z*taul12/std::sqrt(tau);
97 }
98 else
99 {
100 G4double t1 = tau+1.;
101 G4double t2 = tau+2.;
102 G4double tsq = tau*tau;
103 G4double beta2 = tau*t2/(t1*t1);
104 G4double f = 2.*G4Log(tau) - (6.*tau+1.5*tsq-tau*(1.-tsq/3.)/t2
105 -tsq*(0.5-tsq/12.)/(t2*t2))/(t1*t1);
106 dEdx = Z*(G4Log(2.*tau+4.)-2.*ionpotlog+f)/beta2;
107
108 // loss from bremsstrahlung follows
109 G4double cbrem = (cbr1+cbr2*Z)
110 * (cbr3+cbr4*G4Log(kinEnergy/Thigh));
111 dEdx += cbrem*Z*(Z+1.)*bremfactor*tau/beta2;
112 }
113 return dEdx*CLHEP::twopi_mc2_rcl2;
114}
115
116// --------------------------------------------------------------------
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4double G4Log(G4double x)
Definition: G4Log.hh:227
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double logZ(G4int Z) const
Definition: G4Pow.hh:137
G4double ComputeValue(const G4int Z, const G4double kinEnergy) final
const G4ParticleDefinition * theParticle