Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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$
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
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"
60#include "G4LossTableManager.hh"
61#include "G4ModifiedTsai.hh"
62
63#include "G4Physics2DVector.hh"
64
65#include "G4ios.hh"
66#include <fstream>
67#include <iomanip>
68
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
71
72using namespace std;
73
74G4Physics2DVector* G4SeltzerBergerModel::dataSB[101] = {0};
75G4double G4SeltzerBergerModel::ylimit[101] = {0.0};
76G4double G4SeltzerBergerModel::expnumlim = -12.;
77
79 const G4String& nam)
80 : G4eBremsstrahlungRelModel(p,nam),useBicubicInterpolation(false)
81{
83 SetLPMFlag(false);
84 nwarn = 0;
85}
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88
90{
91 for(size_t i=0; i<101; ++i) {
92 if(dataSB[i]) {
93 delete dataSB[i];
94 dataSB[i] = 0;
95 }
96 }
97}
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
100
102 const G4DataVector& cuts)
103{
104 // check environment variable
105 // Build the complete string identifying the file with the data set
106 char* path = getenv("G4LEDATA");
107
108 // Access to elements
109 const G4ElementTable* theElmTable = G4Element::GetElementTable();
110 size_t numOfElm = G4Element::GetNumberOfElements();
111 if(numOfElm > 0) {
112 for(size_t i=0; i<numOfElm; ++i) {
113 G4int Z = G4int(((*theElmTable)[i])->GetZ());
114 if(Z < 1) { Z = 1; }
115 else if(Z > 100) { Z = 100; }
116 //G4cout << "Z= " << Z << G4endl;
117 // Initialisation
118 if(!dataSB[Z]) { ReadData(Z, path); }
119 }
120 }
121
123}
124
125//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
126
127void G4SeltzerBergerModel::ReadData(size_t Z, const char* path)
128{
129 // G4cout << "ReadData Z= " << Z << G4endl;
130 // G4cout << "Status for Z= " << dataSB[Z] << G4endl;
131 //if(path) { G4cout << path << G4endl; }
132 if(dataSB[Z]) { return; }
133 const char* datadir = path;
134
135 if(!datadir) {
136 datadir = getenv("G4LEDATA");
137 if(!datadir) {
138 G4Exception("G4SeltzerBergerModel::ReadData()","em0006",FatalException,
139 "Environment variable G4LEDATA not defined");
140 return;
141 }
142 }
143 std::ostringstream ost;
144 ost << datadir << "/brem_SB/br" << Z;
145 std::ifstream fin(ost.str().c_str());
146 if( !fin.is_open()) {
148 ed << "Bremsstrahlung data file <" << ost.str().c_str()
149 << "> is not opened!" << G4endl;
150 G4Exception("G4SeltzerBergerModel::ReadData()","em0003",FatalException,
151 ed,"G4LEDATA version should be G4EMLOW6.23 or later.");
152 return;
153 }
154 //G4cout << "G4SeltzerBergerModel read from <" << ost.str().c_str()
155 // << ">" << G4endl;
157 const G4double emaxlog = 4*log(10.);
158 if(v->Retrieve(fin)) {
159 if(useBicubicInterpolation) { v->SetBicubicInterpolation(true); }
160 dataSB[Z] = v;
161 ylimit[Z] = v->Value(0.97, emaxlog);
162 } else {
164 ed << "Bremsstrahlung data file <" << ost.str().c_str()
165 << "> is not retrieved!" << G4endl;
166 G4Exception("G4SeltzerBergerModel::ReadData()","em0005",FatalException,
167 ed,"G4LEDATA version should be G4EMLOW6.23 or later.");
168 delete v;
169 }
170 // G4cout << dataSB[Z] << G4endl;
171}
172
173//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
174
176{
177
178 if(gammaEnergy < 0.0 || kinEnergy <= 0.0) { return 0.0; }
179 G4double x = gammaEnergy/kinEnergy;
180 G4double y = log(kinEnergy/MeV);
181 G4int Z = G4int(currentZ);
182 //G4cout << "G4SeltzerBergerModel::ComputeDXSectionPerAtom Z= " << Z
183 // << " x= " << x << " y= " << y << " " << dataSB[Z] << G4endl;
184 if(!dataSB[Z]) { ReadData(Z); }
186 G4double cross = dataSB[Z]->Value(x,y)*invb2*millibarn/bremFactor;
187
188 if(!isElectron) {
189 G4double invbeta1 = sqrt(invb2);
190 G4double e2 = kinEnergy - gammaEnergy;
191 if(e2 > 0.0) {
192 G4double invbeta2 = (e2 + particleMass)/sqrt(e2*(e2 + 2*particleMass));
193 G4double xxx = twopi*fine_structure_const*currentZ*(invbeta1 - invbeta2);
194 if(xxx < expnumlim) { cross = 0.0; }
195 else { cross *= exp(xxx); }
196 } else {
197 cross = 0.0;
198 }
199 }
200
201 return cross;
202}
203
204//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
205
206void
207G4SeltzerBergerModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
208 const G4MaterialCutsCouple* couple,
209 const G4DynamicParticle* dp,
210 G4double cutEnergy,
211 G4double maxEnergy)
212{
213 G4double kineticEnergy = dp->GetKineticEnergy();
214 G4double cut = std::min(cutEnergy, kineticEnergy);
215 G4double emax = std::min(maxEnergy, kineticEnergy);
216 if(cut >= emax) { return; }
217
218 SetupForMaterial(particle, couple->GetMaterial(), kineticEnergy);
219
220 const G4Element* elm =
221 SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
222 SetCurrentElement(elm->GetZ());
223 G4int Z = G4int(currentZ);
224
225 totalEnergy = kineticEnergy + particleMass;
227 G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
228 /*
229 G4cout << "G4SeltzerBergerModel::SampleSecondaries E(MeV)= "
230 << kineticEnergy/MeV
231 << " Z= " << Z << " cut(MeV)= " << cut/MeV
232 << " emax(MeV)= " << emax/MeV << " corr= " << densityCorr << G4endl;
233 */
234 G4double xmin = log(cut*cut + densityCorr);
235 G4double xmax = log(emax*emax + densityCorr);
236 G4double y = log(kineticEnergy/MeV);
237
238 G4double gammaEnergy, v;
239
240 // majoranta
241 G4double x0 = cut/kineticEnergy;
242 G4double vmax = dataSB[Z]->Value(x0, y)*1.02;
243 // G4double invbeta1 = 0;
244
245 const G4double epeaklimit= 300*MeV;
246 const G4double elowlimit = 10*keV;
247
248 // majoranta corrected for e-
249 if(isElectron && x0 < 0.97 &&
250 ((kineticEnergy > epeaklimit) || (kineticEnergy < elowlimit))) {
251 G4double ylim = std::min(ylimit[Z],1.1*dataSB[Z]->Value(0.97, y));
252 if(ylim > vmax) { vmax = ylim; }
253 }
254 if(x0 < 0.05) { vmax *= 1.2; }
255
256 //G4cout<<"y= "<<y<<" xmin= "<<xmin<<" xmax= "<<xmax<<" vmax= "<<vmax<<G4endl;
257 // G4int ncount = 0;
258 do {
259 //++ncount;
260 G4double x = exp(xmin + G4UniformRand()*(xmax - xmin)) - densityCorr;
261 if(x < 0.0) { x = 0.0; }
262 gammaEnergy = sqrt(x);
263 G4double x1 = gammaEnergy/kineticEnergy;
264 v = dataSB[Z]->Value(x1, y);
265
266 // correction for positrons
267 if(!isElectron) {
268 G4double e1 = kineticEnergy - cut;
269 G4double invbeta1 = (e1 + particleMass)/sqrt(e1*(e1 + 2*particleMass));
270 G4double e2 = kineticEnergy - gammaEnergy;
271 G4double invbeta2 = (e2 + particleMass)/sqrt(e2*(e2 + 2*particleMass));
272 G4double xxx = twopi*fine_structure_const*currentZ*(invbeta1 - invbeta2);
273
274 if(xxx < expnumlim) { v = 0.0; }
275 else { v *= exp(xxx); }
276 }
277
278 if (v > 1.05*vmax && nwarn < 20) {
279 ++nwarn;
280 G4cout << "### G4SeltzerBergerModel Warning: Majoranta exceeded! "
281 << v << " > " << vmax << " by " << v/vmax
282 << " Egamma(MeV)= " << gammaEnergy
283 << " Ee(MeV)= " << kineticEnergy
284 << " Z= " << Z << " " << particle->GetParticleName()
285 //<< " ncount= " << ncount
286 << G4endl;
287
288 if ( 20 == nwarn ) {
289 G4cout << "### G4SeltzerBergerModel Warnings will not be printed anymore"
290 << G4endl;
291 }
292 }
293 } while (v < vmax*G4UniformRand());
294
295 //
296 // angles of the emitted gamma. ( Z - axis along the parent particle)
297 // use general interface
298 //
299
300 G4ThreeVector gammaDirection =
302 Z, couple->GetMaterial());
303
304 // create G4DynamicParticle object for the Gamma
305 G4DynamicParticle* gamma =
306 new G4DynamicParticle(theGamma,gammaDirection,gammaEnergy);
307 vdp->push_back(gamma);
308
309 G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
310 - gammaEnergy*gammaDirection).unit();
311
312 /*
313 G4cout << "### G4SBModel: v= "
314 << " Eg(MeV)= " << gammaEnergy
315 << " Ee(MeV)= " << kineticEnergy
316 << " DirE " << direction << " DirG " << gammaDirection
317 << G4endl;
318 */
319 // energy of primary
320 G4double finalE = kineticEnergy - gammaEnergy;
321
322 // stop tracking and create new secondary instead of primary
323 if(gammaEnergy > SecondaryThreshold()) {
326 G4DynamicParticle* el =
328 direction, finalE);
329 vdp->push_back(el);
330
331 // continue tracking
332 } else {
335 }
336}
337
338//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
339
340
std::vector< G4Element * > G4ElementTable
@ FatalException
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
static size_t GetNumberOfElements()
Definition: G4Element.cc:406
static const G4ElementTable * GetElementTable()
Definition: G4Element.cc:399
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
const G4String & GetParticleName() const
G4bool Retrieve(std::ifstream &fIn)
G4double Value(G4double x, G4double y)
void SetBicubicInterpolation(G4bool)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy)
G4SeltzerBergerModel(const G4ParticleDefinition *p=0, const G4String &nam="eBremSB")
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:508
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:634
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:286
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:557
void ProposeTrackStatus(G4TrackStatus status)
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
const G4ParticleDefinition * particle
G4ParticleChangeForLoss * fParticleChange
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76