Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
97 SetMinKinEnergy( 0*CLHEP::MeV );
98 SetMaxKinEnergy( 20*CLHEP::MeV );
99
100 onFlightDB = true;
101 theCrossSections = 0;
102 theProjectile=projectile;
103
104 theHPData = NULL;
105 instanceOfWorker = false;
107 theHPData = new G4ParticleHPData( theProjectile );
108 } else {
109 instanceOfWorker = true;
110 }
111 element_cache = NULL;
112 material_cache = NULL;
113 ke_cache = 0.0;
114 xs_cache = 0.0;
115}
116
118{
119 if ( theCrossSections != NULL && instanceOfWorker != true ) {
120 theCrossSections->clearAndDestroy();
121 delete theCrossSections;
122 theCrossSections = NULL;
123 }
124 if ( theHPData != NULL && instanceOfWorker != true ) {
125 delete theHPData;
126 theHPData = NULL;
127 }
128}
129
131 G4int /*Z*/ , G4int /*A*/ ,
132 const G4Element* /*elm*/ ,
133 const G4Material* /*mat*/ )
134{
135 G4double eKin = dp->GetKineticEnergy();
136 if ( eKin > GetMaxKinEnergy()
137 || eKin < GetMinKinEnergy()
138 || dp->GetDefinition() != theProjectile ) return false;
139
140 return true;
141}
142
144 G4int /*Z*/ , G4int /*A*/ ,
145 const G4Isotope* /*iso*/ ,
146 const G4Element* element ,
147 const G4Material* material )
148{
149 if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && 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
159/*
160G4bool G4ParticleHPInelasticData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
161{
162 G4bool result = true;
163 G4double eKin = aP->GetKineticEnergy();
164 if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
165 return result;
166}
167*/
168
169//void G4ParticleHPInelasticData::BuildPhysicsTableHP(G4ParticleDefinition* projectile,const char* /* dataDirVariable */)
171{
172 // if(&projectile!=G4Neutron::Neutron())
173 // throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
174
175//080428
176 if ( G4ParticleHPManager::GetInstance()->GetNeglectDoppler() )
177 {
178 onFlightDB = false;
179 #ifdef G4VERBOSE
181 G4cout << "Find a flag of \"G4PHP_NEGLECT_DOPPLER\"." << G4endl;
182 G4cout << "On the fly Doppler broadening will be neglect in the cross section calculation of inelastic scattering of neutrons (<20MeV)." << G4endl;
183 }
184 #endif
185 }
186
188 theCrossSections = G4ParticleHPManager::GetInstance()->GetInelasticCrossSections( &projectile );
189 return;
190 } else {
191 if ( theHPData == NULL ) theHPData = G4ParticleHPData::Instance( const_cast<G4ParticleDefinition*> ( &projectile ) );
192 }
193
194
195
196 size_t numberOfElements = G4Element::GetNumberOfElements();
197// theCrossSections = new G4PhysicsTable( numberOfElements );
198// TKDB
199 //if ( theCrossSections == 0 )
200 //{ theCrossSections = new G4PhysicsTable( numberOfElements ); }
201 if ( theCrossSections == NULL )
202 theCrossSections = new G4PhysicsTable( numberOfElements );
203 else
204 theCrossSections->clearAndDestroy();
205
206 // make a PhysicsVector for each element
207
208 //G4ParticleHPData* hpData = new G4ParticleHPData(projectile); //NEW
209 static G4ThreadLocal G4ElementTable *theElementTable = 0 ;
210 if (!theElementTable) theElementTable= G4Element::GetElementTable();
211 for( size_t i=0; i<numberOfElements; ++i )
212 {
213 //NEW G4PhysicsVector* physVec = G4ParticleHPData::
214 //NEW Instance(projectile, dataDirVariable)->MakePhysicsVector((*theElementTable)[i], this);
215 //G4PhysicsVector* physVec = hpData->MakePhysicsVector((*theElementTable)[i], this);
216 G4PhysicsVector* physVec = theHPData->MakePhysicsVector((*theElementTable)[i], this);
217 theCrossSections->push_back(physVec);
218 }
219
220 G4ParticleHPManager::GetInstance()->RegisterInelasticCrossSections( &projectile , theCrossSections );
221}
222
224{
225 if(&projectile!=theProjectile)
226 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use ParticleHP data for a wrong projectile!!!");
227
228 #ifdef G4VERBOSE
229 if ( G4HadronicParameters::Instance()->GetVerboseLevel() == 0 ) return;
230//
231// Dump element based cross section
232// range 10e-5 eV to 20 MeV
233// 10 point per decade
234// in barn
235//
236
237 G4cout << G4endl;
238 G4cout << G4endl;
239 G4cout << "Inelastic Cross Section of Neutron HP"<< G4endl;
240 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
241 G4cout << G4endl;
242 G4cout << "Name of Element" << G4endl;
243 G4cout << "Energy[eV] XS[barn]" << G4endl;
244 G4cout << G4endl;
245
246 size_t numberOfElements = G4Element::GetNumberOfElements();
247 static G4ThreadLocal G4ElementTable *theElementTable = 0 ;
248 if (!theElementTable) theElementTable= G4Element::GetElementTable();
249
250 for ( size_t i = 0 ; i < numberOfElements ; ++i )
251 {
252
253 G4cout << (*theElementTable)[i]->GetName() << G4endl;
254
255 G4int ie = 0;
256
257 for ( ie = 0 ; ie < 130 ; ie++ )
258 {
259 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *CLHEP::eV;
260 G4bool outOfRange = false;
261
262 if ( eKinetic < 20*CLHEP::MeV )
263 {
264 G4cout << eKinetic/CLHEP::eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/CLHEP::barn << G4endl;
265 }
266
267 }
268
269 G4cout << G4endl;
270 }
271
272 //G4cout << "G4ParticleHPInelasticData::DumpPhysicsTable still to be implemented"<<G4endl;
273 #endif
274}
275
276#include "G4NucleiProperties.hh"
277
279GetCrossSection(const G4DynamicParticle* projectile, const G4Element*anE, G4double aT)
280{
281 G4double result = 0;
282 G4bool outOfRange;
283 G4int index = anE->GetIndex();
284
285 // prepare neutron
286 G4double eKinetic = projectile->GetKineticEnergy();
287
288 if ( !onFlightDB )
289 {
290 //NEGLECT_DOPPLER
291 G4double factor = 1.0;
292 if ( eKinetic < aT * CLHEP::k_Boltzmann )
293 {
294 // below 0.1 eV neutrons
295 // Have to do some, but now just igonre.
296 // Will take care after performance check.
297 // factor = factor * targetV;
298 }
299 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
300
301 }
302
303 G4ReactionProduct theNeutron( projectile->GetDefinition() );
304 theNeutron.SetMomentum( projectile->GetMomentum() );
305 theNeutron.SetKineticEnergy( eKinetic );
306
307 // prepare thermal nucleus
308 G4Nucleus aNuc;
309 G4double eps = 0.0001;
310 G4double theA = anE->GetN();
311 G4double theZ = anE->GetZ();
312 G4double eleMass;
313 eleMass = G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theA+eps), static_cast<G4int>(theZ+eps) );
314
315 G4ReactionProduct boosted;
316 G4double aXsection;
317
318 // MC integration loop
319 G4int counter = 0;
320 G4int failCount = 0;
321 G4double buffer = 0; G4int size = G4int(std::max(10., aT/60*CLHEP::kelvin));
322 G4ThreeVector neutronVelocity = 1./theProjectile->GetPDGMass()*theNeutron.GetMomentum();
323 G4double neutronVMag = neutronVelocity.mag();
324
325 // G4cout << " G4ParticleHPInelasticData 2 " << size << G4endl;//GDEB
326#ifndef G4PHP_DOPPLER_LOOP_ONCE
327 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer) // Loop checking, 11.05.2015, T. Koi
328 {
329 if(counter) buffer = result/counter;
330 while (counter<size) // Loop checking, 11.05.2015, T. Koi
331 {
332 counter ++;
333#endif
334 //G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus( eleMass/theProjectile->GetPDGMass(), aT );
335 //G4Nucleus::GetThermalNucleus requests normalization of mass in neutron mass
336 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus( eleMass/G4Neutron::Neutron()->GetPDGMass(), aT );
337 boosted.Lorentz(theNeutron, aThermalNuc);
338 G4double theEkin = boosted.GetKineticEnergy();
339 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
340 // G4cout << " G4ParticleHPInelasticData aXsection " << aXsection << " index " << index << " theEkin " << theEkin << " outOfRange " << outOfRange <<G4endl;//GDEB
341 if(aXsection <0)
342 {
343 if(failCount<1000)
344 {
345 failCount++;
346#ifndef G4PHP_DOPPLER_LOOP_ONCE
347 counter--;
348 continue;
349#endif
350 }
351 else
352 {
353 aXsection = 0;
354 }
355 }
356 // velocity correction.
357 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
358 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
359 result += aXsection;
360 }
361#ifndef G4PHP_DOPPLER_LOOP_ONCE
362 size += size;
363 }
364 result /= counter;
365#endif
366
367/*
368 // Checking impact of G4PHP_NEGLECT_DOPPLER
369 G4cout << " result " << result << " "
370 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
371 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
372*/
373// G4cout << this << " G4ParticleHPInelasticData result " << result << G4endl; //GDEB
374
375 return result;
376}
377
379{
381}
383{
385}
387{
388 outFile << "Extension of High Precision cross section for inelastic reaction of proton, deuteron, triton, He3 and alpha below 20MeV\n";
389}
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:130
static size_t GetNumberOfElements()
Definition: G4Element.cc:404
size_t GetIndex() const
Definition: G4Element.hh:181
G4double GetN() const
Definition: G4Element.hh:134
static G4HadronicParameters * Instance()
static G4He3 * He3()
Definition: G4He3.cc:93
G4double GetTemperature() const
Definition: G4Material.hh:180
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:143
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:94
void SetMaxKinEnergy(G4double value)
void SetName(const G4String &)
void SetMinKinEnergy(G4double value)
G4bool IsWorkerThread()
Definition: G4Threading.cc:123
G4bool IsMasterThread()
Definition: G4Threading.cc:124
#define G4ThreadLocal
Definition: tls.hh:77
#define buffer
Definition: xmlparse.cc:628