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 "G4HadronicException.hh"
00041 #include "G4NeutronElasticXS.hh"
00042 #include "G4Neutron.hh"
00043 #include "G4DynamicParticle.hh"
00044 #include "G4Element.hh"
00045 #include "G4ElementTable.hh"
00046 #include "G4PhysicsLogVector.hh"
00047 #include "G4PhysicsVector.hh"
00048 #include "G4GlauberGribovCrossSection.hh"
00049 #include "G4HadronNucleonXsc.hh"
00050 #include "G4NistManager.hh"
00051 #include "G4Proton.hh"
00052
00053 #include <iostream>
00054 #include <fstream>
00055 #include <sstream>
00056
00057 using namespace std;
00058
00059 G4NeutronElasticXS::G4NeutronElasticXS()
00060 : G4VCrossSectionDataSet("G4NeutronElasticXS"),
00061 proton(G4Proton::Proton()), maxZ(92)
00062 {
00063
00064 if(verboseLevel > 0){
00065 G4cout << "G4NeutronElasticXS::G4NeutronElasticXS Initialise for Z < "
00066 << maxZ + 1 << G4endl;
00067 }
00068 data.resize(maxZ+1, 0);
00069 coeff.resize(maxZ+1, 1.0);
00070 ggXsection = new G4GlauberGribovCrossSection();
00071 fNucleon = new G4HadronNucleonXsc();
00072 isInitialized = false;
00073 }
00074
00075 G4NeutronElasticXS::~G4NeutronElasticXS()
00076 {
00077 delete fNucleon;
00078 for(G4int i=0; i<=maxZ; ++i) {
00079 delete data[i];
00080 }
00081 }
00082
00083 void G4NeutronElasticXS::CrossSectionDescription(std::ostream& outFile) const
00084 {
00085 outFile << "G4NeutronElasticXS calculates the neutron elastic scattering\n"
00086 << "cross section on nuclei using data from the high precision\n"
00087 << "neutron database. These data are simplified and smoothed over\n"
00088 << "the resonance region in order to reduce CPU time.\n"
00089 << "G4NeutronElasticXS is valid for energies up to 20 MeV, for all\n"
00090 << "targets through U.\n";
00091 }
00092
00093 G4bool
00094 G4NeutronElasticXS::IsElementApplicable(const G4DynamicParticle*,
00095 G4int, const G4Material*)
00096 {
00097 return true;
00098 }
00099
00100 G4double
00101 G4NeutronElasticXS::GetElementCrossSection(const G4DynamicParticle* aParticle,
00102 G4int Z, const G4Material*)
00103 {
00104 G4double xs = 0.0;
00105 G4double ekin = aParticle->GetKineticEnergy();
00106
00107 if(Z < 1 || Z > maxZ) { return xs; }
00108
00109 G4int Amean = G4int(G4NistManager::Instance()->GetAtomicMassAmu(Z)+0.5);
00110 G4PhysicsVector* pv = data[Z];
00111
00112
00113
00114 if(!pv) {
00115 Initialise(Z);
00116 pv = data[Z];
00117 if(!pv) { return xs; }
00118 }
00119
00120 G4double e1 = pv->Energy(0);
00121 if(ekin <= e1) { return (*pv)[0]; }
00122
00123 G4int n = pv->GetVectorLength() - 1;
00124 G4double e2 = pv->Energy(n);
00125
00126 if(ekin <= e2) {
00127 xs = pv->Value(ekin);
00128 } else if(1 == Z) {
00129 fNucleon->GetHadronNucleonXscPDG(aParticle, proton);
00130 xs = coeff[1]*fNucleon->GetElasticHadronNucleonXsc();
00131 } else {
00132 ggXsection->GetIsoCrossSection(aParticle, Z, Amean);
00133 xs = coeff[Z]*ggXsection->GetElasticGlauberGribovXsc();
00134 }
00135
00136 if(verboseLevel > 0){
00137 G4cout << "ekin= " << ekin << ", XSinel= " << xs << G4endl;
00138 }
00139 return xs;
00140 }
00141
00142 void
00143 G4NeutronElasticXS::BuildPhysicsTable(const G4ParticleDefinition& p)
00144 {
00145 if(isInitialized) { return; }
00146 if(verboseLevel > 0){
00147 G4cout << "G4NeutronElasticXS::BuildPhysicsTable for "
00148 << p.GetParticleName() << G4endl;
00149 }
00150 if(p.GetParticleName() != "neutron") {
00151 throw G4HadronicException(__FILE__, __LINE__,"Wrong particle type");
00152 return;
00153 }
00154 isInitialized = true;
00155
00156
00157
00158 char* path = getenv("G4NEUTRONXSDATA");
00159 if (!path){
00160 throw G4HadronicException(__FILE__, __LINE__,
00161 "G4NEUTRONXSDATA environment variable not defined");
00162 return;
00163 }
00164
00165 G4DynamicParticle* dynParticle =
00166 new G4DynamicParticle(G4Neutron::Neutron(),G4ThreeVector(1,0,0),1);
00167
00168
00169 const G4ElementTable* theElmTable = G4Element::GetElementTable();
00170 size_t numOfElm = G4Element::GetNumberOfElements();
00171 if(numOfElm > 0) {
00172 for(size_t i=0; i<numOfElm; ++i) {
00173 G4int Z = G4int(((*theElmTable)[i])->GetZ());
00174 if(Z < 1) { Z = 1; }
00175 else if(Z > maxZ) { Z = maxZ; }
00176
00177
00178 if(!data[Z]) { Initialise(Z, dynParticle, path); }
00179 }
00180 }
00181 delete dynParticle;
00182 }
00183
00184 void
00185 G4NeutronElasticXS::Initialise(G4int Z, G4DynamicParticle* dp,
00186 const char* p)
00187 {
00188 if(data[Z]) { return; }
00189 const char* path = p;
00190 if(!p) {
00191
00192
00193 path = getenv("G4NEUTRONXSDATA");
00194 if (!path) {
00195 throw G4HadronicException(__FILE__, __LINE__,
00196 "G4NEUTRONXSDATA environment variable not defined");
00197 return;
00198 }
00199 }
00200 G4DynamicParticle* dynParticle = dp;
00201 if(!dp) {
00202 dynParticle =
00203 new G4DynamicParticle(G4Neutron::Neutron(),G4ThreeVector(1,0,0),1);
00204 }
00205
00206 G4int Amean = G4int(G4NistManager::Instance()->GetAtomicMassAmu(Z)+0.5);
00207
00208
00209 data[Z] = new G4PhysicsLogVector();
00210
00211 std::ostringstream ost;
00212 ost << path << "/elast" << Z ;
00213 std::ifstream filein(ost.str().c_str());
00214 if (!(filein)) {
00215 G4cout << ost.str()
00216 << " is not opened by G4NeutronElasticXS" << G4endl;
00217 throw G4HadronicException(__FILE__, __LINE__,"NO data sets opened");
00218 return;
00219 }else{
00220 if(verboseLevel > 1) {
00221 G4cout << "file " << ost.str()
00222 << " is opened by G4NeutronElasticXS" << G4endl;
00223 }
00224
00225
00226 data[Z]->Retrieve(filein, true);
00227
00228
00229 size_t n = data[Z]->GetVectorLength() - 1;
00230 G4double emax = data[Z]->Energy(n);
00231 G4double sig1 = (*data[Z])[n];
00232 dynParticle->SetKineticEnergy(emax);
00233 G4double sig2 = 0.0;
00234 if(1 == Z) {
00235 fNucleon->GetHadronNucleonXscPDG(dynParticle, proton);
00236 sig2 = fNucleon->GetElasticHadronNucleonXsc();
00237 } else {
00238 ggXsection->GetIsoCrossSection(dynParticle, Z, Amean);
00239 sig2 = ggXsection->GetElasticGlauberGribovXsc();
00240 }
00241 if(sig2 > 0.) { coeff[Z] = sig1/sig2; }
00242 }
00243 if(!dp) { delete dynParticle; }
00244 }