Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4SeltzerBergerModel.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 // $Id: G4SeltzerBergerModel.cc 98737 2016-08-09 12:51:38Z 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 
79  const G4String& nam)
80  : G4eBremsstrahlungRelModel(p,nam),useBicubicInterpolation(false)
81 {
84  SetLPMFlag(false);
85  nwarn = 0;
86  idx = idy = 0;
87 }
88 
89 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
90 
92 {
93  if(IsMaster()) {
94  for(size_t i=0; i<101; ++i) {
95  if(dataSB[i]) {
96  delete dataSB[i];
97  dataSB[i] = nullptr;
98  }
99  }
100  }
101 }
102 
103 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
104 
106  const G4DataVector& cuts)
107 {
108  // Access to elements
109  if(IsMaster()) {
110 
111  // check environment variable
112  // Build the complete string identifying the file with the data set
113  char* path = getenv("G4LEDATA");
114 
115  const G4ElementTable* theElmTable = G4Element::GetElementTable();
116  size_t numOfElm = G4Element::GetNumberOfElements();
117  if(numOfElm > 0) {
118  for(size_t i=0; i<numOfElm; ++i) {
119  G4int Z = G4lrint(((*theElmTable)[i])->GetZ());
120  if(Z < 1) { Z = 1; }
121  else if(Z > 100) { Z = 100; }
122  //G4cout << "Z= " << Z << G4endl;
123  // Initialisation
124  if(nullptr == dataSB[Z]) { ReadData(Z, path); }
125  }
126  }
127  }
128 
130 }
131 
132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
133 
135 {
136  return "/brem_SB/br";
137 }
138 
139 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
140 
141 void G4SeltzerBergerModel::ReadData(G4int Z, const char* path)
142 {
143  // G4cout << "ReadData Z= " << Z << G4endl;
144  // G4cout << "Status for Z= " << dataSB[Z] << G4endl;
145  //if(path) { G4cout << path << G4endl; }
146  if(dataSB[Z]) { return; }
147  const char* datadir = path;
148 
149  if(!datadir) {
150  datadir = getenv("G4LEDATA");
151  if(!datadir) {
152  G4Exception("G4SeltzerBergerModel::ReadData()","em0006",FatalException,
153  "Environment variable G4LEDATA not defined");
154  return;
155  }
156  }
157  std::ostringstream ost;
158  ost << datadir << DirectoryPath() << Z;
159  std::ifstream fin(ost.str().c_str());
160  if( !fin.is_open()) {
162  ed << "Bremsstrahlung data file <" << ost.str().c_str()
163  << "> is not opened!";
164  G4Exception("G4SeltzerBergerModel::ReadData()","em0003",FatalException,
165  ed,"G4LEDATA version should be G4EMLOW6.23 or later.");
166  return;
167  }
168  //G4cout << "G4SeltzerBergerModel read from <" << ost.str().c_str()
169  // << ">" << G4endl;
171  if(v->Retrieve(fin)) {
173  dataSB[Z] = v;
174  static const G4double emaxlog = 4*G4Log(10.);
175  ylimit[Z] = v->Value(0.97, emaxlog, idx, idy);
176  } else {
178  ed << "Bremsstrahlung data file <" << ost.str().c_str()
179  << "> is not retrieved!";
180  G4Exception("G4SeltzerBergerModel::ReadData()","em0005",FatalException,
181  ed,"G4LEDATA version should be G4EMLOW6.23 or later.");
182  delete v;
183  }
184  // G4cout << dataSB[Z] << G4endl;
185 }
186 
187 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
188 
190 {
191 
192  if(gammaEnergy < 0.0 || kinEnergy <= 0.0) { return 0.0; }
193  G4double x = gammaEnergy/kinEnergy;
195  G4int Z = G4lrint(currentZ);
196 
197  //G4cout << "G4SeltzerBergerModel::ComputeDXSectionPerAtom Z= " << Z
198  // << " x= " << x << " y= " << y << " " << dataSB[Z] << G4endl;
199  if(nullptr == dataSB[Z]) { InitialiseForElement(0, Z); }
200  /*
201  G4ExceptionDescription ed;
202  ed << "Bremsstrahlung data for Z= " << Z
203  << " are not initialized!";
204  G4Exception("G4SeltzerBergerModel::ComputeDXSectionPerAtom()","em0005",
205  FatalException, ed,
206  "G4LEDATA version should be G4EMLOW6.23 or later.");
207  }
208  */
209  G4double invb2 =
211  G4double cross = dataSB[Z]->Value(x,y,idx,idy)*invb2*millibarn/bremFactor;
212 
213  if(!isElectron) {
214  G4double invbeta1 = sqrt(invb2);
215  G4double e2 = kinEnergy - gammaEnergy;
216  if(e2 > 0.0) {
217  G4double invbeta2 = (e2 + particleMass)/sqrt(e2*(e2 + 2*particleMass));
219  G4double xxx = alpha*currentZ*(invbeta1 - invbeta2);
220  if(xxx < expnumlim) { cross = 0.0; }
221  else { cross *= G4Exp(xxx); }
222  } else {
223  cross = 0.0;
224  }
225  }
226 
227  return cross;
228 }
229 
230 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
231 
232 void
233 G4SeltzerBergerModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
234  const G4MaterialCutsCouple* couple,
235  const G4DynamicParticle* dp,
236  G4double cutEnergy,
237  G4double maxEnergy)
238 {
240  G4double cut = std::min(cutEnergy, kineticEnergy);
241  G4double emax = std::min(maxEnergy, kineticEnergy);
242  if(cut >= emax) { return; }
243 
245 
246  const G4Element* elm =
247  SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
248  SetCurrentElement(elm->GetZasInt());
249 
250  totalEnergy = kineticEnergy + particleMass;
251  densityCorr = densityFactor*totalEnergy*totalEnergy;
252  G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
253  /*
254  G4cout << "G4SeltzerBergerModel::SampleSecondaries E(MeV)= "
255  << kineticEnergy/MeV
256  << " Z= " << Z << " cut(MeV)= " << cut/MeV
257  << " emax(MeV)= " << emax/MeV << " corr= " << densityCorr << G4endl;
258  */
259  G4double xmin = G4Log(cut*cut + densityCorr);
260  G4double xmax = G4Log(emax*emax + densityCorr);
261  G4double y = G4Log(kineticEnergy/MeV);
262 
263  G4double gammaEnergy, v;
264 
265  // majoranta
266  G4double x0 = cut/kineticEnergy;
267  G4double vmax;
268  if(currentZ <= 92) {
269  vmax = dataSB[currentZ]->Value(x0, y, idx, idy)*1.02;
270  } else {
271  idx = idy = 0;
272  vmax = dataSB[currentZ]->Value(x0, y, idx, idy)*1.2;
273  }
274 
275  static const G4double epeaklimit= 300*CLHEP::MeV;
276  static const G4double elowlimit = 20*CLHEP::keV;
277 
278  // majoranta corrected for e-
279  if(isElectron && x0 < 0.97 &&
280  ((kineticEnergy > epeaklimit) || (kineticEnergy < elowlimit))) {
281  G4double ylim = std::min(ylimit[currentZ],1.1*dataSB[currentZ]->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  static const G4int ncountmax = 100;
289  CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
290  G4double rndm[2];
291 
292  for(G4int nn=0; nn<ncountmax; ++nn) {
293  rndmEngine->flatArray(2, rndm);
294  G4double x = G4Exp(xmin + rndm[0]*(xmax - xmin)) - densityCorr;
295  if(x < 0.0) { x = 0.0; }
296  gammaEnergy = sqrt(x);
297  G4double x1 = gammaEnergy/kineticEnergy;
298  v = dataSB[currentZ]->Value(x1, y, idx, idy);
299 
300  // correction for positrons
301  if(!isElectron) {
302  G4double e1 = kineticEnergy - cut;
303  G4double invbeta1 = (e1 + particleMass)/sqrt(e1*(e1 + 2*particleMass));
304  G4double e2 = kineticEnergy - gammaEnergy;
305  G4double invbeta2 = (e2 + particleMass)/sqrt(e2*(e2 + 2*particleMass));
306  G4double xxx = twopi*fine_structure_const*currentZ*(invbeta1 - invbeta2);
307 
308  if(xxx < expnumlim) { v = 0.0; }
309  else { v *= G4Exp(xxx); }
310  }
311 
312  if (v > 1.05*vmax && nwarn < 5) {
313  ++nwarn;
315  ed << "### G4SeltzerBergerModel Warning: Majoranta exceeded! "
316  << v << " > " << vmax << " by " << v/vmax
317  << " Niter= " << nn
318  << " Egamma(MeV)= " << gammaEnergy
319  << " Ee(MeV)= " << kineticEnergy
320  << " Z= " << currentZ << " " << particle->GetParticleName();
321 
322  if ( 20 == nwarn ) {
323  ed << "\n ### G4SeltzerBergerModel Warnings stopped";
324  }
325  G4Exception("G4SeltzerBergerModel::SampleScattering","em0044",
326  JustWarning, ed,"");
327 
328  }
329  if(v >= vmax*rndm[1]) { break; }
330  }
331 
332  //
333  // angles of the emitted gamma. ( Z - axis along the parent particle)
334  // use general interface
335  //
336 
337  G4ThreeVector gammaDirection =
338  GetAngularDistribution()->SampleDirection(dp, totalEnergy-gammaEnergy,
339  currentZ, couple->GetMaterial());
340 
341  // create G4DynamicParticle object for the Gamma
342  G4DynamicParticle* gamma =
343  new G4DynamicParticle(theGamma,gammaDirection,gammaEnergy);
344  vdp->push_back(gamma);
345 
346  G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
347  - gammaEnergy*gammaDirection).unit();
348 
349  /*
350  G4cout << "### G4SBModel: v= "
351  << " Eg(MeV)= " << gammaEnergy
352  << " Ee(MeV)= " << kineticEnergy
353  << " DirE " << direction << " DirG " << gammaDirection
354  << G4endl;
355  */
356  // energy of primary
357  G4double finalE = kineticEnergy - gammaEnergy;
358 
359  // stop tracking and create new secondary instead of primary
360  if(gammaEnergy > SecondaryThreshold()) {
363  G4DynamicParticle* el =
364  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
365  direction, finalE);
366  vdp->push_back(el);
367 
368  // continue tracking
369  } else {
372  }
373 }
374 
375 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
376 
377 #include "G4AutoLock.hh"
378 namespace { G4Mutex SeltzerBergerModelMutex = G4MUTEX_INITIALIZER; }
380  G4int Z)
381 {
382  G4AutoLock l(&SeltzerBergerModelMutex);
383  // G4cout << "G4SeltzerBergerModel::InitialiseForElement Z= " << Z << G4endl;
384  if(nullptr == dataSB[Z]) { ReadData(Z); }
385 }
386 
387 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
388 
389 
Float_t x
Definition: compare.C:6
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
virtual G4String DirectoryPath() const
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:45
const G4ParticleDefinition * particle
static constexpr double MeV
Definition: G4SIunits.hh:214
static const G4double elowlimit
static G4Physics2DVector * dataSB[101]
G4SeltzerBergerModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="eBremSB")
static constexpr double keV
Definition: G4SIunits.hh:216
Float_t x1[n_points_granero]
Definition: compare.C:5
const G4ThreeVector & GetMomentumDirection() const
Float_t y
Definition: compare.C:6
const char * p
Definition: xmltok.h:285
const G4String & GetParticleName() const
virtual void flatArray(const int size, double *vect)=0
G4bool Retrieve(std::ifstream &fIn)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:552
static const G4double emax
void SetProposedMomentumDirection(const G4ThreeVector &dir)
static const G4double emaxlog
void ReadData(G4int Z, const char *path=0)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static const G4double epeaklimit
fin
Definition: AddMC.C:2
Float_t Z
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:756
double G4double
Definition: G4Types.hh:76
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:88
static constexpr double millibarn
Definition: G4SIunits.hh:106
static constexpr double electron_mass_c2
std::vector< G4Element * > G4ElementTable
static const G4double alpha
static constexpr double MeV
static constexpr double twopi
Definition: G4SIunits.hh:76
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z) override
G4double Value(G4double x, G4double y, size_t &lastidx, size_t &lastidy) const
static constexpr double twopi
Definition: SystemOfUnits.h:55
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
int G4lrint(double ad)
Definition: templates.hh:151
int G4int
Definition: G4Types.hh:78
static G4double ylimit[101]
virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy) override
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:714
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:644
G4double GetKineticEnergy() const
G4bool IsMaster() const
Definition: G4VEmModel.hh:700
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:575
const G4Material * GetMaterial() const
static size_t GetNumberOfElements()
Definition: G4Element.cc:405
G4ParticleChangeForLoss * fParticleChange
static constexpr double fine_structure_const
static constexpr double keV
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
void ProposeTrackStatus(G4TrackStatus status)
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:413
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void SetBicubicInterpolation(G4bool)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:398
std::mutex G4Mutex
Definition: G4Threading.hh:84