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

#include <G4NeutronHPorLEInelasticData.hh>

+ Inheritance diagram for G4NeutronHPorLEInelasticData:

Public Member Functions

 G4NeutronHPorLEInelasticData ()
 
 G4NeutronHPorLEInelasticData (G4NeutronHPChannelList *, std::set< G4String > *)
 
 ~G4NeutronHPorLEInelasticData ()
 
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 56 of file G4NeutronHPorLEInelasticData.hh.

Constructor & Destructor Documentation

◆ G4NeutronHPorLEInelasticData() [1/2]

G4NeutronHPorLEInelasticData::G4NeutronHPorLEInelasticData ( )

Definition at line 42 of file G4NeutronHPorLEInelasticData.cc.

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

◆ G4NeutronHPorLEInelasticData() [2/2]

G4NeutronHPorLEInelasticData::G4NeutronHPorLEInelasticData ( G4NeutronHPChannelList pChannel,
std::set< G4String > *  pSet 
)

Definition at line 90 of file G4NeutronHPorLEInelasticData.cc.

91:G4VCrossSectionDataSet("NeutronHPorLEInelasticXS")
92{
93 theInelasticChannel = pChannel;
94 unavailable_elements = pSet;
95 //BuildPhysicsTable(*G4Neutron::Neutron());
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}

◆ ~G4NeutronHPorLEInelasticData()

G4NeutronHPorLEInelasticData::~G4NeutronHPorLEInelasticData ( )

Definition at line 54 of file G4NeutronHPorLEInelasticData.cc.

55{
56// delete theCrossSections;
57}

Member Function Documentation

◆ BuildPhysicsTable()

void G4NeutronHPorLEInelasticData::BuildPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 121 of file G4NeutronHPorLEInelasticData.cc.

122{
123 if( &aP!=G4Neutron::Neutron() )
124 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
125
126 size_t numberOfElements = G4Element::GetNumberOfElements();
127 theCrossSections = new G4PhysicsTable( numberOfElements );
128
129 static const G4ElementTable *theElementTable = G4Element::GetElementTable();
130 for ( size_t i=0 ; i < numberOfElements; ++i )
131 {
132 G4PhysicsVector* thePhysVec = new G4LPhysicsFreeVector(0, 0, 0);
133
134 if ( unavailable_elements->find( (*theElementTable)[i]->GetName() ) == unavailable_elements->end() )
135 {
136
137 G4NeutronHPElementData* theElementData = new G4NeutronHPElementData();
138 theElementData->Init( (*theElementTable)[i] );
139
140 G4NeutronHPVector* theHPVector = theElementData->GetData( (G4NeutronHPInelasticData*)this );
141
142 G4int len = theHPVector->GetVectorLength();
143
144 if ( len!=0 )
145 {
146 G4double emin = theHPVector->GetX(0);
147 G4double emax = theHPVector->GetX(len-1);
148
149 G4LPhysicsFreeVector* aPhysVector= new G4LPhysicsFreeVector ( len , emin , emax );
150 for ( G4int ii=0; ii<len; ii++ )
151 {
152 aPhysVector->PutValues( ii , theHPVector->GetX(ii) , theHPVector->GetY(ii) );
153 }
154 delete thePhysVec;
155 thePhysVec = aPhysVector;
156 }
157
158 //G4PhysicsVector* physVec = G4NeutronHPData::
159 //Instance()->MakePhysicsVector((*theElementTable)[i], this);
160 //theCrossSections->push_back(physVec);
161 }
162
163 theCrossSections->push_back(thePhysVec);
164 }
165}
std::vector< G4Element * > G4ElementTable
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
static size_t GetNumberOfElements()
Definition: G4Element.cc:406
static const G4ElementTable * GetElementTable()
Definition: G4Element.cc:399
void PutValues(size_t binNumber, G4double binValue, G4double dataValue)
void Init(G4Element *theElement)
G4NeutronHPVector * GetData(G4NeutronHPFissionData *)
G4int GetVectorLength() const
G4double GetX(G4int i) const
G4double GetY(G4double x)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
void push_back(G4PhysicsVector *)

◆ DumpPhysicsTable()

void G4NeutronHPorLEInelasticData::DumpPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 169 of file G4NeutronHPorLEInelasticData.cc.

170{
171 if(&aP!=G4Neutron::Neutron())
172 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
173// G4cout << "G4NeutronHPorLEInelasticData::DumpPhysicsTable still to be implemented"<<G4endl;
174}

◆ GetCrossSection()

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

Definition at line 183 of file G4NeutronHPorLEInelasticData.cc.

185{
186
187 // G4cout << "Choice G4NeutronHPorLEInelasticData for element " << anE->GetName() << G4endl;
188 G4double result = 0;
189 //G4bool outOfRange;
190 G4int index = anE->GetIndex();
191
192 // prepare neutron
193 G4double eKinetic = aP->GetKineticEnergy();
194 G4ReactionProduct theNeutron( aP->GetDefinition() );
195 theNeutron.SetMomentum( aP->GetMomentum() );
196 theNeutron.SetKineticEnergy( eKinetic );
197
198 // prepare thermal nucleus
199 G4Nucleus aNuc;
200 G4double eps = 0.0001;
201 G4double theA = anE->GetN();
202 G4double theZ = anE->GetZ();
203 G4double eleMass;
204 eleMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theA+eps), static_cast<G4int>(theZ+eps))
206
207 G4ReactionProduct boosted;
208 G4double aXsection;
209
210 // MC integration loop
211 G4int counter = 0;
212 G4double buffer = 0;
213 G4int size = G4int(std::max(10., aT/60*kelvin));
214 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
215 G4double neutronVMag = neutronVelocity.mag();
216 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer)
217 {
218 if(counter) buffer = result/counter;
219 while (counter<size)
220 {
221 counter ++;
222 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
223 boosted.Lorentz(theNeutron, aThermalNuc);
224 G4double theEkin = boosted.GetKineticEnergy();
225 //aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
226 aXsection = theInelasticChannel[index].GetXsec( theEkin );
227 // velocity correction.
228 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
229 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
230 result += aXsection;
231 }
232 size += size;
233 }
234 result /= counter;
235 //return result;
236 return result*barn;
237}
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 anEnergy)
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 G4NeutronHPorLEInelasticData::GetIsoCrossSection ( const G4DynamicParticle dp,
G4int  ,
G4int  ,
const G4Isotope ,
const G4Element element,
const G4Material material 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 73 of file G4NeutronHPorLEInelasticData.cc.

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

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 59 of file G4NeutronHPorLEInelasticData.cc.

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

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