| 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.9 2008/07/10 18:06:39 vnivanch Exp $
|
|---|
| 27 | // GEANT4 tag $Name: geant4-09-02-ref-02 $
|
|---|
| 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(G4Vee2hadrons* m, G4int ver,
|
|---|
| 68 | const G4String& nam)
|
|---|
| 69 | : G4VEmModel(nam),
|
|---|
| 70 | model(m),
|
|---|
| 71 | crossPerElectron(0),
|
|---|
| 72 | crossBornPerElectron(0),
|
|---|
| 73 | isInitialised(false),
|
|---|
| 74 | nbins(100),
|
|---|
| 75 | verbose(ver)
|
|---|
| 76 | {
|
|---|
| 77 | theGamma = G4Gamma::Gamma();
|
|---|
| 78 | }
|
|---|
| 79 |
|
|---|
| 80 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 81 |
|
|---|
| 82 | G4eeToHadronsModel::~G4eeToHadronsModel()
|
|---|
| 83 | {
|
|---|
| 84 | delete model;
|
|---|
| 85 | delete crossPerElectron;
|
|---|
| 86 | delete crossBornPerElectron;
|
|---|
| 87 | }
|
|---|
| 88 |
|
|---|
| 89 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 90 |
|
|---|
| 91 | void G4eeToHadronsModel::Initialise(const G4ParticleDefinition*,
|
|---|
| 92 | const G4DataVector&)
|
|---|
| 93 | {
|
|---|
| 94 | if(isInitialised) return;
|
|---|
| 95 | isInitialised = true;
|
|---|
| 96 |
|
|---|
| 97 | // Lab system
|
|---|
| 98 | highKinEnergy = HighEnergyLimit();
|
|---|
| 99 | lowKinEnergy = LowEnergyLimit();
|
|---|
| 100 |
|
|---|
| 101 | // CM system
|
|---|
| 102 | emin = model->LowEnergy();
|
|---|
| 103 | emax = model->HighEnergy();
|
|---|
| 104 |
|
|---|
| 105 | G4double emin0 =
|
|---|
| 106 | 2.0*electron_mass_c2*sqrt(1.0 + 0.5*lowKinEnergy/electron_mass_c2);
|
|---|
| 107 | G4double emax0 =
|
|---|
| 108 | 2.0*electron_mass_c2*sqrt(1.0 + 0.5*highKinEnergy/electron_mass_c2);
|
|---|
| 109 |
|
|---|
| 110 | // recompute low energy
|
|---|
| 111 | if(emin0 > emax) {
|
|---|
| 112 | emin0 = emax;
|
|---|
| 113 | model->SetLowEnergy(emin0);
|
|---|
| 114 | }
|
|---|
| 115 | if(emin > emin0) {
|
|---|
| 116 | emin0 = emin;
|
|---|
| 117 | lowKinEnergy = 0.5*emin*emin/electron_mass_c2 - 2.0*electron_mass_c2;
|
|---|
| 118 | SetLowEnergyLimit(lowKinEnergy);
|
|---|
| 119 | }
|
|---|
| 120 |
|
|---|
| 121 | // recompute high energy
|
|---|
| 122 | if(emax < emax0) {
|
|---|
| 123 | emax0 = emax;
|
|---|
| 124 | highKinEnergy = 0.5*emax*emax/electron_mass_c2 - 2.0*electron_mass_c2;
|
|---|
| 125 | SetHighEnergyLimit(highKinEnergy);
|
|---|
| 126 | }
|
|---|
| 127 |
|
|---|
| 128 | // peak energy
|
|---|
| 129 | epeak = std::min(model->PeakEnergy(), emax);
|
|---|
| 130 | peakKinEnergy = 0.5*epeak*epeak/electron_mass_c2 - 2.0*electron_mass_c2;
|
|---|
| 131 |
|
|---|
| 132 | if(verbose>0) {
|
|---|
| 133 | G4cout << "G4eeToHadronsModel::Initialise: " << G4endl;
|
|---|
| 134 | G4cout << "LabSystem: emin(GeV)= " << lowKinEnergy/GeV
|
|---|
| 135 | << " epeak(GeV)= " << peakKinEnergy/GeV
|
|---|
| 136 | << " emax(GeV)= " << highKinEnergy/GeV
|
|---|
| 137 | << G4endl;
|
|---|
| 138 | G4cout << "SM System: emin(MeV)= " << emin/MeV
|
|---|
| 139 | << " epeak(MeV)= " << epeak/MeV
|
|---|
| 140 | << " emax(MeV)= " << emax/MeV
|
|---|
| 141 | << G4endl;
|
|---|
| 142 | }
|
|---|
| 143 |
|
|---|
| 144 | if(lowKinEnergy < peakKinEnergy) {
|
|---|
| 145 | crossBornPerElectron = model->PhysicsVector(emin, emax);
|
|---|
| 146 | crossPerElectron = model->PhysicsVector(emin, emax);
|
|---|
| 147 | nbins = crossPerElectron->GetVectorLength();
|
|---|
| 148 | for(G4int i=0; i<nbins; i++) {
|
|---|
| 149 | G4double e = crossPerElectron->GetLowEdgeEnergy(i);
|
|---|
| 150 | G4double cs = model->ComputeCrossSection(e);
|
|---|
| 151 | crossBornPerElectron->PutValue(i, cs);
|
|---|
| 152 | }
|
|---|
| 153 | ComputeCMCrossSectionPerElectron();
|
|---|
| 154 | }
|
|---|
| 155 | if(verbose>1) {
|
|---|
| 156 | G4cout << "G4eeToHadronsModel: Cross secsions per electron"
|
|---|
| 157 | << " nbins= " << nbins
|
|---|
| 158 | << " emin(MeV)= " << emin/MeV
|
|---|
| 159 | << " emax(MeV)= " << emax/MeV
|
|---|
| 160 | << G4endl;
|
|---|
| 161 | G4bool b;
|
|---|
| 162 | for(G4int i=0; i<nbins; i++) {
|
|---|
| 163 | G4double e = crossPerElectron->GetLowEdgeEnergy(i);
|
|---|
| 164 | G4double s1 = crossPerElectron->GetValue(e, b);
|
|---|
| 165 | G4double s2 = crossBornPerElectron->GetValue(e, b);
|
|---|
| 166 | G4cout << "E(MeV)= " << e/MeV
|
|---|
| 167 | << " cross(nb)= " << s1/nanobarn
|
|---|
| 168 | << " crossBorn(nb)= " << s2/nanobarn
|
|---|
| 169 | << G4endl;
|
|---|
| 170 | }
|
|---|
| 171 | }
|
|---|
| 172 | }
|
|---|
| 173 |
|
|---|
| 174 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 175 |
|
|---|
| 176 | G4double G4eeToHadronsModel::CrossSectionPerVolume(
|
|---|
| 177 | const G4Material* mat,
|
|---|
| 178 | const G4ParticleDefinition* p,
|
|---|
| 179 | G4double kineticEnergy,
|
|---|
| 180 | G4double, G4double)
|
|---|
| 181 | {
|
|---|
| 182 | return mat->GetElectronDensity()*
|
|---|
| 183 | ComputeCrossSectionPerElectron(p, kineticEnergy);
|
|---|
| 184 | }
|
|---|
| 185 |
|
|---|
| 186 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 187 |
|
|---|
| 188 | G4double G4eeToHadronsModel::ComputeCrossSectionPerAtom(
|
|---|
| 189 | const G4ParticleDefinition* p,
|
|---|
| 190 | G4double kineticEnergy,
|
|---|
| 191 | G4double Z, G4double,
|
|---|
| 192 | G4double, G4double)
|
|---|
| 193 | {
|
|---|
| 194 | return Z*ComputeCrossSectionPerElectron(p, kineticEnergy);
|
|---|
| 195 | }
|
|---|
| 196 |
|
|---|
| 197 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 198 |
|
|---|
| 199 | G4double G4eeToHadronsModel::ComputeCrossSectionPerElectron(
|
|---|
| 200 | const G4ParticleDefinition*,
|
|---|
| 201 | G4double kineticEnergy,
|
|---|
| 202 | G4double, G4double)
|
|---|
| 203 | {
|
|---|
| 204 | G4double cross = 0.0;
|
|---|
| 205 | if(crossPerElectron) {
|
|---|
| 206 | G4bool b;
|
|---|
| 207 | G4double e = 2.0*electron_mass_c2*
|
|---|
| 208 | sqrt(1.0 + 0.5*kineticEnergy/electron_mass_c2);
|
|---|
| 209 | cross = crossPerElectron->GetValue(e, b);
|
|---|
| 210 | }
|
|---|
| 211 | // G4cout << "e= " << kineticEnergy << " cross= " << cross << G4endl;
|
|---|
| 212 | return cross;
|
|---|
| 213 | }
|
|---|
| 214 |
|
|---|
| 215 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 216 |
|
|---|
| 217 | void G4eeToHadronsModel::SampleSecondaries(std::vector<G4DynamicParticle*>* newp,
|
|---|
| 218 | const G4MaterialCutsCouple*,
|
|---|
| 219 | const G4DynamicParticle* dParticle,
|
|---|
| 220 | G4double,
|
|---|
| 221 | G4double)
|
|---|
| 222 | {
|
|---|
| 223 | if(crossPerElectron) {
|
|---|
| 224 | G4double t = dParticle->GetKineticEnergy();
|
|---|
| 225 | G4double e = 2.0*electron_mass_c2*sqrt(1.0 + 0.5*t/electron_mass_c2);
|
|---|
| 226 | G4LorentzVector inlv = dParticle->Get4Momentum();
|
|---|
| 227 | G4ThreeVector inBoost = inlv.boostVector();
|
|---|
| 228 | if(e > emin) {
|
|---|
| 229 | G4DynamicParticle* gamma = GenerateCMPhoton(e);
|
|---|
| 230 | G4LorentzVector gLv = gamma->Get4Momentum();
|
|---|
| 231 | G4LorentzVector lv(0.0,0.0,0.0,e);
|
|---|
| 232 | lv -= gLv;
|
|---|
| 233 | G4double m = lv.m();
|
|---|
| 234 | G4ThreeVector boost = lv.boostVector();
|
|---|
| 235 | const G4ThreeVector dir = gamma->GetMomentumDirection();
|
|---|
| 236 | model->SampleSecondaries(newp, m, dir);
|
|---|
| 237 | G4int np = newp->size();
|
|---|
| 238 | for(G4int j=0; j<np; j++) {
|
|---|
| 239 | G4DynamicParticle* dp = (*newp)[j];
|
|---|
| 240 | G4LorentzVector v = dp->Get4Momentum();
|
|---|
| 241 | v.boost(boost);
|
|---|
| 242 | v.boost(inBoost);
|
|---|
| 243 | dp->Set4Momentum(v);
|
|---|
| 244 | }
|
|---|
| 245 | gLv.boost(inBoost);
|
|---|
| 246 | gamma->Set4Momentum(gLv);
|
|---|
| 247 | newp->push_back(gamma);
|
|---|
| 248 | }
|
|---|
| 249 | }
|
|---|
| 250 | }
|
|---|
| 251 |
|
|---|
| 252 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 253 |
|
|---|
| 254 | void G4eeToHadronsModel::ComputeCMCrossSectionPerElectron()
|
|---|
| 255 | {
|
|---|
| 256 | G4bool b;
|
|---|
| 257 | for(G4int i=0; i<nbins; i++) {
|
|---|
| 258 | G4double e = crossPerElectron->GetLowEdgeEnergy(i);
|
|---|
| 259 | G4double cs = 0.0;
|
|---|
| 260 | if(i > 0) {
|
|---|
| 261 | G4double L = 2.0*log(e/electron_mass_c2);
|
|---|
| 262 | G4double bt = 2.0*fine_structure_const*(L - 1.0)/pi;
|
|---|
| 263 | G4double btm1= bt - 1.0;
|
|---|
| 264 | G4double del = 1. + fine_structure_const*(1.5*L + pi*pi/3. -2.)/pi;
|
|---|
| 265 | G4double s1 = crossBornPerElectron->GetValue(e, b);
|
|---|
| 266 | G4double e1 = crossPerElectron->GetLowEdgeEnergy(i-1);
|
|---|
| 267 | G4double x1 = 1. - e1/e;
|
|---|
| 268 | cs += s1*(del*pow(x1,bt) - bt*(x1 - 0.25*x1*x1));
|
|---|
| 269 | if(i > 1) {
|
|---|
| 270 | G4double e2 = e1;
|
|---|
| 271 | G4double x2 = x1;
|
|---|
| 272 | G4double s2 = crossBornPerElectron->GetValue(e2, b);
|
|---|
| 273 | G4double w2 = bt*(del*pow(x2,btm1) - 1.0 + 0.5*x2);
|
|---|
| 274 |
|
|---|
| 275 | for(G4int j=i-2; j>=0; j--) {
|
|---|
| 276 | e1 = crossPerElectron->GetLowEdgeEnergy(j);
|
|---|
| 277 | x1 = 1. - e1/e;
|
|---|
| 278 | G4double s1 = crossBornPerElectron->GetValue(e1, b);
|
|---|
| 279 | G4double w1 = bt*(del*pow(x1,btm1) - 1.0 + 0.5*x1);
|
|---|
| 280 | cs += 0.5*(x1 - x2)*(w2*s2 + w1*s1);
|
|---|
| 281 | e2 = e1;
|
|---|
| 282 | x2 = x1;
|
|---|
| 283 | s2 = s1;
|
|---|
| 284 | w2 = w1;
|
|---|
| 285 | }
|
|---|
| 286 | }
|
|---|
| 287 | }
|
|---|
| 288 | crossPerElectron->PutValue(i, cs);
|
|---|
| 289 | // G4cout << "e= " << e << " cs= " << cs << G4endl;
|
|---|
| 290 | }
|
|---|
| 291 | }
|
|---|
| 292 |
|
|---|
| 293 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 294 |
|
|---|
| 295 | G4DynamicParticle* G4eeToHadronsModel::GenerateCMPhoton(G4double e)
|
|---|
| 296 | {
|
|---|
| 297 | G4bool b;
|
|---|
| 298 | G4double x;
|
|---|
| 299 | G4DynamicParticle* gamma = 0;
|
|---|
| 300 | G4double L = 2.0*log(e/electron_mass_c2);
|
|---|
| 301 | G4double bt = 2.0*fine_structure_const*(L - 1.)/pi;
|
|---|
| 302 | G4double btm1= bt - 1.0;
|
|---|
| 303 | G4double del = 1. + fine_structure_const*(1.5*L + pi*pi/3. -2.)/pi;
|
|---|
| 304 |
|
|---|
| 305 | G4double s0 = crossBornPerElectron->GetValue(e, b);
|
|---|
| 306 | G4double de = (emax - emin)/(G4double)nbins;
|
|---|
| 307 | G4double x0 = min(de,e - emin)/e;
|
|---|
| 308 | G4double ds = crossBornPerElectron->GetValue(e, b)
|
|---|
| 309 | *(del*pow(x0,bt) - bt*(x0 - 0.25*x0*x0));
|
|---|
| 310 | G4double e1 = e*(1. - x0);
|
|---|
| 311 |
|
|---|
| 312 | if(e1 < emax && s0*G4UniformRand()<ds) {
|
|---|
| 313 | x = x0*pow(G4UniformRand(),1./bt);
|
|---|
| 314 | } else {
|
|---|
| 315 |
|
|---|
| 316 | x = 1. - e1/e;
|
|---|
| 317 | G4double s1 = crossBornPerElectron->GetValue(e1, b);
|
|---|
| 318 | G4double w1 = bt*(del*pow(x,btm1) - 1.0 + 0.5*x);
|
|---|
| 319 | G4double grej = s1*w1;
|
|---|
| 320 | G4double f;
|
|---|
| 321 | // G4cout << "e= " << e/GeV << " epeak= " << epeak/GeV
|
|---|
| 322 | // << " s1= " << s1 << " w1= " << w1
|
|---|
| 323 | // << " grej= " << grej << G4endl;
|
|---|
| 324 | // Above emax cross section is 0
|
|---|
| 325 | if(e1 > emax) {
|
|---|
| 326 | x = 1. - emax/e;
|
|---|
| 327 | G4double s2 = crossBornPerElectron->GetValue(emax, b);
|
|---|
| 328 | G4double w2 = bt*(del*pow(x,btm1) - 1.0 + 0.5*x);
|
|---|
| 329 | grej = s2*w2;
|
|---|
| 330 | // G4cout << "emax= " << emax << " s2= " << s2 << " w2= " << w2
|
|---|
| 331 | // << " grej= " << grej << G4endl;
|
|---|
| 332 | }
|
|---|
| 333 |
|
|---|
| 334 | if(e1 > epeak) {
|
|---|
| 335 | x = 1. - epeak/e;
|
|---|
| 336 | G4double s2 = crossBornPerElectron->GetValue(epeak, b);
|
|---|
| 337 | G4double w2 = bt*(del*pow(x,btm1) - 1.0 + 0.5*x);
|
|---|
| 338 | grej = max(grej,s2*w2);
|
|---|
| 339 | //G4cout << "epeak= " << epeak << " s2= " << s2 << " w2= " << w2
|
|---|
| 340 | // << " grej= " << grej << G4endl;
|
|---|
| 341 | }
|
|---|
| 342 | G4double xmin = 1. - e1/e;
|
|---|
| 343 | if(e1 > emax) xmin = 1. - emax/e;
|
|---|
| 344 | G4double xmax = 1. - emin/e;
|
|---|
| 345 | do {
|
|---|
| 346 | x = xmin + G4UniformRand()*(xmax - xmin);
|
|---|
| 347 | G4double s2 = crossBornPerElectron->GetValue((1.0 - x)*e, b);
|
|---|
| 348 | G4double w2 = bt*(del*pow(x,btm1) - 1.0 + 0.5*x);
|
|---|
| 349 | //G4cout << "x= " << x << " xmin= " << xmin << " xmax= " << xmax
|
|---|
| 350 | // << " s2= " << s2 << " w2= " << w2
|
|---|
| 351 | // << G4endl;
|
|---|
| 352 | f = s2*w2;
|
|---|
| 353 | if(f > grej) {
|
|---|
| 354 | G4cout << "G4DynamicParticle* G4eeToHadronsModel:WARNING "
|
|---|
| 355 | << f << " > " << grej << " majorant is`small!"
|
|---|
| 356 | << G4endl;
|
|---|
| 357 | }
|
|---|
| 358 | } while (f < grej*G4UniformRand());
|
|---|
| 359 | }
|
|---|
| 360 |
|
|---|
| 361 | G4ThreeVector dir(0.0,0.0,1.0);
|
|---|
| 362 | gamma = new G4DynamicParticle(theGamma,dir,x*e);
|
|---|
| 363 | return gamma;
|
|---|
| 364 | }
|
|---|
| 365 |
|
|---|
| 366 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 367 |
|
|---|