Geant4  10.02.p01
G4LivermoreBremsstrahlungModel.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: G4LivermoreBremsstrahlungModel.cc 76220 2013-11-08 10:15:00Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4LivermoreBremsstrahlungModel
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 
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 "G4Generator2BS.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 G4ParticleDefinition* p, const G4String& nam)
85  : G4eBremsstrahlungRelModel(p,nam),useBicubicInterpolation(false)
86 {
87  SetLowEnergyLimit(10.0*eV);
88  SetLPMFlag(false);
89  nwarn = 0;
90  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] = 0;
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 = G4int(((*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(!dataSB[Z]) { ReadData(Z, path); }
130  }
131  }
132  }
133 
135 }
136 
137 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
138 
140 {
141  return "/livermore/brem/br";
142 }
143 
144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
145 
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("G4LivermoreBremsstrahlungModel::ReadData()","em0006",
158  FatalException,"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("G4LivermoreBremsstrahlungModel::ReadData()","em0003",
170  FatalException,ed,
171  "G4LEDATA version should be G4EMLOW6.23 or later.");
172  return;
173  }
174  //G4cout << "G4LivermoreBremsstrahlungModel read from <" << ost.str().c_str()
175  // << ">" << G4endl;
177  if(v->Retrieve(fin)) {
179  dataSB[Z] = v;
180  ylimit[Z] = v->Value(0.97, emaxlog, idx, idy);
181  } else {
183  ed << "Bremsstrahlung data file <" << ost.str().c_str()
184  << "> is not retrieved!";
185  G4Exception("G4LivermoreBremsstrahlungModel::ReadData()","em0005",
186  FatalException,ed,
187  "G4LEDATA version should be G4EMLOW6.23 or later.");
188  delete v;
189  }
190  // G4cout << dataSB[Z] << G4endl;
191 }
192 
193 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
194 
195 G4double
197 {
198 
199  if(gammaEnergy < 0.0 || kinEnergy <= 0.0) { return 0.0; }
200  G4double x = gammaEnergy/kinEnergy;
202  G4int Z = G4lrint(currentZ);
203 
204  //G4cout << "G4LivermoreBremsstrahlungModel::ComputeDXSectionPerAtom Z= " << Z
205  // << " x= " << x << " y= " << y << " " << dataSB[Z] << G4endl;
206  if(!dataSB[Z]) { InitialiseForElement(0, Z); }
207  /*
208  G4ExceptionDescription ed;
209  ed << "Bremsstrahlung data for Z= " << Z
210  << " are not initialized!";
211  G4Exception("G4LivermoreBremsstrahlungModel::ComputeDXSectionPerAtom()",
212  "em0005", FatalException, ed,
213  "G4LEDATA version should be G4EMLOW6.23 or later.");
214  }
215  */
216  G4double invb2 =
218  G4double cross = dataSB[Z]->Value(x,y,idx,idy)*invb2*millibarn/bremFactor;
219 
220  if(!isElectron) {
221  G4double invbeta1 = sqrt(invb2);
222  G4double e2 = kinEnergy - gammaEnergy;
223  if(e2 > 0.0) {
224  G4double invbeta2 = (e2 + particleMass)/sqrt(e2*(e2 + 2*particleMass));
225  G4double xxx = alpha*currentZ*(invbeta1 - invbeta2);
226  if(xxx < expnumlim) { cross = 0.0; }
227  else { cross *= G4Exp(xxx); }
228  } else {
229  cross = 0.0;
230  }
231  }
232 
233  return cross;
234 }
235 
236 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
237 
238 void
240  std::vector<G4DynamicParticle*>* vdp,
241  const G4MaterialCutsCouple* couple,
242  const G4DynamicParticle* dp,
243  G4double cutEnergy,
244  G4double maxEnergy)
245 {
246  G4double kineticEnergy = dp->GetKineticEnergy();
247  G4double cut = std::min(cutEnergy, kineticEnergy);
248  G4double emax = std::min(maxEnergy, kineticEnergy);
249  if(cut >= emax) { return; }
250 
251  SetupForMaterial(particle, couple->GetMaterial(), kineticEnergy);
252 
253  const G4Element* elm =
254  SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
255  SetCurrentElement(elm->GetZ());
256  G4int Z = G4int(currentZ);
257 
258  totalEnergy = kineticEnergy + particleMass;
260  G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
261  /*
262  G4cout << "G4LivermoreBremsstrahlungModel::SampleSecondaries E(MeV)= "
263  << kineticEnergy/MeV
264  << " Z= " << Z << " cut(MeV)= " << cut/MeV
265  << " emax(MeV)= " << emax/MeV << " corr= " << densityCorr << G4endl;
266  */
267  G4double xmin = G4Log(cut*cut + densityCorr);
268  G4double xmax = G4Log(emax*emax + densityCorr);
269  G4double y = G4Log(kineticEnergy/MeV);
270 
271  G4double gammaEnergy, v;
272 
273  // majoranta
274  G4double x0 = cut/kineticEnergy;
275  G4double vmax = dataSB[Z]->Value(x0, y, idx, idy)*1.02;
276  // G4double invbeta1 = 0;
277 
278  // majoranta corrected for e-
279  if(isElectron && x0 < 0.97 &&
280  ((kineticEnergy > epeaklimit) || (kineticEnergy < elowlimit))) {
281  G4double ylim = std::min(ylimit[Z],1.1*dataSB[Z]->Value(0.97,y,idx,idy));
282  if(ylim > vmax) { vmax = ylim; }
283  }
284  if(x0 < 0.05) { vmax *= 1.2; }
285 
286  //G4cout<<"y= "<<y<<" xmin= "<<xmin<<" xmax= "<<xmax
287  //<<" vmax= "<<vmax<<G4endl;
288  // G4int ncount = 0;
289  do {
290  //++ncount;
291  G4double x = G4Exp(xmin + G4UniformRand()*(xmax - xmin)) - densityCorr;
292  if(x < 0.0) { x = 0.0; }
293  gammaEnergy = sqrt(x);
294  G4double x1 = gammaEnergy/kineticEnergy;
295  v = dataSB[Z]->Value(x1, y, idx, idy);
296 
297  // correction for positrons
298  if(!isElectron) {
299  G4double e1 = kineticEnergy - cut;
300  G4double invbeta1 = (e1 + particleMass)/sqrt(e1*(e1 + 2*particleMass));
301  G4double e2 = kineticEnergy - gammaEnergy;
302  G4double invbeta2 = (e2 + particleMass)/sqrt(e2*(e2 + 2*particleMass));
303  G4double xxx = twopi*fine_structure_const*currentZ*(invbeta1 - invbeta2);
304 
305  if(xxx < expnumlim) { v = 0.0; }
306  else { v *= G4Exp(xxx); }
307  }
308 
309  if (v > 1.05*vmax && nwarn < 5) {
310  ++nwarn;
312  ed << "### G4LivermoreBremsstrahlungModel Warning: Majoranta exceeded! "
313  << v << " > " << vmax << " by " << v/vmax
314  << " Egamma(MeV)= " << gammaEnergy
315  << " Ee(MeV)= " << kineticEnergy
316  << " Z= " << Z << " " << particle->GetParticleName();
317 
318  if ( 20 == nwarn ) {
319  ed << "\n ### G4LivermoreBremsstrahlungModel Warnings stopped";
320  }
321  G4Exception("G4LivermoreBremsstrahlungModel::SampleScattering","em0044",
322  JustWarning, ed,"");
323 
324  }
325  } while (v < vmax*G4UniformRand());
326 
327  //
328  // angles of the emitted gamma. ( Z - axis along the parent particle)
329  // use general interface
330  //
331 
332  G4ThreeVector gammaDirection =
333  GetAngularDistribution()->SampleDirection(dp, totalEnergy-gammaEnergy,
334  Z, couple->GetMaterial());
335 
336  // create G4DynamicParticle object for the Gamma
337  G4DynamicParticle* gamma =
338  new G4DynamicParticle(theGamma,gammaDirection,gammaEnergy);
339  vdp->push_back(gamma);
340 
341  G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
342  - gammaEnergy*gammaDirection).unit();
343 
344  /*
345  G4cout << "### G4SBModel: v= "
346  << " Eg(MeV)= " << gammaEnergy
347  << " Ee(MeV)= " << kineticEnergy
348  << " DirE " << direction << " DirG " << gammaDirection
349  << G4endl;
350  */
351  // energy of primary
352  G4double finalE = kineticEnergy - gammaEnergy;
353 
354  // stop tracking and create new secondary instead of primary
355  if(gammaEnergy > SecondaryThreshold()) {
358  G4DynamicParticle* el =
359  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
360  direction, finalE);
361  vdp->push_back(el);
362 
363  // continue tracking
364  } else {
367  }
368 }
369 
370 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
371 
372 #include "G4AutoLock.hh"
373 namespace { G4Mutex LivermoreBremsstrahlungModelMutex = G4MUTEX_INITIALIZER; }
375  const G4ParticleDefinition*,
376  G4int Z)
377 {
378  G4AutoLock l(&LivermoreBremsstrahlungModelMutex);
379  //G4cout << "G4LivermoreBremsstrahlungModel::InitialiseForElement Z= "
380  //<< Z << G4endl;
381  if(!dataSB[Z]) { ReadData(Z); }
382  l.unlock();
383 }
384 
385 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
386 
387 
static const double MeV
Definition: G4SIunits.hh:211
void SetBicubicInterpolation(G4bool)
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:669
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
static G4Physics2DVector * dataSB[101]
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:617
static const G4double e2
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
int G4int
Definition: G4Types.hh:78
static const G4double emaxlog
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:175
const G4String & GetParticleName() const
const G4ParticleDefinition * particle
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
#define G4UniformRand()
Definition: Randomize.hh:97
G4double Value(G4double x, G4double y, size_t &lastidx, size_t &lastidy) const
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
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 const double eV
Definition: G4SIunits.hh:212
G4int G4Mutex
Definition: G4Threading.hh:173
void ReadData(G4int Z, const char *path=0)
int G4lrint(double ad)
Definition: templates.hh:163
G4bool Retrieve(std::ifstream &fIn)
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:774
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:624
const G4double x[NPOINTSGL]
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static const double millibarn
Definition: G4SIunits.hh:105
static const G4double epeaklimit
G4ParticleChangeForLoss * fParticleChange
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static const double keV
Definition: G4SIunits.hh:213
virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy)
std::vector< G4Element * > G4ElementTable
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
static const G4double elowlimit
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
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
G4LivermoreBremsstrahlungModel(const G4ParticleDefinition *p=0, const G4String &nam="LowEnBrem")