Geant4 11.4.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 234 of file G4ParticleHPChannel.cc.

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

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

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

228 {
229 if ( wendtFissionGenerator ) return wendtFissionGenerator;
230 else return nullptr;
231}

◆ 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 xa = std::max(xa, 0.0);
198 xp = std::max(xp, 0.0);
199 if (std::abs(xp - xa) < 0.001*xa) {
200 ++p;
201 }
202 }
203 else {
204 tmp = anActive;
205 t = a;
206 anActive = aPassive;
207 a = p;
208 aPassive = tmp;
209 p = t;
210 }
211 }
212 while (a != anActive->GetVectorLength()) // Loop checking, 11.05.2015, T. Koi
213 {
214 theMerge->SetData(m_tmp++, anActive->GetEnergy(a), anActive->GetXsec(a));
215 ++a;
216 }
217 while (p != aPassive->GetVectorLength()) // Loop checking, 11.05.2015, T. Koi
218 {
219 if (std::abs(theMerge->GetEnergy(std::max(0, m_tmp - 1)) -
220 aPassive->GetEnergy(p)) / aPassive->GetEnergy(p) > 0.001)
221 theMerge->SetData(m_tmp++, aPassive->GetEnergy(p), aPassive->GetXsec(p));
222 ++p;
223 }
224 delete theStore;
225 theStore = theMerge;
226}
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: