Geant4  10.01.p03
G4NeutronHPThermalScattering.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 #ifndef G4NeutronHPThermalScattering_h
27 #define G4NeutronHPThermalScattering_h 1
28 
29 // Thermal Neutron Scattering
30 // Koi, Tatsumi (SLAC/SCCS)
31 //
32 // Class Description
33 // Final State Generators for a high precision (based on evaluated data
34 // libraries) description of themal neutron scattering below 4 eV;
35 // Based on Thermal neutron scattering files
36 // from the evaluated nuclear data files ENDF/B-VI, Release2
37 // To be used in your physics list in case you need this physics.
38 // In this case you want to register an object of this class with
39 // the corresponding process.
40 // Class Description - End
41 
42 #include "globals.hh"
44 #include "G4HadronicInteraction.hh"
45 
47 class G4NeutronHPElastic;
48 
49 struct E_isoAng
50 {
53  std::vector < G4double > isoAngle;
55  energy=0.0;
56  n=0;
57  };
58 };
59 
60 struct E_P_E_isoAng
61 {
64  std::vector < G4double > prob;
65  std::vector < E_isoAng* > vE_isoAngle;
66  G4double sum_of_probXdEs; // should be close to 1
68  energy=0.0;
69  n=0;
70  sum_of_probXdEs=0.0;
71  };
72 };
73 
75 {
76  public:
77 
79 
81 
82  G4HadFinalState * ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& aTargetNucleus);
83 
84  virtual const std::pair<G4double, G4double> GetFatalEnergyCheckLevels() const;
85 
86  //For user prepared thermal files
87  //Name of G4Element , Name of NDL file
89 
91 
92  private:
93 
94  void clearCurrentFSData();
95 
97 
98  // Coherent Elastic
99  // ElementID temp BraggE cumulativeP
100  std::map < G4int , std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >* >* coherentFSs;
101  std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >* readACoherentFSDATA( G4String );
102 
103  // Incoherent Elastic
104  // ElementID temp aFS for this temp (and this element)
105  std::map < G4int , std::map < G4double , std::vector < E_isoAng* >* >* >* incoherentFSs;
106  std::map < G4double , std::vector < E_isoAng* >* >* readAnIncoherentFSDATA( G4String );
107  E_isoAng* readAnE_isoAng ( std::istream* );
108 
109  // Inelastic
110  // ElementID temp aFS for this temp (and this element)
111  std::map < G4int , std::map < G4double , std::vector < E_P_E_isoAng* >* >* >* inelasticFSs;
112  std::map < G4double , std::vector < E_P_E_isoAng* >* >* readAnInelasticFSDATA( G4String );
113  E_P_E_isoAng* readAnE_P_E_isoAng ( std::istream* );
114 
116 
118 
119  G4double getMu ( E_isoAng* );
120 
121  std::pair< G4double , G4double > find_LH ( G4double , std::vector<G4double>* );
122  G4double get_linear_interpolated ( G4double , std::pair < G4double , G4double > , std::pair < G4double , G4double > );
123 
124  E_isoAng create_E_isoAng_from_energy( G4double , std::vector< E_isoAng* >* );
125 
127 
128  std::pair< G4double , E_isoAng > create_sE_and_EPM_from_pE_and_vE_P_E_isoAng ( G4double , G4double , std::vector < E_P_E_isoAng* >* );
129 
130  std::map < std::pair < const G4Material* , const G4Element* > , G4int > dic;
131  void buildPhysicsTable();
132  G4int getTS_ID( const G4Material* , const G4Element* );
133 
134  //size_t sizeOfMaterialTable;
135 
137 
138  //In order to judge whether rebuild of physics table is a necessity
139  size_t nMaterial;
140  size_t nElement;
141 };
142 
143 #endif
std::map< G4int, std::map< G4double, std::vector< E_P_E_isoAng * > * > * > * inelasticFSs
G4int getTS_ID(const G4Material *, const G4Element *)
E_isoAng * readAnE_isoAng(std::istream *)
std::map< G4int, std::map< G4double, std::vector< E_isoAng * > * > * > * incoherentFSs
std::map< G4double, std::vector< E_isoAng * > * > * readAnIncoherentFSDATA(G4String)
std::vector< G4double > isoAngle
E_isoAng create_E_isoAng_from_energy(G4double, std::vector< E_isoAng * > *)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
std::pair< G4double, E_isoAng > create_sE_and_EPM_from_pE_and_vE_P_E_isoAng(G4double, G4double, std::vector< E_P_E_isoAng * > *)
int G4int
Definition: G4Types.hh:78
std::pair< G4double, G4double > find_LH(G4double, std::vector< G4double > *)
std::map< std::pair< const G4Material *, const G4Element * >, G4int > dic
void BuildPhysicsTable(const G4ParticleDefinition &)
std::map< G4double, std::vector< E_P_E_isoAng * > * > * readAnInelasticFSDATA(G4String)
G4double get_linear_interpolated(G4double, std::pair< G4double, G4double >, std::pair< G4double, G4double >)
bool G4bool
Definition: G4Types.hh:79
E_P_E_isoAng * readAnE_P_E_isoAng(std::istream *)
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
std::map< G4double, std::vector< std::pair< G4double, G4double > * > * > * readACoherentFSDATA(G4String)
G4NeutronHPThermalScatteringNames names
double G4double
Definition: G4Types.hh:76
void AddUserThermalScatteringFile(G4String, G4String)
G4NeutronHPThermalScatteringData * theXSection
std::vector< E_isoAng * > vE_isoAngle
std::map< G4int, std::map< G4double, std::vector< std::pair< G4double, G4double > * > * > * > * coherentFSs
G4double get_secondary_energy_from_E_P_E_isoAng(G4double, E_P_E_isoAng *)
std::vector< G4double > prob