Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPInelastic.cc
Go to the documentation of this file.
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//
43
83#include "G4SystemOfUnits.hh"
84#include "G4AutoLock.hh"
85
86G4bool G4ParticleHPInelastic::fLock[] = {false, false, false, false, false, false};
87std::vector<G4ParticleHPChannelList*>*
88G4ParticleHPInelastic::theInelastic[] = {nullptr, nullptr, nullptr, nullptr, nullptr, nullptr};
89
90static std::once_flag applyOnce;
91
92namespace
93{
94 G4Mutex theHPInelastic = G4MUTEX_INITIALIZER;
95}
96
98 : G4HadronicInteraction(name), theProjectile(p)
99{
101 dirName = fManager->GetParticleHPPath(theProjectile) + "/Inelastic";
102 indexP = fManager->GetPHPIndex(theProjectile);
103
104#ifdef G4VERBOSE
105 if (fManager->GetVerboseLevel() > 1)
106 G4cout << "@@@ G4ParticleHPInelastic instantiated for "
107 << p->GetParticleName() << " indexP=" << indexP
108 << "/n data directory " << dirName << G4endl;
109#endif
110}
111
113{
114 // Vector is shared, only one delete
115 if (isFirst) {
116 ClearData();
117 }
118}
119
120void G4ParticleHPInelastic::ClearData()
121{
122 for (G4int i=0; i<6; ++i) {
123 if (theInelastic[i] != nullptr) {
124 for (auto const& p : *(theInelastic[i])) {
125 delete p;
126 }
127 delete theInelastic[i];
128 theInelastic[i] = nullptr;
129 }
130 }
131}
132
134 G4Nucleus& aNucleus)
135{
137 const G4Material* theMaterial = aTrack.GetMaterial();
138 auto n = (G4int)theMaterial->GetNumberOfElements();
139 auto elm = theMaterial->GetElement(0);
140 std::size_t index = elm->GetIndex();
141 G4int it = 0;
142 /*
143 G4cout << "G4ParticleHPInelastic::ApplyYourself n=" << n << " index=" << index
144 << " indexP=" << indexP << " "
145 << aTrack.GetDefinition()->GetParticleName() << G4endl;
146 */
147 if (n != 1) {
148 auto xSec = new G4double[n];
149 G4double sum = 0;
150 G4int i;
151 const G4double* NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
152 G4double rWeight;
153 G4double xs;
154 G4ParticleHPThermalBoost aThermalE;
155 for (i = 0; i < n; ++i) {
156 elm = theMaterial->GetElement(i);
157 index = elm->GetIndex();
158 /*
159 G4cout << "i=" << i << " index=" << index << " " << elm->GetName()
160 << " " << (*(theInelastic[indexP]))[index] << G4endl;
161 */
162 rWeight = NumAtomsPerVolume[i];
163 if (aTrack.GetDefinition() == G4Neutron::Neutron()) {
164 xs = (*(theInelastic[indexP]))[index]->GetXsec(aThermalE.GetThermalEnergy(aTrack, elm,
165 theMaterial->GetTemperature()));
166 }
167 else {
168 xs = (*(theInelastic[indexP]))[index]->GetXsec(aTrack.GetKineticEnergy());
169 }
170 xs *= rWeight;
171 sum += xs;
172 xSec[i] = sum;
173#ifdef G4VERBOSE
174 if (fManager->GetDEBUG())
175 G4cout << " G4ParticleHPInelastic XSEC ELEM " << i << " = " << xSec[i] << G4endl;
176#endif
177 }
178 sum *= G4UniformRand();
179 for (it = 0; it < n; ++it) {
180 elm = theMaterial->GetElement(it);
181 index = elm->GetIndex();
182 if (sum <= xSec[it]) break;
183 }
184 delete[] xSec;
185 }
186
187#ifdef G4VERBOSE
188 if (fManager->GetDEBUG())
189 G4cout << " G4ParticleHPInelastic: Elem it=" << it << " "
190 << elm->GetName() << " index=" << index
191 << " from material " << theMaterial->GetName()
192 << G4endl;
193#endif
194
195 G4HadFinalState* result =
196 (*(theInelastic[indexP]))[index]->ApplyYourself(elm, aTrack);
197
198 aNucleus.SetParameters(fManager->GetReactionWhiteBoard()->GetTargA(),
199 fManager->GetReactionWhiteBoard()->GetTargZ());
200
201 const G4Element* target_element = (*G4Element::GetElementTable())[index];
202 const G4Isotope* target_isotope = nullptr;
203 auto iele = (G4int)target_element->GetNumberOfIsotopes();
204 for (G4int j = 0; j != iele; ++j) {
205 target_isotope = target_element->GetIsotope(j);
206 if (target_isotope->GetN() == fManager->GetReactionWhiteBoard()->GetTargA())
207 break;
208 }
209 aNucleus.SetIsotope(target_isotope);
210
212
213 return result;
214}
215
216const std::pair<G4double, G4double> G4ParticleHPInelastic::GetFatalEnergyCheckLevels() const
217{
218 // max energy non-conservation is mass of heavy nucleus
219 return std::pair<G4double, G4double>(10.0 * perCent, 350.0 * CLHEP::GeV);
220}
221
223{
224 // only one object destroy all list of final state models
225 std::call_once(applyOnce, [this]() { isFirst = true; });
226
228 if (isFirst && nelm != numEle) {
229 G4AutoLock l(&theHPInelastic);
230 if (nelm != numEle) {
231 for (G4int i=0; i<6; ++i) { fLock[i] = false; }
232 }
233 l.unlock();
234 }
235 if (fLock[indexP]) { return; }
236
237 // extra elements should be initialized
238 G4AutoLock l(&theHPInelastic);
239 if (fLock[indexP]) { return; }
240
241 fLock[indexP] = true;
242 G4int n0 = numEle;
243 numEle = nelm;
244
245 if (nullptr == theInelastic[indexP]) {
246 theInelastic[indexP] = new std::vector<G4ParticleHPChannelList*>;
247 }
248
249 if (fManager->GetVerboseLevel() > 0 && isFirst) {
250 fManager->DumpSetting();
251 G4cout << "@@@ G4ParticleHPInelastic instantiated for particle "
252 << theProjectile->GetParticleName() << "/n data directory is "
253 << dirName << G4endl;
254 }
255
256 auto table = G4Element::GetElementTable();
257 for (G4int i = n0; i < nelm; ++i) {
258 auto clist = new G4ParticleHPChannelList(36, theProjectile);
259 theInelastic[indexP]->push_back(clist);
260 clist->Init((*table)[i], dirName, theProjectile);
261 clist->Register(new G4ParticleHPNInelasticFS, "F01/"); // has
262 clist->Register(new G4ParticleHPNXInelasticFS, "F02/");
263 clist->Register(new G4ParticleHP2NDInelasticFS, "F03/");
264 clist->Register(new G4ParticleHP2NInelasticFS, "F04/"); // has, E Done
265 clist->Register(new G4ParticleHP3NInelasticFS, "F05/"); // has, E Done
266 clist->Register(new G4ParticleHPNAInelasticFS, "F06/");
267 clist->Register(new G4ParticleHPN3AInelasticFS, "F07/");
268 clist->Register(new G4ParticleHP2NAInelasticFS, "F08/");
269 clist->Register(new G4ParticleHP3NAInelasticFS, "F09/");
270 clist->Register(new G4ParticleHPNPInelasticFS, "F10/");
271 clist->Register(new G4ParticleHPN2AInelasticFS, "F11/");
272 clist->Register(new G4ParticleHP2N2AInelasticFS, "F12/");
273 clist->Register(new G4ParticleHPNDInelasticFS, "F13/");
274 clist->Register(new G4ParticleHPNTInelasticFS, "F14/");
275 clist->Register(new G4ParticleHPNHe3InelasticFS, "F15/");
276 clist->Register(new G4ParticleHPND2AInelasticFS, "F16/");
277 clist->Register(new G4ParticleHPNT2AInelasticFS, "F17/");
278 clist->Register(new G4ParticleHP4NInelasticFS, "F18/"); // has, E Done
279 clist->Register(new G4ParticleHP2NPInelasticFS, "F19/");
280 clist->Register(new G4ParticleHP3NPInelasticFS, "F20/");
281 clist->Register(new G4ParticleHPN2PInelasticFS, "F21/");
282 clist->Register(new G4ParticleHPNPAInelasticFS, "F22/");
283 clist->Register(new G4ParticleHPPInelasticFS, "F23/");
284 clist->Register(new G4ParticleHPDInelasticFS, "F24/");
285 clist->Register(new G4ParticleHPTInelasticFS, "F25/");
286 clist->Register(new G4ParticleHPHe3InelasticFS, "F26/");
287 clist->Register(new G4ParticleHPAInelasticFS, "F27/");
288 clist->Register(new G4ParticleHP2AInelasticFS, "F28/");
289 clist->Register(new G4ParticleHP3AInelasticFS, "F29/");
290 clist->Register(new G4ParticleHP2PInelasticFS, "F30/");
291 clist->Register(new G4ParticleHPPAInelasticFS, "F31/");
292 clist->Register(new G4ParticleHPD2AInelasticFS, "F32/");
293 clist->Register(new G4ParticleHPT2AInelasticFS, "F33/");
294 clist->Register(new G4ParticleHPPDInelasticFS, "F34/");
295 clist->Register(new G4ParticleHPPTInelasticFS, "F35/");
296 clist->Register(new G4ParticleHPDAInelasticFS, "F36/");
297#ifdef G4VERBOSE
298 if (fManager->GetVerboseLevel() > 1) {
299 G4cout << "ParticleHP::Inelastic for "
300 << theProjectile->GetParticleName() << " off "
301 << (*table)[i]->GetName() << G4endl;
302 }
303#endif
304 }
305 fManager->RegisterInelasticFinalStates( &projectile , theInelastic[indexP] );
306 l.unlock();
307}
308
309void G4ParticleHPInelastic::ModelDescription(std::ostream& outFile) const
310{
311 outFile << "High Precision (HP) model for inelastic reaction of "
312 << theProjectile->GetParticleName() << " below 20MeV\n";
313}
G4TemplateAutoLock< G4Mutex > G4AutoLock
#define G4MUTEX_INITIALIZER
std::mutex G4Mutex
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
static std::size_t GetNumberOfElements()
Definition G4Element.cc:408
std::size_t GetNumberOfIsotopes() const
Definition G4Element.hh:143
std::size_t GetIndex() const
Definition G4Element.hh:159
const G4Isotope * GetIsotope(G4int iso) const
Definition G4Element.hh:151
static const G4ElementTable * GetElementTable()
Definition G4Element.cc:401
const G4Material * GetMaterial() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4HadronicInteraction(const G4String &modelName="HadronicModel")
G4int GetN() const
Definition G4Isotope.hh:83
G4double GetTemperature() const
const G4Element * GetElement(G4int iel) const
const G4double * GetVecNbOfAtomsPerVolume() const
std::size_t GetNumberOfElements() const
const G4String & GetName() const
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
void SetParameters(const G4double A, const G4double Z, const G4int numberOfLambdas=0)
Definition G4Nucleus.cc:319
void SetIsotope(const G4Isotope *iso)
Definition G4Nucleus.hh:93
const G4String & GetParticleName() const
G4ParticleHPManager * fManager
void ModelDescription(std::ostream &outFile) const override
static std::vector< G4ParticleHPChannelList * > * theInelastic[6]
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus) override
const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const override
G4ParticleHPInelastic(G4ParticleDefinition *p=G4Neutron::Neutron(), const char *name="NeutronHPInelastic")
static G4ParticleHPManager * GetInstance()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)