Geant4 11.3.0
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) const
 
G4double GetWeightedXsec (G4double energy, G4int isoNumber) const
 
G4double GetFSCrossSection (G4double energy, G4int isoNumber) const
 
G4bool IsActive (G4int isoNumber) const
 
G4bool HasFSData (G4int isoNumber) const
 
G4bool HasAnyData (G4int isoNumber) const
 
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 () const
 
void DumpInfo () const
 
const G4StringGetFSType ()
 
G4ParticleHPFinalState ** GetFinalStates () const
 
G4WendtFissionFragmentGeneratorGetWendtFissionGenerator () 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{
55 if (fManager->GetUseWendtFissionModel()) {
56 wendtFissionGenerator = G4WendtFissionFragmentGenerator::GetInstance();
57 // Make sure both fission fragment models are not active at same time
58 fManager->SetProduceFissionFragments(false);
59 }
60 theProjectile = (nullptr == p) ? G4Neutron::Neutron() : p;
61 theChannelData = new G4ParticleHPVector;
62}
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
G4ParticleHPManager * fManager
static G4ParticleHPManager * GetInstance()
static G4WendtFissionFragmentGenerator * GetInstance()

Referenced by G4ParticleHPChannel(), and operator=().

◆ ~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 232 of file G4ParticleHPChannel.cc.

234{
235 //G4cout << "G4ParticleHPChannel::ApplyYourself niso=" << niso
236 // << " ni=" << anIsotope << " isElastic=" << isElastic <<G4endl;
237 if (anIsotope != -1 && anIsotope != -2) {
238 // Inelastic Case
239 //G4cout << "G4ParticleHPChannel Inelastic Case"
240 //<< " Z= " << GetZ(anIsotope) << " A = " << GetN(anIsotope) << G4endl;
241 fManager->GetReactionWhiteBoard()->SetTargA((G4int)GetN(anIsotope));
242 fManager->GetReactionWhiteBoard()->SetTargZ((G4int)GetZ(anIsotope));
243 return theFinalStates[anIsotope]->ApplyYourself(theTrack);
244 }
245 G4double sum = 0;
246 G4int it = 0;
247 auto xsec = new G4double[niso];
248 G4ParticleHPThermalBoost aThermalE;
249 for (G4int i = 0; i < niso; i++) {
250 if (theFinalStates[i]->HasAnyData()) {
251 /*
252 G4cout << "FS: " << i << theTrack.GetDefinition()->GetParticleName()
253 << " Z=" << theFinalStates[i]->GetZ()
254 << " A=" << theFinalStates[i]->GetN()
255 << G4endl;
256 */
257 xsec[i] = theIsotopeWiseData[i].GetXsec(
258 aThermalE.GetThermalEnergy(theTrack, theFinalStates[i]->GetN(),
259 theFinalStates[i]->GetZ(),
260 theTrack.GetMaterial()->GetTemperature()));
261 sum += xsec[i];
262 }
263 else {
264 xsec[i] = 0;
265 }
266 }
267 if (sum == 0) {
268 it = G4lrint(niso * G4UniformRand());
269 }
270 else {
271 G4double random = G4UniformRand();
272 G4double running = 0;
273 for (G4int ix = 0; ix < niso; ix++) {
274 running += xsec[ix];
275 if (sum == 0 || random <= running / sum) {
276 it = ix;
277 break;
278 }
279 }
280 if (it == niso) it--;
281 }
282 delete[] xsec;
283 G4HadFinalState* theFinalState = nullptr;
284 const auto A = (G4int)this->GetN(it);
285 const auto Z = (G4int)this->GetZ(it);
286 const auto M = (G4int)this->GetM(it);
287
288 //-2:Marker for Fission
289 if ((wendtFissionGenerator != nullptr) && anIsotope == -2) {
290 theFinalState = wendtFissionGenerator->ApplyYourself(theTrack, Z, A);
291 }
292
293 // Use the standard procedure if the G4FissionFragmentGenerator model fails
294 if (theFinalState == nullptr) {
295 G4int icounter = 0;
296 G4int icounter_max = 1024;
297 while (theFinalState == nullptr) // Loop checking, 11.05.2015, T. Koi
298 {
299 icounter++;
300 if (icounter > icounter_max) {
301 G4cout << "Loop-counter exceeded the threshold value at "
302 << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
303 break;
304 }
305 if (isElastic) {
306 // Register 0 K cross-section for DBRC for Doppler broadened elastic scattering kernel
307 G4ParticleHPVector* xsforFS = theIsotopeWiseData[it].MakeChannelData();
308 // Only G4ParticleHPElasticFS has the RegisterCrossSection method
309 static_cast<G4ParticleHPElasticFS*>(theFinalStates[it])->RegisterCrossSection(xsforFS);
310 }
311 theFinalState = theFinalStates[it]->ApplyYourself(theTrack);
312 }
313 }
314
315 // G4cout <<"THE IMPORTANT RETURN"<<G4endl;
316 // G4cout << "TK G4ParticleHPChannel Elastic, Capture and Fission Cases "
317 //<< " Z= " << this->GetZ(it) << " A = " << this->GetN(it) << G4endl;
318 fManager->GetReactionWhiteBoard()->SetTargA(A);
319 fManager->GetReactionWhiteBoard()->SetTargZ(Z);
320 fManager->GetReactionWhiteBoard()->SetTargM(M);
321
322 return theFinalState;
323}
#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) const
G4double GetN(G4int i) const
G4double GetM(G4int i) const
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
int G4lrint(double ad)
Definition templates.hh:134

Referenced by G4NeutronHPElasticVI::ApplyYourself().

◆ DumpInfo()

void G4ParticleHPChannel::DumpInfo ( ) const

Definition at line 325 of file G4ParticleHPChannel.cc.

326{
327 G4cout << " Element: " << theElement->GetName() << G4endl;
328 G4cout << " Directory name: " << theDir << G4endl;
329 G4cout << " FS name: " << theFSType << G4endl;
330 G4cout << " Number of Isotopes: " << niso << G4endl;
331 G4cout << " Have cross sections: " << G4endl;
332 for (int i = 0; i < niso; i++) {
333 G4cout << theFinalStates[i]->HasXsec() << " ";
334 }
335 G4cout << G4endl;
336 if (theChannelData != nullptr) {
337 G4cout << " Cross Section (total for this channel):" << G4endl;
338 int np = theChannelData->GetVectorLength();
339 G4cout << np << G4endl;
340 for (int i = 0; i < np; i++) {
341 G4cout << theChannelData->GetEnergy(i) / eV << " " << theChannelData->GetXsec(i) << G4endl;
342 }
343 }
344}

◆ GetFinalStates()

G4ParticleHPFinalState ** G4ParticleHPChannel::GetFinalStates ( ) const
inline

Definition at line 129 of file G4ParticleHPChannel.hh.

129{ return theFinalStates; }

◆ GetFSCrossSection()

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

Definition at line 90 of file G4ParticleHPChannel.cc.

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

◆ GetFSType()

const G4String & G4ParticleHPChannel::GetFSType ( )
inline

Definition at line 127 of file G4ParticleHPChannel.hh.

127{ return theFSType; }

◆ GetM()

G4double G4ParticleHPChannel::GetM ( G4int i) const
inline

Definition at line 111 of file G4ParticleHPChannel.hh.

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

Referenced by ApplyYourself().

◆ GetN()

G4double G4ParticleHPChannel::GetN ( G4int i) const
inline

Definition at line 109 of file G4ParticleHPChannel.hh.

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

Referenced by ApplyYourself().

◆ GetNiso()

G4int G4ParticleHPChannel::GetNiso ( ) const
inline

Definition at line 107 of file G4ParticleHPChannel.hh.

107{ return niso; }

◆ GetWeightedXsec()

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

Definition at line 84 of file G4ParticleHPChannel.cc.

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

◆ GetWendtFissionGenerator()

G4WendtFissionFragmentGenerator * G4ParticleHPChannel::GetWendtFissionGenerator ( ) const

Definition at line 226 of file G4ParticleHPChannel.cc.

226 {
227 if ( wendtFissionGenerator ) return wendtFissionGenerator;
228 else return nullptr;
229}

◆ GetXsec()

G4double G4ParticleHPChannel::GetXsec ( G4double energy) const

Definition at line 79 of file G4ParticleHPChannel.cc.

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

◆ GetZ()

G4double G4ParticleHPChannel::GetZ ( G4int i) const
inline

Definition at line 110 of file G4ParticleHPChannel.hh.

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

Referenced by 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;
186 G4ParticleHPVector* tmp;
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}
G4double GetXsec(G4int i)
G4double GetEnergy(G4int i) const
G4int GetVectorLength() const

Referenced by UpdateData().

◆ HasAnyData()

G4bool G4ParticleHPChannel::HasAnyData ( G4int isoNumber) const
inline

Definition at line 80 of file G4ParticleHPChannel.hh.

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

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

◆ HasDataInAnyFinalState()

G4bool G4ParticleHPChannel::HasDataInAnyFinalState ( ) const
inline

Definition at line 113 of file G4ParticleHPChannel.hh.

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

Referenced by Register().

◆ HasFSData()

G4bool G4ParticleHPChannel::HasFSData ( G4int isoNumber) const
inline

Definition at line 75 of file G4ParticleHPChannel.hh.

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

◆ 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 Init().

◆ 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 70 of file G4ParticleHPChannel.hh.

71 {
72 return active[isoNumber];
73 }

◆ 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}
void UpdateData(G4int A, G4int Z, G4int index, G4double abundance, G4ParticleDefinition *projectile)
G4bool HasDataInAnyFinalState() const
virtual G4ParticleHPFinalState * New()=0
void SetProjectile(G4ParticleDefinition *projectile)

◆ UpdateData() [1/2]

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

Definition at line 92 of file G4ParticleHPChannel.hh.

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

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 const 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)

Member Data Documentation

◆ fManager

G4ParticleHPManager* G4ParticleHPChannel::fManager
protected

Definition at line 139 of file G4ParticleHPChannel.hh.

Referenced by ApplyYourself(), and G4ParticleHPChannel().


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