Geant4 10.7.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// -------------------------------------------------------------------
27//
28// GEANT4 Class file
29//
30//
31// File name: G4SeltzerBergerModel
32//
33// Author: Vladimir Ivanchenko use inheritance from Andreas Schaelicke
34// base class implementing ultra relativistic bremsstrahlung
35// model
36//
37// Creation date: 04.10.2011
38//
39// Modifications:
40//
41// 24.07.2018 Introduced possibility to use sampling tables to sample the
42// emitted photon energy (instead of using rejectio) from the
43// Seltzer-Berger scalled DCS for bremsstrahlung photon emission.
44// Using these sampling tables option gives faster(30-70%) final
45// state generation than the original rejection but takes some
46// extra memory (+ ~6MB in the case of the full CMS detector).
47// (M Novak)
48//
49// -------------------------------------------------------------------
50//
51
54#include "G4SystemOfUnits.hh"
55#include "Randomize.hh"
56#include "G4Material.hh"
57#include "G4Element.hh"
58#include "G4ElementVector.hh"
60#include "G4SBBremTable.hh"
61#include "G4ModifiedTsai.hh"
62
63#include "G4EmParameters.hh"
65
66#include "G4Physics2DVector.hh"
67#include "G4Exp.hh"
68#include "G4Log.hh"
69
70#include "G4ios.hh"
71
72#include <fstream>
73#include <iomanip>
74#include <sstream>
75
76G4Physics2DVector* G4SeltzerBergerModel::gSBDCSData[] = { nullptr };
77G4SBBremTable* G4SeltzerBergerModel::gSBSamplingTable = nullptr;
78G4double G4SeltzerBergerModel::gYLimitData[] = { 0.0 };
79G4String G4SeltzerBergerModel::gDataDirectory = "";
80
81#ifdef G4MULTITHREADED
82 G4Mutex G4SeltzerBergerModel::theSBMutex = G4MUTEX_INITIALIZER;
83#endif
84
85static const G4double kMC2 = CLHEP::electron_mass_c2;
86static const G4double kAlpha = CLHEP::twopi*CLHEP::fine_structure_const;
87
89 const G4String& nam)
90: G4eBremsstrahlungRelModel(p,nam), fIsUseBicubicInterpolation(false),
91 fIsUseSamplingTables(true), fNumWarnings(0), fIndx(0), fIndy(0)
92{
93 fLowestKinEnergy = 1.0*keV;
95 SetLPMFlag(false);
97}
98
100{
101 // delete SB-DCS data per Z
102 if (IsMaster()) {
103 for (size_t iz = 0; iz < gMaxZet; ++iz) {
104 if (gSBDCSData[iz]) {
105 delete gSBDCSData[iz];
106 gSBDCSData[iz] = nullptr;
107 }
108 }
109 if (gSBSamplingTable) {
110 delete gSBSamplingTable;
111 gSBSamplingTable = nullptr;
112 }
113 }
114}
115
117 const G4DataVector& cuts)
118{
119 if (p) {
120 SetParticle(p);
121 }
122 fIsUseSamplingTables = G4EmParameters::Instance()->EnableSamplingTable();
123 // Access to elements
124 if (IsMaster()) {
125
126 auto theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
127 size_t numOfCouples = theCoupleTable->GetTableSize();
128 for(size_t j=0; j<numOfCouples; ++j) {
129 auto mat = theCoupleTable->GetMaterialCutsCouple(j)->GetMaterial();
130 auto elmVec = mat->GetElementVector();
131 for (auto & elm : *elmVec) {
132 G4int Z = std::max(1,std::min(elm->GetZasInt(), gMaxZet-1));
133 // load SB-DCS data for this atomic number if it has not been loaded yet
134 if (gSBDCSData[Z] == nullptr) ReadData(Z);
135 }
136 }
137 // elem.selectr. only for master: base class init-local will set for workers
140 }
141 // init sampling tables if it was requested
142 if (fIsUseSamplingTables) {
143 if (!gSBSamplingTable) {
144 gSBSamplingTable = new G4SBBremTable();
145 }
146 gSBSamplingTable->Initialize(std::max(fLowestKinEnergy,LowEnergyLimit()),
148 }
149 }
150 //
152 if (GetTripletModel()) {
153 GetTripletModel()->Initialise(p, cuts);
154 fIsScatOffElectron = true;
155 }
156}
157
158const G4String& G4SeltzerBergerModel::FindDirectoryPath()
159{
160 // check environment variable
161 // build the complete string identifying the file with the data set
162 if(gDataDirectory.empty()) {
163 const char* path = std::getenv("G4LEDATA");
164 if (path) {
165 std::ostringstream ost;
166 ost << path << "/brem_SB/br";
167 gDataDirectory = ost.str();
168 } else {
169 G4Exception("G4SeltzerBergerModel::FindDirectoryPath()","em0006",
171 "Environment variable G4LEDATA not defined");
172 }
173 }
174 return gDataDirectory;
175}
176
177void G4SeltzerBergerModel::ReadData(G4int Z) {
178 // return if it has been already loaded
179 if (gSBDCSData[Z] != nullptr) return;
180
181#ifdef G4MULTITHREADED
182 G4MUTEXLOCK(&theSBMutex);
183 if (gSBDCSData[Z] != nullptr) return;
184#endif
185
186 std::ostringstream ost;
187 ost << FindDirectoryPath() << Z;
188 std::ifstream fin(ost.str().c_str());
189 if (!fin.is_open()) {
191 ed << "Bremsstrahlung data file <" << ost.str().c_str()
192 << "> is not opened!";
193 G4Exception("G4SeltzerBergerModel::ReadData()","em0003",FatalException,
194 ed,"G4LEDATA version should be G4EMLOW6.23 or later.");
195 return;
196 }
197 //G4cout << "G4SeltzerBergerModel read from <" << ost.str().c_str()
198 // << ">" << G4endl;
200 if (v->Retrieve(fin)) {
201 v->SetBicubicInterpolation(fIsUseBicubicInterpolation);
202 static const G4double emaxlog = 4*G4Log(10.);
203 gYLimitData[Z] = v->Value(0.97, emaxlog, fIndx, fIndy);
204 gSBDCSData[Z] = v;
205 } else {
207 ed << "Bremsstrahlung data file <" << ost.str().c_str()
208 << "> is not retrieved!";
209 G4Exception("G4SeltzerBergerModel::ReadData()","em0005",FatalException,
210 ed,"G4LEDATA version should be G4EMLOW6.23 or later.");
211 delete v;
212 }
213#ifdef G4MULTITHREADED
214 G4MUTEXUNLOCK(&theSBMutex);
215#endif
216}
217
219{
220 G4double dxsec = 0.0;
221 if (gammaEnergy < 0.0 || fPrimaryKinEnergy <= 0.0) {
222 return dxsec;
223 }
224 // reduced photon energy
225 const G4double x = gammaEnergy/fPrimaryKinEnergy;
226 // l-kinetic energy of the e-/e+
227 const G4double y = G4Log(fPrimaryKinEnergy/CLHEP::MeV);
228 // make sure that the Z-related SB-DCS are loaded
229 // NOTE: fCurrentIZ should have been set before.
230 fCurrentIZ = std::max(std::min(fCurrentIZ, gMaxZet-1), 1);
231 if (nullptr == gSBDCSData[fCurrentIZ]) {
232 ReadData(fCurrentIZ);
233 }
234 // NOTE: SetupForMaterial should have been called before!
235 const G4double pt2 = fPrimaryKinEnergy*(fPrimaryKinEnergy+2.*kMC2);
237 G4double val = gSBDCSData[fCurrentIZ]->Value(x,y,fIndx,fIndy);
238 dxsec = val*invb2*CLHEP::millibarn/gBremFactor;
239 // e+ correction
240 if (!fIsElectron) {
241 const G4double invbeta1 = std::sqrt(invb2);
242 const G4double e2 = fPrimaryKinEnergy-gammaEnergy;
243 if (e2 > 0.0) {
244 const G4double invbeta2 = (e2+kMC2)/std::sqrt(e2*(e2+2.0*kMC2));
245 const G4double dum0 = kAlpha*fCurrentIZ*(invbeta1-invbeta2);
246 if (dum0 < gExpNumLimit) {
247 dxsec = 0.0;
248 } else {
249 dxsec *= G4Exp(dum0);
250 }
251 } else {
252 dxsec = 0.0;
253 }
254 }
255 return dxsec;
256}
257
258void
259G4SeltzerBergerModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
260 const G4MaterialCutsCouple* couple,
261 const G4DynamicParticle* dp,
262 G4double cutEnergy,
263 G4double maxEnergy)
264{
265 const G4double kinEnergy = dp->GetKineticEnergy();
266 const G4double logKinEnergy = dp->GetLogKineticEnergy();
267 const G4double tmin = std::min(cutEnergy, kinEnergy);
268 const G4double tmax = std::min(maxEnergy, kinEnergy);
269 if (tmin >= tmax) {
270 return;
271 }
272 // set local variables and select target element
273 SetupForMaterial(fPrimaryParticle, couple->GetMaterial(), kinEnergy);
274 const G4Element* elm = SelectTargetAtom(couple, fPrimaryParticle, kinEnergy,
275 logKinEnergy, tmin, tmax);
276 fCurrentIZ = std::max(std::min(elm->GetZasInt(),gMaxZet-1), 1);
277 //
278 const G4double totMomentum = std::sqrt(kinEnergy*(fPrimaryTotalEnergy+kMC2));
279 /*
280 G4cout << "G4SeltzerBergerModel::SampleSecondaries E(MeV)= "
281 << kinEnergy/MeV
282 << " Z= " << fCurrentIZ << " cut(MeV)= " << tmin/MeV
283 << " emax(MeV)= " << tmax/MeV << " corr= " << fDensityCorr << G4endl;
284 */
285 // sample emitted photon energy either by rejection or from samplign tables
286 const G4double gammaEnergy = !fIsUseSamplingTables
287 ? SampleEnergyTransfer(kinEnergy, logKinEnergy, tmin, tmax)
288 : gSBSamplingTable->SampleEnergyTransfer(kinEnergy, logKinEnergy, tmin,
290 // should never happen under normal conditions but protect it
291 if (gammaEnergy <= 0.) {
292 return;
293 }
294 //
295 // angles of the emitted gamma. ( Z - axis along the parent particle) use
296 // general interface
298 fPrimaryTotalEnergy-gammaEnergy, fCurrentIZ, couple->GetMaterial());
299 // create G4DynamicParticle object for the emitted Gamma
301 gammaEnergy);
302 vdp->push_back(gamma);
303 //
304 // compute post-interaction kinematics of the primary e-/e+
305 G4ThreeVector dir =
306 (totMomentum*dp->GetMomentumDirection()-gammaEnergy*gamDir).unit();
307 const G4double finalE = kinEnergy - gammaEnergy;
308 /*
309 G4cout << "### G4SBModel: v= "
310 << " Eg(MeV)= " << gammaEnergy
311 << " Ee(MeV)= " << kineticEnergy
312 << " DirE " << direction << " DirG " << gammaDirection
313 << G4endl;
314 */
315 // if secondary gamma energy is higher than threshold(very high by default)
316 // then stop tracking the primary particle and create new secondary e-/e+
317 // instead of the primary
318 if (gammaEnergy > SecondaryThreshold()) {
322 const_cast<G4ParticleDefinition*>(fPrimaryParticle), dir, finalE);
323 vdp->push_back(el);
324 } else { // continue tracking the primary e-/e+ otherwise
327 }
328}
329
330// sample emitted photon energy by usign rejection
332G4SeltzerBergerModel::SampleEnergyTransfer(const G4double kinEnergy,
333 const G4double logKinEnergy,
334 const G4double tmin,
335 const G4double tmax)
336{
337 // min max of the transformed variable: x(k) = ln(k^2+k_p^2) that is in
338 // [ln(k_c^2+k_p^2), ln(E_k^2+k_p^2)]
339 const G4double xmin = G4Log(tmin*tmin+fDensityCorr);
340 const G4double xrange = G4Log(tmax*tmax+fDensityCorr)-xmin;
341 const G4double y = logKinEnergy;
342 // majoranta
343 const G4double x0 = tmin/kinEnergy;
344 G4double vmax;
345 if (nullptr == gSBDCSData[fCurrentIZ]) {
346 ReadData(fCurrentIZ);
347 }
348 vmax = gSBDCSData[fCurrentIZ]->Value(x0, y, fIndx, fIndy)*1.02;
349 //
350 static const G4double kEPeakLim = 300.*CLHEP::MeV;
351 static const G4double kELowLim = 20.*CLHEP::keV;
352 // majoranta corrected for e-
353 if (fIsElectron && x0 < 0.97 &&
354 ((kinEnergy>kEPeakLim) || (kinEnergy<kELowLim))) {
355 G4double ylim = std::min(gYLimitData[fCurrentIZ],
356 1.1*gSBDCSData[fCurrentIZ]->Value(0.97,y,fIndx,fIndy));
357 vmax = std::max(vmax, ylim);
358 }
359 if (x0 < 0.05) {
360 vmax *= 1.2;
361 }
362 //G4cout<<"y= "<<y<<" xmin= "<<xmin<<" xmax= "<<xmax
363 //<<" vmax= "<<vmax<<G4endl;
364 static const G4int kNCountMax = 100;
365 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
366 G4double rndm[2];
367 G4double gammaEnergy, v;
368 for (G4int nn = 0; nn < kNCountMax; ++nn) {
369 rndmEngine->flatArray(2, rndm);
370 gammaEnergy =
371 std::sqrt(std::max(G4Exp(xmin + rndm[0]*xrange)-fDensityCorr,0.));
372 v = gSBDCSData[fCurrentIZ]->Value(gammaEnergy/kinEnergy, y, fIndx, fIndy);
373 // e+ correction
374 if (!fIsElectron) {
375 const G4double e1 = kinEnergy - tmin;
376 const G4double invbeta1 = (e1+kMC2)/std::sqrt(e1*(e1+2.*kMC2));
377 const G4double e2 = kinEnergy-gammaEnergy;
378 const G4double invbeta2 = (e2+kMC2)/std::sqrt(e2*(e2+2.*kMC2));
379 const G4double dum0 = kAlpha*fCurrentIZ*(invbeta1-invbeta2);
380 if (dum0 < gExpNumLimit) {
381 v = 0.0;
382 } else {
383 v *= G4Exp(dum0);
384 }
385 }
386 if (v > 1.05*vmax && fNumWarnings < 11) {
387 ++fNumWarnings;
389 ed << "### G4SeltzerBergerModel Warning: Majoranta exceeded! "
390 << v << " > " << vmax << " by " << v/vmax
391 << " Niter= " << nn
392 << " Egamma(MeV)= " << gammaEnergy
393 << " Ee(MeV)= " << kinEnergy
394 << " Z= " << fCurrentIZ << " " << fPrimaryParticle->GetParticleName();
395 //
396 if (10 == fNumWarnings) {
397 ed << "\n ### G4SeltzerBergerModel Warnings stopped";
398 }
399 G4Exception("G4SeltzerBergerModel::SampleScattering","em0044",
400 JustWarning, ed,"");
401 }
402 if (v >= vmax*rndm[1]) {
403 break;
404 }
405 }
406 return gammaEnergy;
407}
408
410 const G4Material* mat,
411 G4double kineticEnergy)
412{
414 // calculate threshold for density effect: gamma*k_p = sqrt(fDensityCorr)
415 fPrimaryKinEnergy = kineticEnergy;
416 fPrimaryTotalEnergy = kineticEnergy+CLHEP::electron_mass_c2;
419}
420
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
std::mutex G4Mutex
Definition: G4Threading.hh:81
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
virtual void flatArray(const int size, double *vect)=0
const G4ThreeVector & GetMomentumDirection() const
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
G4int GetZasInt() const
Definition: G4Element.hh:131
static G4EmParameters * Instance()
G4bool EnableSamplingTable() const
const G4Material * GetMaterial() const
G4double GetElectronDensity() const
Definition: G4Material.hh:215
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
const G4String & GetParticleName() const
G4bool Retrieve(std::ifstream &fIn)
G4double Value(G4double x, G4double y, std::size_t &lastidx, std::size_t &lastidy) const
void SetBicubicInterpolation(G4bool)
static G4ProductionCutsTable * GetProductionCutsTable()
double SampleEnergyTransfer(const G4double eekin, const G4double leekin, const G4double gcut, const G4double dielSupConst, const G4int izet, const G4int matCutIndx, const bool iselectron)
void Initialize(const G4double lowe, const G4double highe)
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
G4double ComputeDXSectionPerAtom(G4double gammaEnergy) override
G4SeltzerBergerModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="eBremSB")
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:611
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
G4VEmModel * GetTripletModel()
Definition: G4VEmModel.hh:628
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:806
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:415
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
G4bool LPMFlag() const
Definition: G4VEmModel.hh:687
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:587
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:680
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118
void ProposeTrackStatus(G4TrackStatus status)
void SetParticle(const G4ParticleDefinition *p)
const G4ParticleDefinition * fPrimaryParticle
G4ParticleDefinition * fGammaParticle
G4ParticleChangeForLoss * fParticleChange