Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes
G4RToEConvForPositron Class Reference

#include <G4RToEConvForPositron.hh>

Inheritance diagram for G4RToEConvForPositron:
G4VRangeToEnergyConverter

Public Member Functions

 G4RToEConvForPositron ()
 
virtual ~G4RToEConvForPositron ()
 
- Public Member Functions inherited from G4VRangeToEnergyConverter
 G4VRangeToEnergyConverter ()
 
 G4VRangeToEnergyConverter (const G4VRangeToEnergyConverter &right)
 
G4VRangeToEnergyConverteroperator= (const G4VRangeToEnergyConverter &right)
 
virtual ~G4VRangeToEnergyConverter ()
 
G4int operator== (const G4VRangeToEnergyConverter &right) const
 
G4int operator!= (const G4VRangeToEnergyConverter &right) const
 
virtual G4double Convert (G4double rangeCut, const G4Material *material)
 
const G4ParticleDefinitionGetParticleType () const
 
const G4PhysicsTableGetLossTable () const
 
virtual void Reset ()
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Protected Member Functions

virtual G4double ComputeLoss (G4double AtomicNumber, G4double KineticEnergy)
 
- Protected Member Functions inherited from G4VRangeToEnergyConverter
virtual void BuildLossTable ()
 
virtual void BuildRangeVector (const G4Material *aMaterial, G4RangeVector *rangeVector)
 
G4double ConvertCutToKineticEnergy (G4RangeVector *theRangeVector, G4double theCutInLength, size_t materialIndex) const
 

Protected Attributes

G4double Mass
 
G4double Z
 
G4double taul
 
G4double ionpot
 
G4double ionpotlog
 
G4double bremfactor
 
- Protected Attributes inherited from G4VRangeToEnergyConverter
G4double fMaxEnergyCut
 
const G4ParticleDefinitiontheParticle
 
G4LossTabletheLossTable
 
G4int NumberOfElements
 
const G4int TotBin
 
std::vector< G4RangeVector * > fRangeVectorStore
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VRangeToEnergyConverter
static void SetEnergyRange (G4double lowedge, G4double highedge)
 
static G4double GetLowEdgeEnergy ()
 
static G4double GetHighEdgeEnergy ()
 
static G4double GetMaxEnergyCut ()
 
static void SetMaxEnergyCut (G4double value)
 
- Protected Types inherited from G4VRangeToEnergyConverter
typedef G4PhysicsTable G4LossTable
 
typedef G4PhysicsLogVector G4LossVector
 
typedef G4PhysicsLogVector G4RangeVector
 
- Static Protected Attributes inherited from G4VRangeToEnergyConverter
static G4double LowestEnergy = 0.99e-3*MeV
 
static G4double HighestEnergy = 100.0e6*MeV
 
static G4double MaxEnergyCut = 10.0*GeV
 

Detailed Description

Definition at line 51 of file G4RToEConvForPositron.hh.

Constructor & Destructor Documentation

G4RToEConvForPositron::G4RToEConvForPositron ( )

Definition at line 45 of file G4RToEConvForPositron.cc.

References G4ParticleTable::FindParticle(), G4cout, G4endl, G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetPDGMass(), G4VRangeToEnergyConverter::GetVerboseLevel(), Mass, and G4VRangeToEnergyConverter::theParticle.

47  Mass(0.0),
48  Z(-1.),
49  taul(0.0),
50  ionpot(0.0),
51  ionpotlog(-1.0e-10),
52  bremfactor(0.1)
53 {
55  if (theParticle ==0) {
56 #ifdef G4VERBOSE
57  if (GetVerboseLevel()>0) {
58  G4cout << " G4RToEConvForPositron::G4RToEConvForPositron() ";
59  G4cout << " Positron is not defined !!" << G4endl;
60  }
61 #endif
62  } else {
64  }
65 }
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
G4GLOB_DLL std::ostream G4cout
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
#define G4endl
Definition: G4ios.hh:61
const G4ParticleDefinition * theParticle
G4RToEConvForPositron::~G4RToEConvForPositron ( )
virtual

Definition at line 67 of file G4RToEConvForPositron.cc.

68 {
69 }

Member Function Documentation

G4double G4RToEConvForPositron::ComputeLoss ( G4double  AtomicNumber,
G4double  KineticEnergy 
)
protectedvirtual

Implements G4VRangeToEnergyConverter.

Definition at line 77 of file G4RToEConvForPositron.cc.

References bremfactor, python.hepunit::GeV, ionpot, ionpotlog, python.hepunit::keV, Mass, python.hepunit::MeV, plottest35::t1, taul, python.hepunit::twopi_mc2_rcl2, and Z.

79 {
80  const G4double cbr1=0.02, cbr2=-5.7e-5, cbr3=1., cbr4=0.072;
81  const G4double Tlow=10.*keV, Thigh=1.*GeV;
82 
83  // calculate dE/dx for electrons
84  if( std::fabs(AtomicNumber-Z)>0.1 ) {
85  Z = AtomicNumber;
86  taul = Tlow/Mass;
87  ionpot = 1.6e-5*MeV*std::exp(0.9*std::log(Z))/Mass;
88  ionpotlog = std::log(ionpot);
89  }
90 
91  G4double tau = KineticEnergy/Mass;
92  G4double dEdx;
93 
94  if(tau<taul)
95  {
96  G4double t1 = taul+1.;
97  G4double t2 = taul+2.;
98  G4double tsq = taul*taul;
99  G4double beta2 = taul*t2/(t1*t1);
100  G4double f = 2.*std::log(taul)
101  -(6.*taul+1.5*tsq-taul*(1.-tsq/3.)/t2-tsq*(0.5-tsq/12.)/
102  (t2*t2))/(t1*t1);
103  dEdx = (std::log(2.*taul+4.)-2.*ionpotlog+f)/beta2;
104  dEdx = twopi_mc2_rcl2*Z*dEdx;
105  G4double clow = dEdx*std::sqrt(taul);
106  dEdx = clow/std::sqrt(KineticEnergy/Mass);
107 
108  } else {
109  G4double t1 = tau+1.;
110  G4double t2 = tau+2.;
111  G4double tsq = tau*tau;
112  G4double beta2 = tau*t2/(t1*t1);
113  G4double f = 2.*std::log(tau)
114  - (6.*tau+1.5*tsq-tau*(1.-tsq/3.)/t2-tsq*(0.5-tsq/12.)/
115  (t2*t2))/(t1*t1);
116  dEdx = (std::log(2.*tau+4.)-2.*ionpotlog+f)/beta2;
117  dEdx = twopi_mc2_rcl2*Z*dEdx;
118 
119  // loss from bremsstrahlung follows
120  G4double cbrem = (cbr1+cbr2*Z)
121  *(cbr3+cbr4*std::log(KineticEnergy/Thigh));
122  cbrem = Z*(Z+1.)*cbrem*tau/beta2;
123  cbrem *= bremfactor ;
124  dEdx += twopi_mc2_rcl2*cbrem;
125  }
126  return dEdx;
127 }
tuple t1
Definition: plottest35.py:33
double G4double
Definition: G4Types.hh:76

Field Documentation

G4double G4RToEConvForPositron::bremfactor
protected

Definition at line 72 of file G4RToEConvForPositron.hh.

Referenced by ComputeLoss().

G4double G4RToEConvForPositron::ionpot
protected

Definition at line 70 of file G4RToEConvForPositron.hh.

Referenced by ComputeLoss().

G4double G4RToEConvForPositron::ionpotlog
protected

Definition at line 71 of file G4RToEConvForPositron.hh.

Referenced by ComputeLoss().

G4double G4RToEConvForPositron::Mass
protected

Definition at line 67 of file G4RToEConvForPositron.hh.

Referenced by ComputeLoss(), and G4RToEConvForPositron().

G4double G4RToEConvForPositron::taul
protected

Definition at line 69 of file G4RToEConvForPositron.hh.

Referenced by ComputeLoss().

G4double G4RToEConvForPositron::Z
protected

Definition at line 68 of file G4RToEConvForPositron.hh.

Referenced by ComputeLoss().


The documentation for this class was generated from the following files: