Geant4  10.01.p03
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 88977 2015-03-17 10:03:17Z gcosmo $
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 #include "G4Log.hh"
64 #include "G4Exp.hh"
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67 
68 using namespace std;
69 
71  const G4String& nam)
72  : G4VEmModel(nam),
73  model(mod),
74  crossPerElectron(0),
75  crossBornPerElectron(0),
76  isInitialised(false),
77  nbins(100),
78  verbose(ver)
79 {
86  epeak = emax;
87  //verbose = 1;
88 }
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91 
93 {
94  delete model;
95 }
96 
97 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
98 
100  const G4DataVector&)
101 {
102  if(isInitialised) { return; }
103  isInitialised = true;
104 
105  // CM system
106  emin = model->LowEnergy();
107  emax = model->HighEnergy();
108 
109  // peak energy
111 
112  if(verbose>0) {
113  G4cout << "G4eeToHadronsModel::Initialise: " << G4endl;
114  G4cout << "CM System: emin(MeV)= " << emin/MeV
115  << " epeak(MeV)= " << epeak/MeV
116  << " emax(MeV)= " << emax/MeV
117  << G4endl;
118  }
119 
123  for(G4int i=0; i<nbins; ++i) {
127  }
129 
130  if(verbose>1) {
131  G4cout << "G4eeToHadronsModel: Cross secsions per electron"
132  << " nbins= " << nbins
133  << " emin(MeV)= " << emin/MeV
134  << " emax(MeV)= " << emax/MeV
135  << G4endl;
136  for(G4int i=0; i<nbins; ++i) {
140  G4cout << "E(MeV)= " << e/MeV
141  << " cross(nb)= " << s1/nanobarn
142  << " crossBorn(nb)= " << s2/nanobarn
143  << G4endl;
144  }
145  }
146 }
147 
148 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
149 
151  const G4Material* mat,
152  const G4ParticleDefinition* p,
153  G4double kineticEnergy,
155 {
156  return mat->GetElectronDensity()*
157  ComputeCrossSectionPerElectron(p, kineticEnergy);
158 }
159 
160 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
161 
163  const G4ParticleDefinition* p,
164  G4double kineticEnergy,
165  G4double Z, G4double,
167 {
168  return Z*ComputeCrossSectionPerElectron(p, kineticEnergy);
169 }
170 
171 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
172 
174  const G4ParticleDefinition*,
177 {
178  G4double cross = 0.0;
179  if(crossPerElectron) {
180  cross = crossPerElectron->Value(energy);
181  }
182  return cross;
183 }
184 
185 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
186 
187 void G4eeToHadronsModel::SampleSecondaries(std::vector<G4DynamicParticle*>* newp,
188  const G4MaterialCutsCouple*,
189  const G4DynamicParticle* dParticle,
190  G4double,
191  G4double)
192 {
193  if(crossPerElectron) {
194  G4double t = dParticle->GetKineticEnergy() + 2*electron_mass_c2;
195  G4LorentzVector inlv = dParticle->Get4Momentum() +
196  G4LorentzVector(0.0,0.0,0.0,electron_mass_c2);
197  G4double e = inlv.m();
198  G4ThreeVector inBoost = inlv.boostVector();
199  //G4cout << "G4eeToHadronsModel::SampleSecondaries e= " << e
200  // << " " << inlv << " " << inBoost <<G4endl;
201  if(e > emin) {
203  G4LorentzVector gLv = gamma->Get4Momentum();
204  G4LorentzVector lv(0.0,0.0,0.0,e);
205  lv -= gLv;
206  G4double mass = lv.m();
207  //G4cout << "mass= " << mass << " " << lv << G4endl;
208  G4ThreeVector boost = lv.boostVector();
209  //G4cout << "mass= " << mass << " " << boost << G4endl;
210  const G4ThreeVector dir = gamma->GetMomentumDirection();
211  model->SampleSecondaries(newp, mass, dir);
212  G4int np = newp->size();
213  for(G4int j=0; j<np; ++j) {
214  G4DynamicParticle* dp = (*newp)[j];
215  G4LorentzVector v = dp->Get4Momentum();
216  v.boost(boost);
217  //G4cout << j << ". " << v << G4endl;
218  v.boost(inBoost);
219  //G4cout << " " << v << G4endl;
220  dp->Set4Momentum(v);
221  t -= v.e();
222  }
223  //G4cout << "Gamma " << gLv << G4endl;
224  gLv.boost(inBoost);
225  //G4cout << " " << gLv << G4endl;
226  gamma->Set4Momentum(gLv);
227  t -= gLv.e();
228  newp->push_back(gamma);
229  if(fabs(t) > MeV) {
230  G4cout << "G4eeToHadronsModel::SampleSecondaries: Ebalance(MeV)= "
231  << t/MeV << " primary 4-momentum: " << inlv << G4endl;
232  }
233  }
234  }
235 }
236 
237 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
238 
240 {
241  for(G4int i=0; i<nbins; i++) {
243  G4double cs = 0.0;
244  if(i > 0) {
245  G4double L = 2.0*G4Log(e/electron_mass_c2);
246  G4double bt = 2.0*fine_structure_const*(L - 1.0)/pi;
247  G4double btm1= bt - 1.0;
248  G4double del = 1. + fine_structure_const*(1.5*L + pi*pi/3. -2.)/pi;
251  G4double x1 = 1. - e1/e;
252  cs += s1*(del*G4Exp(G4Log(x1)*bt) - bt*(x1 - 0.25*x1*x1));
253  if(i > 1) {
254  G4double e2 = e1;
255  G4double x2 = x1;
257  G4double w2 = bt*(del*G4Exp(G4Log(x2)*btm1) - 1.0 + 0.5*x2);
258  G4double w1;
259 
260  for(G4int j=i-2; j>=0; --j) {
261  e1 = crossPerElectron->Energy(j);
262  x1 = 1. - e1/e;
263  s1 = crossBornPerElectron->Value(e1);
264  w1 = bt*(del*G4Exp(G4Log(x1)*btm1) - 1.0 + 0.5*x1);
265  cs += 0.5*(x1 - x2)*(w2*s2 + w1*s1);
266  e2 = e1;
267  x2 = x1;
268  s2 = s1;
269  w2 = w1;
270  }
271  }
272  }
273  crossPerElectron->PutValue(i, cs);
274  }
275 }
276 
277 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
278 
280 {
281  G4double x;
282  G4DynamicParticle* gamma = 0;
283  G4double L = 2.0*G4Log(e/electron_mass_c2);
284  G4double bt = 2.0*fine_structure_const*(L - 1.)/pi;
285  G4double btm1= bt - 1.0;
286  G4double del = 1. + fine_structure_const*(1.5*L + pi*pi/3. -2.)/pi;
287 
289  G4double de = (emax - emin)/(G4double)nbins;
290  G4double xmax = 0.5*(1.0 - (emin*emin)/(e*e));
291  G4double xmin = std::min(de/e, xmax);
292  G4double ds = s0*(del*G4Exp(G4Log(xmin)*bt) - bt*(xmin - 0.25*xmin*xmin));
293  G4double e1 = e*(1. - xmin);
294 
295  //G4cout << "e1= " << e1 << G4endl;
296  if(e1 < emax && s0*G4UniformRand()<ds) {
297  x = xmin*G4Exp(G4Log(G4UniformRand())/bt);
298  } else {
299 
300  x = xmin;
302  G4double w1 = bt*(del*G4Exp(G4Log(x)*btm1) - 1.0 + 0.5*x);
303  G4double grej = s1*w1;
304  G4double f;
305  /*
306  G4cout << "e(GeV)= " << e/GeV << " epeak(GeV)= " << epeak/GeV
307  << " s1= " << s1 << " w1= " << w1
308  << " grej= " << grej << G4endl;
309  */
310  // Above emax cross section is const
311  if(e1 > emax) {
312  x = 0.5*(1. - (emax*emax)/(e*e));
314  G4double w2 = bt*(del*G4Exp(G4Log(x)*btm1) - 1.0 + 0.5*x);
315  grej = s2*w2;
316  //G4cout << "emax= " << emax << " s2= " << s2 << " w2= " << w2
317  // << " grej= " << grej << G4endl;
318  }
319 
320  if(e1 > epeak) {
321  x = 0.5*(1.0 - (epeak*epeak)/(e*e));
323  G4double w2 = bt*(del*G4Exp(G4Log(x)*btm1) - 1.0 + 0.5*x);
324  grej = std::max(grej,s2*w2);
325  //G4cout << "epeak= " << epeak << " s2= " << s2 << " w2= " << w2
326  // << " grej= " << grej << G4endl;
327  }
328  do {
329  x = xmin + G4UniformRand()*(xmax - xmin);
330 
331  G4double s2 = crossBornPerElectron->Value(sqrt(1.0 - 2*x)*e);
332  G4double w2 = bt*(del*G4Exp(G4Log(x)*btm1) - 1.0 + 0.5*x);
333  /*
334  G4cout << "x= " << x << " xmin= " << xmin << " xmax= " << xmax
335  << " s2= " << s2 << " w2= " << w2 << G4endl;
336  */
337  f = s2*w2;
338  if(f > grej) {
339  G4cout << "G4DynamicParticle* G4eeToHadronsModel:WARNING "
340  << f << " > " << grej << " majorant is`small!"
341  << G4endl;
342  }
343  } while (f < grej*G4UniformRand());
344  }
345 
346  G4ThreeVector dir(0.0,0.0,1.0);
347  if(G4UniformRand() > 0.5) { dir.set(0.0,0.0,-1.0); }
348  //G4cout << "Egamma(MeV)= " << x*e << " " << dir << G4endl;
349  gamma = new G4DynamicParticle(theGamma,dir,x*e);
350  return gamma;
351 }
352 
353 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
354 
G4double LowEnergy() const
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:623
static const double MeV
Definition: G4SIunits.hh:193
static const double nanobarn
Definition: G4SIunits.hh:98
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:616
G4eeToHadronsModel(G4Vee2hadrons *, G4int ver=0, const G4String &nam="eeToHadrons")
G4PhysicsVector * PhysicsVector() const
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 &)
size_t GetVectorLength() 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)
static const G4int L[nN]
G4PhysicsVector * crossPerElectron
#define G4UniformRand()
Definition: Randomize.hh:93
G4GLOB_DLL std::ostream G4cout
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, G4double, const G4ThreeVector &)=0
G4double GetElectronDensity() const
Definition: G4Material.hh:217
const G4ThreeVector & GetMomentumDirection() const
void PutValue(size_t index, G4double theValue)
G4PhysicsVector * crossBornPerElectron
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4double Energy(size_t index) const
G4double Value(G4double theEnergy, size_t &lastidx) const
G4LorentzVector Get4Momentum() const
static const G4double e1
virtual G4double PeakEnergy() const =0
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
void Set4Momentum(const G4LorentzVector &momentum)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double energy(const ThreeVector &p, const G4double m)
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
CLHEP::HepLorentzVector G4LorentzVector