Geant4-11
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 "G4Pow.hh"
46
49{
50 const char* dataDirVariable;
51 G4String particleName;
52 if( projectile == G4Neutron::Neutron() ) {
53 dataDirVariable = "G4NEUTRONHPDATA";
54 }else if( projectile == G4Proton::Proton() ) {
55 dataDirVariable = "G4PROTONHPDATA";
56 particleName = "Proton";
57 }else if( projectile == G4Deuteron::Deuteron() ) {
58 dataDirVariable = "G4DEUTERONHPDATA";
59 particleName = "Deuteron";
60 }else if( projectile == G4Triton::Triton() ) {
61 dataDirVariable = "G4TRITONHPDATA";
62 particleName = "Triton";
63 }else if( projectile == G4He3::He3() ) {
64 dataDirVariable = "G4HE3HPDATA";
65 particleName = "He3";
66 }else if( projectile == G4Alpha::Alpha() ) {
67 dataDirVariable = "G4ALPHAHPDATA";
68 particleName = "Alpha";
69 } else {
70 G4String message("G4ParticleHPInelasticData may only be called for neutron, proton, deuteron, triton, He3 or alpha, while it is called for " + projectile->GetParticleName());
71 throw G4HadronicException(__FILE__, __LINE__,message.c_str());
72 }
73 // G4cout << this << " G4ParticleHPInelasticData::G4ParticleHPInelasticData " << projectile->GetParticleName() << " DATADIR " << dataDirVariable << G4endl;//GDEB
74 G4String dataName = projectile->GetParticleName()+"HPInelasticXS";
75 dataName.at(0) = toupper(dataName.at(0)) ;
76 SetName( dataName );
77
78 if ( !std::getenv(dataDirVariable) && !std::getenv( "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 ( std::getenv(dataDirVariable) ) {
86 dirName = std::getenv(dataDirVariable);
87 } else {
88 G4String baseName = std::getenv( "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
99
101 theProjectile=projectile;
102
103 theHPData = NULL;
104 instanceOfWorker = false;
107 } else {
108 instanceOfWorker = true;
109 }
110 element_cache = NULL;
111 material_cache = NULL;
112 ke_cache = 0.0;
113 xs_cache = 0.0;
114}
115
117{
118 if ( theCrossSections != NULL && instanceOfWorker != true ) {
120 delete theCrossSections;
121 theCrossSections = NULL;
122 }
123 if ( theHPData != NULL && instanceOfWorker != true ) {
124 delete theHPData;
125 theHPData = NULL;
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 && material == material_cache ) return xs_cache;
149
151 element_cache = element;
153 G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
154 xs_cache = xs;
155 return xs;
156}
157
158/*
159G4bool G4ParticleHPInelasticData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
160{
161 G4bool result = true;
162 G4double eKin = aP->GetKineticEnergy();
163 if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
164 return result;
165}
166*/
167
168//void G4ParticleHPInelasticData::BuildPhysicsTableHP(G4ParticleDefinition* projectile,const char* /* dataDirVariable */)
170{
171 // if(&projectile!=G4Neutron::Neutron())
172 // throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
173
176 return;
177 } else {
178 if ( theHPData == NULL ) theHPData = G4ParticleHPData::Instance( const_cast<G4ParticleDefinition*> ( &projectile ) );
179 }
180
181 size_t numberOfElements = G4Element::GetNumberOfElements();
182// theCrossSections = new G4PhysicsTable( numberOfElements );
183// TKDB
184 //if ( theCrossSections == 0 )
185 //{ theCrossSections = new G4PhysicsTable( numberOfElements ); }
186 if ( theCrossSections == NULL )
187 theCrossSections = new G4PhysicsTable( numberOfElements );
188 else
190
191 // make a PhysicsVector for each element
192
193 //G4ParticleHPData* hpData = new G4ParticleHPData(projectile); //NEW
194 static G4ThreadLocal G4ElementTable *theElementTable = 0 ;
195 if (!theElementTable) theElementTable= G4Element::GetElementTable();
196 for( size_t i=0; i<numberOfElements; ++i )
197 {
198 //NEW G4PhysicsVector* physVec = G4ParticleHPData::
199 //NEW Instance(projectile, dataDirVariable)->MakePhysicsVector((*theElementTable)[i], this);
200 //G4PhysicsVector* physVec = hpData->MakePhysicsVector((*theElementTable)[i], this);
201 G4PhysicsVector* physVec = theHPData->MakePhysicsVector((*theElementTable)[i], this);
202 theCrossSections->push_back(physVec);
203 }
204
206}
207
209{
210 if(&projectile!=theProjectile)
211 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use ParticleHP data for a wrong projectile!!!");
212
213 #ifdef G4VERBOSE
214 if ( G4HadronicParameters::Instance()->GetVerboseLevel() == 0 ) return;
215//
216// Dump element based cross section
217// range 10e-5 eV to 20 MeV
218// 10 point per decade
219// in barn
220//
221
222 G4cout << G4endl;
223 G4cout << G4endl;
224 G4cout << "Inelastic Cross Section of Neutron HP"<< G4endl;
225 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
226 G4cout << G4endl;
227 G4cout << "Name of Element" << G4endl;
228 G4cout << "Energy[eV] XS[barn]" << G4endl;
229 G4cout << G4endl;
230
231 size_t numberOfElements = G4Element::GetNumberOfElements();
232 static G4ThreadLocal G4ElementTable *theElementTable = 0 ;
233 if (!theElementTable) theElementTable= G4Element::GetElementTable();
234
235 for ( size_t i = 0 ; i < numberOfElements ; ++i )
236 {
237
238 G4cout << (*theElementTable)[i]->GetName() << G4endl;
239
240 G4int ie = 0;
241
242 for ( ie = 0 ; ie < 130 ; ie++ )
243 {
244 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *CLHEP::eV;
245 G4bool outOfRange = false;
246
247 if ( eKinetic < 20*CLHEP::MeV )
248 {
249 G4cout << eKinetic/CLHEP::eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/CLHEP::barn << G4endl;
250 }
251
252 }
253
254 G4cout << G4endl;
255 }
256
257 //G4cout << "G4ParticleHPInelasticData::DumpPhysicsTable still to be implemented"<<G4endl;
258 #endif
259}
260
261#include "G4NucleiProperties.hh"
262
264GetCrossSection(const G4DynamicParticle* projectile, const G4Element*anE, G4double aT)
265{
266 G4double result = 0;
267 G4bool outOfRange;
268 G4int index = anE->GetIndex();
269
270 // prepare neutron
271 G4double eKinetic = projectile->GetKineticEnergy();
272
273 if ( G4ParticleHPManager::GetInstance()->GetNeglectDoppler() )
274 {
275 //NEGLECT_DOPPLER
276 G4double factor = 1.0;
277 if ( eKinetic < aT * CLHEP::k_Boltzmann )
278 {
279 // below 0.1 eV neutrons
280 // Have to do some, but now just igonre.
281 // Will take care after performance check.
282 // factor = factor * targetV;
283 }
284 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
285 }
286
287 G4ReactionProduct theNeutron( projectile->GetDefinition() );
288 theNeutron.SetMomentum( projectile->GetMomentum() );
289 theNeutron.SetKineticEnergy( eKinetic );
290
291 // prepare thermal nucleus
292 G4Nucleus aNuc;
293 G4double eps = 0.0001;
294 G4double theA = anE->GetN();
295 G4double theZ = anE->GetZ();
296 G4double eleMass;
297 eleMass = G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theA+eps), static_cast<G4int>(theZ+eps) );
298
299 G4ReactionProduct boosted;
300 G4double aXsection;
301
302 // MC integration loop
303 G4int counter = 0;
304 G4int failCount = 0;
305 G4double buffer = 0; G4int size = G4int(std::max(10., aT/60*CLHEP::kelvin));
306 G4ThreeVector neutronVelocity = 1./theProjectile->GetPDGMass()*theNeutron.GetMomentum();
307 G4double neutronVMag = neutronVelocity.mag();
308
309 // G4cout << " G4ParticleHPInelasticData 2 " << size << G4endl;//GDEB
310#ifndef G4PHP_DOPPLER_LOOP_ONCE
311 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer) // Loop checking, 11.05.2015, T. Koi
312 {
313 if(counter) buffer = result/counter;
314 while (counter<size) // Loop checking, 11.05.2015, T. Koi
315 {
316 counter ++;
317#endif
318 //G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus( eleMass/theProjectile->GetPDGMass(), aT );
319 //G4Nucleus::GetThermalNucleus requests normalization of mass in neutron mass
320 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus( eleMass/G4Neutron::Neutron()->GetPDGMass(), aT );
321 boosted.Lorentz(theNeutron, aThermalNuc);
322 G4double theEkin = boosted.GetKineticEnergy();
323 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
324 // G4cout << " G4ParticleHPInelasticData aXsection " << aXsection << " index " << index << " theEkin " << theEkin << " outOfRange " << outOfRange <<G4endl;//GDEB
325 if(aXsection <0)
326 {
327 if(failCount<1000)
328 {
329 failCount++;
330#ifndef G4PHP_DOPPLER_LOOP_ONCE
331 counter--;
332 continue;
333#endif
334 }
335 else
336 {
337 aXsection = 0;
338 }
339 }
340 // velocity correction.
341 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
342 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
343 result += aXsection;
344 }
345#ifndef G4PHP_DOPPLER_LOOP_ONCE
346 size += size;
347 }
348 result /= counter;
349#endif
350
351/*
352 // Checking impact of G4PHP_NEGLECT_DOPPLER
353 G4cout << " result " << result << " "
354 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
355 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
356*/
357// G4cout << this << " G4ParticleHPInelasticData result " << result << G4endl; //GDEB
358
359 return result;
360}
361
363{
365}
367{
369}
371{
372 outFile << "Extension of High Precision cross section for inelastic reaction of proton, deuteron, triton, He3 and alpha below 20MeV\n";
373}
static const G4double eps
std::vector< G4Element * > G4ElementTable
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:397
G4double GetZ() const
Definition: G4Element.hh:131
static size_t GetNumberOfElements()
Definition: G4Element.cc:404
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
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 &)
G4ParticleDefinition * theProjectile
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)
static constexpr double k_Boltzmann
static constexpr double barn
Definition: SystemOfUnits.h:86
static constexpr double kelvin
static constexpr double MeV
static constexpr double eV
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4bool IsWorkerThread()
Definition: G4Threading.cc:123
G4bool IsMasterThread()
Definition: G4Threading.cc:124
string material
Definition: eplot.py:19
#define G4ThreadLocal
Definition: tls.hh:77
#define buffer
Definition: xmlparse.cc:628