Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPChannel Class Reference

#include <G4ParticleHPChannel.hh>

Public Member Functions

 G4ParticleHPChannel (G4ParticleDefinition *projectile=nullptr)
 
 ~G4ParticleHPChannel ()
 
G4double GetXsec (G4double energy)
 
G4double GetWeightedXsec (G4double energy, G4int isoNumber)
 
G4double GetFSCrossSection (G4double energy, G4int isoNumber)
 
G4bool IsActive (G4int isoNumber) const
 
G4bool HasFSData (G4int isoNumber)
 
G4bool HasAnyData (G4int isoNumber)
 
G4bool Register (G4ParticleHPFinalState *theFS)
 
void Init (G4Element *theElement, const G4String &dirName)
 
void Init (G4Element *theElement, const G4String &dirName, const G4String &fsType)
 
void UpdateData (G4int A, G4int Z, G4int index, G4double abundance, G4ParticleDefinition *projectile)
 
void UpdateData (G4int A, G4int Z, G4int M, G4int index, G4double abundance, G4ParticleDefinition *projectile)
 
void Harmonise (G4ParticleHPVector *&theStore, G4ParticleHPVector *theNew)
 
G4HadFinalStateApplyYourself (const G4HadProjectile &theTrack, G4int isoNumber=-1, G4bool isElastic=false)
 
G4int GetNiso () const
 
G4double GetN (G4int i) const
 
G4double GetZ (G4int i) const
 
G4double GetM (G4int i) const
 
G4bool HasDataInAnyFinalState ()
 
void DumpInfo ()
 
G4StringGetFSType ()
 
G4ParticleHPFinalState ** GetFinalStates () const
 
 G4ParticleHPChannel (G4ParticleHPChannel &)=delete
 
G4ParticleHPChanneloperator= (const G4ParticleHPChannel &right)=delete
 

Protected Attributes

G4ParticleHPManagerfManager
 

Detailed Description

Definition at line 56 of file G4ParticleHPChannel.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPChannel() [1/2]

G4ParticleHPChannel::G4ParticleHPChannel ( G4ParticleDefinition * projectile = nullptr)

Definition at line 52 of file G4ParticleHPChannel.cc.

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}
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
G4ParticleHPManager * fManager
void SetProduceFissionFragments(G4bool val)
G4bool GetUseWendtFissionModel() const
static G4ParticleHPManager * GetInstance()
static G4WendtFissionFragmentGenerator * GetInstance()

◆ ~G4ParticleHPChannel()

G4ParticleHPChannel::~G4ParticleHPChannel ( )

Definition at line 64 of file G4ParticleHPChannel.cc.

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}
int G4int
Definition G4Types.hh:85

◆ G4ParticleHPChannel() [2/2]

G4ParticleHPChannel::G4ParticleHPChannel ( G4ParticleHPChannel & )
delete

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4ParticleHPChannel::ApplyYourself ( const G4HadProjectile & theTrack,
G4int isoNumber = -1,
G4bool isElastic = false )

Definition at line 227 of file G4ParticleHPChannel.cc.

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}
#define M(row, col)
double G4double
Definition G4Types.hh:83
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
const G4Material * GetMaterial() const
G4double GetTemperature() const
G4double GetZ(G4int i) const
G4bool HasAnyData(G4int isoNumber)
G4double GetN(G4int i) const
G4double GetM(G4int i) const
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &)
G4ParticleHPVector * MakeChannelData()
G4double GetXsec(G4double energy)
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
G4HadFinalState * ApplyYourself(const G4HadProjectile &projectile, G4int Z, G4int A)
int G4lrint(double ad)
Definition templates.hh:134

Referenced by G4FissLib::ApplyYourself(), G4NeutronHPElasticVI::ApplyYourself(), and G4ParticleHPChannelList::ApplyYourself().

◆ DumpInfo()

void G4ParticleHPChannel::DumpInfo ( )

Definition at line 320 of file G4ParticleHPChannel.cc.

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}
const G4String & GetName() const
Definition G4Element.hh:115
G4double GetXsec(G4int i)
G4double GetEnergy(G4int i) const
G4int GetVectorLength() const

Referenced by G4ParticleHPChannelList::DumpInfo().

◆ GetFinalStates()

G4ParticleHPFinalState ** G4ParticleHPChannel::GetFinalStates ( ) const
inline

Definition at line 128 of file G4ParticleHPChannel.hh.

128{ return theFinalStates; }

◆ GetFSCrossSection()

G4double G4ParticleHPChannel::GetFSCrossSection ( G4double energy,
G4int isoNumber )

Definition at line 90 of file G4ParticleHPChannel.cc.

92{
93 return theFinalStates[isoNumber]->GetXsec(energy);
94}
virtual G4double GetXsec(G4double)

Referenced by G4ParticleHPChannelList::ApplyYourself().

◆ GetFSType()

G4String & G4ParticleHPChannel::GetFSType ( )
inline

Definition at line 126 of file G4ParticleHPChannel.hh.

126{ return theFSType; }

◆ GetM()

G4double G4ParticleHPChannel::GetM ( G4int i) const
inline

Definition at line 110 of file G4ParticleHPChannel.hh.

110{ return theFinalStates[i]->GetM(); }

Referenced by ApplyYourself().

◆ GetN()

G4double G4ParticleHPChannel::GetN ( G4int i) const
inline

Definition at line 108 of file G4ParticleHPChannel.hh.

108{ return theFinalStates[i]->GetN(); }

Referenced by ApplyYourself(), and G4ParticleHPChannelList::ApplyYourself().

◆ GetNiso()

G4int G4ParticleHPChannel::GetNiso ( ) const
inline

Definition at line 106 of file G4ParticleHPChannel.hh.

106{ return niso; }

Referenced by G4ParticleHPChannelList::ApplyYourself().

◆ GetWeightedXsec()

G4double G4ParticleHPChannel::GetWeightedXsec ( G4double energy,
G4int isoNumber )

Definition at line 84 of file G4ParticleHPChannel.cc.

86{
87 return theIsotopeWiseData[isoNumber].GetXsec(energy);
88}

Referenced by G4ParticleHPChannelList::ApplyYourself().

◆ GetXsec()

G4double G4ParticleHPChannel::GetXsec ( G4double energy)

Definition at line 79 of file G4ParticleHPChannel.cc.

80{
81 return std::max(0., theChannelData->GetXsec(energy));
82}

Referenced by G4FissLib::ApplyYourself().

◆ GetZ()

G4double G4ParticleHPChannel::GetZ ( G4int i) const
inline

Definition at line 109 of file G4ParticleHPChannel.hh.

109{ return theFinalStates[i]->GetZ(); }

Referenced by ApplyYourself(), and G4ParticleHPChannelList::ApplyYourself().

◆ Harmonise()

void G4ParticleHPChannel::Harmonise ( G4ParticleHPVector *& theStore,
G4ParticleHPVector * theNew )

Definition at line 179 of file G4ParticleHPChannel.cc.

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}

Referenced by UpdateData().

◆ HasAnyData()

G4bool G4ParticleHPChannel::HasAnyData ( G4int isoNumber)
inline

Definition at line 79 of file G4ParticleHPChannel.hh.

80 {
81 return theFinalStates[isoNumber]->HasAnyData();
82 }

Referenced by ApplyYourself(), HasDataInAnyFinalState(), and UpdateData().

◆ HasDataInAnyFinalState()

G4bool G4ParticleHPChannel::HasDataInAnyFinalState ( )
inline

Definition at line 112 of file G4ParticleHPChannel.hh.

113 {
114 G4bool result = false;
115 for (G4int i = 0; i < niso; ++i) {
116 if (theFinalStates[i]->HasAnyData()) {
117 result = true;
118 break;
119 }
120 }
121 return result;
122 }
bool G4bool
Definition G4Types.hh:86

Referenced by Register().

◆ HasFSData()

G4bool G4ParticleHPChannel::HasFSData ( G4int isoNumber)
inline

Definition at line 74 of file G4ParticleHPChannel.hh.

75 {
76 return theFinalStates[isoNumber]->HasFSData();
77 }

◆ Init() [1/2]

void G4ParticleHPChannel::Init ( G4Element * theElement,
const G4String & dirName )

Definition at line 103 of file G4ParticleHPChannel.cc.

104{
105 theDir = dirName;
106 theElement = anElement;
107}

Referenced by G4FissLib::G4FissLib(), Init(), and G4ParticleHPChannelList::Register().

◆ Init() [2/2]

void G4ParticleHPChannel::Init ( G4Element * theElement,
const G4String & dirName,
const G4String & fsType )

Definition at line 96 of file G4ParticleHPChannel.cc.

98{
99 theFSType = aFSType;
100 Init(anElement, dirName);
101}
void Init(G4Element *theElement, const G4String &dirName)

◆ IsActive()

G4bool G4ParticleHPChannel::IsActive ( G4int isoNumber) const
inline

Definition at line 69 of file G4ParticleHPChannel.hh.

70 {
71 return active[isoNumber];
72 }

◆ operator=()

G4ParticleHPChannel & G4ParticleHPChannel::operator= ( const G4ParticleHPChannel & right)
delete

◆ Register()

G4bool G4ParticleHPChannel::Register ( G4ParticleHPFinalState * theFS)

Definition at line 109 of file G4ParticleHPChannel.cc.

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}
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
G4int GetZasInt() const
Definition G4Element.hh:120
G4int Getm() const
Definition G4Isotope.hh:89
G4int GetN() const
Definition G4Isotope.hh:83
void UpdateData(G4int A, G4int Z, G4int index, G4double abundance, G4ParticleDefinition *projectile)
void SetA_Z(G4double anA, G4double aZ, G4int aM=0)
virtual G4ParticleHPFinalState * New()=0
void SetProjectile(G4ParticleDefinition *projectile)

Referenced by G4FissLib::G4FissLib(), and G4ParticleHPChannelList::Register().

◆ UpdateData() [1/2]

void G4ParticleHPChannel::UpdateData ( G4int A,
G4int Z,
G4int index,
G4double abundance,
G4ParticleDefinition * projectile )
inline

Definition at line 91 of file G4ParticleHPChannel.hh.

93 {
94 UpdateData(A, Z, 0, index, abundance, projectile);
95 }

Referenced by Register(), and UpdateData().

◆ UpdateData() [2/2]

void G4ParticleHPChannel::UpdateData ( G4int A,
G4int Z,
G4int M,
G4int index,
G4double abundance,
G4ParticleDefinition * projectile )

Definition at line 149 of file G4ParticleHPChannel.cc.

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}
void Harmonise(G4ParticleHPVector *&theStore, G4ParticleHPVector *theNew)
void Init(G4double A, G4double Z, G4String &dirName, G4String &aFSType, G4ParticleDefinition *p)
G4bool Init(G4int A, G4int Z, G4double abun, G4String dirName, G4String aFSType)
void FillChannelData(G4ParticleHPVector *aBuffer)
void Times(G4double factor)
void InitializeANucleus(const G4int A, const G4int Z, const G4int M, const G4String &dataDirectory)

Member Data Documentation

◆ fManager

G4ParticleHPManager* G4ParticleHPChannel::fManager
protected

Definition at line 135 of file G4ParticleHPChannel.hh.

Referenced by ApplyYourself(), and G4ParticleHPChannel().


The documentation for this class was generated from the following files: