Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
47#include "globals.hh"
48
51
53{
56 std::vector<G4double> isoAngle;
58 {
59 energy = 0.0;
60 n = 0;
61 };
62};
63
65{
68 std::vector<G4double> prob;
69 std::vector<E_isoAng*> vE_isoAngle;
70 G4double sum_of_probXdEs; // should be close to 1
71 std::vector<G4double> secondary_energy_cdf;
72 std::vector<G4double> secondary_energy_pdf;
73 std::vector<G4double> secondary_energy_value;
76 {
77 energy = 0.0;
78 n = 0;
79 sum_of_probXdEs = 0.0;
81 };
82};
83
85{
86 public:
88
90
92 G4Nucleus& aTargetNucleus) override;
93
94 const std::pair<G4double, G4double> GetFatalEnergyCheckLevels() const override;
95
96 // For user prepared thermal files
97 // Name of G4Element , Name of NDL file
99
100 void BuildPhysicsTable(const G4ParticleDefinition&) override;
101
102 void ModelDescription(std::ostream& outFile) const override;
103
104 private:
105 void clearCurrentFSData();
106
108
109 // Coherent Elastic
110 // ElementID temp BraggE cumulativeP
111 std::map<G4int, std::map<G4double, std::vector<std::pair<G4double, G4double>*>*>*>* coherentFSs{
112 nullptr};
113 std::map<G4double, std::vector<std::pair<G4double, G4double>*>*>* readACoherentFSDATA(G4String);
114
115 // Incoherent Elastic
116 // ElementID temp aFS for this temp (and this element)
117 std::map<G4int, std::map<G4double, std::vector<E_isoAng*>*>*>* incoherentFSs{nullptr};
118 std::map<G4double, std::vector<E_isoAng*>*>* readAnIncoherentFSDATA(G4String);
119 E_isoAng* readAnE_isoAng(std::istream*);
120
121 // Inelastic
122 // ElementID temp aFS for this temp (and this element)
123 std::map<G4int, std::map<G4double, std::vector<E_P_E_isoAng*>*>*>* inelasticFSs{nullptr};
124 std::map<G4double, std::vector<E_P_E_isoAng*>*>* readAnInelasticFSDATA(G4String);
125 E_P_E_isoAng* readAnE_P_E_isoAng(std::istream*);
126
128
129 G4ParticleHPElastic* theHPElastic;
130
131 G4double getMu(E_isoAng*);
132 G4double getMu(G4double rndm1, G4double rndm2, E_isoAng* anEPM);
133
134 std::pair<G4double, G4double> find_LH(G4double, std::vector<G4double>*);
135 G4double get_linear_interpolated(G4double, std::pair<G4double, G4double>,
136 std::pair<G4double, G4double>);
137
138 E_isoAng create_E_isoAng_from_energy(G4double, std::vector<E_isoAng*>*);
139
140 G4double get_secondary_energy_from_E_P_E_isoAng(G4double random, E_P_E_isoAng* anE_P_E_isoAng);
141
142 std::pair<G4double, G4double> sample_inelastic_E_mu(G4double pE,
143 std::vector<E_P_E_isoAng*>* vNEP_EPM);
144 std::pair<G4double, G4int> sample_inelastic_E(G4double rndm1, G4double rndm2,
145 E_P_E_isoAng* anE_P_E_isoAng);
146
147 std::pair<G4double, E_isoAng>
148 create_sE_and_EPM_from_pE_and_vE_P_E_isoAng(G4double, G4double, std::vector<E_P_E_isoAng*>*);
149
150 std::map<std::pair<const G4Material*, const G4Element*>, G4int> dic;
151 void buildPhysicsTable();
152 G4int getTS_ID(const G4Material*, const G4Element*);
153
154 // size_t sizeOfMaterialTable;
155
156 G4bool check_E_isoAng(E_isoAng*);
157
158 // In order to judge whether the rebuilding of physics table is a necessity or not
159 size_t nMaterial;
160 size_t nElement;
161};
162
163#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void ModelDescription(std::ostream &outFile) const override
const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const override
void AddUserThermalScatteringFile(G4String, G4String)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
std::vector< E_isoAng * > vE_isoAngle
std::vector< G4double > secondary_energy_pdf
std::vector< G4double > secondary_energy_value
std::vector< G4double > secondary_energy_cdf
std::vector< G4double > prob
std::vector< G4double > isoAngle