Geant4-11
G4ParticleHPThermalScattering.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// P. Arce, June-2014 Conversion neutron_hp to particle_hp
28//
29#ifndef G4ParticleHPThermalScattering_h
30#define G4ParticleHPThermalScattering_h 1
31
32// Thermal Neutron Scattering
33// Koi, Tatsumi (SLAC/SCCS)
34//
35// Class Description
36// Final State Generators for a high precision (based on evaluated data
37// libraries) description of themal neutron scattering below 4 eV;
38// Based on Thermal neutron scattering files
39// from the evaluated nuclear data files ENDF/B-VI, Release2
40// To be used in your physics list in case you need this physics.
41// In this case you want to register an object of this class with
42// the corresponding process.
43// Class Description - End
44
45#include "globals.hh"
48
51
52struct E_isoAng
53{
56 std::vector < G4double > isoAngle;
58 energy=0.0;
59 n=0;
60 };
61};
62
64{
67 std::vector < G4double > prob;
68 std::vector < E_isoAng* > vE_isoAngle;
69 G4double sum_of_probXdEs; // should be close to 1
70 std::vector< G4double > secondary_energy_cdf;
71 std::vector< G4double > secondary_energy_pdf;
72 std::vector< G4double > secondary_energy_value;
75 energy=0.0;
76 n=0;
78 };
79};
80
82{
83 public:
84
86
88
89 G4HadFinalState * ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& aTargetNucleus);
90
91 virtual const std::pair<G4double, G4double> GetFatalEnergyCheckLevels() const;
92
93 //For user prepared thermal files
94 //Name of G4Element , Name of NDL file
96
98
99 virtual void ModelDescription(std::ostream& outFile) const;
100
101 private:
102
103 void clearCurrentFSData();
104
106
107 // Coherent Elastic
108 // ElementID temp BraggE cumulativeP
109 std::map < G4int , std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >* >* coherentFSs;
110 std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >* readACoherentFSDATA( G4String );
111
112 // Incoherent Elastic
113 // ElementID temp aFS for this temp (and this element)
114 std::map < G4int , std::map < G4double , std::vector < E_isoAng* >* >* >* incoherentFSs;
115 std::map < G4double , std::vector < E_isoAng* >* >* readAnIncoherentFSDATA( G4String );
116 E_isoAng* readAnE_isoAng ( std::istream* );
117
118 // Inelastic
119 // ElementID temp aFS for this temp (and this element)
120 std::map < G4int , std::map < G4double , std::vector < E_P_E_isoAng* >* >* >* inelasticFSs;
121 std::map < G4double , std::vector < E_P_E_isoAng* >* >* readAnInelasticFSDATA( G4String );
122 E_P_E_isoAng* readAnE_P_E_isoAng ( std::istream* );
123
125
127
129 G4double getMu ( G4double rndm1 , G4double rndm2 , E_isoAng* anEPM );
130
131 std::pair< G4double , G4double > find_LH ( G4double , std::vector<G4double>* );
132 G4double get_linear_interpolated ( G4double , std::pair < G4double , G4double > , std::pair < G4double , G4double > );
133
134 E_isoAng create_E_isoAng_from_energy( G4double , std::vector< E_isoAng* >* );
135
137
138 std::pair< G4double, G4double > sample_inelastic_E_mu( G4double pE , std::vector< E_P_E_isoAng* >* vNEP_EPM );
139 std::pair< G4double, G4int > sample_inelastic_E( G4double rndm1 , G4double rndm2 , E_P_E_isoAng* anE_P_E_isoAng );
140
141 std::pair< G4double , E_isoAng > create_sE_and_EPM_from_pE_and_vE_P_E_isoAng ( G4double , G4double , std::vector < E_P_E_isoAng* >* );
142
143 std::map < std::pair < const G4Material* , const G4Element* > , G4int > dic;
144 void buildPhysicsTable();
145 G4int getTS_ID( const G4Material* , const G4Element* );
146
147 //size_t sizeOfMaterialTable;
148
150
151 //In order to judge whether the rebuilding of physics table is a necessity or not
152 size_t nMaterial;
153 size_t nElement;
154
155};
156
157#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
std::map< G4double, std::vector< std::pair< G4double, G4double > * > * > * readACoherentFSDATA(G4String)
std::pair< G4double, G4double > find_LH(G4double, std::vector< G4double > *)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
G4int getTS_ID(const G4Material *, const G4Element *)
E_isoAng create_E_isoAng_from_energy(G4double, std::vector< E_isoAng * > *)
void AddUserThermalScatteringFile(G4String, G4String)
std::map< G4double, std::vector< E_isoAng * > * > * readAnIncoherentFSDATA(G4String)
std::pair< G4double, E_isoAng > create_sE_and_EPM_from_pE_and_vE_P_E_isoAng(G4double, G4double, std::vector< E_P_E_isoAng * > *)
E_P_E_isoAng * readAnE_P_E_isoAng(std::istream *)
G4ParticleHPThermalScatteringNames names
std::map< G4int, std::map< G4double, std::vector< std::pair< G4double, G4double > * > * > * > * coherentFSs
virtual void ModelDescription(std::ostream &outFile) const
G4double get_linear_interpolated(G4double, std::pair< G4double, G4double >, std::pair< G4double, G4double >)
std::pair< G4double, G4double > sample_inelastic_E_mu(G4double pE, std::vector< E_P_E_isoAng * > *vNEP_EPM)
std::pair< G4double, G4int > sample_inelastic_E(G4double rndm1, G4double rndm2, E_P_E_isoAng *anE_P_E_isoAng)
std::map< G4int, std::map< G4double, std::vector< E_P_E_isoAng * > * > * > * inelasticFSs
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
std::map< G4int, std::map< G4double, std::vector< E_isoAng * > * > * > * incoherentFSs
std::map< G4double, std::vector< E_P_E_isoAng * > * > * readAnInelasticFSDATA(G4String)
G4double get_secondary_energy_from_E_P_E_isoAng(G4double random, E_P_E_isoAng *anE_P_E_isoAng)
void BuildPhysicsTable(const G4ParticleDefinition &)
std::map< std::pair< const G4Material *, const G4Element * >, G4int > dic
G4ParticleHPThermalScatteringData * theXSection
std::vector< G4double > prob
std::vector< E_isoAng * > vE_isoAngle
std::vector< G4double > secondary_energy_pdf
std::vector< G4double > secondary_energy_cdf
std::vector< G4double > secondary_energy_value
std::vector< G4double > isoAngle