Geant4  10.00.p02
G4eeToHadronsModel.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 // $Id: G4eeToHadronsModel.cc 66241 2012-12-13 18:34:42Z gunter $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4eeToHadronsModel
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 12.08.2003
38 //
39 // Modifications:
40 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
41 // 18-05-05 Use optimized interfaces (V.Ivantchenko)
42 //
43 //
44 // -------------------------------------------------------------------
45 //
46 
47 
48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
50 
51 #include "G4eeToHadronsModel.hh"
52 #include "Randomize.hh"
53 #include "G4PhysicalConstants.hh"
54 #include "G4SystemOfUnits.hh"
55 #include "G4Electron.hh"
56 #include "G4Gamma.hh"
57 #include "G4Positron.hh"
58 #include "G4PionPlus.hh"
59 #include "Randomize.hh"
60 #include "G4Vee2hadrons.hh"
61 #include "G4PhysicsVector.hh"
62 #include "G4PhysicsLogVector.hh"
63 
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65 
66 using namespace std;
67 
69  const G4String& nam)
70  : G4VEmModel(nam),
71  model(mod),
72  crossPerElectron(0),
73  crossBornPerElectron(0),
74  isInitialised(false),
75  nbins(100),
76  verbose(ver)
77 {
84  epeak = emax;
85 }
86 
87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88 
90 {
91  delete model;
92  delete crossPerElectron;
93  delete crossBornPerElectron;
94 }
95 
96 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
97 
99  const G4DataVector&)
100 {
101  if(isInitialised) { return; }
102  isInitialised = true;
103 
104  // Lab system
107 
108  // CM system
109  emin = model->LowEnergy();
110  emax = model->HighEnergy();
111 
112  G4double emin0 =
113  2.0*electron_mass_c2*sqrt(1.0 + 0.5*lowKinEnergy/electron_mass_c2);
114  G4double emax0 =
115  2.0*electron_mass_c2*sqrt(1.0 + 0.5*highKinEnergy/electron_mass_c2);
116 
117  // recompute low energy
118  if(emin0 > emax) {
119  emin0 = emax;
120  model->SetLowEnergy(emin0);
121  }
122  if(emin > emin0) {
123  emin0 = emin;
124  lowKinEnergy = 0.5*emin*emin/electron_mass_c2 - 2.0*electron_mass_c2;
126  }
127 
128  // recompute high energy
129  if(emax < emax0) {
130  emax0 = emax;
131  highKinEnergy = 0.5*emax*emax/electron_mass_c2 - 2.0*electron_mass_c2;
133  }
134 
135  // peak energy
137  peakKinEnergy = 0.5*epeak*epeak/electron_mass_c2 - 2.0*electron_mass_c2;
138 
139  if(verbose>0) {
140  G4cout << "G4eeToHadronsModel::Initialise: " << G4endl;
141  G4cout << "LabSystem: emin(GeV)= " << lowKinEnergy/GeV
142  << " epeak(GeV)= " << peakKinEnergy/GeV
143  << " emax(GeV)= " << highKinEnergy/GeV
144  << G4endl;
145  G4cout << "SM System: emin(MeV)= " << emin/MeV
146  << " epeak(MeV)= " << epeak/MeV
147  << " emax(MeV)= " << emax/MeV
148  << G4endl;
149  }
150 
151  if(lowKinEnergy < peakKinEnergy) {
155  for(G4int i=0; i<nbins; i++) {
159  }
161  }
162  if(verbose>1) {
163  G4cout << "G4eeToHadronsModel: Cross secsions per electron"
164  << " nbins= " << nbins
165  << " emin(MeV)= " << emin/MeV
166  << " emax(MeV)= " << emax/MeV
167  << G4endl;
168  G4bool b;
169  for(G4int i=0; i<nbins; i++) {
171  G4double s1 = crossPerElectron->GetValue(e, b);
173  G4cout << "E(MeV)= " << e/MeV
174  << " cross(nb)= " << s1/nanobarn
175  << " crossBorn(nb)= " << s2/nanobarn
176  << G4endl;
177  }
178  }
179 }
180 
181 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
182 
184  const G4Material* mat,
185  const G4ParticleDefinition* p,
186  G4double kineticEnergy,
188 {
189  return mat->GetElectronDensity()*
190  ComputeCrossSectionPerElectron(p, kineticEnergy);
191 }
192 
193 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
194 
196  const G4ParticleDefinition* p,
197  G4double kineticEnergy,
198  G4double Z, G4double,
200 {
201  return Z*ComputeCrossSectionPerElectron(p, kineticEnergy);
202 }
203 
204 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
205 
207  const G4ParticleDefinition*,
208  G4double kineticEnergy,
210 {
211  G4double cross = 0.0;
212  if(crossPerElectron) {
213  G4bool b;
214  G4double e = 2.0*electron_mass_c2*
215  sqrt(1.0 + 0.5*kineticEnergy/electron_mass_c2);
216  cross = crossPerElectron->GetValue(e, b);
217  }
218  // G4cout << "e= " << kineticEnergy << " cross= " << cross << G4endl;
219  return cross;
220 }
221 
222 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
223 
224 void G4eeToHadronsModel::SampleSecondaries(std::vector<G4DynamicParticle*>* newp,
225  const G4MaterialCutsCouple*,
226  const G4DynamicParticle* dParticle,
227  G4double,
228  G4double)
229 {
230  if(crossPerElectron) {
231  G4double t = dParticle->GetKineticEnergy();
232  G4double e = 2.0*electron_mass_c2*sqrt(1.0 + 0.5*t/electron_mass_c2);
233  G4LorentzVector inlv = dParticle->Get4Momentum();
234  G4ThreeVector inBoost = inlv.boostVector();
235  if(e > emin) {
237  G4LorentzVector gLv = gamma->Get4Momentum();
238  G4LorentzVector lv(0.0,0.0,0.0,e);
239  lv -= gLv;
240  G4double mass = lv.m();
241  G4ThreeVector boost = lv.boostVector();
242  const G4ThreeVector dir = gamma->GetMomentumDirection();
243  model->SampleSecondaries(newp, mass, dir);
244  G4int np = newp->size();
245  for(G4int j=0; j<np; j++) {
246  G4DynamicParticle* dp = (*newp)[j];
247  G4LorentzVector v = dp->Get4Momentum();
248  v.boost(boost);
249  v.boost(inBoost);
250  dp->Set4Momentum(v);
251  }
252  gLv.boost(inBoost);
253  gamma->Set4Momentum(gLv);
254  newp->push_back(gamma);
255  }
256  }
257 }
258 
259 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
260 
262 {
263  G4bool b;
264  for(G4int i=0; i<nbins; i++) {
266  G4double cs = 0.0;
267  if(i > 0) {
268  G4double L = 2.0*log(e/electron_mass_c2);
269  G4double bt = 2.0*fine_structure_const*(L - 1.0)/pi;
270  G4double btm1= bt - 1.0;
271  G4double del = 1. + fine_structure_const*(1.5*L + pi*pi/3. -2.)/pi;
274  G4double x1 = 1. - e1/e;
275  cs += s1*(del*pow(x1,bt) - bt*(x1 - 0.25*x1*x1));
276  if(i > 1) {
277  G4double e2 = e1;
278  G4double x2 = x1;
280  G4double w2 = bt*(del*pow(x2,btm1) - 1.0 + 0.5*x2);
281  G4double w1;
282 
283  for(G4int j=i-2; j>=0; j--) {
285  x1 = 1. - e1/e;
286  s1 = crossBornPerElectron->GetValue(e1, b);
287  w1 = bt*(del*pow(x1,btm1) - 1.0 + 0.5*x1);
288  cs += 0.5*(x1 - x2)*(w2*s2 + w1*s1);
289  e2 = e1;
290  x2 = x1;
291  s2 = s1;
292  w2 = w1;
293  }
294  }
295  }
296  crossPerElectron->PutValue(i, cs);
297  // G4cout << "e= " << e << " cs= " << cs << G4endl;
298  }
299 }
300 
301 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
302 
304 {
305  G4bool b;
306  G4double x;
307  G4DynamicParticle* gamma = 0;
308  G4double L = 2.0*log(e/electron_mass_c2);
309  G4double bt = 2.0*fine_structure_const*(L - 1.)/pi;
310  G4double btm1= bt - 1.0;
311  G4double del = 1. + fine_structure_const*(1.5*L + pi*pi/3. -2.)/pi;
312 
314  G4double de = (emax - emin)/(G4double)nbins;
315  G4double x0 = min(de,e - emin)/e;
317  *(del*pow(x0,bt) - bt*(x0 - 0.25*x0*x0));
318  G4double e1 = e*(1. - x0);
319 
320  if(e1 < emax && s0*G4UniformRand()<ds) {
321  x = x0*pow(G4UniformRand(),1./bt);
322  } else {
323 
324  x = 1. - e1/e;
326  G4double w1 = bt*(del*pow(x,btm1) - 1.0 + 0.5*x);
327  G4double grej = s1*w1;
328  G4double f;
329  // G4cout << "e= " << e/GeV << " epeak= " << epeak/GeV
330  // << " s1= " << s1 << " w1= " << w1
331  // << " grej= " << grej << G4endl;
332  // Above emax cross section is 0
333  if(e1 > emax) {
334  x = 1. - emax/e;
336  G4double w2 = bt*(del*pow(x,btm1) - 1.0 + 0.5*x);
337  grej = s2*w2;
338  // G4cout << "emax= " << emax << " s2= " << s2 << " w2= " << w2
339  // << " grej= " << grej << G4endl;
340  }
341 
342  if(e1 > epeak) {
343  x = 1. - epeak/e;
345  G4double w2 = bt*(del*pow(x,btm1) - 1.0 + 0.5*x);
346  grej = max(grej,s2*w2);
347  //G4cout << "epeak= " << epeak << " s2= " << s2 << " w2= " << w2
348  // << " grej= " << grej << G4endl;
349  }
350  G4double xmin = 1. - e1/e;
351  if(e1 > emax) xmin = 1. - emax/e;
352  G4double xmax = 1. - emin/e;
353  do {
354  x = xmin + G4UniformRand()*(xmax - xmin);
355  G4double s2 = crossBornPerElectron->GetValue((1.0 - x)*e, b);
356  G4double w2 = bt*(del*pow(x,btm1) - 1.0 + 0.5*x);
357  //G4cout << "x= " << x << " xmin= " << xmin << " xmax= " << xmax
358  // << " s2= " << s2 << " w2= " << w2
359  // << G4endl;
360  f = s2*w2;
361  if(f > grej) {
362  G4cout << "G4DynamicParticle* G4eeToHadronsModel:WARNING "
363  << f << " > " << grej << " majorant is`small!"
364  << G4endl;
365  }
366  } while (f < grej*G4UniformRand());
367  }
368 
369  G4ThreeVector dir(0.0,0.0,1.0);
370  gamma = new G4DynamicParticle(theGamma,dir,x*e);
371  return gamma;
372 }
373 
374 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
375 
G4double LowEnergy() const
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
static const double MeV
Definition: G4SIunits.hh:193
G4double GetValue(G4double theEnergy, G4bool &isOutRange) const
static const double nanobarn
Definition: G4SIunits.hh:98
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
void SetLowEnergy(G4double val)
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
G4eeToHadronsModel(G4Vee2hadrons *, G4int ver=0, const G4String &nam="eeToHadrons")
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double maxEnergy=DBL_MAX)
const G4double pi
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
static const G4double e2
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4PhysicsVector * PhysicsVector(G4double, G4double) const =0
size_t GetVectorLength() const
G4double GetLowEdgeEnergy(size_t binNumber) const
int G4int
Definition: G4Types.hh:78
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
static const G4int L[nN]
G4PhysicsVector * crossPerElectron
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, G4double, const G4ThreeVector &)=0
G4double GetElectronDensity() const
Definition: G4Material.hh:215
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
void PutValue(size_t index, G4double theValue)
G4PhysicsVector * crossBornPerElectron
static const double GeV
Definition: G4SIunits.hh:196
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4LorentzVector Get4Momentum() const
static const G4double e1
virtual G4double PeakEnergy() const =0
void Set4Momentum(const G4LorentzVector &momentum)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4DynamicParticle * GenerateCMPhoton(G4double)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
#define G4endl
Definition: G4ios.hh:61
G4double HighEnergy() const
G4ParticleDefinition * theGamma
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
double G4double
Definition: G4Types.hh:76
virtual G4double ComputeCrossSection(G4double) const =0
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
CLHEP::HepLorentzVector G4LorentzVector