Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PenelopeBremsstrahlungAngular.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// -------------------------------------------------------------------
28// $Id$
29//
30// Author: L.Pandola
31//
32// History:
33// -----------
34// 23 Nov 2010 L. Pandola 1st implementation
35// 24 May 2011 L. Pandola Renamed (make default Penelope)
36// 13 Mar 2012 L. Pandola Made a derived class of G4VEmAngularDistribution
37// and update the interface accordingly
38// 18 Jul 2012 L. Pandola Migrate to the new interface of G4VEmAngularDistribution
39//
40// Class description:
41// Calculation of angular distribution for Penelope Bremsstrahlung
42// version 2008
43// --------------------------------------------------------------
44
45
46#ifndef G4PENELOPEBREMSSTRAHLUNGANGULAR_HH
47#define G4PENELOPEBREMSSTRAHLUNGANGULAR_HH 1
48#include "globals.hh"
49#include <map>
51
52class G4PhysicsTable;
53class G4Material;
54
56{
57
58public:
61
62 //! The Initialize() method forces the cleaning of tables
63 void Initialize();
64
65 //! Old interface, backwards compatibility. Will not work in this case
66 //! it will produce a G4Exception().
67 G4double PolarAngle(const G4double initial_energy,
68 const G4double final_energy,
69 const G4int Z);
70
71 //! Samples the direction of the outgoing photon (in global coordinates).
72 //! Forces the calculation of tables, if they are not available
74 G4double out_energy,
75 G4int Z,
76 const G4Material* mat = 0);
77
78 //! Set/Get Verbosity level
79 void SetVerbosityLevel(G4int vl){verbosityLevel = vl;};
80 G4int GetVerbosityLevel(){return verbosityLevel;};
81
82private:
83 void PrepareInterpolationTables(G4double Zeq);
84 void ClearTables();
85
86 G4double GetEffectiveZ(const G4Material* material);
87 std::map<const G4Material*,G4double> *theEffectiveZSq;
88
89 //Tables containing the Lorentz sampling coefficients
90 //The key is the effective Z of the material
91 std::map<G4double,G4PhysicsTable*> *theLorentzTables1;
92 std::map<G4double,G4PhysicsTable*> *theLorentzTables2;
93
94 void ReadDataFile();
95 G4bool dataRead;
96
97 static const G4int NumberofZPoints=6;
98 static const G4int NumberofEPoints=6;
99 static const G4int NumberofKPoints=4;
100
101 G4double QQ1[NumberofZPoints][NumberofEPoints][NumberofKPoints];
102 G4double QQ2[NumberofZPoints][NumberofEPoints][NumberofKPoints];
103
104 G4int verbosityLevel;
105
106};
107
108
109
110#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=0)
void SetVerbosityLevel(G4int vl)
Set/Get Verbosity level.
G4double PolarAngle(const G4double initial_energy, const G4double final_energy, const G4int Z)
void Initialize()
The Initialize() method forces the cleaning of tables.