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

#include <G4ParticleHPCaptureData.hh>

+ Inheritance diagram for G4ParticleHPCaptureData:

Public Member Functions

 G4ParticleHPCaptureData ()
 
 ~G4ParticleHPCaptureData ()
 
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 &)
 
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 51 of file G4ParticleHPCaptureData.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPCaptureData()

G4ParticleHPCaptureData::G4ParticleHPCaptureData ( )

Definition at line 52 of file G4ParticleHPCaptureData.cc.

53:G4VCrossSectionDataSet("NeutronHPCaptureXS")
54{
55 SetMinKinEnergy( 0*MeV );
56 SetMaxKinEnergy( 20*MeV );
57
58 theCrossSections = 0;
59
60 instanceOfWorker = false;
62 instanceOfWorker = true;
63 }
64
65 element_cache = nullptr;
66 material_cache = nullptr;
67 ke_cache = 0.0;
68 xs_cache = 0.0;
69
70 //BuildPhysicsTable(*G4Neutron::Neutron());
71}
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
G4bool IsWorkerThread()
Definition: G4Threading.cc:123

◆ ~G4ParticleHPCaptureData()

G4ParticleHPCaptureData::~G4ParticleHPCaptureData ( )

Definition at line 73 of file G4ParticleHPCaptureData.cc.

74{
75 if ( theCrossSections != nullptr && instanceOfWorker != true ) {
76 theCrossSections->clearAndDestroy();
77 delete theCrossSections;
78 theCrossSections = nullptr;
79 }
80}
void clearAndDestroy()

Member Function Documentation

◆ BuildPhysicsTable()

void G4ParticleHPCaptureData::BuildPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 121 of file G4ParticleHPCaptureData.cc.

122{
123 if(&aP!=G4Neutron::Neutron())
124 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
125
128 return;
129 }
130
131 std::size_t numberOfElements = G4Element::GetNumberOfElements();
132 // G4cout << "CALLED G4ParticleHPCaptureData::BuildPhysicsTable "<<numberOfElements<<G4endl;
133 // TKDB
134 //if ( theCrossSections == 0 ) theCrossSections = new G4PhysicsTable( numberOfElements );
135 if ( theCrossSections == nullptr )
136 theCrossSections = new G4PhysicsTable( numberOfElements );
137 else
138 theCrossSections->clearAndDestroy();
139
140 // make a PhysicsVector for each element
141
142 static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
143 for( std::size_t i=0; i<numberOfElements; ++i )
144 {
145 #ifdef G4VERBOSE
146 if(std::getenv("CaptureDataIndexDebug"))
147 {
148 std::size_t index_debug = ((*theElementTable)[i])->GetIndex();
149 if ( G4HadronicParameters::Instance()->GetVerboseLevel() > 0 ) G4cout << "IndexDebug "<< i <<" "<<index_debug<<G4endl;
150 }
151 #endif
153 Instance(G4Neutron::Neutron())->MakePhysicsVector((*theElementTable)[i], this);
154 theCrossSections->push_back(physVec);
155 }
156
158}
std::vector< G4Element * > G4ElementTable
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
static G4HadronicParameters * Instance()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4ParticleHPData * Instance(G4ParticleDefinition *projectile)
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4ParticleHPFissionData *theP)
G4PhysicsTable * GetCaptureCrossSections()
void RegisterCaptureCrossSections(G4PhysicsTable *val)
static G4ParticleHPManager * GetInstance()
void push_back(G4PhysicsVector *)
#define G4ThreadLocal
Definition: tls.hh:77

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 297 of file G4ParticleHPCaptureData.cc.

298{
299 outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF) for radiative capture reaction of neutrons below 20MeV\n" ;
300}

◆ DumpPhysicsTable()

void G4ParticleHPCaptureData::DumpPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 160 of file G4ParticleHPCaptureData.cc.

161{
162 if(&aP!=G4Neutron::Neutron())
163 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
164
165 #ifdef G4VERBOSE
166 if ( G4HadronicParameters::Instance()->GetVerboseLevel() == 0 ) return;
167
168//
169// Dump element based cross section
170// range 10e-5 eV to 20 MeV
171// 10 point per decade
172// in barn
173//
174
175 G4cout << G4endl;
176 G4cout << G4endl;
177 G4cout << "Capture Cross Section of Neutron HP"<< G4endl;
178 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
179 G4cout << G4endl;
180 G4cout << "Name of Element" << G4endl;
181 G4cout << "Energy[eV] XS[barn]" << G4endl;
182 G4cout << G4endl;
183
184 std::size_t numberOfElements = G4Element::GetNumberOfElements();
185 static G4ThreadLocal G4ElementTable *theElementTable = 0 ;
186 if (!theElementTable) theElementTable= G4Element::GetElementTable();
187
188 for ( std::size_t i = 0 ; i < numberOfElements ; ++i )
189 {
190
191 G4cout << (*theElementTable)[i]->GetName() << G4endl;
192
193 G4int ie = 0;
194
195 for ( ie = 0 ; ie < 130 ; ++ie )
196 {
197 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *eV;
198 G4bool outOfRange = false;
199
200 if ( eKinetic < 20*MeV )
201 {
202 G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
203 }
204
205 }
206
207 G4cout << G4endl;
208 }
209
210 //G4cout << "G4ParticleHPCaptureData::DumpPhysicsTable still to be implemented"<<G4endl;
211 #endif
212}
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 G4ParticleHPCaptureData::GetCrossSection ( const G4DynamicParticle aP,
const G4Element anE,
G4double  aT 
)

Definition at line 214 of file G4ParticleHPCaptureData.cc.

216{
217 G4double result = 0;
218 G4bool outOfRange;
219 G4int index = (G4int)anE->GetIndex();
220
221 // prepare neutron
222 G4double eKinetic = aP->GetKineticEnergy();
223
225 {
226 //NEGLECT_DOPPLER
227 G4double factor = 1.0;
228 if ( eKinetic < aT * k_Boltzmann )
229 {
230 // below 0.1 eV neutrons
231 // Have to do some, but now just igonre.
232 // Will take care after performance check.
233 // factor = factor * targetV;
234 }
235 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
236 }
237
238 G4ReactionProduct theNeutron( aP->GetDefinition() );
239 theNeutron.SetMomentum( aP->GetMomentum() );
240 theNeutron.SetKineticEnergy( eKinetic );
241
242 // prepare thermal nucleus
243 G4Nucleus aNuc;
244 G4double eps = 0.0001;
245 G4double theA = anE->GetN();
246 G4double theZ = anE->GetZ();
247 G4double eleMass;
248 eleMass = G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) ) / G4Neutron::Neutron()->GetPDGMass();
249
250 G4ReactionProduct boosted;
251 G4double aXsection;
252
253 // MC integration loop
254 G4int counter = 0;
255 G4double buffer = 0;
256 G4int size = G4int(std::max(10., aT/60*kelvin));
257 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
258 G4double neutronVMag = neutronVelocity.mag();
259
260 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer) // Loop checking, 11.05.2015, T. Koi
261 {
262 if(counter) buffer = result/counter;
263 while (counter<size) // Loop checking, 11.05.2015, T. Koi
264 {
265 ++counter;
266 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
267 boosted.Lorentz(theNeutron, aThermalNuc);
268 G4double theEkin = boosted.GetKineticEnergy();
269 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
270 // velocity correction, or luminosity factor...
271 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
272 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
273 result += aXsection;
274 }
275 size += size;
276 }
277 result /= counter;
278/*
279 // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
280 G4cout << " result " << result << " "
281 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
282 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
283*/
284 return result;
285}
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 G4ParticleHPCaptureData::GetIsoCrossSection ( const G4DynamicParticle dp,
G4int  ,
G4int  ,
const G4Isotope ,
const G4Element element,
const G4Material material 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 95 of file G4ParticleHPCaptureData.cc.

100{
101 if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
102
103 ke_cache = dp->GetKineticEnergy();
104 element_cache = element;
105 material_cache = material;
106 G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
107 xs_cache = xs;
108 return xs;
109}
G4double GetTemperature() const
Definition: G4Material.hh:177
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)

◆ GetVerboseLevel()

G4int G4ParticleHPCaptureData::GetVerboseLevel ( ) const

Definition at line 287 of file G4ParticleHPCaptureData.cc.

Referenced by DumpPhysicsTable().

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 82 of file G4ParticleHPCaptureData.cc.

86{
87 G4double eKin = dp->GetKineticEnergy();
88 if ( eKin > GetMaxKinEnergy()
89 || eKin < GetMinKinEnergy()
90 || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
91
92 return true;
93}

◆ SetVerboseLevel()

void G4ParticleHPCaptureData::SetVerboseLevel ( G4int  newValue)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 292 of file G4ParticleHPCaptureData.cc.

293{
295}

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