Geant4 11.2.2
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[] = {true, true, true, true, true, true};
87std::vector<G4ParticleHPChannelList*>*
88G4ParticleHPInelastic::theInelastic[] = {nullptr, nullptr, nullptr, nullptr, nullptr, nullptr};
89
90namespace
91{
92 G4Mutex theHPInelastic = G4MUTEX_INITIALIZER;
93}
94
96 : G4HadronicInteraction(name), theProjectile(p)
97{
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
111{
112 // Vector is shared, only one delete
113 if (isFirst) {
114 ClearData();
115 }
116}
117
118void 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
130 G4Nucleus& aNucleus)
131{
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
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
208
209 return result;
210}
211
212const 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
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
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) {
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 l.unlock();
298}
299
300void G4ParticleHPInelastic::ModelDescription(std::ostream& outFile) const
301{
302 outFile << "High Precision (HP) model for inelastic reaction of "
303 << theProjectile->GetParticleName() << " below 20MeV\n";
304}
#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 G4ElementTable * GetElementTable()
Definition G4Element.cc:389
static size_t GetNumberOfElements()
Definition G4Element.cc:393
const G4Isotope * GetIsotope(G4int iso) const
Definition G4Element.hh:151
size_t GetIndex() const
Definition G4Element.hh:159
size_t GetNumberOfIsotopes() const
Definition G4Element.hh:143
const G4Material * GetMaterial() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
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:114
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")
const G4String & GetParticleHPPath(const G4ParticleDefinition *) const
G4int GetPHPIndex(const G4ParticleDefinition *) const
static G4ParticleHPManager * GetInstance()
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)