Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/particle_hp/src/G4ParticleHPElastic.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 // neutron_hp -- source file
 27 // J.P. Wellisch, Nov-1996
 28 // A prototype of the low energy neutron transport model.
 29 //
 30 // 070523 bug fix for G4FPE_DEBUG on by A. Howard ( and T. Koi)
 31 // 080319 Compilation warnings - gcc-4.3.0 fix by T. Koi
 32 //
 33 // P. Arce, June-2014 Conversion neutron_hp to particle_hp
 34 //
 35 #include "G4ParticleHPElastic.hh"
 36 
 37 #include "G4ParticleHPElasticFS.hh"
 38 #include "G4ParticleHPManager.hh"
 39 #include "G4ParticleHPThermalBoost.hh"
 40 #include "G4SystemOfUnits.hh"
 41 #include "G4Threading.hh"
 42 
 43 G4ParticleHPElastic::G4ParticleHPElastic() : G4HadronicInteraction("NeutronHPElastic")
 44 {
 45   overrideSuspension = false;
 46   SetMinEnergy(0. * eV);
 47   SetMaxEnergy(20. * MeV);
 48 }
 49 
 50 G4ParticleHPElastic::~G4ParticleHPElastic()
 51 {
 52   // the vectror is shared among threads, only master deletes
 53   if (!G4Threading::IsWorkerThread()) {
 54     if (theElastic != nullptr) {
 55       for (auto it = theElastic->cbegin(); it != theElastic->cend(); ++it) {
 56         delete *it;
 57       }
 58       theElastic->clear();
 59     }
 60   }
 61 }
 62 
 63 G4HadFinalState* G4ParticleHPElastic::ApplyYourself(const G4HadProjectile& aTrack,
 64                                                     G4Nucleus& aNucleus)
 65 {
 66   return this->ApplyYourself(aTrack, aNucleus, false);
 67 }
 68 
 69 //--------------------------------------------------------
 70 // New method added by L. Thulliez (CEA-Saclay) 2021/05/04
 71 //--------------------------------------------------------
 72 G4HadFinalState* G4ParticleHPElastic::ApplyYourself(const G4HadProjectile& aTrack,
 73                                                     G4Nucleus& aNucleus, G4bool isFromTSL)
 74 {
 75   G4ParticleHPManager::GetInstance()->OpenReactionWhiteBoard();
 76   const G4Material* theMaterial = aTrack.GetMaterial();
 77   auto n = (G4int)theMaterial->GetNumberOfElements();
 78   std::size_t index = theMaterial->GetElement(0)->GetIndex();
 79 
 80   if (!isFromTSL) {
 81     if (n != 1) {
 82       G4int i;
 83       auto xSec = new G4double[n];
 84       G4double sum = 0;
 85       const G4double* NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
 86       G4double rWeight;
 87       G4ParticleHPThermalBoost aThermalE;
 88       for (i = 0; i < n; ++i) {
 89         index = theMaterial->GetElement(i)->GetIndex();
 90         rWeight = NumAtomsPerVolume[i];
 91         xSec[i] = ((*theElastic)[index])
 92                     ->GetXsec(aThermalE.GetThermalEnergy(aTrack, theMaterial->GetElement(i),
 93                                                          theMaterial->GetTemperature()));
 94         xSec[i] *= rWeight;
 95         sum += xSec[i];
 96       }
 97       G4double random = G4UniformRand();
 98       G4double running = 0;
 99       for (i = 0; i < n; ++i) {
100         running += xSec[i];
101         index = theMaterial->GetElement(i)->GetIndex();
102         if (sum == 0 || random <= running / sum) break;
103       }
104       delete[] xSec;
105     }
106   }
107   else {
108     G4int i;
109     if (n != 1) {
110       for (i = 0; i < n; ++i) {
111         if (aNucleus.GetZ_asInt() == (G4int)(theMaterial->GetElement(i)->GetZ())) {
112           index = theMaterial->GetElement(i)->GetIndex();
113         }
114       }
115     }
116   }
117 
118   // The boolean "true", as last argument, specifies to G4ParticleHPChannel::ApplyYourself
119   // that it is an elastic channel: this is needed for the special DBRC treatment.
120   G4HadFinalState* finalState = ((*theElastic)[index])->ApplyYourself(aTrack, -1, true);
121 
122   if (overrideSuspension) finalState->SetStatusChange(isAlive);
123 
124   // Overwrite target parameters
125   aNucleus.SetParameters(G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),
126                          G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
127   const G4Element* target_element = (*G4Element::GetElementTable())[index];
128   const G4Isotope* target_isotope = nullptr;
129   auto iele = (G4int)target_element->GetNumberOfIsotopes();
130   for (G4int j = 0; j != iele; ++j) {
131     target_isotope = target_element->GetIsotope(j);
132     if (target_isotope->GetN()
133         == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA())
134       break;
135   }
136   aNucleus.SetIsotope(target_isotope);
137 
138   G4ParticleHPManager::GetInstance()->CloseReactionWhiteBoard();
139   return finalState;
140 }
141 
142 const std::pair<G4double, G4double> G4ParticleHPElastic::GetFatalEnergyCheckLevels() const
143 {
144   // max energy non-conservation is mass of heavy nucleus
145   return std::pair<G4double, G4double>(10.0 * perCent, 350.0 * CLHEP::GeV);
146 }
147 
148 G4int G4ParticleHPElastic::GetVerboseLevel() const
149 {
150   return G4ParticleHPManager::GetInstance()->GetVerboseLevel();
151 }
152 
153 void G4ParticleHPElastic::SetVerboseLevel(G4int newValue)
154 {
155   G4ParticleHPManager::GetInstance()->SetVerboseLevel(newValue);
156 }
157 
158 void G4ParticleHPElastic::BuildPhysicsTable(const G4ParticleDefinition&)
159 {
160   G4ParticleHPManager* hpmanager = G4ParticleHPManager::GetInstance();
161 
162   theElastic = hpmanager->GetElasticFinalStates();
163 
164   if (G4Threading::IsMasterThread()) {
165     if (theElastic == nullptr) theElastic = new std::vector<G4ParticleHPChannel*>;
166 
167     if (numEle == (G4int)G4Element::GetNumberOfElements()) return;
168 
169     if (theElastic->size() == G4Element::GetNumberOfElements()) {
170       numEle = (G4int)G4Element::GetNumberOfElements();
171       return;
172     }
173 
174     auto theFS = new G4ParticleHPElasticFS;
175     if (G4FindDataDir("G4NEUTRONHPDATA") == nullptr)
176       throw G4HadronicException(
177         __FILE__, __LINE__,
178         "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
179     dirName = G4FindDataDir("G4NEUTRONHPDATA");
180     G4String tString = "/Elastic";
181     dirName = dirName + tString;
182     for (G4int i = numEle; i < (G4int)G4Element::GetNumberOfElements(); ++i) {
183       theElastic->push_back(new G4ParticleHPChannel);
184       ((*theElastic)[i])->Init((*(G4Element::GetElementTable()))[i], dirName);
185       // while(!((*theElastic)[i])->Register(theFS)) ;
186       ((*theElastic)[i])->Register(theFS);
187     }
188     delete theFS;
189     hpmanager->RegisterElasticFinalStates(theElastic);
190   }
191   numEle = (G4int)G4Element::GetNumberOfElements();
192 }
193 
194 void G4ParticleHPElastic::ModelDescription(std::ostream& outFile) const
195 {
196   outFile << "High Precision model based on Evaluated Nuclear Data Files (ENDF) for inelastic "
197              "reaction of neutrons below 20MeV\n";
198 }
199