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
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"
47
48struct E_isoAng
49{
52 std::vector < G4double > isoAngle;
53};
54
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 //G4int GetNiso() {return theElastic[0].GetNiso();}
76
77 private:
78
80
81 // Coherent Elastic
82 // ElementID temp BraggE cumulativeP
83 std::map < G4int , std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >* > coherentFSs;
84 std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >* readACoherentFSDATA( G4String );
85
86 // Incoherent Elastic
87 // ElementID temp aFS for this temp (and this element)
88 std::map < G4int , std::map < G4double , std::vector < E_isoAng* >* >* > incoherentFSs;
89 std::map < G4double , std::vector < E_isoAng* >* >* readAnIncoherentFSDATA( G4String );
90 E_isoAng* readAnE_isoAng ( std::ifstream* );
91
92 // Inelastic
93 // ElementID temp aFS for this temp (and this element)
94 std::map < G4int , std::map < G4double , std::vector < E_P_E_isoAng* >* >* > inelasticFSs;
95 std::map < G4double , std::vector < E_P_E_isoAng* >* >* readAnInelasticFSDATA( G4String );
96 E_P_E_isoAng* readAnE_P_E_isoAng ( std::ifstream* );
97
98
100 G4double * xSec;
101 //G4String dirName;
102 G4int numEle;
103
104
105
106 G4NeutronHPElastic* theHPElastic;
107
108 G4double getMu ( E_isoAng* );
109
110 std::pair< G4double , G4double > find_LH ( G4double , std::vector<G4double>* );
111 G4double get_linear_interpolated ( G4double , std::pair < G4double , G4double > , std::pair < G4double , G4double > );
112
113 E_isoAng create_E_isoAng_from_energy( G4double , std::vector< E_isoAng* >* );
114
115 G4double get_secondary_energy_from_E_P_E_isoAng ( G4double , E_P_E_isoAng* );
116
117 std::pair< G4double , E_isoAng > create_sE_and_EPM_from_pE_and_vE_P_E_isoAng ( G4double , G4double , std::vector < E_P_E_isoAng* >* );
118
119 std::map < std::pair < const G4Material* , const G4Element* > , G4int > dic;
120 void buildPhysicsTable();
121 G4int getTS_ID( const G4Material* , const G4Element* );
122
123};
124
125#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
std::vector< G4double > prob
std::vector< E_isoAng * > vE_isoAngle
std::vector< G4double > isoAngle