source: trunk/source/processes/cuts/src/G4RToEConvForPositron.cc@ 1351

Last change on this file since 1351 was 1337, checked in by garnier, 15 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 4.5 KB
Line 
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: G4RToEConvForPositron.cc,v 1.7 2009/09/11 15:21:39 kurasige Exp $
28// GEANT4 tag $Name: geant4-09-04-beta-01 $
29//
30//
31// --------------------------------------------------------------
32// GEANT 4 class implementation file/ History:
33// 5 Oct. 2002, H.Kuirashige : Structure created based on object model
34// --------------------------------------------------------------
35
36#include "G4RToEConvForPositron.hh"
37#include "G4ParticleDefinition.hh"
38#include "G4ParticleTable.hh"
39#include "G4Material.hh"
40#include "G4PhysicsLogVector.hh"
41
42#include "G4ios.hh"
43
44G4RToEConvForPositron::G4RToEConvForPositron() : G4VRangeToEnergyConverter()
45{
46 theParticle = G4ParticleTable::GetParticleTable()->FindParticle("e+");
47 if (theParticle ==0) {
48#ifdef G4VERBOSE
49 if (GetVerboseLevel()>0) {
50 G4cout << " G4RToEConvForPositron::G4RToEConvForPositron() ";
51 G4cout << " Positron is not defined !!" << G4endl;
52 }
53#endif
54 }
55}
56
57G4RToEConvForPositron::~G4RToEConvForPositron()
58{
59}
60
61
62
63
64// **********************************************************************
65// ************************* ComputeLoss ********************************
66// **********************************************************************
67G4double G4RToEConvForPositron::ComputeLoss(G4double AtomicNumber,
68 G4double KineticEnergy) const
69{
70 static G4double Z;
71 static G4double taul, ionpot, ionpotlog;
72 const G4double cbr1=0.02, cbr2=-5.7e-5, cbr3=1., cbr4=0.072;
73 const G4double Tlow=10.*keV, Thigh=1.*GeV;
74 static G4double bremfactor = 0.1 ;
75
76 G4double Mass = theParticle->GetPDGMass();
77 // calculate dE/dx for electrons
78 if( std::fabs(AtomicNumber-Z)>0.1 ) {
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-tsq*(0.5-tsq/12.)/
96 (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 G4double t1 = tau+1.;
104 G4double t2 = tau+2.;
105 G4double tsq = tau*tau;
106 G4double beta2 = tau*t2/(t1*t1);
107 G4double f = 2.*std::log(tau)
108 - (6.*tau+1.5*tsq-tau*(1.-tsq/3.)/t2-tsq*(0.5-tsq/12.)/
109 (t2*t2))/(t1*t1);
110 dEdx = (std::log(2.*tau+4.)-2.*ionpotlog+f)/beta2;
111 dEdx = twopi_mc2_rcl2*Z*dEdx;
112
113 // loss from bremsstrahlung follows
114 G4double cbrem = (cbr1+cbr2*Z)
115 *(cbr3+cbr4*std::log(KineticEnergy/Thigh));
116 cbrem = Z*(Z+1.)*cbrem*tau/beta2;
117 cbrem *= bremfactor ;
118 dEdx += twopi_mc2_rcl2*cbrem;
119 }
120 return dEdx;
121}
122
123
Note: See TracBrowser for help on using the repository browser.