Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ParticleHPInelasticData.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// particle_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// P. Arce, June-2014 Conversion neutron_hp to particle_hp
38//
41#include "G4Neutron.hh"
42#include "G4ElementTable.hh"
43#include "G4ParticleHPData.hh"
45#include "G4NucleiProperties.hh"
46#include "G4Pow.hh"
47
50{
51 const char* dataDirVariable;
52 G4String particleName;
53 if( projectile == G4Neutron::Neutron() ) {
54 dataDirVariable = "G4NEUTRONHPDATA";
55 }else if( projectile == G4Proton::Proton() ) {
56 dataDirVariable = "G4PROTONHPDATA";
57 particleName = "Proton";
58 }else if( projectile == G4Deuteron::Deuteron() ) {
59 dataDirVariable = "G4DEUTERONHPDATA";
60 particleName = "Deuteron";
61 }else if( projectile == G4Triton::Triton() ) {
62 dataDirVariable = "G4TRITONHPDATA";
63 particleName = "Triton";
64 }else if( projectile == G4He3::He3() ) {
65 dataDirVariable = "G4HE3HPDATA";
66 particleName = "He3";
67 }else if( projectile == G4Alpha::Alpha() ) {
68 dataDirVariable = "G4ALPHAHPDATA";
69 particleName = "Alpha";
70 } else {
71 G4String message("G4ParticleHPInelasticData may only be called for neutron, proton, deuteron, triton, He3 or alpha, while it is called for " + projectile->GetParticleName());
72 throw G4HadronicException(__FILE__, __LINE__,message.c_str());
73 }
74 G4String dataName = projectile->GetParticleName()+"HPInelasticXS";
75 dataName.at(0) = (char)std::toupper(dataName.at(0)) ;
76 SetName( dataName );
77
78 if ( !G4FindDataDir(dataDirVariable) && !G4FindDataDir( "G4PARTICLEHPDATA" ) ){
79 G4String message("Please setenv G4PARTICLEHPDATA (recommended) or, at least setenv " +
80 G4String(dataDirVariable) + " to point to the " + projectile->GetParticleName() + " cross-section files.");
81 throw G4HadronicException(__FILE__, __LINE__,message.c_str());
82 }
83
84 G4String dirName;
85 if ( G4FindDataDir(dataDirVariable) ) {
86 dirName = G4FindDataDir(dataDirVariable);
87 } else {
88 G4String baseName = G4FindDataDir( "G4PARTICLEHPDATA" );
89 dirName = baseName + "/" + particleName;
90 }
91 #ifdef G4VERBOSE
93 G4cout << "@@@ G4ParticleHPInelasticData instantiated for particle " << projectile->GetParticleName() << " data directory variable is " << dataDirVariable << " pointing to " << dirName << G4endl;
94 }
95 #endif
96
97 SetMinKinEnergy( 0*CLHEP::MeV );
98 SetMaxKinEnergy( 20*CLHEP::MeV );
99
100 theCrossSections = 0;
101 theProjectile=projectile;
102
103 theHPData = nullptr;
104 instanceOfWorker = false;
106 theHPData = new G4ParticleHPData( theProjectile );
107 } else {
108 instanceOfWorker = true;
109 }
110 element_cache = nullptr;
111 material_cache = nullptr;
112 ke_cache = 0.0;
113 xs_cache = 0.0;
114}
115
117{
118 if ( theCrossSections != nullptr && instanceOfWorker != true ) {
119 theCrossSections->clearAndDestroy();
120 delete theCrossSections;
121 theCrossSections = nullptr;
122 }
123 if ( theHPData != nullptr && instanceOfWorker != true ) {
124 delete theHPData;
125 theHPData = nullptr;
126 }
127}
128
130 G4int /*Z*/ , G4int /*A*/ ,
131 const G4Element* /*elm*/ ,
132 const G4Material* /*mat*/ )
133{
134 G4double eKin = dp->GetKineticEnergy();
135 if ( eKin > GetMaxKinEnergy()
136 || eKin < GetMinKinEnergy()
137 || dp->GetDefinition() != theProjectile ) return false;
138
139 return true;
140}
141
143 G4int /*Z*/ , G4int /*A*/ ,
144 const G4Isotope* /*iso*/ ,
145 const G4Element* element ,
146 const G4Material* material )
147{
148 if ( dp->GetKineticEnergy() == ke_cache && element == element_cache
149 && material == material_cache ) return xs_cache;
150
151 ke_cache = dp->GetKineticEnergy();
152 element_cache = element;
153 material_cache = material;
154 G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
155 xs_cache = xs;
156 return xs;
157}
158
160{
162 theCrossSections = G4ParticleHPManager::GetInstance()->GetInelasticCrossSections( &projectile );
163 return;
164 } else {
165 if ( theHPData == nullptr ) theHPData = G4ParticleHPData::Instance( const_cast<G4ParticleDefinition*> ( &projectile ) );
166 }
167
168 std::size_t numberOfElements = G4Element::GetNumberOfElements();
169 if ( theCrossSections == nullptr )
170 theCrossSections = new G4PhysicsTable( numberOfElements );
171 else
172 theCrossSections->clearAndDestroy();
173
174 // make a PhysicsVector for each element
175
176 static G4ThreadLocal G4ElementTable *theElementTable = nullptr ;
177 if (!theElementTable) theElementTable= G4Element::GetElementTable();
178 for( std::size_t i=0; i<numberOfElements; ++i )
179 {
180 G4PhysicsVector* physVec = theHPData->MakePhysicsVector((*theElementTable)[i], this);
181 theCrossSections->push_back(physVec);
182 }
183 G4ParticleHPManager::GetInstance()->RegisterInelasticCrossSections( &projectile , theCrossSections );
184}
185
187{
188 if(&projectile!=theProjectile)
189 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use ParticleHP data for a wrong projectile!!!");
190
191#ifdef G4VERBOSE
192 if ( G4HadronicParameters::Instance()->GetVerboseLevel() == 0 ) return;
193
194 // Dump element based cross section
195 // range 10e-5 eV to 20 MeV
196 // 10 point per decade
197 // in barn
198
199 G4cout << G4endl;
200 G4cout << G4endl;
201 G4cout << "Inelastic Cross Section of Neutron HP"<< G4endl;
202 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
203 G4cout << G4endl;
204 G4cout << "Name of Element" << G4endl;
205 G4cout << "Energy[eV] XS[barn]" << G4endl;
206 G4cout << G4endl;
207
208 std::size_t numberOfElements = G4Element::GetNumberOfElements();
209 static G4ThreadLocal G4ElementTable *theElementTable = nullptr ;
210 if (!theElementTable) theElementTable= G4Element::GetElementTable();
211
212 for ( std::size_t i = 0 ; i < numberOfElements ; ++i )
213 {
214 G4cout << (*theElementTable)[i]->GetName() << G4endl;
215
216 G4int ie = 0;
217
218 for ( ie = 0 ; ie < 130 ; ie++ )
219 {
220 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *CLHEP::eV;
221 G4bool outOfRange = false;
222
223 if ( eKinetic < 20*CLHEP::MeV )
224 {
225 G4cout << eKinetic/CLHEP::eV << " "
226 << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/CLHEP::barn
227 << G4endl;
228 }
229 }
230 G4cout << G4endl;
231 }
232#endif
233}
234
236GetCrossSection(const G4DynamicParticle* projectile, const G4Element*anE, G4double aT)
237{
238 G4double result = 0;
239 G4bool outOfRange;
240 std::size_t index = anE->GetIndex();
241
242 // prepare neutron
243 G4double eKinetic = projectile->GetKineticEnergy();
244
245 if ( G4ParticleHPManager::GetInstance()->GetNeglectDoppler() )
246 {
247 //NEGLECT_DOPPLER
248 G4double factor = 1.0;
249 if ( eKinetic < aT * CLHEP::k_Boltzmann )
250 {
251 // below 0.1 eV neutrons
252 // Have to do some, but now just igonre.
253 // Will take care after performance check.
254 // factor = factor * targetV;
255 }
256 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
257 }
258
259 G4ReactionProduct theNeutron( projectile->GetDefinition() );
260 theNeutron.SetMomentum( projectile->GetMomentum() );
261 theNeutron.SetKineticEnergy( eKinetic );
262
263 // prepare thermal nucleus
264 G4Nucleus aNuc;
265 G4double eps = 0.0001;
266 G4double theA = anE->GetN();
267 G4double theZ = anE->GetZ();
268 G4double eleMass;
269 eleMass = G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theA+eps), static_cast<G4int>(theZ+eps) );
270
271 G4ReactionProduct boosted;
272 G4double aXsection;
273
274 // MC integration loop
275 G4int counter = 0;
276 G4int failCount = 0;
277 G4double buffer = 0; G4int size = G4int(std::max(10., aT/60*CLHEP::kelvin));
278 G4ThreeVector neutronVelocity = 1./theProjectile->GetPDGMass()*theNeutron.GetMomentum();
279 G4double neutronVMag = neutronVelocity.mag();
280
281#ifndef G4PHP_DOPPLER_LOOP_ONCE
282 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer) // Loop checking, 11.05.2015, T. Koi
283 {
284 if(counter) buffer = result/counter;
285 while (counter<size) // Loop checking, 11.05.2015, T. Koi
286 {
287 ++counter;
288#endif
289 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus( eleMass/G4Neutron::Neutron()->GetPDGMass(), aT );
290 boosted.Lorentz(theNeutron, aThermalNuc);
291 G4double theEkin = boosted.GetKineticEnergy();
292 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
293 if(aXsection <0)
294 {
295 if(failCount<1000)
296 {
297 ++failCount;
298#ifndef G4PHP_DOPPLER_LOOP_ONCE
299 --counter;
300 continue;
301#endif
302 }
303 else
304 {
305 aXsection = 0;
306 }
307 }
308 // velocity correction.
309 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
310 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
311 result += aXsection;
312 }
313#ifndef G4PHP_DOPPLER_LOOP_ONCE
314 size += size;
315 }
316 result /= counter;
317#endif
318
319 return result;
320}
321
323{
325}
326
328{
330}
331
333{
334 outFile << "Extension of High Precision cross section for inelastic reaction of proton, deuteron, triton, He3 and alpha below 20MeV\n";
335}
std::vector< G4Element * > G4ElementTable
const char * G4FindDataDir(const char *)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double mag() const
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:93
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
G4double GetZ() const
Definition: G4Element.hh:131
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
size_t GetIndex() const
Definition: G4Element.hh:182
G4double GetN() const
Definition: G4Element.hh:135
static G4HadronicParameters * Instance()
static G4He3 * He3()
Definition: G4He3.cc:93
G4double GetTemperature() const
Definition: G4Material.hh:177
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:236
const G4String & GetParticleName() const
static G4ParticleHPData * Instance(G4ParticleDefinition *projectile)
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4ParticleHPFissionData *theP)
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
virtual void CrossSectionDescription(std::ostream &) const
void DumpPhysicsTable(const G4ParticleDefinition &)
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)
G4ParticleHPInelasticData(G4ParticleDefinition *projectile=G4Neutron::Neutron())
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4PhysicsTable * GetInelasticCrossSections(const G4ParticleDefinition *)
void RegisterInelasticCrossSections(const G4ParticleDefinition *, G4PhysicsTable *)
static G4ParticleHPManager * GetInstance()
void push_back(G4PhysicsVector *)
void clearAndDestroy()
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:230
static G4Proton * Proton()
Definition: G4Proton.cc:92
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
static G4Triton * Triton()
Definition: G4Triton.cc:93
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
void SetName(const G4String &nam)
G4bool IsWorkerThread()
Definition: G4Threading.cc:123
G4bool IsMasterThread()
Definition: G4Threading.cc:124
#define G4ThreadLocal
Definition: tls.hh:77