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
G4PenelopeBremsstrahlungFS.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// Author: Luciano Pandola
28//
29// History:
30// -----------
31// 20 Oct 2010 L. Pandola 1st implementation.
32// 02 May 2011 L. Pandola Remove dependency on CLHEP::HepMatrix
33// 24 May 2011 L. Pandola Renamed (make default Penelope)
34// 03 Oct 2013 L. Pandola Migration to MT
35// 07 Oct 2013 L. Pandola Add verbosity and ismaster flag for the
36// master-only methods
37// 30 Oct 2013 L. Pandola Add a G4Cache of temporary vectors as
38// private member.
39//
40// -------------------------------------------------------------------
41//
42// Class description:
43// Helper class for the calculation of final state (energy and angular
44// distribution) for bremsstrahlung, Penelope Model, version 2008
45// -------------------------------------------------------------------
46
47#ifndef G4PENELOPEBREMSSTRAHLUNGFS_HH
48#define G4PENELOPEBREMSSTRAHLUNGFS_HH 1
49
50#include "globals.hh"
51#include "G4DataVector.hh"
52#include "G4Cache.hh"
53#include <map>
54
57class G4Material;
58class G4PhysicsTable;
59
61{
62public:
63 //! Only master models are supposed to create instances
64 explicit G4PenelopeBremsstrahlungFS(G4int verbosity=0);
66
67 //!
68 //! Master and workers (do not touch tables)
69 //! All of them are const
70 //!
72 size_t GetNBinsX() const {return fNBinsX;};
74 G4double up,G4int momOrder) const;
76 const G4double cut) const;
78 const G4Material*, const G4double cut) const;
79
80 //! Reserved for the master model: they build and handle tables
81 void ClearTables(G4bool isMaster=true);
82 void BuildScaledXSTable(const G4Material* material,G4double cut,
83 G4bool isMaster);
84
85 void SetVerbosity(G4int ver){fVerbosity=ver;};
86 G4int GetVerbosity(){return fVerbosity;};
87
90
91private:
92 void ReadDataFile(G4int Z);
93 //Tables for energy sampling
94 void InitializeEnergySampling(const G4Material*,G4double cut);
95
96 //Differential cross section tables
97 //Table contains G4PhysicsVectors of log(XS(E,x)) vs. log(E)
98 //for a grid of 32 values in x (= reduced photon energy)
99 std::map< std::pair<const G4Material*,G4double> ,
100 G4PhysicsTable*> *fReducedXSTable;
101
102 std::map<const G4Material*,G4double> *fEffectiveZSq;
103
104 //Map of element data vs. Z.
105 //This is conceptually an array [57][33] with 57 energy
106 //points and 32 points in x. The 33-th column gives the total XS vs. E.
107 //It is implemented as a one-dimensional array of dimension
108 //57*33=1881 elements. data[e][x] --> theElementData[e*(nBinsX+1)+x]
109 std::map<G4int,G4DataVector*> *fElementData;
110
111 //Table contains G4PhysicsVectors of integral_XS(E,x) vs. x for a grid of
112 //57 values in energy
113 std::map< std::pair<const G4Material*,G4double> ,
114 G4PhysicsTable*> *fSamplingTable;
115 std::map< std::pair<const G4Material*,G4double> ,
116 G4PhysicsFreeVector* > *fPBcut;
117
118 //temporary vector. Used as member variable to avoid to book/release the
119 //memory on the fly. This vector is over-written at every call of
120 //SampleGammaEnergy(). It is thread-local (each thread has its own)
121 //and managed by G4Cache
123
124 //Element data table
125 static const size_t fNBinsE = 57;
126 static const size_t fNBinsX = 32;
127 //x and E grids used in the data tables
128 G4double theXGrid[fNBinsX];
129 G4double theEGrid[fNBinsE];
130
131 G4int fVerbosity;
132
133};
134
135#endif
136
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
void BuildScaledXSTable(const G4Material *material, G4double cut, G4bool isMaster)
G4double GetEffectiveZSquared(const G4Material *mat) const
const G4PhysicsTable * GetScaledXSTable(const G4Material *, const G4double cut) const
G4PenelopeBremsstrahlungFS & operator=(const G4PenelopeBremsstrahlungFS &right)=delete
void ClearTables(G4bool isMaster=true)
Reserved for the master model: they build and handle tables.
G4PenelopeBremsstrahlungFS(const G4PenelopeBremsstrahlungFS &)=delete
G4double SampleGammaEnergy(G4double energy, const G4Material *, const G4double cut) const
G4double GetMomentumIntegral(G4double *y, G4double up, G4int momOrder) const