Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 //
27 // $Id$
28 //
29 //
30 // --------------------------------------------------------------
31 // GEANT 4 class implementation file/ History:
32 // 5 Oct. 2002, H.Kuirashige : Structure created based on object model
33 // --------------------------------------------------------------
34 
35 #include "G4RToEConvForPositron.hh"
36 #include "G4ParticleDefinition.hh"
37 #include "G4ParticleTable.hh"
38 #include "G4Material.hh"
39 #include "G4PhysicsLogVector.hh"
40 
41 #include "G4ios.hh"
42 #include "G4PhysicalConstants.hh"
43 #include "G4SystemOfUnits.hh"
44 
46 {
48  if (theParticle ==0) {
49 #ifdef G4VERBOSE
50  if (GetVerboseLevel()>0) {
51  G4cout << " G4RToEConvForPositron::G4RToEConvForPositron() ";
52  G4cout << " Positron is not defined !!" << G4endl;
53  }
54 #endif
55  }
56 }
57 
59 {
60 }
61 
62 
63 
64 
65 // **********************************************************************
66 // ************************* ComputeLoss ********************************
67 // **********************************************************************
69  G4double KineticEnergy) const
70 {
71  static G4double Z;
72  static G4double taul, ionpot, ionpotlog;
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  static G4double bremfactor = 0.1 ;
76 
77  G4double Mass = theParticle->GetPDGMass();
78  // calculate dE/dx for electrons
79  if( std::fabs(AtomicNumber-Z)>0.1 ) {
80  Z = AtomicNumber;
81  taul = Tlow/Mass;
82  ionpot = 1.6e-5*MeV*std::exp(0.9*std::log(Z))/Mass;
83  ionpotlog = std::log(ionpot);
84  }
85 
86  G4double tau = KineticEnergy/Mass;
87  G4double dEdx;
88 
89  if(tau<taul)
90  {
91  G4double t1 = taul+1.;
92  G4double t2 = taul+2.;
93  G4double tsq = taul*taul;
94  G4double beta2 = taul*t2/(t1*t1);
95  G4double f = 2.*std::log(taul)
96  -(6.*taul+1.5*tsq-taul*(1.-tsq/3.)/t2-tsq*(0.5-tsq/12.)/
97  (t2*t2))/(t1*t1);
98  dEdx = (std::log(2.*taul+4.)-2.*ionpotlog+f)/beta2;
99  dEdx = twopi_mc2_rcl2*Z*dEdx;
100  G4double clow = dEdx*std::sqrt(taul);
101  dEdx = clow/std::sqrt(KineticEnergy/Mass);
102 
103  } else {
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-tsq*(0.5-tsq/12.)/
110  (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 }
123 
124