Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PenelopeBremsstrahlungModel.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// 23 Aug 2010 L. Pandola 1st implementation.
32// 24 May 2011 L. Pandola Renamed (make v2008 as default Penelope)
33// 02 Oct 2013 L. Pandola Migrated to multi-thread
34//
35// -------------------------------------------------------------------
36//
37// Class description:
38// Low Energy Electromagnetic Physics, e+ and e- bremsstrahlung
39// with Penelope Model, version 2008
40// -------------------------------------------------------------------
41
42#ifndef G4PENELOPEBREMSSTRAHLUNGMODEL_HH
43#define G4PENELOPEBREMSSTRAHLUNGMODEL_HH 1
44
45#include "globals.hh"
46#include "G4VEmModel.hh"
47#include "G4DataVector.hh"
49
50
56class G4Material;
61
63{
64
65public:
66
68 const G4String& processName ="PenBrem");
69
71
72 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
73 virtual void InitialiseLocal(const G4ParticleDefinition*,
74 G4VEmModel*);
75
76 //DUMMY METHOD
78 G4double kinEnergy,
79 G4double Z,
80 G4double A=0,
81 G4double cut=0,
82 G4double emax=DBL_MAX);
83
84
85 virtual G4double CrossSectionPerVolume(const G4Material* material,
86 const G4ParticleDefinition* theParticle,
87 G4double kineticEnergy,
88 G4double cutEnergy,
89 G4double maxEnergy = DBL_MAX);
90
91
92 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
94 const G4DynamicParticle*,
95 G4double tmin,
96 G4double maxEnergy);
97
100 G4double kineticEnergy,
101 G4double cutEnergy);
102
103 // Min cut in kinetic energy allowed by the model
105 const G4MaterialCutsCouple*);
106
107 void SetVerbosityLevel(G4int lev){verboseLevel = lev;};
108 G4int GetVerbosityLevel(){return verboseLevel;};
109
110protected:
113
114private:
115 void ClearTables();
116
119
120 G4PenelopeCrossSection* GetCrossSectionTableForCouple(const G4ParticleDefinition*,
121 const G4Material*,G4double cut);
122 void SetParticle(const G4ParticleDefinition*);
123
124 //Intrinsic energy limits of the model: cannot be extended by the parent process
125 G4double fIntrinsicLowEnergyLimit;
126 G4double fIntrinsicHighEnergyLimit;
127
128 G4int verboseLevel;
129
130 G4bool isInitialised;
131
132 G4PenelopeOscillatorManager* oscManager;
133
134 //
135 //Members to handle and store cross section tables
136 void BuildXSTable(const G4Material* material,G4double cut);
137
138 //This is the main energy grid
139 G4PhysicsLogVector* energyGrid;
140 size_t nBins;
141 //G4PenelopeCrossSection takes care of the logs
142 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*> *XSTableElectron;
143 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*> *XSTablePositron;
144 G4double GetPositronXSCorrection(const G4Material*,G4double energy);
145
146 //Helpers
147 G4PenelopeBremsstrahlungFS* fPenelopeFSHelper;
148 G4PenelopeBremsstrahlungAngular* fPenelopeAngular;
149
150 //Used only for G4EmCalculator and Unit Tests
151 G4bool fLocalTable;
152
153};
154
155#endif
156
double A(double temperature)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *theParticle, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy=DBL_MAX)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *theParticle, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *)
#define DBL_MAX
Definition: templates.hh:62