ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4LivermoreBremsstrahlungModel.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4LivermoreBremsstrahlungModel.cc
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 file
30 //
31 //
32 // File name: G4LivermoreBremsstrahlungModel
33 //
34 // Author: Vladimir Ivanchenko use inheritance from Andreas Schaelicke
35 // base class implementing ultra relativistic bremsstrahlung
36 // model
37 //
38 // Creation date: 04.10.2011
39 //
40 // Modifications:
41 //
42 // -------------------------------------------------------------------
43 //
44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46 
48 #include "G4PhysicalConstants.hh"
49 #include "G4SystemOfUnits.hh"
50 #include "G4Electron.hh"
51 #include "G4Positron.hh"
52 #include "G4Gamma.hh"
53 #include "Randomize.hh"
54 #include "G4Material.hh"
55 #include "G4Element.hh"
56 #include "G4ElementVector.hh"
57 #include "G4ProductionCutsTable.hh"
59 #include "G4Generator2BS.hh"
60 
61 #include "G4Physics2DVector.hh"
62 #include "G4Exp.hh"
63 #include "G4Log.hh"
64 
65 #include "G4ios.hh"
66 #include <fstream>
67 #include <iomanip>
68 
69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70 
71 using namespace std;
72 
76 
77 static const G4double emaxlog = 4*G4Log(10.);
79 static const G4double epeaklimit= 300*CLHEP::MeV;
80 static const G4double elowlimit = 20*CLHEP::keV;
81 
83  const G4ParticleDefinition* p, const G4String& nam)
84  : G4eBremsstrahlungRelModel(p,nam),useBicubicInterpolation(false)
85 {
86  SetLowEnergyLimit(10.0*eV);
87  SetLPMFlag(false);
88  nwarn = 0;
89  idx = idy = 0;
91 }
92 
93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
94 
96 {
97  if(IsMaster()) {
98  for(size_t i=0; i<101; ++i) {
99  if(dataSB[i]) {
100  delete dataSB[i];
101  dataSB[i] = 0;
102  }
103  }
104  }
105 }
106 
107 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
108 
110  const G4DataVector& cuts)
111 {
112  // Access to elements
113  if(IsMaster()) {
114 
115  // check environment variable
116  // Build the complete string identifying the file with the data set
117  char* path = std::getenv("G4LEDATA");
118 
119  const G4ElementTable* theElmTable = G4Element::GetElementTable();
120  size_t numOfElm = G4Element::GetNumberOfElements();
121  if(numOfElm > 0) {
122  for(size_t i=0; i<numOfElm; ++i) {
123  G4int Z = (*theElmTable)[i]->GetZasInt();
124  if(Z < 1) { Z = 1; }
125  else if(Z > 100) { Z = 100; }
126  //G4cout << "Z= " << Z << G4endl;
127  // Initialisation
128  if(!dataSB[Z]) { ReadData(Z, path); }
129  }
130  }
131  }
132 
134 }
135 
136 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
137 
139 {
140  return "/livermore/brem/br";
141 }
142 
143 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
144 
146 {
147  // G4cout << "ReadData Z= " << Z << G4endl;
148  // G4cout << "Status for Z= " << dataSB[Z] << G4endl;
149  //if(path) { G4cout << path << G4endl; }
150  if(dataSB[Z]) { return; }
151  const char* datadir = path;
152 
153  if(!datadir) {
154  datadir = std::getenv("G4LEDATA");
155  if(!datadir) {
156  G4Exception("G4LivermoreBremsstrahlungModel::ReadData()","em0006",
157  FatalException,"Environment variable G4LEDATA not defined");
158  return;
159  }
160  }
161  std::ostringstream ost;
162  ost << datadir << DirectoryPath() << Z;
163  std::ifstream fin(ost.str().c_str());
164  if( !fin.is_open()) {
166  ed << "Bremsstrahlung data file <" << ost.str().c_str()
167  << "> is not opened!";
168  G4Exception("G4LivermoreBremsstrahlungModel::ReadData()","em0003",
169  FatalException,ed,
170  "G4LEDATA version should be G4EMLOW6.23 or later.");
171  return;
172  }
173  //G4cout << "G4LivermoreBremsstrahlungModel 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("G4LivermoreBremsstrahlungModel::ReadData()","em0005",
185  FatalException,ed,
186  "G4LEDATA version should be G4EMLOW6.23 or later.");
187  delete v;
188  }
189  // G4cout << dataSB[Z] << G4endl;
190 }
191 
192 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
193 
194 G4double
196 {
197 
198  if(gammaEnergy < 0.0 || fPrimaryKinEnergy <= 0.0) { return 0.0; }
199  G4double x = gammaEnergy/fPrimaryKinEnergy;
201  G4int Z = fCurrentIZ;
202 
203  //G4cout << "G4LivermoreBremsstrahlungModel::ComputeDXSectionPerAtom Z= " << Z
204  // << " x= " << x << " y= " << y << " " << dataSB[Z] << G4endl;
205  if(!dataSB[Z]) { InitialiseForElement(0, Z); }
206  /*
207  G4ExceptionDescription ed;
208  ed << "Bremsstrahlung data for Z= " << Z
209  << " are not initialized!";
210  G4Exception("G4LivermoreBremsstrahlungModel::ComputeDXSectionPerAtom()",
211  "em0005", FatalException, ed,
212  "G4LEDATA version should be G4EMLOW6.23 or later.");
213  }
214  */
218 
219  if(!fIsElectron) {
220  G4double invbeta1 = sqrt(invb2);
221  G4double e2 = fPrimaryKinEnergy - gammaEnergy;
222  if(e2 > 0.0) {
223  G4double invbeta2 = (e2 + fPrimaryParticleMass)
224  /sqrt(e2*(e2 + 2.*fPrimaryParticleMass));
225  G4double xxx = alpha*fCurrentIZ*(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  // sets total energy, kinetic energy and density correction
251  SetupForMaterial(fPrimaryParticle, couple->GetMaterial(), kineticEnergy);
252 
253  const G4Element* elm =
254  SelectRandomAtom(couple,fPrimaryParticle,kineticEnergy,cut,emax);
255  fCurrentIZ = elm->GetZasInt();
256  G4int Z = fCurrentIZ;
257 
258  G4double totMomentum = sqrt(kineticEnergy*(fPrimaryTotalEnergy+electron_mass_c2));
259  /*
260  G4cout << "G4LivermoreBremsstrahlungModel::SampleSecondaries E(MeV)= "
261  << kineticEnergy/MeV
262  << " Z= " << Z << " cut(MeV)= " << cut/MeV
263  << " emax(MeV)= " << emax/MeV << " corr= " << fDensityCorr << G4endl;
264  */
265  G4double xmin = G4Log(cut*cut + fDensityCorr);
266  G4double xmax = G4Log(emax*emax + fDensityCorr);
267  G4double y = G4Log(kineticEnergy/MeV);
268 
269  G4double gammaEnergy, v;
270 
271  // majoranta
272  G4double x0 = cut/kineticEnergy;
273  G4double vmax = dataSB[Z]->Value(x0, y, idx, idy)*1.02;
274  // G4double invbeta1 = 0;
275 
276  // majoranta corrected for e-
277  if(fIsElectron && x0 < 0.97 &&
278  ((kineticEnergy > epeaklimit) || (kineticEnergy < elowlimit))) {
279  G4double ylim = std::min(ylimit[Z],1.1*dataSB[Z]->Value(0.97,y,idx,idy));
280  if(ylim > vmax) { vmax = ylim; }
281  }
282  if(x0 < 0.05) { vmax *= 1.2; }
283 
284  //G4cout<<"y= "<<y<<" xmin= "<<xmin<<" xmax= "<<xmax
285  //<<" vmax= "<<vmax<<G4endl;
286  // G4int ncount = 0;
287  do {
288  //++ncount;
289  G4double x = G4Exp(xmin + G4UniformRand()*(xmax - xmin)) - fDensityCorr;
290  if(x < 0.0) { x = 0.0; }
291  gammaEnergy = sqrt(x);
292  G4double x1 = gammaEnergy/kineticEnergy;
293  v = dataSB[Z]->Value(x1, y, idx, idy);
294 
295  // correction for positrons
296  if(!fIsElectron) {
297  G4double e1 = kineticEnergy - cut;
298  G4double invbeta1 = (e1 + fPrimaryParticleMass)
299  /sqrt(e1*(e1 + 2*fPrimaryParticleMass));
300  G4double e2 = kineticEnergy - gammaEnergy;
301  G4double invbeta2 = (e2 + fPrimaryParticleMass)
302  /sqrt(e2*(e2 + 2*fPrimaryParticleMass));
303  G4double xxx = twopi*fine_structure_const*fCurrentIZ*(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 << " " << fPrimaryParticle->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 =
334  Z, couple->GetMaterial());
335 
336  // create G4DynamicParticle object for the Gamma
337  G4DynamicParticle* gamma =
338  new G4DynamicParticle(fGammaParticle,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*>(fPrimaryParticle),
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