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

#include <G4NeutronHPorLFissionData.hh>

+ Inheritance diagram for G4NeutronHPorLFissionData:

Public Member Functions

 G4NeutronHPorLFissionData ()
 
 G4NeutronHPorLFissionData (G4NeutronHPChannel *, std::set< G4String > *)
 
 ~G4NeutronHPorLFissionData ()
 
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 &)
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
virtual G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
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 55 of file G4NeutronHPorLFissionData.hh.

Constructor & Destructor Documentation

◆ G4NeutronHPorLFissionData() [1/2]

G4NeutronHPorLFissionData::G4NeutronHPorLFissionData ( )

Definition at line 43 of file G4NeutronHPorLFissionData.cc.

44{
45 SetMinKinEnergy( 0*MeV );
46 SetMaxKinEnergy( 20*MeV );
47
48 ke_cache = 0.0;
49 xs_cache = 0.0;
50 element_cache = NULL;
51 material_cache = NULL;
52// BuildPhysicsTable(*G4Neutron::Neutron());
53}
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)

◆ G4NeutronHPorLFissionData() [2/2]

G4NeutronHPorLFissionData::G4NeutronHPorLFissionData ( G4NeutronHPChannel pChannel,
std::set< G4String > *  pSet 
)

Definition at line 91 of file G4NeutronHPorLFissionData.cc.

92:G4VCrossSectionDataSet("NeutronHPorLFissionXS")
93{
94 theFissionChannel = pChannel;
95 unavailable_elements = pSet;
96
97 SetMinKinEnergy( 0*MeV );
98 SetMaxKinEnergy( 20*MeV );
99
100 ke_cache = 0.0;
101 xs_cache = 0.0;
102 element_cache = NULL;
103 material_cache = NULL;
104}

◆ ~G4NeutronHPorLFissionData()

G4NeutronHPorLFissionData::~G4NeutronHPorLFissionData ( )

Definition at line 55 of file G4NeutronHPorLFissionData.cc.

56{
57// delete theCrossSections;
58}

Member Function Documentation

◆ BuildPhysicsTable()

void G4NeutronHPorLFissionData::BuildPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 119 of file G4NeutronHPorLFissionData.cc.

120{
121 if( &aP!=G4Neutron::Neutron() )
122 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
123}
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104

◆ DumpPhysicsTable()

void G4NeutronHPorLFissionData::DumpPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 127 of file G4NeutronHPorLFissionData.cc.

128{
129 if(&aP!=G4Neutron::Neutron())
130 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
131// G4cout << "G4NeutronHPorLFissionData::DumpPhysicsTable still to be implemented"<<G4endl;
132}

◆ GetCrossSection()

G4double G4NeutronHPorLFissionData::GetCrossSection ( const G4DynamicParticle aP,
const G4Element anE,
G4double  aT 
)

Definition at line 141 of file G4NeutronHPorLFissionData.cc.

143{
144 G4double result = 0;
145 if ( anE->GetZ() < 90 ) return result;
146
147// G4bool outOfRange;
148 G4int index = anE->GetIndex();
149
150 // prepare neutron
151 G4double eKinetic = aP->GetKineticEnergy();
152 G4ReactionProduct theNeutron( aP->GetDefinition() );
153 theNeutron.SetMomentum( aP->GetMomentum() );
154 theNeutron.SetKineticEnergy( eKinetic );
155
156 // prepare thermal nucleus
157 G4Nucleus aNuc;
158 G4double eps = 0.0001;
159 G4double theA = anE->GetN();
160 G4double theZ = anE->GetZ();
161 G4double eleMass;
162 eleMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theA+eps), static_cast<G4int>(theZ+eps))
164
165 G4ReactionProduct boosted;
166 G4double aXsection;
167
168 // MC integration loop
169 G4int counter = 0;
170 G4double buffer = 0;
171 G4int size = G4int(std::max(10., aT/60*kelvin));
172 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
173 G4double neutronVMag = neutronVelocity.mag();
174 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer)
175 {
176 if(counter) buffer = result/counter;
177 while (counter<size)
178 {
179 counter ++;
180 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
181 boosted.Lorentz(theNeutron, aThermalNuc);
182 G4double theEkin = boosted.GetKineticEnergy();
183 //aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
184 aXsection = theFissionChannel[index].GetXsec( theEkin );
185 // velocity correction.
186 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
187 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
188 result += aXsection;
189 }
190 size += size;
191 }
192 result /= counter;
193 return result;
194}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
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:134
G4double GetXsec(G4double energy)
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:130
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4double GetMass() const
#define buffer
Definition: xmlparse.cc:611

Referenced by GetIsoCrossSection().

◆ GetIsoCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 74 of file G4NeutronHPorLFissionData.cc.

79{
80 if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
81
82 ke_cache = dp->GetKineticEnergy();
83 element_cache = element;
84 material_cache = material;
85 G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
86 xs_cache = xs;
87 return xs;
88 //return GetCrossSection( dp , element , material->GetTemperature() );
89}
G4double GetTemperature() const
Definition: G4Material.hh:181
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 60 of file G4NeutronHPorLFissionData.cc.

64{
65 G4double eKin = dp->GetKineticEnergy();
66 if ( eKin > GetMaxKinEnergy()
67 || eKin < GetMinKinEnergy()
68 || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
69 if ( unavailable_elements->find( element->GetName() ) != unavailable_elements->end() ) return false;
70
71 return true;
72}
const G4String & GetName() const
Definition: G4Element.hh:127

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