[819] | 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 | // $Id: G4eeToHadronsModel.cc,v 1.8 2007/05/22 17:37:30 vnivanch Exp $ |
---|
| 27 | // GEANT4 tag $Name: $ |
---|
| 28 | // |
---|
| 29 | // ------------------------------------------------------------------- |
---|
| 30 | // |
---|
| 31 | // GEANT4 Class header file |
---|
| 32 | // |
---|
| 33 | // |
---|
| 34 | // File name: G4eeToHadronsModel |
---|
| 35 | // |
---|
| 36 | // Author: Vladimir Ivanchenko |
---|
| 37 | // |
---|
| 38 | // Creation date: 12.08.2003 |
---|
| 39 | // |
---|
| 40 | // Modifications: |
---|
| 41 | // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko) |
---|
| 42 | // 18-05-05 Use optimized interfaces (V.Ivantchenko) |
---|
| 43 | // |
---|
| 44 | // |
---|
| 45 | // ------------------------------------------------------------------- |
---|
| 46 | // |
---|
| 47 | |
---|
| 48 | |
---|
| 49 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 50 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 51 | |
---|
| 52 | #include "G4eeToHadronsModel.hh" |
---|
| 53 | #include "Randomize.hh" |
---|
| 54 | #include "G4Electron.hh" |
---|
| 55 | #include "G4Gamma.hh" |
---|
| 56 | #include "G4Positron.hh" |
---|
| 57 | #include "G4PionPlus.hh" |
---|
| 58 | #include "Randomize.hh" |
---|
| 59 | #include "G4Vee2hadrons.hh" |
---|
| 60 | #include "G4PhysicsVector.hh" |
---|
| 61 | #include "G4PhysicsLogVector.hh" |
---|
| 62 | |
---|
| 63 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 64 | |
---|
| 65 | using namespace std; |
---|
| 66 | |
---|
| 67 | G4eeToHadronsModel::G4eeToHadronsModel(const G4Vee2hadrons* m, |
---|
| 68 | G4int ver, |
---|
| 69 | const G4String& nam) |
---|
| 70 | : G4VEmModel(nam), |
---|
| 71 | model(m), |
---|
| 72 | crossPerElectron(0), |
---|
| 73 | crossBornPerElectron(0), |
---|
| 74 | isInitialised(false), |
---|
| 75 | nbins(100), |
---|
| 76 | verbose(ver) |
---|
| 77 | { |
---|
| 78 | theGamma = G4Gamma::Gamma(); |
---|
| 79 | } |
---|
| 80 | |
---|
| 81 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 82 | |
---|
| 83 | G4eeToHadronsModel::~G4eeToHadronsModel() |
---|
| 84 | { |
---|
| 85 | delete model; |
---|
| 86 | delete crossPerElectron; |
---|
| 87 | delete crossBornPerElectron; |
---|
| 88 | } |
---|
| 89 | |
---|
| 90 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 91 | |
---|
| 92 | void G4eeToHadronsModel::Initialise(const G4ParticleDefinition*, |
---|
| 93 | const G4DataVector&) |
---|
| 94 | { |
---|
| 95 | if(isInitialised) return; |
---|
| 96 | isInitialised = true; |
---|
| 97 | |
---|
| 98 | highKinEnergy = HighEnergyLimit(); |
---|
| 99 | lowKinEnergy = LowEnergyLimit(); |
---|
| 100 | |
---|
| 101 | emin = model->ThresholdEnergy(); |
---|
| 102 | emax = 2.0*electron_mass_c2*sqrt(1.0 + 0.5*highKinEnergy/electron_mass_c2); |
---|
| 103 | if(emin > emax) emin = emax; |
---|
| 104 | |
---|
| 105 | lowKinEnergy = 0.5*emin*emin/electron_mass_c2 - 2.0*electron_mass_c2; |
---|
| 106 | |
---|
| 107 | epeak = min(model->PeakEnergy(), emax); |
---|
| 108 | peakKinEnergy = 0.5*epeak*epeak/electron_mass_c2 - 2.0*electron_mass_c2; |
---|
| 109 | |
---|
| 110 | if(verbose>0) { |
---|
| 111 | G4cout << "G4eeToHadronsModel::Initialise: " << G4endl; |
---|
| 112 | G4cout << "LabSystem: emin(GeV)= " << lowKinEnergy/GeV |
---|
| 113 | << " epeak(GeV)= " << peakKinEnergy/GeV |
---|
| 114 | << " emax(GeV)= " << highKinEnergy/GeV |
---|
| 115 | << G4endl; |
---|
| 116 | G4cout << "SM System: emin(MeV)= " << emin/MeV |
---|
| 117 | << " epeak(MeV)= " << epeak/MeV |
---|
| 118 | << " emax(MeV)= " << emax/MeV |
---|
| 119 | << G4endl; |
---|
| 120 | } |
---|
| 121 | |
---|
| 122 | if(lowKinEnergy < peakKinEnergy) { |
---|
| 123 | crossBornPerElectron = model->PhysicsVector(emin, emax); |
---|
| 124 | crossPerElectron = model->PhysicsVector(emin, emax); |
---|
| 125 | nbins = crossPerElectron->GetVectorLength(); |
---|
| 126 | for(G4int i=0; i<nbins; i++) { |
---|
| 127 | G4double e = crossPerElectron->GetLowEdgeEnergy(i); |
---|
| 128 | G4double cs = model->ComputeCrossSection(e); |
---|
| 129 | crossBornPerElectron->PutValue(i, cs); |
---|
| 130 | } |
---|
| 131 | ComputeCMCrossSectionPerElectron(); |
---|
| 132 | } |
---|
| 133 | if(verbose>1) { |
---|
| 134 | G4cout << "G4eeToHadronsModel: Cross secsions per electron" |
---|
| 135 | << " nbins= " << nbins |
---|
| 136 | << " emin(MeV)= " << emin/MeV |
---|
| 137 | << " emax(MeV)= " << emax/MeV |
---|
| 138 | << G4endl; |
---|
| 139 | G4bool b; |
---|
| 140 | for(G4int i=0; i<nbins; i++) { |
---|
| 141 | G4double e = crossPerElectron->GetLowEdgeEnergy(i); |
---|
| 142 | G4double s1 = crossPerElectron->GetValue(e, b); |
---|
| 143 | G4double s2 = crossBornPerElectron->GetValue(e, b); |
---|
| 144 | G4cout << "E(MeV)= " << e/MeV |
---|
| 145 | << " cross(nb)= " << s1/nanobarn |
---|
| 146 | << " crossBorn(nb)= " << s2/nanobarn |
---|
| 147 | << G4endl; |
---|
| 148 | } |
---|
| 149 | } |
---|
| 150 | } |
---|
| 151 | |
---|
| 152 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 153 | |
---|
| 154 | G4double G4eeToHadronsModel::ComputeCrossSectionPerElectron( |
---|
| 155 | const G4ParticleDefinition*, |
---|
| 156 | G4double kineticEnergy, |
---|
| 157 | G4double, G4double) |
---|
| 158 | { |
---|
| 159 | G4double cross = 0.0; |
---|
| 160 | if(crossPerElectron) { |
---|
| 161 | G4bool b; |
---|
| 162 | G4double e = 2.0*electron_mass_c2* |
---|
| 163 | sqrt(1.0 + 0.5*kineticEnergy/electron_mass_c2); |
---|
| 164 | cross = crossPerElectron->GetValue(e, b); |
---|
| 165 | } |
---|
| 166 | // G4cout << "e= " << kineticEnergy << " cross= " << cross << G4endl; |
---|
| 167 | return cross; |
---|
| 168 | } |
---|
| 169 | |
---|
| 170 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 171 | |
---|
| 172 | void G4eeToHadronsModel::SampleSecondaries(std::vector<G4DynamicParticle*>* newp, |
---|
| 173 | const G4MaterialCutsCouple*, |
---|
| 174 | const G4DynamicParticle* dParticle, |
---|
| 175 | G4double, |
---|
| 176 | G4double) |
---|
| 177 | { |
---|
| 178 | if(crossPerElectron) { |
---|
| 179 | G4double t = dParticle->GetKineticEnergy(); |
---|
| 180 | G4double e = 2.0*electron_mass_c2*sqrt(1.0 + 0.5*t/electron_mass_c2); |
---|
| 181 | G4LorentzVector inlv = dParticle->Get4Momentum(); |
---|
| 182 | G4ThreeVector inBoost = inlv.boostVector(); |
---|
| 183 | if(e > emin) { |
---|
| 184 | G4DynamicParticle* gamma = GenerateCMPhoton(e); |
---|
| 185 | G4LorentzVector gLv = gamma->Get4Momentum(); |
---|
| 186 | G4LorentzVector lv(0.0,0.0,0.0,e); |
---|
| 187 | lv -= gLv; |
---|
| 188 | G4double m = lv.m(); |
---|
| 189 | G4ThreeVector boost = lv.boostVector(); |
---|
| 190 | const G4ThreeVector dir = gamma->GetMomentumDirection(); |
---|
| 191 | model->SampleSecondaries(newp, m, dir); |
---|
| 192 | G4int np = newp->size(); |
---|
| 193 | for(G4int j=0; j<np; j++) { |
---|
| 194 | G4DynamicParticle* dp = (*newp)[j]; |
---|
| 195 | G4LorentzVector v = dp->Get4Momentum(); |
---|
| 196 | v.boost(boost); |
---|
| 197 | v.boost(inBoost); |
---|
| 198 | dp->Set4Momentum(v); |
---|
| 199 | } |
---|
| 200 | gLv.boost(inBoost); |
---|
| 201 | gamma->Set4Momentum(gLv); |
---|
| 202 | newp->push_back(gamma); |
---|
| 203 | } |
---|
| 204 | } |
---|
| 205 | } |
---|
| 206 | |
---|
| 207 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 208 | |
---|
| 209 | void G4eeToHadronsModel::ComputeCMCrossSectionPerElectron() |
---|
| 210 | { |
---|
| 211 | G4bool b; |
---|
| 212 | for(G4int i=0; i<nbins; i++) { |
---|
| 213 | G4double e = crossPerElectron->GetLowEdgeEnergy(i); |
---|
| 214 | G4double cs = 0.0; |
---|
| 215 | if(i > 0) { |
---|
| 216 | G4double L = 2.0*log(e/electron_mass_c2); |
---|
| 217 | G4double bt = 2.0*fine_structure_const*(L - 1.0)/pi; |
---|
| 218 | G4double btm1= bt - 1.0; |
---|
| 219 | G4double del = 1. + fine_structure_const*(1.5*L + pi*pi/3. -2.)/pi; |
---|
| 220 | G4double s1 = crossBornPerElectron->GetValue(e, b); |
---|
| 221 | G4double e1 = crossPerElectron->GetLowEdgeEnergy(i-1); |
---|
| 222 | G4double x1 = 1. - e1/e; |
---|
| 223 | cs += s1*(del*pow(x1,bt) - bt*(x1 - 0.25*x1*x1)); |
---|
| 224 | if(i > 1) { |
---|
| 225 | G4double e2 = e1; |
---|
| 226 | G4double x2 = x1; |
---|
| 227 | G4double s2 = crossBornPerElectron->GetValue(e2, b); |
---|
| 228 | G4double w2 = bt*(del*pow(x2,btm1) - 1.0 + 0.5*x2); |
---|
| 229 | |
---|
| 230 | for(G4int j=i-2; j>=0; j--) { |
---|
| 231 | e1 = crossPerElectron->GetLowEdgeEnergy(j); |
---|
| 232 | x1 = 1. - e1/e; |
---|
| 233 | G4double s1 = crossBornPerElectron->GetValue(e1, b); |
---|
| 234 | G4double w1 = bt*(del*pow(x1,btm1) - 1.0 + 0.5*x1); |
---|
| 235 | cs += 0.5*(x1 - x2)*(w2*s2 + w1*s1); |
---|
| 236 | e2 = e1; |
---|
| 237 | x2 = x1; |
---|
| 238 | s2 = s1; |
---|
| 239 | w2 = w1; |
---|
| 240 | } |
---|
| 241 | } |
---|
| 242 | } |
---|
| 243 | crossPerElectron->PutValue(i, cs); |
---|
| 244 | // G4cout << "e= " << e << " cs= " << cs << G4endl; |
---|
| 245 | } |
---|
| 246 | } |
---|
| 247 | |
---|
| 248 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 249 | |
---|
| 250 | G4DynamicParticle* G4eeToHadronsModel::GenerateCMPhoton(G4double e) |
---|
| 251 | { |
---|
| 252 | G4bool b; |
---|
| 253 | G4double x; |
---|
| 254 | G4DynamicParticle* gamma = 0; |
---|
| 255 | G4double L = 2.0*log(e/electron_mass_c2); |
---|
| 256 | G4double bt = 2.0*fine_structure_const*(L - 1.)/pi; |
---|
| 257 | G4double btm1= bt - 1.0; |
---|
| 258 | G4double del = 1. + fine_structure_const*(1.5*L + pi*pi/3. -2.)/pi; |
---|
| 259 | |
---|
| 260 | G4double s0 = crossBornPerElectron->GetValue(e, b); |
---|
| 261 | G4double de = (emax - emin)/(G4double)nbins; |
---|
| 262 | G4double x0 = min(de,e - emin)/e; |
---|
| 263 | G4double ds = crossBornPerElectron->GetValue(e, b) |
---|
| 264 | *(del*pow(x0,bt) - bt*(x0 - 0.25*x0*x0)); |
---|
| 265 | G4double e1 = e*(1. - x0); |
---|
| 266 | |
---|
| 267 | if(e1 < emax && s0*G4UniformRand()<ds) { |
---|
| 268 | x = x0*pow(G4UniformRand(),1./bt); |
---|
| 269 | } else { |
---|
| 270 | |
---|
| 271 | x = 1. - e1/e; |
---|
| 272 | G4double s1 = crossBornPerElectron->GetValue(e1, b); |
---|
| 273 | G4double w1 = bt*(del*pow(x,btm1) - 1.0 + 0.5*x); |
---|
| 274 | G4double grej = s1*w1; |
---|
| 275 | G4double f; |
---|
| 276 | // G4cout << "e= " << e/GeV << " epeak= " << epeak/GeV |
---|
| 277 | // << " s1= " << s1 << " w1= " << w1 |
---|
| 278 | // << " grej= " << grej << G4endl; |
---|
| 279 | // Above emax cross section is 0 |
---|
| 280 | if(e1 > emax) { |
---|
| 281 | x = 1. - emax/e; |
---|
| 282 | G4double s2 = crossBornPerElectron->GetValue(emax, b); |
---|
| 283 | G4double w2 = bt*(del*pow(x,btm1) - 1.0 + 0.5*x); |
---|
| 284 | grej = s2*w2; |
---|
| 285 | // G4cout << "emax= " << emax << " s2= " << s2 << " w2= " << w2 |
---|
| 286 | // << " grej= " << grej << G4endl; |
---|
| 287 | } |
---|
| 288 | |
---|
| 289 | if(e1 > epeak) { |
---|
| 290 | x = 1. - epeak/e; |
---|
| 291 | G4double s2 = crossBornPerElectron->GetValue(epeak, b); |
---|
| 292 | G4double w2 = bt*(del*pow(x,btm1) - 1.0 + 0.5*x); |
---|
| 293 | grej = max(grej,s2*w2); |
---|
| 294 | //G4cout << "epeak= " << epeak << " s2= " << s2 << " w2= " << w2 |
---|
| 295 | // << " grej= " << grej << G4endl; |
---|
| 296 | } |
---|
| 297 | G4double xmin = 1. - e1/e; |
---|
| 298 | if(e1 > emax) xmin = 1. - emax/e; |
---|
| 299 | G4double xmax = 1. - emin/e; |
---|
| 300 | do { |
---|
| 301 | x = xmin + G4UniformRand()*(xmax - xmin); |
---|
| 302 | G4double s2 = crossBornPerElectron->GetValue((1.0 - x)*e, b); |
---|
| 303 | G4double w2 = bt*(del*pow(x,btm1) - 1.0 + 0.5*x); |
---|
| 304 | //G4cout << "x= " << x << " xmin= " << xmin << " xmax= " << xmax |
---|
| 305 | // << " s2= " << s2 << " w2= " << w2 |
---|
| 306 | // << G4endl; |
---|
| 307 | f = s2*w2; |
---|
| 308 | if(f > grej) { |
---|
| 309 | G4cout << "G4DynamicParticle* G4eeToHadronsModel:WARNING " |
---|
| 310 | << f << " > " << grej << " majorant is`small!" |
---|
| 311 | << G4endl; |
---|
| 312 | } |
---|
| 313 | } while (f < grej*G4UniformRand()); |
---|
| 314 | } |
---|
| 315 | |
---|
| 316 | G4ThreeVector dir(0.0,0.0,1.0); |
---|
| 317 | gamma = new G4DynamicParticle(theGamma,dir,x*e); |
---|
| 318 | return gamma; |
---|
| 319 | } |
---|
| 320 | |
---|
| 321 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 322 | |
---|