Geant4  10.02.p02
G4SeltzerBergerModel.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: G4SeltzerBergerModel.cc 93567 2015-10-26 14:51:41Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4SeltzerBergerModel
34 //
35 // Author: Vladimir Ivanchenko use inheritance from Andreas Schaelicke
36 // base class implementing ultra relativistic bremsstrahlung
37 // model
38 //
39 // Creation date: 04.10.2011
40 //
41 // Modifications:
42 //
43 // -------------------------------------------------------------------
44 //
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47 
48 #include "G4SeltzerBergerModel.hh"
49 #include "G4PhysicalConstants.hh"
50 #include "G4SystemOfUnits.hh"
51 #include "G4Electron.hh"
52 #include "G4Positron.hh"
53 #include "G4Gamma.hh"
54 #include "Randomize.hh"
55 #include "G4Material.hh"
56 #include "G4Element.hh"
57 #include "G4ElementVector.hh"
58 #include "G4ProductionCutsTable.hh"
60 #include "G4ModifiedTsai.hh"
61 
62 #include "G4Physics2DVector.hh"
63 #include "G4Exp.hh"
64 #include "G4Log.hh"
65 
66 #include "G4ios.hh"
67 #include <fstream>
68 #include <iomanip>
69 
70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
71 
72 using namespace std;
73 
77 
78 static const G4double emaxlog = 4*G4Log(10.);
79 static const G4double alpha = CLHEP::twopi*CLHEP::fine_structure_const;
80 static const G4double epeaklimit= 300*CLHEP::MeV;
81 static const G4double elowlimit = 20*CLHEP::keV;
82 
84  const G4String& nam)
85  : G4eBremsstrahlungRelModel(p,nam),useBicubicInterpolation(false)
86 {
89  SetLPMFlag(false);
90  nwarn = 0;
91  idx = idy = 0;
92 }
93 
94 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
95 
97 {
98  if(IsMaster()) {
99  for(size_t i=0; i<101; ++i) {
100  if(dataSB[i]) {
101  delete dataSB[i];
102  dataSB[i] = nullptr;
103  }
104  }
105  }
106 }
107 
108 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
109 
111  const G4DataVector& cuts)
112 {
113  // Access to elements
114  if(IsMaster()) {
115 
116  // check environment variable
117  // Build the complete string identifying the file with the data set
118  char* path = getenv("G4LEDATA");
119 
120  const G4ElementTable* theElmTable = G4Element::GetElementTable();
121  size_t numOfElm = G4Element::GetNumberOfElements();
122  if(numOfElm > 0) {
123  for(size_t i=0; i<numOfElm; ++i) {
124  G4int Z = G4lrint(((*theElmTable)[i])->GetZ());
125  if(Z < 1) { Z = 1; }
126  else if(Z > 100) { Z = 100; }
127  //G4cout << "Z= " << Z << G4endl;
128  // Initialisation
129  if(nullptr == dataSB[Z]) { ReadData(Z, path); }
130  }
131  }
132  }
133 
135 }
136 
137 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
138 
140 {
141  return "/brem_SB/br";
142 }
143 
144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
145 
146 void G4SeltzerBergerModel::ReadData(G4int Z, const char* path)
147 {
148  // G4cout << "ReadData Z= " << Z << G4endl;
149  // G4cout << "Status for Z= " << dataSB[Z] << G4endl;
150  //if(path) { G4cout << path << G4endl; }
151  if(dataSB[Z]) { return; }
152  const char* datadir = path;
153 
154  if(!datadir) {
155  datadir = getenv("G4LEDATA");
156  if(!datadir) {
157  G4Exception("G4SeltzerBergerModel::ReadData()","em0006",FatalException,
158  "Environment variable G4LEDATA not defined");
159  return;
160  }
161  }
162  std::ostringstream ost;
163  ost << datadir << DirectoryPath() << Z;
164  std::ifstream fin(ost.str().c_str());
165  if( !fin.is_open()) {
167  ed << "Bremsstrahlung data file <" << ost.str().c_str()
168  << "> is not opened!";
169  G4Exception("G4SeltzerBergerModel::ReadData()","em0003",FatalException,
170  ed,"G4LEDATA version should be G4EMLOW6.23 or later.");
171  return;
172  }
173  //G4cout << "G4SeltzerBergerModel read from <" << ost.str().c_str()
174  // << ">" << G4endl;
176  if(v->Retrieve(fin)) {
178  dataSB[Z] = v;
179  ylimit[Z] = v->Value(0.97, emaxlog, idx, idy);
180  } else {
182  ed << "Bremsstrahlung data file <" << ost.str().c_str()
183  << "> is not retrieved!";
184  G4Exception("G4SeltzerBergerModel::ReadData()","em0005",FatalException,
185  ed,"G4LEDATA version should be G4EMLOW6.23 or later.");
186  delete v;
187  }
188  // G4cout << dataSB[Z] << G4endl;
189 }
190 
191 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
192 
194 {
195 
196  if(gammaEnergy < 0.0 || kinEnergy <= 0.0) { return 0.0; }
197  G4double x = gammaEnergy/kinEnergy;
199  G4int Z = G4lrint(currentZ);
200 
201  //G4cout << "G4SeltzerBergerModel::ComputeDXSectionPerAtom Z= " << Z
202  // << " x= " << x << " y= " << y << " " << dataSB[Z] << G4endl;
203  if(nullptr == dataSB[Z]) { InitialiseForElement(0, Z); }
204  /*
205  G4ExceptionDescription ed;
206  ed << "Bremsstrahlung data for Z= " << Z
207  << " are not initialized!";
208  G4Exception("G4SeltzerBergerModel::ComputeDXSectionPerAtom()","em0005",
209  FatalException, ed,
210  "G4LEDATA version should be G4EMLOW6.23 or later.");
211  }
212  */
213  G4double invb2 =
215  G4double cross = dataSB[Z]->Value(x,y,idx,idy)*invb2*millibarn/bremFactor;
216 
217  if(!isElectron) {
218  G4double invbeta1 = sqrt(invb2);
219  G4double e2 = kinEnergy - gammaEnergy;
220  if(e2 > 0.0) {
221  G4double invbeta2 = (e2 + particleMass)/sqrt(e2*(e2 + 2*particleMass));
222  G4double xxx = alpha*currentZ*(invbeta1 - invbeta2);
223  if(xxx < expnumlim) { cross = 0.0; }
224  else { cross *= G4Exp(xxx); }
225  } else {
226  cross = 0.0;
227  }
228  }
229 
230  return cross;
231 }
232 
233 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
234 
235 void
236 G4SeltzerBergerModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
237  const G4MaterialCutsCouple* couple,
238  const G4DynamicParticle* dp,
239  G4double cutEnergy,
240  G4double maxEnergy)
241 {
242  G4double kineticEnergy = dp->GetKineticEnergy();
243  G4double cut = std::min(cutEnergy, kineticEnergy);
244  G4double emax = std::min(maxEnergy, kineticEnergy);
245  if(cut >= emax) { return; }
246 
247  SetupForMaterial(particle, couple->GetMaterial(), kineticEnergy);
248 
249  const G4Element* elm =
250  SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
251  SetCurrentElement(elm->GetZ());
252  G4int Z = G4int(currentZ);
253 
254  totalEnergy = kineticEnergy + particleMass;
256  G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
257  /*
258  G4cout << "G4SeltzerBergerModel::SampleSecondaries E(MeV)= "
259  << kineticEnergy/MeV
260  << " Z= " << Z << " cut(MeV)= " << cut/MeV
261  << " emax(MeV)= " << emax/MeV << " corr= " << densityCorr << G4endl;
262  */
263  G4double xmin = G4Log(cut*cut + densityCorr);
264  G4double xmax = G4Log(emax*emax + densityCorr);
265  G4double y = G4Log(kineticEnergy/MeV);
266 
267  G4double gammaEnergy, v;
268 
269  // majoranta
270  G4double x0 = cut/kineticEnergy;
271  G4double vmax = dataSB[Z]->Value(x0, y, idx, idy)*1.02;
272  // G4double invbeta1 = 0;
273 
274  // majoranta corrected for e-
275  if(isElectron && x0 < 0.97 &&
276  ((kineticEnergy > epeaklimit) || (kineticEnergy < elowlimit))) {
277  G4double ylim = std::min(ylimit[Z],1.1*dataSB[Z]->Value(0.97,y,idx,idy));
278  if(ylim > vmax) { vmax = ylim; }
279  }
280  if(x0 < 0.05) { vmax *= 1.2; }
281 
282  //G4cout<<"y= "<<y<<" xmin= "<<xmin<<" xmax= "<<xmax
283  //<<" vmax= "<<vmax<<G4endl;
284  static const G4int ncountmax = 100;
285  CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
286  G4double rndm[2];
287 
288  for(G4int nn=0; nn<ncountmax; ++nn) {
289  rndmEngine->flatArray(2, rndm);
290  G4double x = G4Exp(xmin + rndm[0]*(xmax - xmin)) - densityCorr;
291  if(x < 0.0) { x = 0.0; }
292  gammaEnergy = sqrt(x);
293  G4double x1 = gammaEnergy/kineticEnergy;
294  v = dataSB[Z]->Value(x1, y, idx, idy);
295 
296  // correction for positrons
297  if(!isElectron) {
298  G4double e1 = kineticEnergy - cut;
299  G4double invbeta1 = (e1 + particleMass)/sqrt(e1*(e1 + 2*particleMass));
300  G4double e2 = kineticEnergy - gammaEnergy;
301  G4double invbeta2 = (e2 + particleMass)/sqrt(e2*(e2 + 2*particleMass));
302  G4double xxx = twopi*fine_structure_const*currentZ*(invbeta1 - invbeta2);
303 
304  if(xxx < expnumlim) { v = 0.0; }
305  else { v *= G4Exp(xxx); }
306  }
307 
308  if (v > 1.05*vmax && nwarn < 5) {
309  ++nwarn;
311  ed << "### G4SeltzerBergerModel Warning: Majoranta exceeded! "
312  << v << " > " << vmax << " by " << v/vmax
313  << " Niter= " << nn
314  << " Egamma(MeV)= " << gammaEnergy
315  << " Ee(MeV)= " << kineticEnergy
316  << " Z= " << Z << " " << particle->GetParticleName();
317 
318  if ( 20 == nwarn ) {
319  ed << "\n ### G4SeltzerBergerModel Warnings stopped";
320  }
321  G4Exception("G4SeltzerBergerModel::SampleScattering","em0044",
322  JustWarning, ed,"");
323 
324  }
325  if(v >= vmax*rndm[1]) { break; }
326  }
327 
328  //
329  // angles of the emitted gamma. ( Z - axis along the parent particle)
330  // use general interface
331  //
332 
333  G4ThreeVector gammaDirection =
334  GetAngularDistribution()->SampleDirection(dp, totalEnergy-gammaEnergy,
335  Z, couple->GetMaterial());
336 
337  // create G4DynamicParticle object for the Gamma
338  G4DynamicParticle* gamma =
339  new G4DynamicParticle(theGamma,gammaDirection,gammaEnergy);
340  vdp->push_back(gamma);
341 
342  G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
343  - gammaEnergy*gammaDirection).unit();
344 
345  /*
346  G4cout << "### G4SBModel: v= "
347  << " Eg(MeV)= " << gammaEnergy
348  << " Ee(MeV)= " << kineticEnergy
349  << " DirE " << direction << " DirG " << gammaDirection
350  << G4endl;
351  */
352  // energy of primary
353  G4double finalE = kineticEnergy - gammaEnergy;
354 
355  // stop tracking and create new secondary instead of primary
356  if(gammaEnergy > SecondaryThreshold()) {
359  G4DynamicParticle* el =
360  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
361  direction, finalE);
362  vdp->push_back(el);
363 
364  // continue tracking
365  } else {
368  }
369 }
370 
371 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
372 
373 #include "G4AutoLock.hh"
374 namespace { G4Mutex SeltzerBergerModelMutex = G4MUTEX_INITIALIZER; }
376  G4int Z)
377 {
378  G4AutoLock l(&SeltzerBergerModelMutex);
379  // G4cout << "G4SeltzerBergerModel::InitialiseForElement Z= " << Z << G4endl;
380  if(nullptr == dataSB[Z]) { ReadData(Z); }
381 }
382 
383 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
384 
385 
static const double MeV
Definition: G4SIunits.hh:211
void SetBicubicInterpolation(G4bool)
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:669
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
void ReadData(G4int Z, const char *path=0)
G4SeltzerBergerModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="eBremSB")
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:617
static const G4double e2
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
int G4int
Definition: G4Types.hh:78
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:175
const G4String & GetParticleName() const
const G4ParticleDefinition * particle
G4double Value(G4double x, G4double y, size_t &lastidx, size_t &lastidy) const
virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy)
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
static const G4double elowlimit
const G4ThreeVector & GetMomentumDirection() const
static const double twopi
Definition: G4SIunits.hh:75
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
static const G4double e1
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const G4double emax
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static G4double ylimit[101]
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
G4int G4Mutex
Definition: G4Threading.hh:173
int G4lrint(double ad)
Definition: templates.hh:163
G4bool Retrieve(std::ifstream &fIn)
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:774
const G4double x[NPOINTSGL]
static G4Physics2DVector * dataSB[101]
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static const double millibarn
Definition: G4SIunits.hh:105
G4ParticleChangeForLoss * fParticleChange
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static const double keV
Definition: G4SIunits.hh:213
static const G4double epeaklimit
std::vector< G4Element * > G4ElementTable
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
virtual G4String DirectoryPath() const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy)
static const G4double alpha
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:544
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:369
const G4Material * GetMaterial() const
static const G4double emaxlog