Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPChannel.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// 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// 071031 bug fix T. Koi on behalf of A. Howard
32// 081203 bug fix in Register method by T. Koi
33//
34// P. Arce, June-2014 Conversion neutron_hp to particle_hp
35//
36// June-2019 - E. Mendoza --> Modification to allow using an incomplete
37// data library if the G4NEUTRONHP_SKIP_MISSING_ISOTOPES environmental
38// flag is defined. The missing XS are set to 0.
39
41
42#include "G4HadTmpUtil.hh"
47#include "G4SystemOfUnits.hh"
48#include "globals.hh"
49
50#include <cstdlib>
51
53{
56 wendtFissionGenerator = G4WendtFissionFragmentGenerator::GetInstance();
57 // Make sure both fission fragment models are not active at same time
59 }
60 theProjectile = (nullptr == p) ? G4Neutron::Neutron() : p;
61 theChannelData = new G4ParticleHPVector;
62}
63
65{
66 delete theChannelData;
67 // Following statement disabled to avoid SEGV
68 // theBuffer is also deleted as "theChannelData" in
69 delete[] theIsotopeWiseData;
70 if (theFinalStates != nullptr) {
71 for (G4int i = 0; i < niso; i++) {
72 delete theFinalStates[i];
73 }
74 delete[] theFinalStates;
75 }
76 delete[] active;
77}
78
80{
81 return std::max(0., theChannelData->GetXsec(energy));
82}
83
85 G4int isoNumber)
86{
87 return theIsotopeWiseData[isoNumber].GetXsec(energy);
88}
89
91 G4int isoNumber)
92{
93 return theFinalStates[isoNumber]->GetXsec(energy);
94}
95
97 const G4String& dirName, const G4String& aFSType)
98{
99 theFSType = aFSType;
100 Init(anElement, dirName);
101}
102
103void G4ParticleHPChannel::Init(G4Element* anElement, const G4String& dirName)
104{
105 theDir = dirName;
106 theElement = anElement;
107}
108
110{
111 ++registerCount;
112 G4int Z = theElement->GetZasInt();
113
114 niso = (G4int)theElement->GetNumberOfIsotopes();
115 const std::size_t nsize = niso > 0 ? niso : 1;
116
117 delete[] theIsotopeWiseData;
118 theIsotopeWiseData = new G4ParticleHPIsoData[nsize];
119 delete[] active;
120 active = new G4bool[nsize];
121
122 delete[] theFinalStates;
123 theFinalStates = new G4ParticleHPFinalState*[nsize];
124 delete theChannelData;
125 theChannelData = new G4ParticleHPVector;
126 for (G4int i = 0; i < niso; ++i) {
127 theFinalStates[i] = theFS->New();
128 theFinalStates[i]->SetProjectile(theProjectile);
129 }
130 if (niso != 0 && registerCount == 0) {
131 for (G4int i1 = 0; i1 < niso; ++i1) {
132 G4int A = theElement->GetIsotope(i1)->GetN();
133 G4int M = theElement->GetIsotope(i1)->Getm();
134 //G4cout <<" Init: normal case i=" << i1
135 // << " Z=" << Z << " A=" << A << G4endl;
136 G4double frac = theElement->GetRelativeAbundanceVector()[i1] / perCent;
137 theFinalStates[i1]->SetA_Z(A, Z, M);
138 UpdateData(A, Z, M, i1, frac, theProjectile);
139 }
140 }
142
143 // To avoid issuing hash by worker threads
144 if (result) theChannelData->Hash();
145
146 return result;
147}
148
150 G4double abundance,
151 G4ParticleDefinition* projectile)
152{
153 // Initialze the G4FissionFragment generator for this isomer if needed
154 if (wendtFissionGenerator != nullptr) {
155 wendtFissionGenerator->InitializeANucleus(A, Z, M, theDir);
156 }
157
158 theFinalStates[index]->Init(A, Z, M, theDir, theFSType, projectile);
159 if (!theFinalStates[index]->HasAnyData()) return;
160 // nothing there for exactly this isotope.
161
162 // the above has put the X-sec into the FS
163 theBuffer = nullptr;
164 if (theFinalStates[index]->HasXsec()) {
165 theBuffer = theFinalStates[index]->GetXsec();
166 theBuffer->Times(abundance / 100.);
167 theIsotopeWiseData[index].FillChannelData(theBuffer);
168 }
169 else // get data from CrossSection directory
170 {
171 G4String tString = "/CrossSection";
172 active[index] = theIsotopeWiseData[index].Init(A, Z, M, abundance,
173 theDir, tString);
174 if (active[index]) theBuffer = theIsotopeWiseData[index].MakeChannelData();
175 }
176 if (theBuffer != nullptr) Harmonise(theChannelData, theBuffer);
177}
178
180 G4ParticleHPVector* theNew)
181{
182 G4int s_tmp = 0, n = 0, m_tmp = 0;
183 auto theMerge = new G4ParticleHPVector;
184 G4ParticleHPVector* anActive = theStore;
185 G4ParticleHPVector* aPassive = theNew;
187 G4int a = s_tmp, p = n, t;
188 while (a < anActive->GetVectorLength() && p < aPassive->GetVectorLength())
189 // Loop checking, 11.05.2015, T. Koi
190 {
191 if (anActive->GetEnergy(a) <= aPassive->GetEnergy(p)) {
192 G4double xa = anActive->GetEnergy(a);
193 theMerge->SetData(m_tmp, xa, anActive->GetXsec(a) + std::max(0., aPassive->GetXsec(xa)));
194 m_tmp++;
195 a++;
196 G4double xp = aPassive->GetEnergy(p);
197 if (std::abs(std::abs(xp - xa) / xa) < 0.001) {
198 ++p;
199 }
200 }
201 else {
202 tmp = anActive;
203 t = a;
204 anActive = aPassive;
205 a = p;
206 aPassive = tmp;
207 p = t;
208 }
209 }
210 while (a != anActive->GetVectorLength()) // Loop checking, 11.05.2015, T. Koi
211 {
212 theMerge->SetData(m_tmp++, anActive->GetEnergy(a), anActive->GetXsec(a));
213 ++a;
214 }
215 while (p != aPassive->GetVectorLength()) // Loop checking, 11.05.2015, T. Koi
216 {
217 if (std::abs(theMerge->GetEnergy(std::max(0, m_tmp - 1)) -
218 aPassive->GetEnergy(p)) / aPassive->GetEnergy(p) > 0.001)
219 theMerge->SetData(m_tmp++, aPassive->GetEnergy(p), aPassive->GetXsec(p));
220 ++p;
221 }
222 delete theStore;
223 theStore = theMerge;
224}
225
228 G4int anIsotope, G4bool isElastic)
229{
230 //G4cout << "G4ParticleHPChannel::ApplyYourself niso=" << niso
231 // << " ni=" << anIsotope << " isElastic=" << isElastic <<G4endl;
232 if (anIsotope != -1 && anIsotope != -2) {
233 // Inelastic Case
234 //G4cout << "G4ParticleHPChannel Inelastic Case"
235 //<< " Z= " << GetZ(anIsotope) << " A = " << GetN(anIsotope) << G4endl;
238 return theFinalStates[anIsotope]->ApplyYourself(theTrack);
239 }
240 G4double sum = 0;
241 G4int it = 0;
242 auto xsec = new G4double[niso];
243 G4ParticleHPThermalBoost aThermalE;
244 for (G4int i = 0; i < niso; i++) {
245 if (theFinalStates[i]->HasAnyData()) {
246 /*
247 G4cout << "FS: " << i << theTrack.GetDefinition()->GetParticleName()
248 << " Z=" << theFinalStates[i]->GetZ()
249 << " A=" << theFinalStates[i]->GetN()
250 << G4endl;
251 */
252 xsec[i] = theIsotopeWiseData[i].GetXsec(
253 aThermalE.GetThermalEnergy(theTrack, theFinalStates[i]->GetN(),
254 theFinalStates[i]->GetZ(),
255 theTrack.GetMaterial()->GetTemperature()));
256 sum += xsec[i];
257 }
258 else {
259 xsec[i] = 0;
260 }
261 }
262 if (sum == 0) {
263 it = G4lrint(niso * G4UniformRand());
264 }
265 else {
266 G4double random = G4UniformRand();
267 G4double running = 0;
268 for (G4int ix = 0; ix < niso; ix++) {
269 running += xsec[ix];
270 if (sum == 0 || random <= running / sum) {
271 it = ix;
272 break;
273 }
274 }
275 if (it == niso) it--;
276 }
277 delete[] xsec;
278 G4HadFinalState* theFinalState = nullptr;
279 const auto A = (G4int)this->GetN(it);
280 const auto Z = (G4int)this->GetZ(it);
281 const auto M = (G4int)this->GetM(it);
282
283 //-2:Marker for Fission
284 if ((wendtFissionGenerator != nullptr) && anIsotope == -2) {
285 theFinalState = wendtFissionGenerator->ApplyYourself(theTrack, Z, A);
286 }
287
288 // Use the standard procedure if the G4FissionFragmentGenerator model fails
289 if (theFinalState == nullptr) {
290 G4int icounter = 0;
291 G4int icounter_max = 1024;
292 while (theFinalState == nullptr) // Loop checking, 11.05.2015, T. Koi
293 {
294 icounter++;
295 if (icounter > icounter_max) {
296 G4cout << "Loop-counter exceeded the threshold value at "
297 << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
298 break;
299 }
300 if (isElastic) {
301 // Register 0 K cross-section for DBRC for Doppler broadened elastic scattering kernel
302 G4ParticleHPVector* xsforFS = theIsotopeWiseData[it].MakeChannelData();
303 // Only G4ParticleHPElasticFS has the RegisterCrossSection method
304 static_cast<G4ParticleHPElasticFS*>(theFinalStates[it])->RegisterCrossSection(xsforFS);
305 }
306 theFinalState = theFinalStates[it]->ApplyYourself(theTrack);
307 }
308 }
309
310 // G4cout <<"THE IMPORTANT RETURN"<<G4endl;
311 // G4cout << "TK G4ParticleHPChannel Elastic, Capture and Fission Cases "
312 //<< " Z= " << this->GetZ(it) << " A = " << this->GetN(it) << G4endl;
316
317 return theFinalState;
318}
319
321{
322 G4cout << " Element: " << theElement->GetName() << G4endl;
323 G4cout << " Directory name: " << theDir << G4endl;
324 G4cout << " FS name: " << theFSType << G4endl;
325 G4cout << " Number of Isotopes: " << niso << G4endl;
326 G4cout << " Have cross sections: " << G4endl;
327 for (int i = 0; i < niso; i++) {
328 G4cout << theFinalStates[i]->HasXsec() << " ";
329 }
330 G4cout << G4endl;
331 if (theChannelData != nullptr) {
332 G4cout << " Cross Section (total for this channel):" << G4endl;
333 int np = theChannelData->GetVectorLength();
334 G4cout << np << G4endl;
335 for (int i = 0; i < np; i++) {
336 G4cout << theChannelData->GetEnergy(i) / eV << " " << theChannelData->GetXsec(i) << G4endl;
337 }
338 }
339}
#define M(row, col)
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
G4double * GetRelativeAbundanceVector() const
Definition G4Element.hh:149
const G4Isotope * GetIsotope(G4int iso) const
Definition G4Element.hh:151
size_t GetNumberOfIsotopes() const
Definition G4Element.hh:143
const G4String & GetName() const
Definition G4Element.hh:115
G4int GetZasInt() const
Definition G4Element.hh:120
const G4Material * GetMaterial() const
G4int Getm() const
Definition G4Isotope.hh:89
G4int GetN() const
Definition G4Isotope.hh:83
G4double GetTemperature() const
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
G4double GetZ(G4int i) const
G4bool HasAnyData(G4int isoNumber)
void UpdateData(G4int A, G4int Z, G4int index, G4double abundance, G4ParticleDefinition *projectile)
G4ParticleHPChannel(G4ParticleDefinition *projectile=nullptr)
G4double GetN(G4int i) const
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack, G4int isoNumber=-1, G4bool isElastic=false)
void Harmonise(G4ParticleHPVector *&theStore, G4ParticleHPVector *theNew)
G4double GetWeightedXsec(G4double energy, G4int isoNumber)
G4double GetM(G4int i) const
G4bool Register(G4ParticleHPFinalState *theFS)
G4double GetXsec(G4double energy)
void Init(G4Element *theElement, const G4String &dirName)
G4ParticleHPManager * fManager
G4double GetFSCrossSection(G4double energy, G4int isoNumber)
virtual G4double GetXsec(G4double)
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &)
void SetA_Z(G4double anA, G4double aZ, G4int aM=0)
virtual G4ParticleHPFinalState * New()=0
void SetProjectile(G4ParticleDefinition *projectile)
void Init(G4double A, G4double Z, G4String &dirName, G4String &aFSType, G4ParticleDefinition *p)
G4ParticleHPVector * MakeChannelData()
G4bool Init(G4int A, G4int Z, G4double abun, G4String dirName, G4String aFSType)
void FillChannelData(G4ParticleHPVector *aBuffer)
G4double GetXsec(G4double energy)
void SetProduceFissionFragments(G4bool val)
G4bool GetUseWendtFissionModel() const
static G4ParticleHPManager * GetInstance()
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
void Times(G4double factor)
G4double GetXsec(G4int i)
G4double GetEnergy(G4int i) const
G4int GetVectorLength() const
void InitializeANucleus(const G4int A, const G4int Z, const G4int M, const G4String &dataDirectory)
G4HadFinalState * ApplyYourself(const G4HadProjectile &projectile, G4int Z, G4int A)
static G4WendtFissionFragmentGenerator * GetInstance()
int G4lrint(double ad)
Definition templates.hh:134