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
G4NeutronHPorLCaptureData.cc
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// 05-11-21 NeutronHP or Low Energy Parameterization Models
28// Implemented by T. Koi (SLAC/SCCS)
29// If NeutronHP data do not available for an element, then Low Energy
30// Parameterization models handle the interactions of the element.
31// 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
32//
33
35#include "G4SystemOfUnits.hh"
36#include "G4Neutron.hh"
37#include "G4ElementTable.hh"
38#include "G4NeutronHPData.hh"
39
40#include "G4PhysicsVector.hh"
41
42
44{
45 SetMinKinEnergy( 0*MeV );
46 SetMaxKinEnergy( 20*MeV );
47
48 ke_cache = 0.0;
49 xs_cache = 0.0;
50 element_cache = NULL;
51 material_cache = NULL;
52// BuildPhysicsTable(*G4Neutron::Neutron());
53}
54
56{
57// delete theCrossSections;
58}
59
61 G4int /*Z*/ , G4int /*A*/ ,
62 const G4Element* element ,
63 const G4Material* /*mat*/ )
64{
65 G4double eKin = dp->GetKineticEnergy();
66 if ( eKin > GetMaxKinEnergy()
67 || eKin < GetMinKinEnergy()
68 || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
69 if ( unavailable_elements->find( element->GetName() ) != unavailable_elements->end() ) return false;
70
71 return true;
72}
73
75 G4int /*Z*/ , G4int /*A*/ ,
76 const G4Isotope* /*iso*/ ,
77 const G4Element* element ,
78 const G4Material* material )
79{
80 if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
81
82 ke_cache = dp->GetKineticEnergy();
83 element_cache = element;
84 material_cache = material;
85 G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
86 xs_cache = xs;
87 return xs;
88 //return GetCrossSection( dp , element , material->GetTemperature() );
89}
90
91
93:G4VCrossSectionDataSet("NeutronHPorLCaptureXS")
94{
95 theCaptureChannel = pChannel;
96 unavailable_elements = pSet;
97
98 SetMinKinEnergy( 0*MeV );
99 SetMaxKinEnergy( 20*MeV );
100
101 ke_cache = 0.0;
102 xs_cache = 0.0;
103 element_cache = NULL;
104 material_cache = NULL;
105}
106
107/*
108G4bool G4NeutronHPorLCaptureData::IsApplicable(const G4DynamicParticle*aP, const G4Element* anElement)
109{
110 G4bool result = true;
111 G4double eKin = aP->GetKineticEnergy();
112 if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
113 if ( unavailable_elements->find( anElement->GetName() ) != unavailable_elements->end() ) result = false;
114 return result;
115}
116*/
117
119{
120 if( &aP!=G4Neutron::Neutron() )
121 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
122}
123
124
125
127{
128 if(&aP!=G4Neutron::Neutron())
129 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
130// G4cout << "G4NeutronHPorLCaptureData::DumpPhysicsTable still to be implemented"<<G4endl;
131}
132
133
134
135#include "G4Nucleus.hh"
136#include "G4NucleiProperties.hh"
137#include "G4Neutron.hh"
138#include "G4Electron.hh"
139
141GetCrossSection(const G4DynamicParticle* aP, const G4Element*anE, G4double aT)
142{
143
144 // G4cout << "Choice G4NeutronHPorLCaptureData for element " << anE->GetName() << G4endl;
145 G4double result = 0;
146// G4bool outOfRange;
147 G4int index = anE->GetIndex();
148
149 // prepare neutron
150 G4double eKinetic = aP->GetKineticEnergy();
151 G4ReactionProduct theNeutron( aP->GetDefinition() );
152 theNeutron.SetMomentum( aP->GetMomentum() );
153 theNeutron.SetKineticEnergy( eKinetic );
154
155 // prepare thermal nucleus
156 G4Nucleus aNuc;
157 G4double eps = 0.0001;
158 G4double theA = anE->GetN();
159 G4double theZ = anE->GetZ();
160 G4double eleMass;
161 eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps), static_cast<G4int>(theZ+eps))
163
164 G4ReactionProduct boosted;
165 G4double aXsection;
166
167 // MC integration loop
168 G4int counter = 0;
169 G4double buffer = 0;
170 G4int size = G4int(std::max(10., aT/60*kelvin));
171 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
172 G4double neutronVMag = neutronVelocity.mag();
173 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer)
174 {
175 if(counter) buffer = result/counter;
176 while (counter<size)
177 {
178 counter ++;
179 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
180 boosted.Lorentz(theNeutron, aThermalNuc);
181 G4double theEkin = boosted.GetKineticEnergy();
182 //aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
183 aXsection = theCaptureChannel[index].GetXsec( theEkin );
184 // velocity correction.
185 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
186 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
187 result += aXsection;
188 }
189 size += size;
190 }
191 result /= counter;
192 return result;
193}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
double mag() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetZ() const
Definition: G4Element.hh:131
size_t GetIndex() const
Definition: G4Element.hh:182
const G4String & GetName() const
Definition: G4Element.hh:127
G4double GetN() const
Definition: G4Element.hh:134
G4double GetTemperature() const
Definition: G4Material.hh:181
G4double GetXsec(G4double energy)
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
void DumpPhysicsTable(const G4ParticleDefinition &)
void BuildPhysicsTable(const G4ParticleDefinition &)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:130
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
void SetKineticEnergy(const G4double en)
G4double GetMass() const
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
#define buffer
Definition: xmlparse.cc:611