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
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{
64public:
66 const G4String& processName ="PenBrem");
68
69 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
71 G4VEmModel*) override;
72
73 //DUMMY METHOD
75 G4double kinEnergy,
76 G4double Z,
77 G4double A=0,
78 G4double cut=0,
79 G4double emax=DBL_MAX) override;
80
81
83 const G4ParticleDefinition* theParticle,
84 G4double kineticEnergy,
85 G4double cutEnergy,
86 G4double maxEnergy = DBL_MAX) override;
87
88 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
90 const G4DynamicParticle*,
91 G4double tmin,
92 G4double maxEnergy) override;
93
96 G4double kineticEnergy,
97 G4double cutEnergy) override;
98
99 // Min cut in kinetic energy allowed by the model
101 const G4MaterialCutsCouple*) override;
102
103 void SetVerbosityLevel(G4int lev){fVerboseLevel = lev;};
104 G4int GetVerbosityLevel(){return fVerboseLevel;};
105
107 (const G4PenelopeBremsstrahlungModel &right) = delete;
109
110protected:
113
114private:
115 void ClearTables();
116 G4PenelopeCrossSection* GetCrossSectionTableForCouple(const G4ParticleDefinition*,
117 const G4Material*,G4double cut);
118 void SetParticle(const G4ParticleDefinition*);
119 //
120 //Members to handle and store cross section tables
121 void BuildXSTable(const G4Material* material,G4double cut);
122 G4double GetPositronXSCorrection(const G4Material*,G4double energy);
123
124 //Helpers
125 G4PenelopeOscillatorManager* fOscManager;
126 G4PenelopeBremsstrahlungFS* fPenelopeFSHelper;
127 G4PenelopeBremsstrahlungAngular* fPenelopeAngular;
128
129 //This is the main energy grid
130 G4PhysicsLogVector* fEnergyGrid;
131 size_t nBins;
132 //G4PenelopeCrossSection takes care of the logs
133 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*> *fXSTableElectron;
134 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*> *fXSTablePositron;
135
136 //Intrinsic energy limits of the model: cannot be extended by the parent process
137 G4double fIntrinsicLowEnergyLimit;
138 G4double fIntrinsicHighEnergyLimit;
139
140 G4int fVerboseLevel;
141 G4bool fIsInitialised;
142 //Used only for G4EmCalculator and Unit Tests
143 G4bool fLocalTable;
144};
145
146#endif
147
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *theParticle, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy=DBL_MAX) override
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *theParticle, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4PenelopeBremsstrahlungModel(const G4PenelopeBremsstrahlungModel &)=delete
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
#define DBL_MAX
Definition: templates.hh:62