Geant4_10
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"
45 #include "G4NeutronHPElastic.hh"
46 #include "G4HadronicInteraction.hh"
47 
48 struct E_isoAng
49 {
52  std::vector < G4double > isoAngle;
53 };
54 
55 struct E_P_E_isoAng
56 {
59  std::vector < G4double > prob;
60  std::vector < E_isoAng* > vE_isoAngle;
61  G4double sum_of_probXdEs; // should be close to 1
62 };
63 
65 {
66  public:
67 
69 
71 
72  G4HadFinalState * ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& aTargetNucleus);
73 
74  virtual const std::pair<G4double, G4double> GetFatalEnergyCheckLevels() const;
75 
76  //For user prepared thermal files
77  //Name of G4Element , Name of NDL file
79 
80  private:
81 
83 
84  // Coherent Elastic
85  // ElementID temp BraggE cumulativeP
86  std::map < G4int , std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >* > coherentFSs;
87  std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >* readACoherentFSDATA( G4String );
88 
89  // Incoherent Elastic
90  // ElementID temp aFS for this temp (and this element)
91  std::map < G4int , std::map < G4double , std::vector < E_isoAng* >* >* > incoherentFSs;
92  std::map < G4double , std::vector < E_isoAng* >* >* readAnIncoherentFSDATA( G4String );
93  E_isoAng* readAnE_isoAng ( std::istream* );
94 
95  // Inelastic
96  // ElementID temp aFS for this temp (and this element)
97  std::map < G4int , std::map < G4double , std::vector < E_P_E_isoAng* >* >* > inelasticFSs;
98  std::map < G4double , std::vector < E_P_E_isoAng* >* >* readAnInelasticFSDATA( G4String );
99  E_P_E_isoAng* readAnE_P_E_isoAng ( std::istream* );
100 
102 
103  G4NeutronHPElastic* theHPElastic;
104 
105  G4double getMu ( E_isoAng* );
106 
107  std::pair< G4double , G4double > find_LH ( G4double , std::vector<G4double>* );
108  G4double get_linear_interpolated ( G4double , std::pair < G4double , G4double > , std::pair < G4double , G4double > );
109 
110  E_isoAng create_E_isoAng_from_energy( G4double , std::vector< E_isoAng* >* );
111 
112  G4double get_secondary_energy_from_E_P_E_isoAng ( G4double , E_P_E_isoAng* );
113 
114  std::pair< G4double , E_isoAng > create_sE_and_EPM_from_pE_and_vE_P_E_isoAng ( G4double , G4double , std::vector < E_P_E_isoAng* >* );
115 
116  std::map < std::pair < const G4Material* , const G4Element* > , G4int > dic;
117  void buildPhysicsTable();
118  G4int getTS_ID( const G4Material* , const G4Element* );
119 
120  size_t sizeOfMaterialTable;
121 
122 };
123 
124 #endif
std::vector< G4double > isoAngle
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
int G4int
Definition: G4Types.hh:78
std::vector< E_isoAng * > vE_isoAngle
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
double G4double
Definition: G4Types.hh:76
void AddUserThermalScatteringFile(G4String, G4String)
std::vector< G4double > prob