Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4RToEConvForProton.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: G4RToEConvForProton.cc 70745 2013-06-05 10:54:00Z gcosmo $
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 "G4RToEConvForProton.hh"
36 #include "G4ParticleTable.hh"
37 #include "G4Material.hh"
38 #include "G4PhysicsLogVector.hh"
39 
40 #include "G4ios.hh"
41 #include "G4PhysicalConstants.hh"
42 #include "G4SystemOfUnits.hh"
43 
46  Mass(0.0),
47  Z(-1.),
48  tau0(0.0), taul(0.0), taum(0.0),
49  ionpot(0.0),
50  ca(0.0), cba(0.0), cc(0.0)
51 {
53  if (theParticle ==0) {
54 #ifdef G4VERBOSE
55  if (GetVerboseLevel()>0) {
56  G4cout << " G4RToEConvForProton::G4RToEConvForProton() ";
57  G4cout << " proton is not defined !!" << G4endl;
58  }
59 #endif
60  } else {
62  }
63 }
64 
66 {
67 }
68 
69 
71 {
72  // Simple formula
73  // range = Ekin/(100*keV)*(1*mm);
74  return (rangeCut/(1.0*mm)) * (100.0*keV);
75 }
76 
77 
78 // **********************************************************************
79 // ************************* ComputeLoss ********************************
80 // **********************************************************************
82  G4double KineticEnergy)
83 {
84  // calculate dE/dx
85  const G4double z2Particle = 1.0;
86 
87  if( std::fabs(AtomicNumber-Z)>0.1 ){
88  // recalculate constants
89  Z = AtomicNumber;
90  G4double Z13 = std::exp(std::log(Z)/3.);
91  tau0 = 0.1*Z13*MeV/proton_mass_c2;
92  taum = 0.035*Z13*MeV/proton_mass_c2;
93  taul = 2.*MeV/proton_mass_c2;
94  ionpot = 1.6e-5*MeV*std::exp(0.9*std::log(Z));
95  cc = (taul+1.)*(taul+1.)*std::log(2.*electron_mass_c2*taul*(taul+2.)/ionpot)/(taul*(taul+2.))-1.;
96  cc = 2.*twopi_mc2_rcl2*Z*cc*std::sqrt(taul);
97  ca = cc/((1.-0.5*std::sqrt(tau0/taum))*tau0);
98  cba = -0.5/std::sqrt(taum);
99  }
100 
101  G4double tau = KineticEnergy/Mass;
102  G4double dEdx;
103  if ( tau <= tau0 ) {
104  dEdx = ca*(std::sqrt(tau)+cba*tau);
105  } else {
106  if( tau <= taul ) {
107  dEdx = cc/std::sqrt(tau);
108  } else {
109  dEdx = (tau+1.)*(tau+1.)*
110  std::log(2.*electron_mass_c2*tau*(tau+2.)/ionpot)/(tau*(tau+2.))-1.;
111  dEdx = 2.*twopi_mc2_rcl2*Z*dEdx;
112  }
113  }
114  return dEdx*z2Particle ;
115 }
116 
117 
118 // **********************************************************************
119 // ************************* Reset ********************************
120 // **********************************************************************
122 {
123  // do nothing because loss tables and range vectors are not used
124  return;
125 }
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
virtual G4double Convert(G4double rangeCut, const G4Material *material)
G4GLOB_DLL std::ostream G4cout
float proton_mass_c2
Definition: hepunit.py:275
float electron_mass_c2
Definition: hepunit.py:274
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
virtual G4double ComputeLoss(G4double AtomicNumber, G4double KineticEnergy)
const G4ParticleDefinition * theParticle