00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040 #include <fstream>
00041 #include <sstream>
00042
00043 #include "G4HadronicException.hh"
00044 #include "G4SystemOfUnits.hh"
00045 #include "G4NeutronCaptureXS.hh"
00046 #include "G4Element.hh"
00047 #include "G4ElementTable.hh"
00048 #include "G4PhysicsLogVector.hh"
00049 #include "G4PhysicsVector.hh"
00050 #include "G4DynamicParticle.hh"
00051 #include "Randomize.hh"
00052
00053 using namespace std;
00054
00055 const G4int G4NeutronCaptureXS::amin[] = {0,
00056 0, 0, 6, 0,10,12,14,16, 0, 0,
00057 0, 0, 0,28, 0, 0, 0,36, 0,40,
00058 0, 0, 0, 0, 0,54, 0,58,63,64,
00059 0,70, 0, 0, 0, 0, 0, 0, 0,90,
00060 0, 0, 0, 0, 0, 0,107,106, 0,112,
00061 0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
00062 0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
00063 0, 0, 0,180, 0, 0, 0, 0, 0, 0,
00064 0,204, 0, 0, 0, 0, 0, 0, 0, 0,
00065 0,235};
00066 const G4int G4NeutronCaptureXS::amax[] = {0,
00067 0, 0, 7, 0,11,13,15,18, 0, 0,
00068 0, 0, 0,30, 0, 0, 0,40, 0,48,
00069 0, 0, 0, 0, 0,58, 0,64,65,70,
00070 0,76, 0, 0, 0, 0, 0, 0, 0,96,
00071 0, 0, 0, 0, 0, 0,109,116, 0,124,
00072 0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
00073 0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
00074 0, 0, 0,186, 0, 0, 0, 0, 0, 0,
00075 0,208, 0, 0, 0, 0, 0, 0, 0, 0,
00076 0,238};
00077
00078 G4NeutronCaptureXS::G4NeutronCaptureXS()
00079 : G4VCrossSectionDataSet("G4NeutronCaptureXS"),
00080 emax(20*MeV),maxZ(92)
00081 {
00082
00083 if(verboseLevel > 0){
00084 G4cout << "G4NeutronCaptureXS::G4NeutronCaptureXS: Initialise for Z < "
00085 << maxZ + 1 << G4endl;
00086 }
00087
00088 data.SetName("NeutronCapture");
00089 work.resize(13,0);
00090 temp.resize(13,0.0);
00091 isInitialized = false;
00092 }
00093
00094 G4NeutronCaptureXS::~G4NeutronCaptureXS()
00095 {
00096
00097
00098
00099
00100
00101 }
00102
00103 void G4NeutronCaptureXS::CrossSectionDescription(std::ostream& outFile) const
00104 {
00105 outFile << "G4NeutronCaptureXS calculates the neutron capture cross sections\n"
00106 << "on nuclei using data from the high precision neutron database.\n"
00107 << "These data are simplified and smoothed over the resonance region\n"
00108 << "in order to reduce CPU time. G4NeutronCaptureXS is valid up to\n"
00109 << "20 MeV for all targets through U.\n";
00110 }
00111
00112 G4bool
00113 G4NeutronCaptureXS::IsElementApplicable(const G4DynamicParticle*,
00114 G4int, const G4Material*)
00115 {
00116 return true;
00117 }
00118
00119 G4bool
00120 G4NeutronCaptureXS::IsIsoApplicable(const G4DynamicParticle*,
00121 G4int , G4int ,
00122 const G4Element*, const G4Material*)
00123 {
00124 return false;
00125 }
00126
00127 G4double
00128 G4NeutronCaptureXS::GetElementCrossSection(const G4DynamicParticle* aParticle,
00129 G4int Z, const G4Material*)
00130 {
00131 G4double xs = 0.0;
00132 G4double ekin = aParticle->GetKineticEnergy();
00133 if(ekin > emax || Z < 1 || Z > maxZ) { return xs; }
00134 const G4double elimit = 1.0e-10*eV;
00135 if(ekin < elimit) { ekin = elimit; }
00136
00137
00138 G4PhysicsVector* pv = data.GetElementData(Z);
00139
00140
00141 if(!pv) {
00142 Initialise(Z);
00143
00144 pv = data.GetElementData(Z);
00145 if(!pv) { return xs; }
00146 }
00147
00148 G4double e1 = pv->Energy(0);
00149 if(ekin < e1) { xs = (*pv)[0]*std::sqrt(e1/ekin); }
00150 else { xs = pv->Value(ekin); }
00151
00152 if(verboseLevel > 0){
00153 G4cout << "ekin= " << ekin << ", xs= " << xs << G4endl;
00154 }
00155 return xs;
00156 }
00157
00158 G4double
00159 G4NeutronCaptureXS::GetIsoCrossSection(const G4DynamicParticle* aParticle,
00160 G4int Z, G4int A,
00161 const G4Isotope*, const G4Element*,
00162 const G4Material*)
00163 {
00164 G4double xs = 0.0;
00165 G4double ekin = aParticle->GetKineticEnergy();
00166 if(ekin > emax || Z < 1 || Z > maxZ) { return xs; }
00167 const G4double elimit = 1.0e-10*eV;
00168 if(ekin < elimit) { ekin = elimit; }
00169
00170
00171 G4PhysicsVector* pv = data.GetElementData(Z);
00172
00173
00174 if(!pv) {
00175 Initialise(Z);
00176
00177 pv = data.GetElementData(Z);
00178 if(!pv) { return xs; }
00179 }
00180 pv = data.GetComponentDataByID(Z, A);
00181 if(!pv) { return xs; }
00182
00183 G4double e1 = pv->Energy(0);
00184 if(ekin < e1) { xs = (*pv)[0]*std::sqrt(e1/ekin); }
00185 else { xs = pv->Value(ekin); }
00186
00187 if(verboseLevel > 0){
00188 G4cout << "ekin= " << ekin << ", xs= " << xs << G4endl;
00189 }
00190 return xs;
00191 }
00192
00193 G4Isotope* G4NeutronCaptureXS::SelectIsotope(const G4Element* anElement,
00194 G4double kinEnergy)
00195 {
00196 G4int nIso = anElement->GetNumberOfIsotopes();
00197 G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
00198 G4Isotope* iso = (*isoVector)[0];
00199
00200
00201 if(1 < nIso) {
00202 G4int Z = G4lrint(anElement->GetZ());
00203 if(Z > maxZ) { Z = maxZ; }
00204 G4double* abundVector = anElement->GetRelativeAbundanceVector();
00205 G4double q = G4UniformRand();
00206 G4double sum = 0.0;
00207
00208
00209 if(0 == amin[Z]) {
00210 for (G4int j = 0; j<nIso; ++j) {
00211 sum += abundVector[j];
00212 if(q <= sum) {
00213 iso = (*isoVector)[j];
00214 break;
00215 }
00216 }
00217 } else {
00218 size_t nmax = data.GetNumberOfComponents(Z);
00219 if(temp.size() < nmax) { temp.resize(nmax,0.0); }
00220 for (size_t i=0; i<nmax; ++i) {
00221 G4int A = (*isoVector)[i]->GetN();
00222 G4PhysicsVector* v = data.GetComponentDataByID(Z, A);
00223 if(v) { sum += abundVector[i]*v->Value(kinEnergy); }
00224 temp[i] = sum;
00225 }
00226 sum *= q;
00227 for (size_t j = 0; j<nmax; ++j) {
00228 if(temp[j] >= sum) {
00229 iso = (*isoVector)[j];
00230 break;
00231 }
00232 }
00233 }
00234 }
00235 return iso;
00236 }
00237
00238 void
00239 G4NeutronCaptureXS::BuildPhysicsTable(const G4ParticleDefinition& p)
00240 {
00241 if(isInitialized) { return; }
00242 if(verboseLevel > 0){
00243 G4cout << "G4NeutronCaptureXS::BuildPhysicsTable for "
00244 << p.GetParticleName() << G4endl;
00245 }
00246 if(p.GetParticleName() != "neutron") {
00247 throw G4HadronicException(__FILE__, __LINE__,"Wrong particle type");
00248 return;
00249 }
00250 isInitialized = true;
00251
00252
00253
00254 char* path = getenv("G4NEUTRONXSDATA");
00255 if (!path){
00256 throw G4HadronicException(__FILE__, __LINE__,
00257 "G4NEUTRONXSDATA environment variable not defined");
00258 return;
00259 }
00260
00261
00262 const G4ElementTable* theElmTable = G4Element::GetElementTable();
00263 size_t numOfElm = G4Element::GetNumberOfElements();
00264 if(numOfElm > 0) {
00265 for(size_t i=0; i<numOfElm; ++i) {
00266 G4int Z = G4int(((*theElmTable)[i])->GetZ());
00267 if(Z < 1) { Z = 1; }
00268 else if(Z > maxZ) { Z = maxZ; }
00269
00270
00271
00272 if(!data.GetElementData(Z)) { Initialise(Z, path); }
00273 }
00274 }
00275 }
00276
00277 void
00278 G4NeutronCaptureXS::Initialise(G4int Z, const char* p)
00279 {
00280 if(data.GetElementData(Z)) { return; }
00281 const char* path = p;
00282
00283
00284 if(!p) {
00285 path = getenv("G4NEUTRONXSDATA");
00286 if (!path) {
00287 throw G4HadronicException(__FILE__, __LINE__,
00288 "G4NEUTRONXSDATA environment variable not defined");
00289 return;
00290 }
00291 }
00292
00293
00294 std::ostringstream ost;
00295 ost << path << "/cap" << Z ;
00296 G4PhysicsVector* v = RetrieveVector(ost, true);
00297 data.InitialiseForElement(Z, v);
00298
00299
00300 if(amin[Z] > 0) {
00301 size_t n = 0;
00302 size_t i = 0;
00303 size_t nmax = (size_t)(amax[Z]-amin[Z]+1);
00304 if(work.size() < nmax) { work.resize(nmax,0); }
00305 for(G4int A=amin[Z]; A<=amax[Z]; ++A) {
00306 std::ostringstream ost1;
00307 ost1 << path << "/cap" << Z << "_" << A;
00308 v = RetrieveVector(ost1, false);
00309 if(v) { ++n; }
00310 work[i] = v;
00311 ++i;
00312 }
00313 data.InitialiseForComponent(Z, n);
00314 for(size_t j=0; j<i; ++j) {
00315 if(work[j]) { data.AddComponent(Z, amin[Z]+j, work[j]); }
00316 }
00317 }
00318 }
00319
00320 G4PhysicsVector*
00321 G4NeutronCaptureXS::RetrieveVector(std::ostringstream& ost, G4bool warn)
00322 {
00323 G4PhysicsLogVector* v = 0;
00324 std::ifstream filein(ost.str().c_str());
00325 if (!(filein)) {
00326 if(!warn) { return v; }
00327 G4cout << ost.str() << " is not opened by G4NeutronCaptureXS" << G4endl;
00328 throw G4HadronicException(__FILE__, __LINE__,"NO data sets opened");
00329 }else{
00330 if(verboseLevel > 1) {
00331 G4cout << "File " << ost.str()
00332 << " is opened by G4NeutronCaptureXS" << G4endl;
00333 }
00334
00335 v = new G4PhysicsLogVector();
00336 if(!v->Retrieve(filein, true)) {
00337 G4cout << ost.str() << " is not retrieved in G4NeutronCaptureXS" << G4endl;
00338 throw G4HadronicException(__FILE__, __LINE__,"ERROR: retrieve data fail");
00339 }
00340 }
00341 return v;
00342 }