Geant4 10.7.0
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#include "G4Material.hh"
35#include "G4PhysicsLogVector.hh"
36
37#include "G4ios.hh"
39#include "G4SystemOfUnits.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
67// **********************************************************************
68// ************************* ComputeLoss ********************************
69// **********************************************************************
71 G4double KineticEnergy)
72{
73 const G4double cbr1=0.02, cbr2=-5.7e-5, cbr3=1., cbr4=0.072;
74 const G4double Tlow=10.*keV, Thigh=1.*GeV;
75
76 // calculate dE/dx for electrons
77 if( std::fabs(AtomicNumber-Z)>0.1 )
78 {
79 Z = AtomicNumber;
80 taul = Tlow/Mass;
81 ionpot = 1.6e-5*MeV*std::exp(0.9*std::log(Z))/Mass;
82 ionpotlog = std::log(ionpot);
83 }
84
85 G4double tau = KineticEnergy/Mass;
86 G4double dEdx;
87
88 if(tau<taul)
89 {
90 G4double t1 = taul+1.;
91 G4double t2 = taul+2.;
92 G4double tsq = taul*taul;
93 G4double beta2 = taul*t2/(t1*t1);
94 G4double f = 2.*std::log(taul)
95 -(6.*taul+1.5*tsq-taul*(1.-tsq/3.)/t2
96 -tsq*(0.5-tsq/12.)/(t2*t2))/(t1*t1);
97 dEdx = (std::log(2.*taul+4.)-2.*ionpotlog+f)/beta2;
98 dEdx = twopi_mc2_rcl2*Z*dEdx;
99 G4double clow = dEdx*std::sqrt(taul);
100 dEdx = clow/std::sqrt(KineticEnergy/Mass);
101 }
102 else
103 {
104 G4double t1 = tau+1.;
105 G4double t2 = tau+2.;
106 G4double tsq = tau*tau;
107 G4double beta2 = tau*t2/(t1*t1);
108 G4double f = 2.*std::log(tau)
109 - (6.*tau+1.5*tsq-tau*(1.-tsq/3.)/t2
110 -tsq*(0.5-tsq/12.)/(t2*t2))/(t1*t1);
111 dEdx = (std::log(2.*tau+4.)-2.*ionpotlog+f)/beta2;
112 dEdx = twopi_mc2_rcl2*Z*dEdx;
113
114 // loss from bremsstrahlung follows
115 G4double cbrem = (cbr1+cbr2*Z)
116 * (cbr3+cbr4*std::log(KineticEnergy/Thigh));
117 cbrem = Z*(Z+1.)*cbrem*tau/beta2;
118 cbrem *= bremfactor;
119 dEdx += twopi_mc2_rcl2*cbrem;
120 }
121 return dEdx;
122}
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
virtual G4double ComputeLoss(G4double AtomicNumber, G4double KineticEnergy)
const G4ParticleDefinition * theParticle