Geant4 11.1.1
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#include "G4AutoLock.hh"
70
71#include "G4ios.hh"
72
73#include <fstream>
74#include <iomanip>
75#include <sstream>
76
77G4Physics2DVector* G4SeltzerBergerModel::gSBDCSData[] = { nullptr };
78G4SBBremTable* G4SeltzerBergerModel::gSBSamplingTable = nullptr;
79G4double G4SeltzerBergerModel::gYLimitData[] = { 0.0 };
80G4String G4SeltzerBergerModel::gDataDirectory = "";
81
82namespace
83{
84 G4Mutex theSBMutex = G4MUTEX_INITIALIZER;
85}
86
87static const G4double kMC2 = CLHEP::electron_mass_c2;
88static const G4double kAlpha = CLHEP::twopi*CLHEP::fine_structure_const;
89
91 const G4String& nam)
92: G4eBremsstrahlungRelModel(p,nam), fIsUseBicubicInterpolation(false),
93 fIsUseSamplingTables(true), fNumWarnings(0), fIndx(0), fIndy(0)
94{
95 fLowestKinEnergy = 1.0*keV;
97 SetLPMFlag(false);
99}
100
102{
103 // delete SB-DCS data per Z
104 if (IsMaster()) {
105 for (std::size_t iz = 0; iz < gMaxZet; ++iz) {
106 if (gSBDCSData[iz]) {
107 delete gSBDCSData[iz];
108 gSBDCSData[iz] = nullptr;
109 }
110 }
111 if (gSBSamplingTable) {
112 delete gSBSamplingTable;
113 gSBSamplingTable = nullptr;
114 }
115 }
116}
117
119 const G4DataVector& cuts)
120{
121 if (p) {
122 SetParticle(p);
123 }
124 fIsUseSamplingTables = G4EmParameters::Instance()->EnableSamplingTable();
125 // Access to elements
126 if (IsMaster()) {
127
128 auto theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
129 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
130 for(G4int j=0; j<numOfCouples; ++j) {
131 auto mat = theCoupleTable->GetMaterialCutsCouple(j)->GetMaterial();
132 auto elmVec = mat->GetElementVector();
133 for (auto & elm : *elmVec) {
134 G4int Z = std::max(1,std::min(elm->GetZasInt(), gMaxZet-1));
135 // load SB-DCS data for this atomic number if it has not been loaded yet
136 if (gSBDCSData[Z] == nullptr) ReadData(Z);
137 }
138 }
139 // elem.selectr. only for master: base class init-local will set for workers
142 }
143 // init sampling tables if it was requested
144 if (fIsUseSamplingTables) {
145 if (!gSBSamplingTable) {
146 gSBSamplingTable = new G4SBBremTable();
147 }
148 gSBSamplingTable->Initialize(std::max(fLowestKinEnergy,LowEnergyLimit()),
150 }
151 }
152 //
154 if (GetTripletModel()) {
155 GetTripletModel()->Initialise(p, cuts);
156 fIsScatOffElectron = true;
157 }
158}
159
160const G4String& G4SeltzerBergerModel::FindDirectoryPath()
161{
162 // check environment variable
163 // build the complete string identifying the file with the data set
164 if(gDataDirectory.empty()) {
165 const char* path = G4FindDataDir("G4LEDATA");
166 if (path) {
167 std::ostringstream ost;
168 ost << path << "/brem_SB/br";
169 gDataDirectory = ost.str();
170 } else {
171 G4Exception("G4SeltzerBergerModel::FindDirectoryPath()","em0006",
173 "Environment variable G4LEDATA not defined");
174 }
175 }
176 return gDataDirectory;
177}
178
179void G4SeltzerBergerModel::ReadData(G4int Z) {
180 // return if it has been already loaded
181 if (gSBDCSData[Z] != nullptr) return;
182
183 G4AutoLock l(&theSBMutex);
184 if (gSBDCSData[Z] == nullptr) {
185 std::ostringstream ost;
186 ost << FindDirectoryPath() << Z;
187 std::ifstream fin(ost.str().c_str());
188 if (!fin.is_open()) {
190 ed << "Bremsstrahlung data file <" << ost.str().c_str()
191 << "> is not opened!";
192 G4Exception("G4SeltzerBergerModel::ReadData()","em0003",FatalException,
193 ed,"G4LEDATA version should be G4EMLOW6.23 or later.");
194 return;
195 }
196 //G4cout << "G4SeltzerBergerModel read from <" << ost.str().c_str()
197 // << ">" << G4endl;
198 auto v = new G4Physics2DVector();
199 if (v->Retrieve(fin)) {
200 v->SetBicubicInterpolation(fIsUseBicubicInterpolation);
201 static const G4double emaxlog = 4*G4Log(10.);
202 gYLimitData[Z] = v->Value(0.97, emaxlog, fIndx, fIndy);
203 gSBDCSData[Z] = v;
204 } else {
206 ed << "Bremsstrahlung data file <" << ost.str().c_str()
207 << "> is not retrieved!";
208 G4Exception("G4SeltzerBergerModel::ReadData()","em0005",FatalException,
209 ed,"G4LEDATA version should be G4EMLOW6.23 or later.");
210 delete v;
211 }
212 }
213 l.unlock();
214}
215
217{
218 G4double dxsec = 0.0;
219 if (gammaEnergy < 0.0 || fPrimaryKinEnergy <= 0.0) {
220 return dxsec;
221 }
222 // reduced photon energy
223 const G4double x = gammaEnergy/fPrimaryKinEnergy;
224 // l-kinetic energy of the e-/e+
225 const G4double y = G4Log(fPrimaryKinEnergy/CLHEP::MeV);
226 // make sure that the Z-related SB-DCS are loaded
227 // NOTE: fCurrentIZ should have been set before.
228 fCurrentIZ = std::max(std::min(fCurrentIZ, gMaxZet-1), 1);
229 if (nullptr == gSBDCSData[fCurrentIZ]) {
230 ReadData(fCurrentIZ);
231 }
232 // NOTE: SetupForMaterial should have been called before!
233 const G4double pt2 = fPrimaryKinEnergy*(fPrimaryKinEnergy+2.*kMC2);
235 G4double val = gSBDCSData[fCurrentIZ]->Value(x,y,fIndx,fIndy);
236 dxsec = val*invb2*CLHEP::millibarn/gBremFactor;
237 // e+ correction
238 if (!fIsElectron) {
239 const G4double invbeta1 = std::sqrt(invb2);
240 const G4double e2 = fPrimaryKinEnergy-gammaEnergy;
241 if (e2 > 0.0) {
242 const G4double invbeta2 = (e2+kMC2)/std::sqrt(e2*(e2+2.0*kMC2));
243 const G4double dum0 = kAlpha*fCurrentIZ*(invbeta1-invbeta2);
244 if (dum0 < gExpNumLimit) {
245 dxsec = 0.0;
246 } else {
247 dxsec *= G4Exp(dum0);
248 }
249 } else {
250 dxsec = 0.0;
251 }
252 }
253 return dxsec;
254}
255
256void
257G4SeltzerBergerModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
258 const G4MaterialCutsCouple* couple,
259 const G4DynamicParticle* dp,
260 G4double cutEnergy,
261 G4double maxEnergy)
262{
263 const G4double kinEnergy = dp->GetKineticEnergy();
264 const G4double logKinEnergy = dp->GetLogKineticEnergy();
265 const G4double tmin = std::min(cutEnergy, kinEnergy);
266 const G4double tmax = std::min(maxEnergy, kinEnergy);
267 if (tmin >= tmax) {
268 return;
269 }
270 // set local variables and select target element
271 SetupForMaterial(fPrimaryParticle, couple->GetMaterial(), kinEnergy);
272 const G4Element* elm = SelectTargetAtom(couple, fPrimaryParticle, kinEnergy,
273 logKinEnergy, tmin, tmax);
274 fCurrentIZ = std::max(std::min(elm->GetZasInt(),gMaxZet-1), 1);
275 //
276 const G4double totMomentum = std::sqrt(kinEnergy*(fPrimaryTotalEnergy+kMC2));
277 /*
278 G4cout << "G4SeltzerBergerModel::SampleSecondaries E(MeV)= "
279 << kinEnergy/MeV
280 << " Z= " << fCurrentIZ << " cut(MeV)= " << tmin/MeV
281 << " emax(MeV)= " << tmax/MeV << " corr= " << fDensityCorr << G4endl;
282 */
283 // sample emitted photon energy either by rejection or from samplign tables
284 const G4double gammaEnergy = !fIsUseSamplingTables
285 ? SampleEnergyTransfer(kinEnergy, logKinEnergy, tmin, tmax)
286 : gSBSamplingTable->SampleEnergyTransfer(kinEnergy, logKinEnergy, tmin,
288 // should never happen under normal conditions but protect it
289 if (gammaEnergy <= 0.) {
290 return;
291 }
292 //
293 // angles of the emitted gamma. ( Z - axis along the parent particle) use
294 // general interface
296 fPrimaryTotalEnergy-gammaEnergy, fCurrentIZ, couple->GetMaterial());
297 // create G4DynamicParticle object for the emitted Gamma
298 auto gamma = new G4DynamicParticle(fGammaParticle, gamDir, gammaEnergy);
299 vdp->push_back(gamma);
300 //
301 // compute post-interaction kinematics of the primary e-/e+
302 G4ThreeVector dir =
303 (totMomentum*dp->GetMomentumDirection()-gammaEnergy*gamDir).unit();
304 const G4double finalE = kinEnergy - gammaEnergy;
305 /*
306 G4cout << "### G4SBModel: v= "
307 << " Eg(MeV)= " << gammaEnergy
308 << " Ee(MeV)= " << kineticEnergy
309 << " DirE " << direction << " DirG " << gammaDirection
310 << G4endl;
311 */
312 // if secondary gamma energy is higher than threshold(very high by default)
313 // then stop tracking the primary particle and create new secondary e-/e+
314 // instead of the primary
315 if (gammaEnergy > SecondaryThreshold()) {
318 auto el = new G4DynamicParticle(
319 const_cast<G4ParticleDefinition*>(fPrimaryParticle), dir, finalE);
320 vdp->push_back(el);
321 } else { // continue tracking the primary e-/e+ otherwise
324 }
325}
326
327// sample emitted photon energy by usign rejection
329G4SeltzerBergerModel::SampleEnergyTransfer(const G4double kinEnergy,
330 const G4double logKinEnergy,
331 const G4double tmin,
332 const G4double tmax)
333{
334 // min max of the transformed variable: x(k) = ln(k^2+k_p^2) that is in
335 // [ln(k_c^2+k_p^2), ln(E_k^2+k_p^2)]
336 const G4double xmin = G4Log(tmin*tmin+fDensityCorr);
337 const G4double xrange = G4Log(tmax*tmax+fDensityCorr)-xmin;
338 const G4double y = logKinEnergy;
339 // majoranta
340 const G4double x0 = tmin/kinEnergy;
341 G4double vmax;
342 if (nullptr == gSBDCSData[fCurrentIZ]) {
343 ReadData(fCurrentIZ);
344 }
345 vmax = gSBDCSData[fCurrentIZ]->Value(x0, y, fIndx, fIndy)*1.02;
346 //
347 static const G4double kEPeakLim = 300.*CLHEP::MeV;
348 static const G4double kELowLim = 20.*CLHEP::keV;
349 // majoranta corrected for e-
350 if (fIsElectron && x0 < 0.97 &&
351 ((kinEnergy>kEPeakLim) || (kinEnergy<kELowLim))) {
352 G4double ylim = std::min(gYLimitData[fCurrentIZ],
353 1.1*gSBDCSData[fCurrentIZ]->Value(0.97,y,fIndx,fIndy));
354 vmax = std::max(vmax, ylim);
355 }
356 if (x0 < 0.05) {
357 vmax *= 1.2;
358 }
359 //G4cout<<"y= "<<y<<" xmin= "<<xmin<<" xmax= "<<xmax
360 //<<" vmax= "<<vmax<<G4endl;
361 static const G4int kNCountMax = 100;
362 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
363 G4double rndm[2];
364 G4double gammaEnergy, v;
365 for (G4int nn = 0; nn < kNCountMax; ++nn) {
366 rndmEngine->flatArray(2, rndm);
367 gammaEnergy =
368 std::sqrt(std::max(G4Exp(xmin + rndm[0]*xrange)-fDensityCorr,0.));
369 v = gSBDCSData[fCurrentIZ]->Value(gammaEnergy/kinEnergy, y, fIndx, fIndy);
370 // e+ correction
371 if (!fIsElectron) {
372 const G4double e1 = kinEnergy - tmin;
373 const G4double invbeta1 = (e1+kMC2)/std::sqrt(e1*(e1+2.*kMC2));
374 const G4double e2 = kinEnergy-gammaEnergy;
375 const G4double invbeta2 = (e2+kMC2)/std::sqrt(e2*(e2+2.*kMC2));
376 const G4double dum0 = kAlpha*fCurrentIZ*(invbeta1-invbeta2);
377 if (dum0 < gExpNumLimit) {
378 v = 0.0;
379 } else {
380 v *= G4Exp(dum0);
381 }
382 }
383 if (v > 1.05*vmax && fNumWarnings < 11) {
384 ++fNumWarnings;
386 ed << "### G4SeltzerBergerModel Warning: Majoranta exceeded! "
387 << v << " > " << vmax << " by " << v/vmax
388 << " Niter= " << nn
389 << " Egamma(MeV)= " << gammaEnergy
390 << " Ee(MeV)= " << kinEnergy
391 << " Z= " << fCurrentIZ << " " << fPrimaryParticle->GetParticleName();
392 //
393 if (10 == fNumWarnings) {
394 ed << "\n ### G4SeltzerBergerModel Warnings stopped";
395 }
396 G4Exception("G4SeltzerBergerModel::SampleScattering","em0044",
397 JustWarning, ed,"");
398 }
399 if (v >= vmax*rndm[1]) {
400 break;
401 }
402 }
403 return gammaEnergy;
404}
405
407 const G4Material* mat,
408 G4double kineticEnergy)
409{
411 // calculate threshold for density effect: gamma*k_p = sqrt(fDensityCorr)
412 fPrimaryKinEnergy = kineticEnergy;
413 fPrimaryTotalEnergy = kineticEnergy+CLHEP::electron_mass_c2;
416}
417
const char * G4FindDataDir(const char *)
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4double G4Log(G4double x)
Definition: G4Log.hh:227
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
std::mutex G4Mutex
Definition: G4Threading.hh:81
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
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:132
static G4EmParameters * Instance()
G4bool EnableSamplingTable() const
const G4Material * GetMaterial() const
G4double GetElectronDensity() const
Definition: G4Material.hh:212
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
const G4String & GetParticleName() const
G4double Value(G4double x, G4double y, std::size_t &lastidx, std::size_t &lastidy) const
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:600
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
G4VEmModel * GetTripletModel()
Definition: G4VEmModel.hh:617
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:795
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:349
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
G4bool LPMFlag() const
Definition: G4VEmModel.hh:676
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:577
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:139
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:669
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:109
void ProposeTrackStatus(G4TrackStatus status)
void SetParticle(const G4ParticleDefinition *p)
const G4ParticleDefinition * fPrimaryParticle
G4ParticleDefinition * fGammaParticle
G4ParticleChangeForLoss * fParticleChange