#include <G4NeutronHPCaptureData.hh>
Inheritance diagram for G4NeutronHPCaptureData:
Public Member Functions | |
G4NeutronHPCaptureData () | |
~G4NeutronHPCaptureData () | |
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 &) |
void | IgnoreOnFlightDopplerBroadening () |
void | EnableOnFlightDopplerBroadening () |
Definition at line 50 of file G4NeutronHPCaptureData.hh.
G4NeutronHPCaptureData::G4NeutronHPCaptureData | ( | ) |
Definition at line 44 of file G4NeutronHPCaptureData.cc.
References BuildPhysicsTable(), G4Neutron::Neutron(), G4VCrossSectionDataSet::SetMaxKinEnergy(), and G4VCrossSectionDataSet::SetMinKinEnergy().
00045 :G4VCrossSectionDataSet("NeutronHPCaptureXS") 00046 { 00047 SetMinKinEnergy( 0*MeV ); 00048 SetMaxKinEnergy( 20*MeV ); 00049 00050 ke_cache = 0.0; 00051 xs_cache = 0.0; 00052 element_cache = NULL; 00053 material_cache = NULL; 00054 00055 theCrossSections = 0; 00056 onFlightDB = true; 00057 00058 BuildPhysicsTable(*G4Neutron::Neutron()); 00059 }
G4NeutronHPCaptureData::~G4NeutronHPCaptureData | ( | ) |
Definition at line 61 of file G4NeutronHPCaptureData.cc.
References G4PhysicsTable::clearAndDestroy().
00062 { 00063 if ( theCrossSections != 0 ) theCrossSections->clearAndDestroy(); 00064 00065 delete theCrossSections; 00066 }
void G4NeutronHPCaptureData::BuildPhysicsTable | ( | const G4ParticleDefinition & | ) | [virtual] |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 107 of file G4NeutronHPCaptureData.cc.
References G4PhysicsTable::clearAndDestroy(), G4cout, G4endl, G4Element::GetElementTable(), G4Element::GetNumberOfElements(), G4NeutronHPData::Instance(), G4NeutronHPData::MakePhysicsVector(), G4Neutron::Neutron(), and G4PhysicsTable::push_back().
Referenced by G4NeutronHPCaptureData().
00108 { 00109 if(&aP!=G4Neutron::Neutron()) 00110 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!"); 00111 00112 //080428 00113 if ( getenv( "G4NEUTRONHP_NEGLECT_DOPPLER" ) ) 00114 { 00115 G4cout << "Find environment variable of \"G4NEUTRONHP_NEGLECT_DOPPLER\"." << G4endl; 00116 G4cout << "On the fly Doppler broadening will be neglect in the cross section calculation of capture reaction of neutrons (<20MeV)." << G4endl; 00117 onFlightDB = false; 00118 } 00119 00120 size_t numberOfElements = G4Element::GetNumberOfElements(); 00121 // G4cout << "CALLED G4NeutronHPCaptureData::BuildPhysicsTable "<<numberOfElements<<G4endl; 00122 // TKDB 00123 //if ( theCrossSections == 0 ) theCrossSections = new G4PhysicsTable( numberOfElements ); 00124 if ( theCrossSections == NULL ) 00125 theCrossSections = new G4PhysicsTable( numberOfElements ); 00126 else 00127 theCrossSections->clearAndDestroy(); 00128 00129 // make a PhysicsVector for each element 00130 00131 static const G4ElementTable *theElementTable = G4Element::GetElementTable(); 00132 for( size_t i=0; i<numberOfElements; ++i ) 00133 { 00134 if(getenv("CaptureDataIndexDebug")) 00135 { 00136 G4int index_debug = ((*theElementTable)[i])->GetIndex(); 00137 G4cout << "IndexDebug "<< i <<" "<<index_debug<<G4endl; 00138 } 00139 G4PhysicsVector* physVec = G4NeutronHPData:: 00140 Instance()->MakePhysicsVector((*theElementTable)[i], this); 00141 theCrossSections->push_back(physVec); 00142 } 00143 }
void G4NeutronHPCaptureData::DumpPhysicsTable | ( | const G4ParticleDefinition & | ) | [virtual] |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 145 of file G4NeutronHPCaptureData.cc.
References G4cout, G4endl, G4Element::GetElementTable(), G4VCrossSectionDataSet::GetName(), G4Element::GetNumberOfElements(), and G4Neutron::Neutron().
00146 { 00147 if(&aP!=G4Neutron::Neutron()) 00148 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!"); 00149 00150 // 00151 // Dump element based cross section 00152 // range 10e-5 eV to 20 MeV 00153 // 10 point per decade 00154 // in barn 00155 // 00156 00157 G4cout << G4endl; 00158 G4cout << G4endl; 00159 G4cout << "Capture Cross Section of Neutron HP"<< G4endl; 00160 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl; 00161 G4cout << G4endl; 00162 G4cout << "Name of Element" << G4endl; 00163 G4cout << "Energy[eV] XS[barn]" << G4endl; 00164 G4cout << G4endl; 00165 00166 size_t numberOfElements = G4Element::GetNumberOfElements(); 00167 static const G4ElementTable *theElementTable = G4Element::GetElementTable(); 00168 00169 for ( size_t i = 0 ; i < numberOfElements ; ++i ) 00170 { 00171 00172 G4cout << (*theElementTable)[i]->GetName() << G4endl; 00173 00174 G4int ie = 0; 00175 00176 for ( ie = 0 ; ie < 130 ; ie++ ) 00177 { 00178 G4double eKinetic = 1.0e-5 * std::pow ( 10.0 , ie/10.0 ) *eV; 00179 G4bool outOfRange = false; 00180 00181 if ( eKinetic < 20*MeV ) 00182 { 00183 G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl; 00184 } 00185 00186 } 00187 00188 G4cout << G4endl; 00189 } 00190 00191 00192 // G4cout << "G4NeutronHPCaptureData::DumpPhysicsTable still to be implemented"<<G4endl; 00193 }
void G4NeutronHPCaptureData::EnableOnFlightDopplerBroadening | ( | ) | [inline] |
G4double G4NeutronHPCaptureData::GetCrossSection | ( | const G4DynamicParticle * | , | |
const G4Element * | , | |||
G4double | aT | |||
) |
Definition at line 198 of file G4NeutronHPCaptureData.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().
00199 { 00200 G4double result = 0; 00201 G4bool outOfRange; 00202 G4int index = anE->GetIndex(); 00203 00204 // prepare neutron 00205 G4double eKinetic = aP->GetKineticEnergy(); 00206 00207 //if ( getenv( "G4NEUTRONHP_NEGLECT_DOPPLER" ) ) 00208 //080428 00209 if ( !onFlightDB ) 00210 { 00211 G4double factor = 1.0; 00212 if ( eKinetic < aT * k_Boltzmann ) 00213 { 00214 // below 0.1 eV neutrons 00215 // Have to do some, but now just igonre. 00216 // Will take care after performance check. 00217 // factor = factor * targetV; 00218 } 00219 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor; 00220 } 00221 00222 G4ReactionProduct theNeutron( aP->GetDefinition() ); 00223 theNeutron.SetMomentum( aP->GetMomentum() ); 00224 theNeutron.SetKineticEnergy( eKinetic ); 00225 00226 // prepare thermal nucleus 00227 G4Nucleus aNuc; 00228 G4double eps = 0.0001; 00229 G4double theA = anE->GetN(); 00230 G4double theZ = anE->GetZ(); 00231 G4double eleMass; 00232 eleMass = G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) ) / G4Neutron::Neutron()->GetPDGMass(); 00233 00234 G4ReactionProduct boosted; 00235 G4double aXsection; 00236 00237 // MC integration loop 00238 G4int counter = 0; 00239 G4double buffer = 0; 00240 G4int size = G4int(std::max(10., aT/60*kelvin)); 00241 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum(); 00242 G4double neutronVMag = neutronVelocity.mag(); 00243 00244 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer) 00245 { 00246 if(counter) buffer = result/counter; 00247 while (counter<size) 00248 { 00249 counter ++; 00250 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT); 00251 boosted.Lorentz(theNeutron, aThermalNuc); 00252 G4double theEkin = boosted.GetKineticEnergy(); 00253 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange); 00254 // velocity correction, or luminosity factor... 00255 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum(); 00256 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag; 00257 result += aXsection; 00258 } 00259 size += size; 00260 } 00261 result /= counter; 00262 /* 00263 // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER 00264 G4cout << " result " << result << " " 00265 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " " 00266 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl; 00267 */ 00268 return result; 00269 }
G4double G4NeutronHPCaptureData::GetIsoCrossSection | ( | const G4DynamicParticle * | , | |
G4int | , | |||
G4int | , | |||
const G4Isotope * | , | |||
const G4Element * | , | |||
const G4Material * | ||||
) | [virtual] |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 81 of file G4NeutronHPCaptureData.cc.
References GetCrossSection(), and G4DynamicParticle::GetKineticEnergy().
00086 { 00087 if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache; 00088 00089 ke_cache = dp->GetKineticEnergy(); 00090 element_cache = element; 00091 material_cache = material; 00092 G4double xs = GetCrossSection( dp , element , material->GetTemperature() ); 00093 xs_cache = xs; 00094 return xs; 00095 }
void G4NeutronHPCaptureData::IgnoreOnFlightDopplerBroadening | ( | ) | [inline] |
G4bool G4NeutronHPCaptureData::IsIsoApplicable | ( | const G4DynamicParticle * | , | |
G4int | , | |||
G4int | , | |||
const G4Element * | , | |||
const G4Material * | ||||
) | [virtual] |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 68 of file G4NeutronHPCaptureData.cc.
References G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4VCrossSectionDataSet::GetMaxKinEnergy(), G4VCrossSectionDataSet::GetMinKinEnergy(), and G4Neutron::Neutron().
00072 { 00073 G4double eKin = dp->GetKineticEnergy(); 00074 if ( eKin > GetMaxKinEnergy() 00075 || eKin < GetMinKinEnergy() 00076 || dp->GetDefinition() != G4Neutron::Neutron() ) return false; 00077 00078 return true; 00079 }