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
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// $Id$
27//
28// Author: Luciano Pandola
29//
30// History:
31// -----------
32// 23 Aug 2010 L. Pandola 1st implementation.
33// 24 May 2011 L. Pandola Renamed (make v2008 as default Penelope)
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;
57class G4PhysicsTable;
62
64{
65
66public:
67
69 const G4String& processName ="PenBrem");
70
72
73 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
74
75 //DUMMY METHOD
77 G4double kinEnergy,
78 G4double Z,
79 G4double A=0,
80 G4double cut=0,
81 G4double emax=DBL_MAX);
82
83
84 virtual G4double CrossSectionPerVolume(const G4Material* material,
85 const G4ParticleDefinition* theParticle,
86 G4double kineticEnergy,
87 G4double cutEnergy,
88 G4double maxEnergy = DBL_MAX);
89
90
91 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
93 const G4DynamicParticle*,
94 G4double tmin,
95 G4double maxEnergy);
96
99 G4double kineticEnergy,
100 G4double cutEnergy);
101
102 // Min cut in kinetic energy allowed by the model
104 const G4MaterialCutsCouple*);
105
106 void SetVerbosityLevel(G4int lev){verboseLevel = lev;};
107 G4int GetVerbosityLevel(){return verboseLevel;};
108
109protected:
111
112private:
113 void ClearTables();
114
117
118 G4PenelopeCrossSection* GetCrossSectionTableForCouple(const G4ParticleDefinition*,
119 const G4Material*,G4double cut);
120
121 //Intrinsic energy limits of the model: cannot be extended by the parent process
122 G4double fIntrinsicLowEnergyLimit;
123 G4double fIntrinsicHighEnergyLimit;
124
125 G4int verboseLevel;
126
127 G4bool isInitialised;
128
129 G4PenelopeOscillatorManager* oscManager;
130
131 //
132 //Members to handle and store cross section tables
133 void BuildXSTable(const G4Material*,G4double cut);
134 //This is the main energy grid
135 G4PhysicsLogVector* energyGrid;
136 size_t nBins;
137 //G4PenelopeCrossSection takes care of the logs
138 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*> *XSTableElectron;
139 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*> *XSTablePositron;
140 G4double GetPositronXSCorrection(const G4Material*,G4double energy);
141
142 //Helpers
143 G4PenelopeBremsstrahlungFS* fPenelopeFSHelper;
144 G4PenelopeBremsstrahlungAngular* fPenelopeAngular;
145
146};
147
148#endif
149
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
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)
#define DBL_MAX
Definition: templates.hh:83