Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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// $Id$
27//
28// Author: Luciano Pandola
29//
30// History:
31// -----------
32// 20 Oct 2010 L. Pandola 1st implementation.
33// 02 May 2011 L. Pandola Remove dependency on CLHEP::HepMatrix
34// 24 May 2011 L. Pandola Renamed (make default Penelope)
35//
36// -------------------------------------------------------------------
37//
38// Class description:
39// Helper class for the calculation of final state (energy and angular
40// distribution) for bremsstrahlung, Penelope Model, version 2008
41// -------------------------------------------------------------------
42
43#ifndef G4PENELOPEBREMSSTRAHLUNGFS_HH
44#define G4PENELOPEBREMSSTRAHLUNGFS_HH 1
45
46#include "globals.hh"
47#include "G4DataVector.hh"
48#include <map>
49
52class G4Material;
53class G4PhysicsTable;
54
56{
57public:
58
61
63 void ClearTables();
64 size_t GetNBinsX(){return nBinsX;};
65
66
67 //utilities
69 G4double up,G4int momOrder);
70
73
74private:
75 //assignment operator
77 //copy constructor
79
80 //Differential cross section tables
81 //Table contains G4PhysicsVectors of log(XS(E,x)) vs. log(E) for a grid of 32 values in
82 //x (= reduced photon energy)
83 std::map< std::pair<const G4Material*,G4double> ,
84 G4PhysicsTable*> *theReducedXSTable;
85
86 void BuildScaledXSTable(const G4Material* material,G4double cut);
87 std::map<const G4Material*,G4double> *theEffectiveZSq;
88
89
90 //Element data table
91 static const size_t nBinsE = 57;
92 static const size_t nBinsX = 32;
93 //x and E grids used in the data tables
94 G4double theXGrid[nBinsX];
95 G4double theEGrid[nBinsE];
96 void ReadDataFile(G4int Z);
97
98 //Map of element data vs. Z.
99 //This is conceptually an array [57][33] with 57 energy
100 //points and 32 points in x. The 33-th column gives the total XS vs. E.
101 //It is implemented as a one-dimensional array of dimension
102 //57*33=1881 elements. data[e][x] --> theElementData[e*(nBinsX+1)+x]
103 std::map<G4int,G4DataVector*> *theElementData;
104
105 //Tables for energy sampling
106 void InitializeEnergySampling(const G4Material*,G4double cut);
107
108 //Table contains G4PhysicsVectors of integral_XS(E,x) vs. x for a grid of
109 //57 values in energy
110 std::map< std::pair<const G4Material*,G4double> ,
111 G4PhysicsTable*> *theSamplingTable;
112 std::map< std::pair<const G4Material*,G4double> ,
113 G4PhysicsFreeVector* > *thePBcut;
114
115 //temporary vector. Used as member variable to avoid to book/release the
116 //memory on the fly. This vector is over-written at every call of
117 //SampleGammaEnergy()
118 G4PhysicsFreeVector* theTempVec;
119
120};
121
122#endif
123
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4double GetMomentumIntegral(G4double *y, G4double up, G4int momOrder)
G4double SampleGammaEnergy(G4double energy, const G4Material *, G4double cut)
G4double GetEffectiveZSquared(const G4Material *mat)
G4PhysicsTable * GetScaledXSTable(const G4Material *, G4double cut)