Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPChannelList.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//
32// P. Arce, June-2014 Conversion neutron_hp to particle_hp
33// V. Ivanchenko, July-2023 Basic revision of particle HP classes
34//
35
37
38#include "G4Element.hh"
39#include "G4HadFinalState.hh"
40#include "G4HadProjectile.hh"
44#include "G4Neutron.hh"
45
47{
48 nChannels = n;
49 theChannels = new G4ParticleHPChannel*[n];
50 theProjectile = (nullptr == p) ? G4Neutron::Neutron() : p;
51}
52
54{
55 if (theChannels != nullptr) {
56 for (G4int i = 0; i < nChannels; ++i) {
57 delete theChannels[i];
58 }
59 delete[] theChannels;
60 }
61}
62
64 const G4HadProjectile& aTrack)
65{
67 G4int i, ii;
68
69 // decide on the isotope
70 G4int numberOfIsos(0);
71 for (ii = 0; ii < nChannels; ++ii) {
72 numberOfIsos = theChannels[ii]->GetNiso();
73 if (numberOfIsos != 0) break;
74 }
75 auto running = new G4double[numberOfIsos];
76 running[0] = 0;
77 for (i = 0; i < numberOfIsos; i++) {
78 if (i != 0) running[i] = running[i - 1];
79 for (ii = 0; ii < nChannels; ii++) {
80 if (theChannels[ii]->HasAnyData(i)) {
81 running[i] += theChannels[ii]->GetWeightedXsec(
82 aThermalE.GetThermalEnergy(aTrack, theChannels[ii]->GetN(i), theChannels[ii]->GetZ(i),
83 aTrack.GetMaterial()->GetTemperature()), i);
84 }
85 }
86 }
87 G4int isotope = nChannels - 1;
88 G4double random = G4UniformRand();
89 for (i = 0; i < numberOfIsos; ++i) {
90 isotope = i;
91 if (running[numberOfIsos - 1] != 0)
92 if (random < running[i] / running[numberOfIsos - 1]) break;
93 }
94 delete[] running;
95
96 // decide on the channel
97 running = new G4double[nChannels];
98 running[0] = 0;
99 G4int targA = -1; // For production of unChanged
100 G4int targZ = -1;
101 for (i = 0; i < nChannels; i++) {
102 if (i != 0) running[i] = running[i - 1];
103 if (theChannels[i]->HasAnyData(isotope)) {
104 targA = (G4int)theChannels[i]->GetN(isotope); // Will be simply used the last valid value
105 targZ = (G4int)theChannels[i]->GetZ(isotope);
106 running[i] += theChannels[i]->GetFSCrossSection(
107 aThermalE.GetThermalEnergy(aTrack, targA, targZ, aTrack.GetMaterial()->GetTemperature()),
108 isotope);
109 targA = (G4int)theChannels[i]->GetN(isotope); // Will be simply used the last valid value
110 targZ = (G4int)theChannels[i]->GetZ(isotope);
111 // G4cout << " G4ParticleHPChannelList " << i << " targA " << targA << " targZ " << targZ <<
112 //" running " << running[i] << G4endl;
113 }
114 }
115
116 // TK120607
117 if (running[nChannels - 1] == 0) {
118 // This happened usually by the miss match between the cross section data and model
119 if (targA == -1 && targZ == -1) {
121 __FILE__, __LINE__,
122 "ParticleHP model encounter lethal discrepancy with cross section data");
123 }
124
125 // TK121106
126 G4cout << "Warning from NeutronHP: could not find proper reaction channel. This may cause by "
127 "inconsistency between cross section and model. Unchanged final states are returned."
128 << G4endl;
129 unChanged.Clear();
130
131 // For Ep Check create unchanged final state including rest target
132 G4ParticleDefinition* targ_pd = G4IonTable::GetIonTable()->GetIon(targZ, targA, 0.0);
133 auto targ_dp = new G4DynamicParticle(targ_pd, G4ThreeVector(1, 0, 0), 0.0);
134 unChanged.SetEnergyChange(aTrack.GetKineticEnergy());
135 unChanged.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
136 unChanged.AddSecondary(targ_dp);
137 // TK121106
140 delete[] running;
141 return &unChanged;
142 }
143 // TK120607
144
145 G4int lChan = 0;
146 random = G4UniformRand();
147 for (i = 0; i < nChannels; i++) {
148 lChan = i;
149 if (running[nChannels - 1] != 0)
150 if (random < running[i] / running[nChannels - 1]) break;
151 }
152 delete[] running;
153#ifdef G4PHPDEBUG
154 if (G4ParticleHPManager::GetInstance()->GetDEBUG())
155 G4cout << " G4ParticleHPChannelList SELECTED ISOTOPE " << isotope << " SELECTED CHANNEL "
156 << lChan << G4endl;
157#endif
158 return theChannels[lChan]->ApplyYourself(aTrack, isotope);
159}
160
161void G4ParticleHPChannelList::Init(G4Element* anElement, const G4String& dirName,
162 G4ParticleDefinition* projectile)
163{
164 theDir = dirName;
165 // G4cout << theDir << G4endl;
166 theElement = anElement;
167 // G4cout << theElement->GetName() << G4endl;
168 theProjectile = projectile;
169}
170
172{
173 theChannels[idx] = new G4ParticleHPChannel(theProjectile);
174 theChannels[idx]->Init(theElement, theDir, aName);
175 theChannels[idx]->Register(theFS);
176 ++idx;
177}
178
180{
181 G4cout << "================================================================" << G4endl;
182 G4cout << " Element: " << theElement->GetName() << G4endl;
183 G4cout << " Number of channels: " << nChannels << G4endl;
184 G4cout << " Projectile: " << theProjectile->GetParticleName() << G4endl;
185 G4cout << " Directory name: " << theDir << G4endl;
186 for (G4int i = 0; i < nChannels; ++i) {
187 if (theChannels[i]->HasDataInAnyFinalState()) {
188 G4cout << "----------------------------------------------------------------" << G4endl;
189 theChannels[i]->DumpInfo();
190 G4cout << "----------------------------------------------------------------" << G4endl;
191 }
192 }
193 G4cout << "================================================================" << G4endl;
194}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector unit() const
Hep3Vector vect() const
const G4String & GetName() const
Definition G4Element.hh:115
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4Material * GetMaterial() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
static G4IonTable * GetIonTable()
G4double GetTemperature() const
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
const G4String & GetParticleName() const
G4ParticleHPChannelList(G4int n=0, G4ParticleDefinition *p=nullptr)
G4HadFinalState * ApplyYourself(const G4Element *theElement, const G4HadProjectile &aTrack)
void Register(G4ParticleHPFinalState *theFS, const G4String &aName)
G4double GetZ(G4int i) const
G4double GetN(G4int i) const
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack, G4int isoNumber=-1, G4bool isElastic=false)
G4double GetWeightedXsec(G4double energy, G4int isoNumber)
G4bool Register(G4ParticleHPFinalState *theFS)
void Init(G4Element *theElement, const G4String &dirName)
G4double GetFSCrossSection(G4double energy, G4int isoNumber)
static G4ParticleHPManager * GetInstance()
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)