ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MyMollerBhabhaModel.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file MyMollerBhabhaModel.cc
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 //
28 //
29 //
30 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 
33 #include "MyMollerBhabhaModel.hh"
34 
35 #include "G4PhysicalConstants.hh"
36 #include "G4SystemOfUnits.hh"
37 
38 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
39 
40 using namespace std;
41 
43  const G4String& nam)
44  : G4MollerBhabhaModel(p,nam)
45 {}
46 
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48 
50 {}
51 
52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53 
55  const G4Material* material,
56  const G4ParticleDefinition* p,
57  G4double kineticEnergy,
58  G4double cutEnergy)
59 {
60  if(!particle) SetParticle(p);
61  // calculate the dE/dx due to the ionization by Seltzer-Berger formula
62 
63  G4double electronDensity = material->GetElectronDensity();
64  G4double Zeff = electronDensity/material->GetTotNbOfAtomsPerVolume();
65  G4double th = 0.25*sqrt(Zeff)*keV;
66  G4double tkin = kineticEnergy;
67  G4bool lowEnergy = false;
68  if (kineticEnergy < th) {
69  tkin = th;
70  lowEnergy = true;
71  }
72  G4double tau = tkin/electron_mass_c2;
73  G4double gam = tau + 1.0;
74  G4double gamma2= gam*gam;
75  G4double beta2 = 1. - 1./gamma2;
76  //G4double bg2 = beta2*gamma2;
77 
78  G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
79  eexc /= electron_mass_c2;
80  G4double eexc2 = eexc*eexc;
81 
82  G4double d = min(cutEnergy, MaxSecondaryEnergy(p, tkin))/electron_mass_c2;
83  G4double dedx;
84 
85  // electron
86  if (isElectron) {
87 
88  dedx = log(2.0*(tau + 2.0)/eexc2) - 1.0 - beta2
89  + log((tau-d)*d) + tau/(tau-d)
90  + (0.5*d*d + (2.0*tau + 1.)*log(1. - d/tau))/gamma2;
91 
92  //positron
93  } else {
94 
95  G4double d2 = d*d*0.5;
96  G4double d3 = d2*d/1.5;
97  G4double d4 = d3*d*3.75;
98  G4double y = 1.0/(1.0 + gam);
99  dedx = log(2.0*(tau + 2.0)/eexc2) + log(tau*d)
100  - beta2*(tau + 2.0*d - y*(3.0*d2
101  + y*(d - d3 + y*(d2 - tau*d3 + d4))))/tau;
102  }
103 
104  //do not apply density correction
105  //G4double cden = material->GetIonisation()->GetCdensity();
106  //G4double mden = material->GetIonisation()->GetMdensity();
107  //G4double aden = material->GetIonisation()->GetAdensity();
108  //G4double x0den = material->GetIonisation()->GetX0density();
109  //G4double x1den = material->GetIonisation()->GetX1density();
110  //G4double x = log(bg2)/twoln10;
111 
112  //if (x >= x0den) {
113  // dedx -= twoln10*x - cden;
114  // if (x < x1den) dedx -= aden*pow(x1den-x, mden);
115  //}
116 
117  // now you can compute the total ionization loss
118  dedx *= twopi_mc2_rcl2*electronDensity/beta2;
119  if (dedx < 0.0) dedx = 0.0;
120 
121  // lowenergy extrapolation
122 
123  if (lowEnergy) {
124 
125  if (kineticEnergy >= lowLimit) dedx *= sqrt(tkin/kineticEnergy);
126  else dedx *= sqrt(tkin*kineticEnergy)/lowLimit;
127 
128  }
129  return dedx;
130 }
131 
132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......