Geant4 11.1.1
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 &)
 
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 ComputeCrossSectionPerElement (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, 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 G4double ComputeIsoCrossSection (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, 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 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
 
void SetName (const G4String &nam)
 
G4VCrossSectionDataSetoperator= (const G4VCrossSectionDataSet &right)=delete
 
 G4VCrossSectionDataSet (const G4VCrossSectionDataSet &)=delete
 

Additional Inherited Members

- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 
G4String name
 

Detailed Description

Definition at line 53 of file G4ParticleHPInelasticData.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPInelasticData()

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

Definition at line 48 of file G4ParticleHPInelasticData.cc.

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}
const char * G4FindDataDir(const char *)
#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:93
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
void SetName(const G4String &nam)
G4bool IsMasterThread()
Definition: G4Threading.cc:124

◆ ~G4ParticleHPInelasticData()

G4ParticleHPInelasticData::~G4ParticleHPInelasticData ( )

Definition at line 116 of file G4ParticleHPInelasticData.cc.

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

Member Function Documentation

◆ BuildPhysicsTable()

void G4ParticleHPInelasticData::BuildPhysicsTable ( const G4ParticleDefinition projectile)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 159 of file G4ParticleHPInelasticData.cc.

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}
std::vector< G4Element * > G4ElementTable
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
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 332 of file G4ParticleHPInelasticData.cc.

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

◆ DumpPhysicsTable()

void G4ParticleHPInelasticData::DumpPhysicsTable ( const G4ParticleDefinition projectile)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 186 of file G4ParticleHPInelasticData.cc.

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}
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

◆ GetCrossSection()

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

Definition at line 235 of file G4ParticleHPInelasticData.cc.

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}
double mag() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetZ() const
Definition: G4Element.hh:131
size_t GetIndex() const
Definition: G4Element.hh:182
G4double GetN() const
Definition: G4Element.hh:135
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:236
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4double GetMass() const

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 142 of file G4ParticleHPInelasticData.cc.

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

◆ GetProjectile()

G4ParticleDefinition * G4ParticleHPInelasticData::GetProjectile ( )
inline

Definition at line 85 of file G4ParticleHPInelasticData.hh.

85{return theProjectile;}

◆ GetVerboseLevel()

G4int G4ParticleHPInelasticData::GetVerboseLevel ( ) const

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 129 of file G4ParticleHPInelasticData.cc.

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

◆ SetVerboseLevel()

void G4ParticleHPInelasticData::SetVerboseLevel ( G4int  newValue)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 327 of file G4ParticleHPInelasticData.cc.

328{
330}

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