Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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$
28//
29//
30// --------------------------------------------------------------
31// GEANT 4 class implementation file/ History:
32// 5 Oct. 2002, H.Kuirashige : Structure created based on object model
33// --------------------------------------------------------------
34
36#include "G4ParticleTable.hh"
37#include "G4Material.hh"
38#include "G4PhysicsLogVector.hh"
39
40#include "G4ios.hh"
42#include "G4SystemOfUnits.hh"
43
45{
47 if (theParticle ==0) {
48#ifdef G4VERBOSE
49 if (GetVerboseLevel()>0) {
50 G4cout << " G4RToEConvForProton::G4RToEConvForProton() ";
51 G4cout << " proton is not defined !!" << G4endl;
52 }
53#endif
54 }
55}
56
58{
59}
60
61
63{
64 // Simple formula
65 // range = Ekin/(100*keV)*(1*mm);
66 return (rangeCut/(1.0*mm)) * (100.0*keV);
67}
68
69
70// **********************************************************************
71// ************************* ComputeLoss ********************************
72// **********************************************************************
74 G4double KineticEnergy) const
75{
76 // calculate dE/dx
77
78 static G4double Z;
79 static G4double ionpot, tau0, taum, taul, ca, cba, cc;
80
81 G4double z2Particle = theParticle->GetPDGCharge()/eplus;
82 z2Particle *= z2Particle;
83 if (z2Particle < 0.1) return 0.0;
84
85 if( std::fabs(AtomicNumber-Z)>0.1 ){
86 // recalculate constants
87 Z = AtomicNumber;
88 G4double Z13 = std::exp(std::log(Z)/3.);
89 tau0 = 0.1*Z13*MeV/proton_mass_c2;
90 taum = 0.035*Z13*MeV/proton_mass_c2;
91 taul = 2.*MeV/proton_mass_c2;
92 ionpot = 1.6e-5*MeV*std::exp(0.9*std::log(Z));
93 cc = (taul+1.)*(taul+1.)*std::log(2.*electron_mass_c2*taul*(taul+2.)/ionpot)/(taul*(taul+2.))-1.;
94 cc = 2.*twopi_mc2_rcl2*Z*cc*std::sqrt(taul);
95 ca = cc/((1.-0.5*std::sqrt(tau0/taum))*tau0);
96 cba = -0.5/std::sqrt(taum);
97 }
98
99 G4double tau = KineticEnergy/theParticle->GetPDGMass();
100 G4double dEdx;
101 if ( tau <= tau0 ) {
102 dEdx = ca*(std::sqrt(tau)+cba*tau);
103 } else {
104 if( tau <= taul ) {
105 dEdx = cc/std::sqrt(tau);
106 } else {
107 dEdx = (tau+1.)*(tau+1.)*
108 std::log(2.*electron_mass_c2*tau*(tau+2.)/ionpot)/(tau*(tau+2.))-1.;
109 dEdx = 2.*twopi_mc2_rcl2*Z*dEdx;
110 }
111 }
112 return dEdx*z2Particle ;
113}
114
115
116// **********************************************************************
117// ************************* Reset ********************************
118// **********************************************************************
120{
121 // do nothing because loss tables and range vectors are not used
122 return;
123}
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4double GetPDGCharge() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
virtual G4double Convert(G4double rangeCut, const G4Material *material)
virtual G4double ComputeLoss(G4double AtomicNumber, G4double KineticEnergy) const
const G4ParticleDefinition * theParticle