Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4SeltzerBergerModel.hh
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//
29// GEANT4 Class header file
30//
31//
32// File name: G4SeltzerBergerModel
33//
34// Author: Andreas Schaelicke & Vladimir Ivantchenko
35//
36// Creation date: 04.10.2011
37//
38// Modifications:
39//
40// 24.07.2018 Introduced possibility to use sampling tables to sample the
41// emitted photon energy (instead of using rejectio) from the
42// Seltzer-Berger scalled DCS for bremsstrahlung photon emission.
43// Using these sampling tables option gives faster(30-70%) final
44// state generation than the original rejection but takes some
45// extra memory (+ ~6MB in the case of the full CMS detector).
46// (M Novak)
47//
48// Class Description:
49//
50// Implementation of the bremssrahlung energy spectrum using
51// 1. S.M. Seltzer and M.J. Berger Nucl. Instr. Meth. B12 (1985) 95
52// 2. S.M. Seltzer and M.J. Berger Atomic data and Nuclear Data
53// Tables 35 (1986) 345
54// Cross section computation in the base class G4eBremsstrahlungRelModel
55
56// -------------------------------------------------------------------
57//
58
59#ifndef G4SeltzerBergerModel_h
60#define G4SeltzerBergerModel_h 1
61
63#include "globals.hh"
64
66class G4SBBremTable;
67
69{
70
71public:
72
73 explicit G4SeltzerBergerModel(const G4ParticleDefinition* p = nullptr,
74 const G4String& nam = "eBremSB");
75
76 ~G4SeltzerBergerModel() override;
77
78 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
79
80 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
82 const G4DynamicParticle*,
83 G4double cutEnergy,
84 G4double maxEnergy) override;
85
87 const G4Material*, G4double) override;
88
90 { fIsUseBicubicInterpolation = val; };
91
92 // hide assignment operator and cctr
95
96protected:
97
98 G4double ComputeDXSectionPerAtom(G4double gammaEnergy) override;
99
100private:
101
102 void ReadData(G4int Z);
103
104 const G4String& FindDirectoryPath();
105
106 G4double SampleEnergyTransfer(const G4double kineticEnergy,
107 const G4double logKineticEnergy,
108 const G4double cut,
109 const G4double emax);
110
111 static constexpr G4int gMaxZet = 101;
112 static constexpr G4double gExpNumLimit = -12.;
113 static G4double gYLimitData[gMaxZet];
114 static G4Physics2DVector* gSBDCSData[gMaxZet];
115 static G4SBBremTable* gSBSamplingTable;
116 static G4String gDataDirectory;
117
118 G4bool fIsUseBicubicInterpolation;
119 G4bool fIsUseSamplingTables;
120
121 G4int fNumWarnings;
122
123 size_t fIndx;
124 size_t fIndy;
125};
126
127#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
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
G4SeltzerBergerModel(const G4SeltzerBergerModel &)=delete
G4double ComputeDXSectionPerAtom(G4double gammaEnergy) override
void SetBicubicInterpolationFlag(G4bool val)
G4SeltzerBergerModel & operator=(const G4SeltzerBergerModel &right)=delete