Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4eeToHadronsModel
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 12.08.2003
37//
38// Modifications:
39// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
40// 18-05-05 Use optimized interfaces (V.Ivantchenko)
41//
42//
43// -------------------------------------------------------------------
44//
45
46
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49
50#include "G4eeToHadronsModel.hh"
51#include "Randomize.hh"
53#include "G4SystemOfUnits.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#include "G4Log.hh"
63#include "G4Exp.hh"
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66
67using namespace std;
68
70 const G4String& nam)
71 : G4VEmModel(nam),
72 model(mod),
73 verbose(ver)
74{
75 theGamma = G4Gamma::Gamma();
76 highKinEnergy = HighEnergyLimit();
77 lowKinEnergy = LowEnergyLimit();
78 emin = lowKinEnergy;
79 emax = highKinEnergy;
80 peakKinEnergy = highKinEnergy;
81 epeak = emax;
82 //verbose = 1;
83}
84
85//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86
88{
89 delete model;
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93
95 const G4DataVector&)
96{
97 if(isInitialised) { return; }
98 isInitialised = true;
99
100 // CM system
101 emin = model->LowEnergy();
102 emax = model->HighEnergy();
103
104 // peak energy
105 epeak = std::min(model->PeakEnergy(), emax);
106
107 if(verbose>0) {
108 G4cout << "G4eeToHadronsModel::Initialise: " << G4endl;
109 G4cout << "CM System: emin(MeV)= " << emin/MeV
110 << " epeak(MeV)= " << epeak/MeV
111 << " emax(MeV)= " << emax/MeV
112 << G4endl;
113 }
114
115 crossBornPerElectron = model->PhysicsVector();
116 crossPerElectron = model->PhysicsVector();
117 nbins = (G4int)crossPerElectron->GetVectorLength();
118 for(G4int i=0; i<nbins; ++i) {
119 G4double e = crossPerElectron->Energy(i);
120 G4double cs = model->ComputeCrossSection(e);
121 crossBornPerElectron->PutValue(i, cs);
122 }
123 ComputeCMCrossSectionPerElectron();
124
125 if(verbose>1) {
126 G4cout << "G4eeToHadronsModel: Cross sections per electron"
127 << " nbins= " << nbins
128 << " emin(MeV)= " << emin/MeV
129 << " emax(MeV)= " << emax/MeV
130 << G4endl;
131 for(G4int i=0; i<nbins; ++i) {
132 G4double e = crossPerElectron->Energy(i);
133 G4double s1 = crossPerElectron->Value(e);
134 G4double s2 = crossBornPerElectron->Value(e);
135 G4cout << "E(MeV)= " << e/MeV
136 << " cross(nb)= " << s1/nanobarn
137 << " crossBorn(nb)= " << s2/nanobarn
138 << G4endl;
139 }
140 }
141}
142
143//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
144
146 const G4Material* mat,
147 const G4ParticleDefinition* p,
148 G4double kineticEnergy,
150{
151 return mat->GetElectronDensity()*
152 ComputeCrossSectionPerElectron(p, kineticEnergy);
153}
154
155//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
156
158 const G4ParticleDefinition* p,
159 G4double kineticEnergy,
162{
163 return Z*ComputeCrossSectionPerElectron(p, kineticEnergy);
164}
165
166//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
167
170 G4double energy,
172{
173 return crossPerElectron->Value(energy);
174}
175
176//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
177
178void G4eeToHadronsModel::SampleSecondaries(std::vector<G4DynamicParticle*>* newp,
180 const G4DynamicParticle* dParticle,
181 G4double,
182 G4double)
183{
184 G4double t = dParticle->GetKineticEnergy() + 2*electron_mass_c2;
185 G4LorentzVector inlv = dParticle->Get4Momentum() +
186 G4LorentzVector(0.0,0.0,0.0,electron_mass_c2);
187 G4double e = inlv.m();
188 G4ThreeVector inBoost = inlv.boostVector();
189 //G4cout << "G4eeToHadronsModel::SampleSecondaries e= " << e
190 // << " " << inlv << " " << inBoost <<G4endl;
191 if(e > emin) {
193 G4LorentzVector gLv = gamma->Get4Momentum();
194 G4LorentzVector lv(0.0,0.0,0.0,e);
195 lv -= gLv;
196 G4double mass = lv.m();
197 //G4cout << "mass= " << mass << " " << lv << G4endl;
198 G4ThreeVector boost = lv.boostVector();
199 //G4cout << "mass= " << mass << " " << boost << G4endl;
200 const G4ThreeVector dir = gamma->GetMomentumDirection();
201 model->SampleSecondaries(newp, mass, dir);
202 std::size_t np = newp->size();
203 for(std::size_t j=0; j<np; ++j) {
204 G4DynamicParticle* dp = (*newp)[j];
206 v.boost(boost);
207 //G4cout << j << ". " << v << G4endl;
208 v.boost(inBoost);
209 //G4cout << " " << v << G4endl;
210 dp->Set4Momentum(v);
211 t -= v.e();
212 }
213 //G4cout << "Gamma " << gLv << G4endl;
214 gLv.boost(inBoost);
215 //G4cout << " " << gLv << G4endl;
216 gamma->Set4Momentum(gLv);
217 t -= gLv.e();
218 newp->push_back(gamma);
219 if(std::abs(t) > CLHEP::MeV) {
220 G4cout << "G4eeToHadronsModel::SampleSecondaries: Ebalance(MeV)= "
221 << t/MeV << " primary 4-momentum: " << inlv << G4endl;
222 }
223 }
224}
225
226//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
227
228void G4eeToHadronsModel::ComputeCMCrossSectionPerElectron()
229{
230 for(G4int i=0; i<nbins; i++) {
231 G4double e = crossPerElectron->Energy(i);
232 G4double cs = 0.0;
233 if(i > 0) {
234 G4double LL = 2.0*G4Log(e/electron_mass_c2);
235 G4double bt = 2.0*fine_structure_const*(LL - 1.0)/pi;
236 G4double btm1= bt - 1.0;
237 G4double del = 1. + fine_structure_const*(1.5*LL + pi*pi/3. -2.)/pi;
238 G4double s1 = crossBornPerElectron->Value(e);
239 G4double e1 = crossPerElectron->Energy(i-1);
240 G4double x1 = 1. - e1/e;
241 cs += s1*(del*G4Exp(G4Log(x1)*bt) - bt*(x1 - 0.25*x1*x1));
242 if(i > 1) {
243 G4double e2 = e1;
244 G4double x2 = x1;
245 G4double s2 = crossBornPerElectron->Value(e2);
246 G4double w2 = bt*(del*G4Exp(G4Log(x2)*btm1) - 1.0 + 0.5*x2);
247 G4double w1;
248
249 for(G4int j=i-2; j>=0; --j) {
250 e1 = crossPerElectron->Energy(j);
251 x1 = 1. - e1/e;
252 s1 = crossBornPerElectron->Value(e1);
253 w1 = bt*(del*G4Exp(G4Log(x1)*btm1) - 1.0 + 0.5*x1);
254 cs += 0.5*(x1 - x2)*(w2*s2 + w1*s1);
255 e2 = e1;
256 x2 = x1;
257 s2 = s1;
258 w2 = w1;
259 }
260 }
261 }
262 crossPerElectron->PutValue(i, cs);
263 }
264}
265
266//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
267
269{
270 G4double x;
271 G4DynamicParticle* gamma = nullptr;
272 G4double LL = 2.0*G4Log(e/electron_mass_c2);
273 G4double bt = 2.0*fine_structure_const*(LL - 1.)/pi;
274 G4double btm1= bt - 1.0;
275 G4double del = 1. + fine_structure_const*(1.5*LL + pi*pi/3. -2.)/pi;
276
277 G4double s0 = crossBornPerElectron->Value(e);
278 G4double de = (emax - emin)/(G4double)nbins;
279 G4double xmax = 0.5*(1.0 - (emin*emin)/(e*e));
280 G4double xmin = std::min(de/e, xmax);
281 G4double ds = s0*(del*G4Exp(G4Log(xmin)*bt) - bt*(xmin - 0.25*xmin*xmin));
282 G4double e1 = e*(1. - xmin);
283
284 //G4cout << "e1= " << e1 << G4endl;
285 if(e1 < emax && s0*G4UniformRand()<ds) {
286 x = xmin*G4Exp(G4Log(G4UniformRand())/bt);
287 } else {
288
289 x = xmin;
290 G4double s1 = crossBornPerElectron->Value(e1);
291 G4double w1 = bt*(del*G4Exp(G4Log(x)*btm1) - 1.0 + 0.5*x);
292 G4double grej = s1*w1;
293 G4double f;
294 /*
295 G4cout << "e(GeV)= " << e/GeV << " epeak(GeV)= " << epeak/GeV
296 << " s1= " << s1 << " w1= " << w1
297 << " grej= " << grej << G4endl;
298 */
299 // Above emax cross section is const
300 if(e1 > emax) {
301 x = 0.5*(1. - (emax*emax)/(e*e));
302 G4double s2 = crossBornPerElectron->Value(emax);
303 G4double w2 = bt*(del*G4Exp(G4Log(x)*btm1) - 1.0 + 0.5*x);
304 grej = s2*w2;
305 //G4cout << "emax= " << emax << " s2= " << s2 << " w2= " << w2
306 // << " grej= " << grej << G4endl;
307 }
308
309 if(e1 > epeak) {
310 x = 0.5*(1.0 - (epeak*epeak)/(e*e));
311 G4double s2 = crossBornPerElectron->Value(epeak);
312 G4double w2 = bt*(del*G4Exp(G4Log(x)*btm1) - 1.0 + 0.5*x);
313 grej = std::max(grej,s2*w2);
314 //G4cout << "epeak= " << epeak << " s2= " << s2 << " w2= " << w2
315 // << " grej= " << grej << G4endl;
316 }
317 G4int ii = 0;
318 const G4int iimax = 1000;
319 do {
320 x = xmin + G4UniformRand()*(xmax - xmin);
321
322 G4double s2 = crossBornPerElectron->Value(sqrt(1.0 - 2*x)*e);
323 G4double w2 = bt*(del*G4Exp(G4Log(x)*btm1) - 1.0 + 0.5*x);
324 /*
325 G4cout << "x= " << x << " xmin= " << xmin << " xmax= " << xmax
326 << " s2= " << s2 << " w2= " << w2 << G4endl;
327 */
328 f = s2*w2;
329 if(f > grej) {
330 G4cout << "G4DynamicParticle* G4eeToHadronsModel:WARNING "
331 << f << " > " << grej << " majorant is`small!"
332 << G4endl;
333 }
334 if(++ii >= iimax) { break; }
335 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
336 } while (f < grej*G4UniformRand());
337 }
338
339 G4ThreeVector dir(0.0,0.0,1.0);
340 if(G4UniformRand() > 0.5) { dir.set(0.0,0.0,-1.0); }
341 //G4cout << "Egamma(MeV)= " << x*e << " " << dir << G4endl;
342 gamma = new G4DynamicParticle(theGamma,dir,x*e);
343 return gamma;
344}
345
346//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4double G4Log(G4double x)
Definition: G4Log.hh:227
CLHEP::HepLorentzVector G4LorentzVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
void set(double x, double y, double z)
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
const G4ThreeVector & GetMomentumDirection() const
G4LorentzVector Get4Momentum() const
G4double GetKineticEnergy() const
void Set4Momentum(const G4LorentzVector &momentum)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
G4double GetElectronDensity() const
Definition: G4Material.hh:212
void PutValue(const std::size_t index, const G4double value)
G4double Energy(const std::size_t index) const
G4double Value(const G4double energy, std::size_t &lastidx) const
std::size_t GetVectorLength() const
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
G4PhysicsVector * PhysicsVector() const
G4double LowEnergy() const
virtual G4double ComputeCrossSection(G4double) const =0
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, G4double, const G4ThreeVector &)=0
virtual G4double PeakEnergy() const =0
G4double HighEnergy() const
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double maxEnergy=DBL_MAX) override
G4DynamicParticle * GenerateCMPhoton(G4double)
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) override
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
~G4eeToHadronsModel() override
G4eeToHadronsModel(G4Vee2hadrons *, G4int ver=0, const G4String &nam="eeToHadrons")