source: trunk/source/processes/hadronic/models/neutron_hp/src/G4NeutronHPInelasticData.cc@ 1201

Last change on this file since 1201 was 962, checked in by garnier, 17 years ago

update processes

File size: 7.6 KB
Line 
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// neutron_hp -- source file
27// J.P. Wellisch, Nov-1996
28// A prototype of the low energy neutron transport model.
29//
30// 070523 add neglecting doppler broadening on the fly. T. Koi
31// 070613 fix memory leaking by T. Koi
32// 071002 enable cross section dump by T. Koi
33// 080428 change checking point of "neglecting doppler broadening" flag
34// from GetCrossSection to BuildPhysicsTable by T. Koi
35// 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
36//
37#include "G4NeutronHPInelasticData.hh"
38#include "G4Neutron.hh"
39#include "G4ElementTable.hh"
40#include "G4NeutronHPData.hh"
41
42G4bool G4NeutronHPInelasticData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
43{
44 G4bool result = true;
45 G4double eKin = aP->GetKineticEnergy();
46 if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
47 return result;
48}
49
50G4NeutronHPInelasticData::G4NeutronHPInelasticData()
51{
52// TKDB
53 onFlightDB = true;
54 theCrossSections = 0;
55 BuildPhysicsTable(*G4Neutron::Neutron());
56}
57
58G4NeutronHPInelasticData::~G4NeutronHPInelasticData()
59{
60// TKDB
61 if ( theCrossSections != 0 )
62 { theCrossSections->clearAndDestroy(); }
63 delete theCrossSections;
64}
65
66void G4NeutronHPInelasticData::BuildPhysicsTable(const G4ParticleDefinition& aP)
67{
68 if(&aP!=G4Neutron::Neutron())
69 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
70
71//080428
72 if ( getenv( "G4NEUTRONHP_NEGLECT_DOPPLER" ) ) onFlightDB = false;
73
74 size_t numberOfElements = G4Element::GetNumberOfElements();
75// theCrossSections = new G4PhysicsTable( numberOfElements );
76// TKDB
77 if ( theCrossSections == 0 )
78 { theCrossSections = new G4PhysicsTable( numberOfElements ); }
79
80 // make a PhysicsVector for each element
81
82 static const G4ElementTable *theElementTable = G4Element::GetElementTable();
83 for( size_t i=0; i<numberOfElements; ++i )
84 {
85 G4PhysicsVector* physVec = G4NeutronHPData::
86 Instance()->MakePhysicsVector((*theElementTable)[i], this);
87 theCrossSections->push_back(physVec);
88 }
89}
90
91void G4NeutronHPInelasticData::DumpPhysicsTable(const G4ParticleDefinition& aP)
92{
93 if(&aP!=G4Neutron::Neutron())
94 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
95
96//
97// Dump element based cross section
98// range 10e-5 eV to 20 MeV
99// 10 point per decade
100// in barn
101//
102
103 G4cout << G4endl;
104 G4cout << G4endl;
105 G4cout << "Inelastic Cross Section of Neutron HP"<< G4endl;
106 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
107 G4cout << G4endl;
108 G4cout << "Name of Element" << G4endl;
109 G4cout << "Energy[eV] XS[barn]" << G4endl;
110 G4cout << G4endl;
111
112 size_t numberOfElements = G4Element::GetNumberOfElements();
113 static const G4ElementTable *theElementTable = G4Element::GetElementTable();
114
115 for ( size_t i = 0 ; i < numberOfElements ; ++i )
116 {
117
118 G4cout << (*theElementTable)[i]->GetName() << G4endl;
119
120 G4int ie = 0;
121
122 for ( ie = 0 ; ie < 130 ; ie++ )
123 {
124 G4double eKinetic = 1.0e-5 * std::pow ( 10.0 , ie/10.0 ) *eV;
125 G4bool outOfRange = false;
126
127 if ( eKinetic < 20*MeV )
128 {
129 G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
130 }
131
132 }
133
134 G4cout << G4endl;
135 }
136
137 //G4cout << "G4NeutronHPInelasticData::DumpPhysicsTable still to be implemented"<<G4endl;
138}
139
140#include "G4NucleiProperties.hh"
141
142G4double G4NeutronHPInelasticData::
143GetCrossSection(const G4DynamicParticle* aP, const G4Element*anE, G4double aT)
144{
145 G4double result = 0;
146 G4bool outOfRange;
147 G4int index = anE->GetIndex();
148
149 // prepare neutron
150 G4double eKinetic = aP->GetKineticEnergy();
151
152 // T. K.
153//if ( getenv( "G4NEUTRONHP_NEGLECT_DOPPLER" ) )
154//080428
155 if ( !onFlightDB )
156 {
157 G4double factor = 1.0;
158 if ( eKinetic < aT * k_Boltzmann )
159 {
160 // below 0.1 eV neutrons
161 // Have to do some, but now just igonre.
162 // Will take care after performance check.
163 // factor = factor * targetV;
164 }
165 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
166 }
167
168 G4ReactionProduct theNeutron( aP->GetDefinition() );
169 theNeutron.SetMomentum( aP->GetMomentum() );
170 theNeutron.SetKineticEnergy( eKinetic );
171
172 // prepare thermal nucleus
173 G4Nucleus aNuc;
174 G4double eps = 0.0001;
175 G4double theA = anE->GetN();
176 G4double theZ = anE->GetZ();
177 G4double eleMass;
178 eleMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theA+eps), static_cast<G4int>(theZ+eps))
179 ) / G4Neutron::Neutron()->GetPDGMass();
180
181 G4ReactionProduct boosted;
182 G4double aXsection;
183
184 // MC integration loop
185 G4int counter = 0;
186 G4int failCount = 0;
187 G4double buffer = 0;
188 G4int size = G4int(std::max(10., aT/60*kelvin));
189 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
190 G4double neutronVMag = neutronVelocity.mag();
191
192 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer)
193 {
194 if(counter) buffer = result/counter;
195 while (counter<size)
196 {
197 counter ++;
198 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
199 boosted.Lorentz(theNeutron, aThermalNuc);
200 G4double theEkin = boosted.GetKineticEnergy();
201 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
202 if(aXsection <0)
203 {
204 if(failCount<1000)
205 {
206 failCount++;
207 counter--;
208 continue;
209 }
210 else
211 {
212 aXsection = 0;
213 }
214 }
215 // velocity correction.
216 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
217 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
218 result += aXsection;
219 }
220 size += size;
221 }
222 result /= counter;
223/*
224 // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
225 G4cout << " result " << result << " "
226 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
227 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
228*/
229 return result;
230}
Note: See TracBrowser for help on using the repository browser.