Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPInelasticData Class Reference

#include <G4ParticleHPInelasticData.hh>

+ Inheritance diagram for G4ParticleHPInelasticData:

Public Member Functions

 G4ParticleHPInelasticData (G4ParticleDefinition *projectile=G4Neutron::Neutron())
 
 ~G4ParticleHPInelasticData ()
 
G4bool IsIsoApplicable (const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
 
G4double GetIsoCrossSection (const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, G4double aT)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
void DumpPhysicsTable (const G4ParticleDefinition &)
 
void IgnoreOnFlightDopplerBroadening ()
 
void EnableOnFlightDopplerBroadening ()
 
G4ParticleDefinitionGetProjectile ()
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int)
 
virtual void CrossSectionDescription (std::ostream &) const
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
virtual const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
virtual G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
bool ForAllAtomsAndEnergies () const
 
void SetForAllAtomsAndEnergies (G4bool val)
 
const G4StringGetName () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VCrossSectionDataSet
void SetName (const G4String &)
 
- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 

Detailed Description

Definition at line 53 of file G4ParticleHPInelasticData.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPInelasticData()

G4ParticleHPInelasticData::G4ParticleHPInelasticData ( G4ParticleDefinition projectile = G4Neutron::Neutron())

Definition at line 47 of file G4ParticleHPInelasticData.cc.

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}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:93
static G4HadronicParameters * Instance()
static G4He3 * He3()
Definition: G4He3.cc:93
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
const G4String & GetParticleName() const
static G4Proton * Proton()
Definition: G4Proton.cc:92
static G4Triton * Triton()
Definition: G4Triton.cc:94
void SetMaxKinEnergy(G4double value)
void SetName(const G4String &)
void SetMinKinEnergy(G4double value)
G4bool IsMasterThread()
Definition: G4Threading.cc:124

◆ ~G4ParticleHPInelasticData()

G4ParticleHPInelasticData::~G4ParticleHPInelasticData ( )

Definition at line 117 of file G4ParticleHPInelasticData.cc.

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}
void clearAndDestroy()

Member Function Documentation

◆ BuildPhysicsTable()

void G4ParticleHPInelasticData::BuildPhysicsTable ( const G4ParticleDefinition projectile)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 170 of file G4ParticleHPInelasticData.cc.

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}
std::vector< G4Element * > G4ElementTable
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:397
static size_t GetNumberOfElements()
Definition: G4Element.cc:404
static G4ParticleHPData * Instance(G4ParticleDefinition *projectile)
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4ParticleHPFissionData *theP)
G4PhysicsTable * GetInelasticCrossSections(const G4ParticleDefinition *)
void RegisterInelasticCrossSections(const G4ParticleDefinition *, G4PhysicsTable *)
static G4ParticleHPManager * GetInstance()
void push_back(G4PhysicsVector *)
G4bool IsWorkerThread()
Definition: G4Threading.cc:123
#define G4ThreadLocal
Definition: tls.hh:77

◆ CrossSectionDescription()

void G4ParticleHPInelasticData::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 386 of file G4ParticleHPInelasticData.cc.

387{
388 outFile << "Extension of High Precision cross section for inelastic reaction of proton, deuteron, triton, He3 and alpha below 20MeV\n";
389}

◆ DumpPhysicsTable()

void G4ParticleHPInelasticData::DumpPhysicsTable ( const G4ParticleDefinition projectile)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 223 of file G4ParticleHPInelasticData.cc.

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}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:230

◆ EnableOnFlightDopplerBroadening()

void G4ParticleHPInelasticData::EnableOnFlightDopplerBroadening ( )
inline

Definition at line 86 of file G4ParticleHPInelasticData.hh.

86{ onFlightDB = true; };

◆ GetCrossSection()

G4double G4ParticleHPInelasticData::GetCrossSection ( const G4DynamicParticle projectile,
const G4Element anE,
G4double  aT 
)

Definition at line 278 of file G4ParticleHPInelasticData.cc.

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}
double mag() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetZ() const
Definition: G4Element.hh:130
size_t GetIndex() const
Definition: G4Element.hh:181
G4double GetN() const
Definition: G4Element.hh:134
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:143
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4double GetMass() const
#define buffer
Definition: xmlparse.cc:628

Referenced by GetIsoCrossSection().

◆ GetIsoCrossSection()

G4double G4ParticleHPInelasticData::GetIsoCrossSection ( const G4DynamicParticle dp,
G4int  ,
G4int  ,
const G4Isotope ,
const G4Element element,
const G4Material material 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 143 of file G4ParticleHPInelasticData.cc.

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}
G4double GetTemperature() const
Definition: G4Material.hh:180
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)

◆ GetProjectile()

G4ParticleDefinition * G4ParticleHPInelasticData::GetProjectile ( )
inline

Definition at line 88 of file G4ParticleHPInelasticData.hh.

88{return theProjectile;}

◆ GetVerboseLevel()

G4int G4ParticleHPInelasticData::GetVerboseLevel ( ) const
virtual

◆ IgnoreOnFlightDopplerBroadening()

void G4ParticleHPInelasticData::IgnoreOnFlightDopplerBroadening ( )
inline

Definition at line 85 of file G4ParticleHPInelasticData.hh.

85{ onFlightDB = false; };

◆ IsIsoApplicable()

G4bool G4ParticleHPInelasticData::IsIsoApplicable ( const G4DynamicParticle dp,
G4int  ,
G4int  ,
const G4Element ,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 130 of file G4ParticleHPInelasticData.cc.

134{
135 G4double eKin = dp->GetKineticEnergy();
136 if ( eKin > GetMaxKinEnergy()
137 || eKin < GetMinKinEnergy()
138 || dp->GetDefinition() != theProjectile ) return false;
139
140 return true;
141}

◆ SetVerboseLevel()

void G4ParticleHPInelasticData::SetVerboseLevel ( G4int  newValue)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 382 of file G4ParticleHPInelasticData.cc.

383{
385}

The documentation for this class was generated from the following files: