Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/particle_hp/src/G4ParticleHPInelastic.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 // this code implementation is the intellectual property of
 27 // neutron_hp -- source file
 28 // J.P. Wellisch, Nov-1996
 29 // A prototype of the low energy neutron transport model.
 30 //
 31 // By copying, distributing or modifying the Program (or any work
 32 // based on the Program) you indicate your acceptance of this statement,
 33 // and all its terms.
 34 //
 35 //
 36 // 070523 bug fix for G4FPE_DEBUG on by A. Howard (and T. Koi)
 37 // 081203 limit maximum trial for creating final states add protection for 1H isotope case by T. Koi
 38 //
 39 // P. Arce, June-2014 Conversion neutron_hp to particle_hp
 40 // V. Ivanchenko, July-2023 Basic revision of particle HP classes
 41 //
 42 #include "G4ParticleHPInelastic.hh"
 43 
 44 #include "G4HadronicParameters.hh"
 45 #include "G4ParticleHP2AInelasticFS.hh"
 46 #include "G4ParticleHP2N2AInelasticFS.hh"
 47 #include "G4ParticleHP2NAInelasticFS.hh"
 48 #include "G4ParticleHP2NDInelasticFS.hh"
 49 #include "G4ParticleHP2NInelasticFS.hh"
 50 #include "G4ParticleHP2NPInelasticFS.hh"
 51 #include "G4ParticleHP2PInelasticFS.hh"
 52 #include "G4ParticleHP3AInelasticFS.hh"
 53 #include "G4ParticleHP3NAInelasticFS.hh"
 54 #include "G4ParticleHP3NInelasticFS.hh"
 55 #include "G4ParticleHP3NPInelasticFS.hh"
 56 #include "G4ParticleHP4NInelasticFS.hh"
 57 #include "G4ParticleHPAInelasticFS.hh"
 58 #include "G4ParticleHPD2AInelasticFS.hh"
 59 #include "G4ParticleHPDAInelasticFS.hh"
 60 #include "G4ParticleHPDInelasticFS.hh"
 61 #include "G4ParticleHPHe3InelasticFS.hh"
 62 #include "G4ParticleHPManager.hh"
 63 #include "G4ParticleHPN2AInelasticFS.hh"
 64 #include "G4ParticleHPN2PInelasticFS.hh"
 65 #include "G4ParticleHPN3AInelasticFS.hh"
 66 #include "G4ParticleHPNAInelasticFS.hh"
 67 #include "G4ParticleHPND2AInelasticFS.hh"
 68 #include "G4ParticleHPNDInelasticFS.hh"
 69 #include "G4ParticleHPNHe3InelasticFS.hh"
 70 #include "G4ParticleHPNInelasticFS.hh"
 71 #include "G4ParticleHPNPAInelasticFS.hh"
 72 #include "G4ParticleHPNPInelasticFS.hh"
 73 #include "G4ParticleHPNT2AInelasticFS.hh"
 74 #include "G4ParticleHPNTInelasticFS.hh"
 75 #include "G4ParticleHPNXInelasticFS.hh"
 76 #include "G4ParticleHPPAInelasticFS.hh"
 77 #include "G4ParticleHPPDInelasticFS.hh"
 78 #include "G4ParticleHPPInelasticFS.hh"
 79 #include "G4ParticleHPPTInelasticFS.hh"
 80 #include "G4ParticleHPT2AInelasticFS.hh"
 81 #include "G4ParticleHPTInelasticFS.hh"
 82 #include "G4ParticleHPThermalBoost.hh"
 83 #include "G4SystemOfUnits.hh"
 84 #include "G4AutoLock.hh"
 85 
 86 G4bool G4ParticleHPInelastic::fLock[] = {true, true, true, true, true, true};
 87 std::vector<G4ParticleHPChannelList*>*
 88 G4ParticleHPInelastic::theInelastic[] = {nullptr, nullptr, nullptr, nullptr, nullptr, nullptr};
 89 
 90 namespace
 91 {
 92   G4Mutex theHPInelastic = G4MUTEX_INITIALIZER;
 93 }
 94 
 95 G4ParticleHPInelastic::G4ParticleHPInelastic(G4ParticleDefinition* p, const char* name)
 96   : G4HadronicInteraction(name), theProjectile(p)
 97 {
 98   fManager = G4ParticleHPManager::GetInstance();
 99   dirName = fManager->GetParticleHPPath(theProjectile) + "/Inelastic";
100   indexP = fManager->GetPHPIndex(theProjectile);
101 
102 #ifdef G4VERBOSE
103   if (fManager->GetVerboseLevel() > 1)
104     G4cout << "@@@ G4ParticleHPInelastic instantiated for "
105            << p->GetParticleName() << " indexP=" << indexP
106      << "/n    data directory " << dirName << G4endl;
107 #endif
108 }
109 
110 G4ParticleHPInelastic::~G4ParticleHPInelastic()
111 {
112   // Vector is shared, only one delete
113   if (isFirst) {
114     ClearData();
115   }
116 }
117 
118 void G4ParticleHPInelastic::ClearData()
119 {
120   if (theInelastic[indexP] != nullptr) {
121     for (auto const& p : *(theInelastic[indexP])) {
122       delete p;
123     }
124     delete theInelastic[indexP];
125     theInelastic[indexP] = nullptr;
126   }
127 }
128 
129 G4HadFinalState* G4ParticleHPInelastic::ApplyYourself(const G4HadProjectile& aTrack,
130                                                       G4Nucleus& aNucleus)
131 {
132   G4ParticleHPManager::GetInstance()->OpenReactionWhiteBoard();
133   const G4Material* theMaterial = aTrack.GetMaterial();
134   auto n = (G4int)theMaterial->GetNumberOfElements();
135   auto elm = theMaterial->GetElement(0);
136   std::size_t index = elm->GetIndex();
137   G4int it = 0;
138   /*
139   G4cout << "G4ParticleHPInelastic::ApplyYourself n=" << n << " index=" << index
140    << " indexP=" << indexP << " "
141          << aTrack.GetDefinition()->GetParticleName() << G4endl;
142   */
143   if (n != 1) {
144     auto xSec = new G4double[n];
145     G4double sum = 0;
146     G4int i;
147     const G4double* NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
148     G4double rWeight;
149     G4double xs; 
150     G4ParticleHPThermalBoost aThermalE;
151     for (i = 0; i < n; ++i) {
152       elm = theMaterial->GetElement(i);
153       index = elm->GetIndex();
154       /*
155       G4cout << "i=" << i << "  index=" << index << "  " << elm->GetName() 
156        << "  " << (*(theInelastic[indexP]))[index] << G4endl;
157       */
158       rWeight = NumAtomsPerVolume[i];
159       if (aTrack.GetDefinition() == G4Neutron::Neutron()) {
160         xs = (*(theInelastic[indexP]))[index]->GetXsec(aThermalE.GetThermalEnergy(aTrack, elm,
161                    theMaterial->GetTemperature()));
162       }
163       else {
164         xs = (*(theInelastic[indexP]))[index]->GetXsec(aTrack.GetKineticEnergy());
165       }
166       xs *= rWeight;
167       sum += xs;
168       xSec[i] = sum;
169 #ifdef G4VERBOSE
170       if (fManager->GetDEBUG())
171         G4cout << " G4ParticleHPInelastic XSEC ELEM " << i << " = " << xSec[i] << G4endl;
172 #endif
173     }
174     sum *= G4UniformRand();
175     for (it = 0; it < n; ++it) {
176       elm = theMaterial->GetElement(it);
177       index = elm->GetIndex();
178       if (sum <= xSec[it]) break;
179     }
180     delete[] xSec;
181   }
182 
183 #ifdef G4VERBOSE
184   if (fManager->GetDEBUG())
185     G4cout << " G4ParticleHPInelastic: Elem it=" << it << "  "
186            << elm->GetName() << " index=" << index
187      << " from material " << theMaterial->GetName()
188            << G4endl;
189 #endif
190 
191   G4HadFinalState* result = 
192     (*(theInelastic[indexP]))[index]->ApplyYourself(elm, aTrack);
193 
194   aNucleus.SetParameters(fManager->GetReactionWhiteBoard()->GetTargA(),
195                          fManager->GetReactionWhiteBoard()->GetTargZ());
196   
197   const G4Element* target_element = (*G4Element::GetElementTable())[index];
198   const G4Isotope* target_isotope = nullptr;
199   auto iele = (G4int)target_element->GetNumberOfIsotopes();
200   for (G4int j = 0; j != iele; ++j) {
201     target_isotope = target_element->GetIsotope(j);
202     if (target_isotope->GetN() == fManager->GetReactionWhiteBoard()->GetTargA())
203       break;
204   }
205   aNucleus.SetIsotope(target_isotope);
206 
207   G4ParticleHPManager::GetInstance()->CloseReactionWhiteBoard();
208 
209   return result;
210 }
211 
212 const std::pair<G4double, G4double> G4ParticleHPInelastic::GetFatalEnergyCheckLevels() const
213 {
214   // max energy non-conservation is mass of heavy nucleus
215   return std::pair<G4double, G4double>(10.0 * perCent, 350.0 * CLHEP::GeV);
216 }
217 
218 void G4ParticleHPInelastic::BuildPhysicsTable(const G4ParticleDefinition& projectile)
219 {
220   if (fLock[indexP]) {
221     G4AutoLock l(&theHPInelastic);
222     if (fLock[indexP]) {
223       isFirst = true;
224       fLock[indexP] = false;
225     }
226     l.unlock();
227   }
228 
229   G4int nelm = (G4int)G4Element::GetNumberOfElements();
230   G4int n0 = numEle;
231   numEle = nelm;
232   if (!isFirst || nelm == n0) { return; }
233 
234   // extra elements should be initialized
235   G4AutoLock l(&theHPInelastic);
236 
237   if (nullptr == theInelastic[indexP]) {
238     theInelastic[indexP] = new std::vector<G4ParticleHPChannelList*>;
239   }
240 
241   if (fManager->GetVerboseLevel() > 0 && isFirst) {
242     fManager->DumpSetting();
243     G4cout << "@@@ G4ParticleHPInelastic instantiated for particle "
244      << theProjectile->GetParticleName() << "/n    data directory is "
245      << dirName << G4endl;
246   }
247 
248   auto table = G4Element::GetElementTable();
249   for (G4int i = n0; i < nelm; ++i) {
250     auto clist = new G4ParticleHPChannelList(36, theProjectile);
251     theInelastic[indexP]->push_back(clist);
252     clist->Init((*table)[i], dirName, theProjectile);
253     clist->Register(new G4ParticleHPNInelasticFS, "F01/");  // has
254     clist->Register(new G4ParticleHPNXInelasticFS, "F02/");
255     clist->Register(new G4ParticleHP2NDInelasticFS, "F03/");
256     clist->Register(new G4ParticleHP2NInelasticFS, "F04/");  // has, E Done
257     clist->Register(new G4ParticleHP3NInelasticFS, "F05/");  // has, E Done
258     clist->Register(new G4ParticleHPNAInelasticFS, "F06/");
259     clist->Register(new G4ParticleHPN3AInelasticFS, "F07/");
260     clist->Register(new G4ParticleHP2NAInelasticFS, "F08/");
261     clist->Register(new G4ParticleHP3NAInelasticFS, "F09/");
262     clist->Register(new G4ParticleHPNPInelasticFS, "F10/");
263     clist->Register(new G4ParticleHPN2AInelasticFS, "F11/");
264     clist->Register(new G4ParticleHP2N2AInelasticFS, "F12/");
265     clist->Register(new G4ParticleHPNDInelasticFS, "F13/");
266     clist->Register(new G4ParticleHPNTInelasticFS, "F14/");
267     clist->Register(new G4ParticleHPNHe3InelasticFS, "F15/");
268     clist->Register(new G4ParticleHPND2AInelasticFS, "F16/");
269     clist->Register(new G4ParticleHPNT2AInelasticFS, "F17/");
270     clist->Register(new G4ParticleHP4NInelasticFS, "F18/");  // has, E Done
271     clist->Register(new G4ParticleHP2NPInelasticFS, "F19/");
272     clist->Register(new G4ParticleHP3NPInelasticFS, "F20/");
273     clist->Register(new G4ParticleHPN2PInelasticFS, "F21/");
274     clist->Register(new G4ParticleHPNPAInelasticFS, "F22/");
275     clist->Register(new G4ParticleHPPInelasticFS, "F23/");
276     clist->Register(new G4ParticleHPDInelasticFS, "F24/");
277     clist->Register(new G4ParticleHPTInelasticFS, "F25/");
278     clist->Register(new G4ParticleHPHe3InelasticFS, "F26/");
279     clist->Register(new G4ParticleHPAInelasticFS, "F27/");
280     clist->Register(new G4ParticleHP2AInelasticFS, "F28/");
281     clist->Register(new G4ParticleHP3AInelasticFS, "F29/");
282     clist->Register(new G4ParticleHP2PInelasticFS, "F30/");
283     clist->Register(new G4ParticleHPPAInelasticFS, "F31/");
284     clist->Register(new G4ParticleHPD2AInelasticFS, "F32/");
285     clist->Register(new G4ParticleHPT2AInelasticFS, "F33/");
286     clist->Register(new G4ParticleHPPDInelasticFS, "F34/");
287     clist->Register(new G4ParticleHPPTInelasticFS, "F35/");
288     clist->Register(new G4ParticleHPDAInelasticFS, "F36/");
289 #ifdef G4VERBOSE
290     if (fManager->GetVerboseLevel() > 1) {
291       G4cout << "ParticleHP::Inelastic for " 
292        << theProjectile->GetParticleName() << " off " 
293        << (*table)[i]->GetName() << G4endl;
294     }
295 #endif
296   }
297   fManager->RegisterInelasticFinalStates( &projectile , theInelastic[indexP] );
298   l.unlock();
299 }
300 
301 void G4ParticleHPInelastic::ModelDescription(std::ostream& outFile) const
302 {
303   outFile << "High Precision (HP) model for inelastic reaction of "
304     << theProjectile->GetParticleName() << " below 20MeV\n";
305 }
306