Geant4 10.7.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// G4RToEConvForProton class implementation
27//
28// Author: H.Kurashige, 05 October 2002 - First implementation
29// --------------------------------------------------------------------
30
32#include "G4ParticleTable.hh"
33#include "G4Material.hh"
34#include "G4PhysicsLogVector.hh"
35
36#include "G4ios.hh"
38#include "G4SystemOfUnits.hh"
39
40// --------------------------------------------------------------------
43{
45 if (theParticle == nullptr)
46 {
47#ifdef G4VERBOSE
48 if (GetVerboseLevel()>0)
49 {
50 G4cout << "G4RToEConvForProton::G4RToEConvForProton() - ";
51 G4cout << "Proton is not defined !!" << G4endl;
52 }
53#endif
54 }
55 else
56 {
58 }
59}
60
61// --------------------------------------------------------------------
63{
64}
65
66// --------------------------------------------------------------------
68{
69 // Simple formula - range = Ekin/(100*keV)*(1*mm);
70
71 return (rangeCut/(1.0*mm)) * (100.0*keV);
72}
73
74// **********************************************************************
75// ************************* ComputeLoss ********************************
76// **********************************************************************
78 G4double KineticEnergy)
79{
80 // calculate dE/dx
81 const G4double z2Particle = 1.0;
82
83 if( std::fabs(AtomicNumber-Z)>0.1 )
84 {
85 // recalculate constants
86 Z = AtomicNumber;
87 G4double Z13 = std::exp(std::log(Z)/3.);
88 tau0 = 0.1*Z13*MeV/proton_mass_c2;
89 taum = 0.035*Z13*MeV/proton_mass_c2;
90 taul = 2.*MeV/proton_mass_c2;
91 ionpot = 1.6e-5*MeV*std::exp(0.9*std::log(Z));
92 cc = (taul+1.)*(taul+1.)*std::log(2.*electron_mass_c2*taul*(taul+2.)/ionpot)
93 / (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/Mass;
100 G4double dEdx;
101 if ( tau <= tau0 )
102 {
103 dEdx = ca*(std::sqrt(tau)+cba*tau);
104 }
105 else
106 {
107 if( tau <= taul )
108 {
109 dEdx = cc/std::sqrt(tau);
110 }
111 else
112 {
113 dEdx = (tau+1.)*(tau+1.)
114 * std::log(2.*electron_mass_c2*tau*(tau+2.)/ionpot)
115 / (tau*(tau+2.))-1.;
116 dEdx = 2.*twopi_mc2_rcl2*Z*dEdx;
117 }
118 }
119 return dEdx*z2Particle ;
120}
121
122// **********************************************************************
123// ************************* Reset ********************************
124// **********************************************************************
126{
127 // do nothing because loss tables and range vectors are not used
128
129 return;
130}
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 Convert(G4double rangeCut, const G4Material *material)
virtual G4double ComputeLoss(G4double AtomicNumber, G4double KineticEnergy)
const G4ParticleDefinition * theParticle