Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/cross_sections/src/G4NeutronElasticXS.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  1 //
  2 // ********************************************************************
  3 // * License and Disclaimer                                           *
  4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.                             *
 10 // *                                                                  *
 11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                                                  *
 18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // ********************************************************************
 25 //
 26 // -------------------------------------------------------------------
 27 //
 28 // GEANT4 Class file
 29 //
 30 //
 31 // File name:    G4NeutronElasticXS
 32 //
 33 // Author  Ivantchenko, Geant4, 3-Aug-09
 34 //
 35 // Modifications:
 36 //
 37 
 38 #include "G4NeutronElasticXS.hh"
 39 #include "G4Neutron.hh"
 40 #include "G4DynamicParticle.hh"
 41 #include "G4ElementTable.hh"
 42 #include "G4Material.hh"
 43 #include "G4Element.hh"
 44 #include "G4PhysicsLogVector.hh"
 45 #include "G4CrossSectionDataSetRegistry.hh"
 46 #include "G4ComponentGGHadronNucleusXsc.hh"
 47 #include "G4HadronicParameters.hh"
 48 #include "Randomize.hh"
 49 #include "G4SystemOfUnits.hh"
 50 #include "G4IsotopeList.hh"
 51 #include "G4AutoLock.hh"
 52 
 53 #include <fstream>
 54 #include <sstream>
 55 
 56 G4PhysicsVector* G4NeutronElasticXS::data[] = {nullptr};
 57 G4double G4NeutronElasticXS::coeff[] = {0.0};
 58 G4String G4NeutronElasticXS::gDataDirectory = "";
 59 G4bool G4NeutronElasticXS::fLock = true;
 60 
 61 namespace
 62 {
 63   G4Mutex nElasticXSMutex = G4MUTEX_INITIALIZER;
 64 }
 65 
 66 G4NeutronElasticXS::G4NeutronElasticXS() 
 67  : G4VCrossSectionDataSet(Default_Name()),
 68    neutron(G4Neutron::Neutron())
 69 {
 70   //  verboseLevel = 0;
 71   if (verboseLevel > 0){
 72     G4cout  << "G4NeutronElasticXS::G4NeutronElasticXS Initialise for Z < " 
 73       << MAXZEL << G4endl;
 74   }
 75   ggXsection = 
 76     G4CrossSectionDataSetRegistry::Instance()->GetComponentCrossSection("Glauber-Gribov");
 77   if (ggXsection == nullptr)
 78     ggXsection = new G4ComponentGGHadronNucleusXsc();
 79   SetForAllAtomsAndEnergies(true);
 80   FindDirectoryPath();
 81 }
 82 
 83 G4NeutronElasticXS::~G4NeutronElasticXS()
 84 {
 85   if (isFirst) {
 86     for(G4int i=0; i<MAXZEL; ++i) {
 87       delete data[i];
 88       data[i] = nullptr;
 89     }
 90   }
 91 }
 92 
 93 void G4NeutronElasticXS::CrossSectionDescription(std::ostream& outFile) const
 94 {
 95   outFile << "G4NeutronElasticXS calculates the neutron elastic scattering\n"
 96           << "cross section on nuclei using data from the high precision\n"
 97           << "neutron database.  These data are simplified and smoothed over\n"
 98           << "the resonance region in order to reduce CPU time.\n"
 99           << "For high energies Glauber-Gribiv cross section is used.\n";
100 }
101 
102 G4bool 
103 G4NeutronElasticXS::IsElementApplicable(const G4DynamicParticle*, 
104           G4int, const G4Material*)
105 {
106   return true;
107 }
108 
109 G4bool G4NeutronElasticXS::IsIsoApplicable(const G4DynamicParticle*,
110                                            G4int, G4int,
111                                            const G4Element*, const G4Material*)
112 {
113   return false;
114 }
115 
116 G4double 
117 G4NeutronElasticXS::GetElementCrossSection(const G4DynamicParticle* aParticle,
118              G4int Z, const G4Material*)
119 {
120   return ElementCrossSection(aParticle->GetKineticEnergy(),
121            aParticle->GetLogKineticEnergy(), Z);
122 }
123 
124 G4double
125 G4NeutronElasticXS::ComputeCrossSectionPerElement(G4double ekin, G4double loge,
126               const G4ParticleDefinition*,
127               const G4Element* elm,
128               const G4Material*)
129 {
130   return ElementCrossSection(ekin, loge, elm->GetZasInt());
131 }
132 
133 G4double G4NeutronElasticXS::ElementCrossSection(G4double ekin, G4double loge, G4int ZZ)
134 {
135   G4int Z = (ZZ >= MAXZEL) ? MAXZEL - 1 : ZZ;
136   auto pv = GetPhysicsVector(Z);
137 
138   G4double xs = (ekin <= pv->GetMaxEnergy()) ? pv->LogVectorValue(ekin, loge) 
139     : coeff[Z]*ggXsection->GetElasticElementCrossSection(neutron, ekin,
140                                                          Z, aeff[Z]);
141 
142 #ifdef G4VERBOSE
143   if(verboseLevel > 1) {
144     G4cout  << "Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV 
145       << ",  nElmXSel(b)= " << xs/CLHEP::barn 
146       << G4endl;
147   }
148 #endif
149   return xs;
150 }
151 
152 G4double
153 G4NeutronElasticXS::ComputeIsoCrossSection(G4double ekin, G4double loge,
154                    const G4ParticleDefinition*,
155                    G4int Z, G4int A,
156                    const G4Isotope*, const G4Element*,
157                    const G4Material*)
158 {
159   return ElementCrossSection(ekin, loge, Z)*A/aeff[Z];
160 }
161 
162 G4double
163 G4NeutronElasticXS::GetIsoCrossSection(const G4DynamicParticle* aParticle, 
164                G4int Z, G4int A,
165                const G4Isotope*, const G4Element*,
166                const G4Material*)
167 {
168   return ElementCrossSection(aParticle->GetKineticEnergy(),
169            aParticle->GetLogKineticEnergy(), Z)*A/aeff[Z];
170 
171 }
172 
173 const G4Isotope* G4NeutronElasticXS::SelectIsotope(
174       const G4Element* anElement, G4double, G4double)
175 {
176   G4int nIso = (G4int)anElement->GetNumberOfIsotopes();
177   const G4Isotope* iso = anElement->GetIsotope(0);
178 
179   //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
180   if(1 == nIso) { return iso; }
181 
182   const G4double* abundVector = anElement->GetRelativeAbundanceVector();
183   G4double q = G4UniformRand();
184   G4double sum = 0.0;
185 
186   // isotope wise cross section not used
187   for (G4int j=0; j<nIso; ++j) {
188     sum += abundVector[j];
189     if(q <= sum) {
190       iso = anElement->GetIsotope(j);
191       break;
192     }
193   }
194   return iso;
195 }
196 
197 void 
198 G4NeutronElasticXS::BuildPhysicsTable(const G4ParticleDefinition& p)
199 {
200   if(verboseLevel > 0){
201     G4cout << "G4NeutronElasticXS::BuildPhysicsTable for " 
202      << p.GetParticleName() << G4endl;
203   }
204   if(p.GetParticleName() != "neutron") { 
205     G4ExceptionDescription ed;
206     ed << p.GetParticleName() << " is a wrong particle type -"
207        << " only neutron is allowed";
208     G4Exception("G4NeutronElasticXS::BuildPhysicsTable(..)","had012",
209     FatalException, ed, "");
210     return; 
211   }
212   if (fLock || isFirst) { 
213     G4AutoLock l(&nElasticXSMutex);
214     if (fLock) { 
215       isFirst = true;
216       fLock = false;
217       FindDirectoryPath();
218     }
219 
220     // Access to elements
221     const G4ElementTable* table = G4Element::GetElementTable();
222     for ( auto & elm : *table ) {
223       G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZEL-1) );
224       if ( nullptr == data[Z] ) { Initialise(Z); }
225     }
226     l.unlock();
227   }
228 }
229 
230 const G4String& G4NeutronElasticXS::FindDirectoryPath()
231 {
232   // build the complete string identifying the file with the data set
233   if (gDataDirectory.empty()) {
234     std::ostringstream ost;
235     ost << G4HadronicParameters::Instance()->GetDirPARTICLEXS() << "/neutron/el";
236     gDataDirectory = ost.str();
237   }
238   return gDataDirectory;
239 }
240 
241 void G4NeutronElasticXS::InitialiseOnFly(G4int Z)
242 {
243   G4AutoLock l(&nElasticXSMutex);
244   Initialise(Z);
245   l.unlock();
246 }
247 
248 void G4NeutronElasticXS::Initialise(G4int Z)
249 {
250   if(data[Z] != nullptr) { return; }
251 
252   // upload data from file
253   data[Z] = new G4PhysicsLogVector();
254 
255   std::ostringstream ost;
256   ost << FindDirectoryPath() << Z ;
257   std::ifstream filein(ost.str().c_str());
258   if (!filein.is_open()) {
259     G4ExceptionDescription ed;
260     ed << "Data file <" << ost.str().c_str()
261        << "> is not opened!";
262     G4Exception("G4NeutronElasticXS::Initialise(..)","had014",
263                 FatalException, ed, "Check G4PARTICLEXSDATA");
264     return;
265   }
266   if(verboseLevel > 1) {
267     G4cout << "file " << ost.str() 
268      << " is opened by G4NeutronElasticXS" << G4endl;
269   }
270     
271   // retrieve data from DB
272   if(!data[Z]->Retrieve(filein, true)) {
273     G4ExceptionDescription ed;
274     ed << "Data file <" << ost.str().c_str()
275        << "> is not retrieved!";
276     G4Exception("G4NeutronElasticXS::Initialise(..)","had015",
277     FatalException, ed, "Check G4PARTICLEXSDATA");
278     return;
279   }
280   // smooth transition 
281   G4double sig1  = (*(data[Z]))[data[Z]->GetVectorLength()-1];
282   G4double ehigh = data[Z]->GetMaxEnergy();
283   G4double sig2  = ggXsection->GetElasticElementCrossSection(neutron, 
284                                ehigh, Z, aeff[Z]);
285   coeff[Z] = (sig2 > 0.) ? sig1/sig2 : 1.0;  
286 }
287