#include <G4NeutronHPFissionData.hh>
Inheritance diagram for G4NeutronHPFissionData:
Public Member Functions | |
G4NeutronHPFissionData () | |
~G4NeutronHPFissionData () | |
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 &) |
Definition at line 48 of file G4NeutronHPFissionData.hh.
G4NeutronHPFissionData::G4NeutronHPFissionData | ( | ) |
Definition at line 42 of file G4NeutronHPFissionData.cc.
References BuildPhysicsTable(), G4Neutron::Neutron(), G4VCrossSectionDataSet::SetMaxKinEnergy(), and G4VCrossSectionDataSet::SetMinKinEnergy().
00043 :G4VCrossSectionDataSet("NeutronHPFissionXS") 00044 { 00045 SetMinKinEnergy( 0*MeV ); 00046 SetMaxKinEnergy( 20*MeV ); 00047 00048 ke_cache = 0.0; 00049 xs_cache = 0.0; 00050 element_cache = NULL; 00051 material_cache = NULL; 00052 00053 theCrossSections = 0; 00054 BuildPhysicsTable(*G4Neutron::Neutron()); 00055 }
G4NeutronHPFissionData::~G4NeutronHPFissionData | ( | ) |
Definition at line 57 of file G4NeutronHPFissionData.cc.
References G4PhysicsTable::clearAndDestroy().
00058 { 00059 if ( theCrossSections != NULL ) theCrossSections->clearAndDestroy(); 00060 delete theCrossSections; 00061 }
void G4NeutronHPFissionData::BuildPhysicsTable | ( | const G4ParticleDefinition & | ) | [virtual] |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 102 of file G4NeutronHPFissionData.cc.
References G4PhysicsTable::clearAndDestroy(), G4Element::GetElementTable(), G4Element::GetNumberOfElements(), G4Neutron::Neutron(), and G4PhysicsTable::push_back().
Referenced by G4NeutronHPFissionData().
00103 { 00104 if(&aP!=G4Neutron::Neutron()) 00105 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!"); 00106 size_t numberOfElements = G4Element::GetNumberOfElements(); 00107 //theCrossSections = new G4PhysicsTable( numberOfElements ); 00108 // TKDB 00109 //if ( theCrossSections == NULL ) theCrossSections = new G4PhysicsTable( numberOfElements ); 00110 if ( theCrossSections == NULL ) 00111 theCrossSections = new G4PhysicsTable( numberOfElements ); 00112 else 00113 theCrossSections->clearAndDestroy(); 00114 00115 // make a PhysicsVector for each element 00116 00117 static const G4ElementTable *theElementTable = G4Element::GetElementTable(); 00118 for( size_t i=0; i<numberOfElements; ++i ) 00119 { 00120 G4PhysicsVector* physVec = G4NeutronHPData:: 00121 Instance()->MakePhysicsVector((*theElementTable)[i], this); 00122 theCrossSections->push_back(physVec); 00123 } 00124 }
void G4NeutronHPFissionData::DumpPhysicsTable | ( | const G4ParticleDefinition & | ) | [virtual] |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 126 of file G4NeutronHPFissionData.cc.
References G4cout, G4endl, G4Element::GetElementTable(), G4VCrossSectionDataSet::GetName(), G4Element::GetNumberOfElements(), and G4Neutron::Neutron().
00127 { 00128 if(&aP!=G4Neutron::Neutron()) 00129 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!"); 00130 00131 // 00132 // Dump element based cross section 00133 // range 10e-5 eV to 20 MeV 00134 // 10 point per decade 00135 // in barn 00136 // 00137 00138 G4cout << G4endl; 00139 G4cout << G4endl; 00140 G4cout << "Fission Cross Section of Neutron HP"<< G4endl; 00141 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl; 00142 G4cout << G4endl; 00143 G4cout << "Name of Element" << G4endl; 00144 G4cout << "Energy[eV] XS[barn]" << G4endl; 00145 G4cout << G4endl; 00146 00147 size_t numberOfElements = G4Element::GetNumberOfElements(); 00148 static const G4ElementTable *theElementTable = G4Element::GetElementTable(); 00149 00150 for ( size_t i = 0 ; i < numberOfElements ; ++i ) 00151 { 00152 00153 G4cout << (*theElementTable)[i]->GetName() << G4endl; 00154 00155 if ( (*((*theCrossSections)(i))).GetVectorLength() == 0 ) 00156 { 00157 G4cout << "The cross-section data of the fission of this element is not available." << G4endl; 00158 G4cout << G4endl; 00159 continue; 00160 } 00161 00162 G4int ie = 0; 00163 00164 for ( ie = 0 ; ie < 130 ; ie++ ) 00165 { 00166 G4double eKinetic = 1.0e-5 * std::pow ( 10.0 , ie/10.0 ) *eV; 00167 G4bool outOfRange = false; 00168 00169 if ( eKinetic < 20*MeV ) 00170 { 00171 G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl; 00172 } 00173 00174 } 00175 00176 G4cout << G4endl; 00177 } 00178 00179 //G4cout << "G4NeutronHPFissionData::DumpPhysicsTable still to be implemented"<<G4endl; 00180 }
G4double G4NeutronHPFissionData::GetCrossSection | ( | const G4DynamicParticle * | , | |
const G4Element * | , | |||
G4double | aT | |||
) |
Definition at line 185 of file G4NeutronHPFissionData.cc.
References buffer, G4DynamicParticle::GetDefinition(), G4Element::GetIndex(), G4ReactionProduct::GetKineticEnergy(), G4DynamicParticle::GetKineticEnergy(), G4ReactionProduct::GetMass(), G4ReactionProduct::GetMomentum(), G4DynamicParticle::GetMomentum(), G4Element::GetN(), G4NucleiProperties::GetNuclearMass(), G4ParticleDefinition::GetPDGMass(), G4Nucleus::GetThermalNucleus(), G4Element::GetZ(), G4ReactionProduct::Lorentz(), and G4Neutron::Neutron().
Referenced by GetIsoCrossSection().
00186 { 00187 G4double result = 0; 00188 if(anE->GetZ()<90) return result; 00189 G4bool outOfRange; 00190 G4int index = anE->GetIndex(); 00191 00192 // 100729 TK add safety 00193 if ( ( ( *theCrossSections )( index ) )->GetVectorLength() == 0 ) return result; 00194 00195 // prepare neutron 00196 G4double eKinetic = aP->GetKineticEnergy(); 00197 G4ReactionProduct theNeutron( aP->GetDefinition() ); 00198 theNeutron.SetMomentum( aP->GetMomentum() ); 00199 theNeutron.SetKineticEnergy( eKinetic ); 00200 00201 // prepare thermal nucleus 00202 G4Nucleus aNuc; 00203 G4double eps = 0.0001; 00204 G4double theA = anE->GetN(); 00205 G4double theZ = anE->GetZ(); 00206 G4double eleMass; 00207 eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) ) 00208 ) / G4Neutron::Neutron()->GetPDGMass(); 00209 00210 G4ReactionProduct boosted; 00211 G4double aXsection; 00212 00213 // MC integration loop 00214 G4int counter = 0; 00215 G4double buffer = 0; 00216 G4int size = G4int(std::max(10., aT/60*kelvin)); 00217 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum(); 00218 G4double neutronVMag = neutronVelocity.mag(); 00219 00220 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer) 00221 { 00222 if(counter) buffer = result/counter; 00223 while (counter<size) 00224 { 00225 counter ++; 00226 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT); 00227 boosted.Lorentz(theNeutron, aThermalNuc); 00228 G4double theEkin = boosted.GetKineticEnergy(); 00229 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange); 00230 // velocity correction. 00231 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum(); 00232 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag; 00233 result += aXsection; 00234 } 00235 size += size; 00236 } 00237 result /= counter; 00238 return result; 00239 }
G4double G4NeutronHPFissionData::GetIsoCrossSection | ( | const G4DynamicParticle * | , | |
G4int | , | |||
G4int | , | |||
const G4Isotope * | , | |||
const G4Element * | , | |||
const G4Material * | ||||
) | [virtual] |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 76 of file G4NeutronHPFissionData.cc.
References GetCrossSection(), and G4DynamicParticle::GetKineticEnergy().
00081 { 00082 if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache; 00083 00084 ke_cache = dp->GetKineticEnergy(); 00085 element_cache = element; 00086 material_cache = material; 00087 G4double xs = GetCrossSection( dp , element , material->GetTemperature() ); 00088 xs_cache = xs; 00089 return xs; 00090 }
G4bool G4NeutronHPFissionData::IsIsoApplicable | ( | const G4DynamicParticle * | , | |
G4int | , | |||
G4int | , | |||
const G4Element * | , | |||
const G4Material * | ||||
) | [virtual] |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 63 of file G4NeutronHPFissionData.cc.
References G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4VCrossSectionDataSet::GetMaxKinEnergy(), G4VCrossSectionDataSet::GetMinKinEnergy(), and G4Neutron::Neutron().
00067 { 00068 G4double eKin = dp->GetKineticEnergy(); 00069 if ( eKin > GetMaxKinEnergy() 00070 || eKin < GetMinKinEnergy() 00071 || dp->GetDefinition() != G4Neutron::Neutron() ) return false; 00072 00073 return true; 00074 }