Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPPhotonDist.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 Try to limit sum of secondary photon energy while keeping distribution shape
31// in the of nDiscrete = 1 an nPartial = 1. Most case are satisfied.
32// T. Koi
33// 070606 Add Partial case by T. Koi
34// 070618 fix memory leaking by T. Koi
35// 080801 fix memory leaking by T. Koi
36// 080801 Correcting data disorder which happened when both InitPartial
37// and InitAnglurar methods was called in a same instance by T. Koi
38// 090514 Fix bug in IC electron emission case
39// Contribution from Chao Zhang ([email protected]) and Dongming Mei([email protected])
40// But it looks like never cause real effect in G4NDL3.13 (at least Natural elements) TK
41// 101111 Change warning message for "repFlag == 2 && isoFlag != 1" case
42//
43// there is a lot of unused (and undebugged) code in this file. Kept for the moment just in case. @@
44// P. Arce, June-2014 Conversion neutron_hp to particle_hp
45//
47
48#include "G4Electron.hh"
51#include "G4Poisson.hh"
52#include "G4SystemOfUnits.hh"
53
54#include <numeric>
55
57{
58 G4bool result = true;
59 if (aDataFile >> repFlag) {
60 aDataFile >> targetMass;
61 if (repFlag == 1) {
62 // multiplicities
63 aDataFile >> nDiscrete;
64 const std::size_t msize = nDiscrete > 0 ? nDiscrete : 1;
65 disType = new G4int[msize];
66 energy = new G4double[msize];
67 // actualMult = new G4int[msize];
68 theYield = new G4ParticleHPVector[msize];
69 for (std::size_t i = 0; i < msize; ++i) {
70 aDataFile >> disType[i] >> energy[i];
71 energy[i] *= eV;
72 theYield[i].Init(aDataFile, eV);
73 }
74 }
75 else if (repFlag == 2) {
76 aDataFile >> theInternalConversionFlag;
77 aDataFile >> theBaseEnergy;
78 theBaseEnergy *= eV;
79 aDataFile >> theInternalConversionFlag;
80 aDataFile >> nGammaEnergies;
81 const std::size_t esize = nGammaEnergies > 0 ? nGammaEnergies : 1;
82 theLevelEnergies = new G4double[esize];
83 theTransitionProbabilities = new G4double[esize];
84 if (theInternalConversionFlag == 2)
85 thePhotonTransitionFraction = new G4double[esize];
86 for (std::size_t ii = 0; ii < esize; ++ii) {
87 if (theInternalConversionFlag == 1) {
88 aDataFile >> theLevelEnergies[ii] >> theTransitionProbabilities[ii];
89 theLevelEnergies[ii] *= eV;
90 }
91 else if (theInternalConversionFlag == 2) {
92 aDataFile >> theLevelEnergies[ii] >> theTransitionProbabilities[ii]
93 >> thePhotonTransitionFraction[ii];
94 theLevelEnergies[ii] *= eV;
95 }
96 else {
97 throw G4HadronicException(__FILE__, __LINE__,
98 "G4ParticleHPPhotonDist: Unknown conversion flag");
99 }
100 }
101 }
102 else {
103 G4cout << "Data representation in G4ParticleHPPhotonDist: " << repFlag << G4endl;
105 __FILE__, __LINE__, "G4ParticleHPPhotonDist: This data representation is not implemented.");
106 }
107 }
108 else {
109 result = false;
110 }
111 return result;
112}
113
114void G4ParticleHPPhotonDist::InitAngular(std::istream& aDataFile)
115{
116 G4int i, ii;
117 // angular distributions
118 aDataFile >> isoFlag;
119 if (isoFlag != 1) {
120 if (repFlag == 2)
121 G4cout << "G4ParticleHPPhotonDist: repFlag == 2 && isoFlag != 1 is unexpected! If you use "
122 "G4ND3.x, then please report to Geant4 HyperNews. "
123 << G4endl;
124 aDataFile >> tabulationType >> nDiscrete2 >> nIso;
125 if (theGammas != nullptr && nDiscrete2 != nDiscrete)
126 G4cout << "080731c G4ParticleHPPhotonDist nDiscrete2 != nDiscrete, It looks like something "
127 "wrong in your NDL files. Please update the latest. If you still have this "
128 "messages after the update, then please report to Geant4 Hyper News."
129 << G4endl;
130
131 // The order of cross section (InitPartials) and distribution
132 // (InitAngular here) data are different, we have to re-coordinate
133 // consistent data order.
134 std::vector<G4double> vct_gammas_par;
135 std::vector<G4double> vct_shells_par;
136 std::vector<G4int> vct_primary_par;
137 std::vector<G4int> vct_distype_par;
138 std::vector<G4ParticleHPVector*> vct_pXS_par;
139 if (theGammas != nullptr && theShells != nullptr) {
140 // copy the cross section data
141 for (i = 0; i < nDiscrete; ++i) {
142 vct_gammas_par.push_back(theGammas[i]);
143 vct_shells_par.push_back(theShells[i]);
144 vct_primary_par.push_back(isPrimary[i]);
145 vct_distype_par.push_back(disType[i]);
146 auto hpv = new G4ParticleHPVector;
147 *hpv = thePartialXsec[i];
148 vct_pXS_par.push_back(hpv);
149 }
150 }
151 const std::size_t psize = nDiscrete2 > 0 ? nDiscrete2 : 1;
152 if (theGammas == nullptr) theGammas = new G4double[psize];
153 if (theShells == nullptr) theShells = new G4double[psize];
154
155 for (i = 0; i < nIso; ++i) // isotropic photons
156 {
157 aDataFile >> theGammas[i] >> theShells[i];
158 theGammas[i] *= eV;
159 theShells[i] *= eV;
160 }
161 const std::size_t tsize = nDiscrete2 - nIso > 0 ? nDiscrete2 - nIso : 1;
162 nNeu = new G4int[tsize];
163 if (tabulationType == 1) theLegendre = new G4ParticleHPLegendreTable*[tsize];
164 if (tabulationType == 2) theAngular = new G4ParticleHPAngularP*[tsize];
165 for (i = nIso; i < nDiscrete2; ++i) {
166 if (tabulationType == 1) {
167 aDataFile >> theGammas[i] >> theShells[i] >> nNeu[i - nIso];
168 theGammas[i] *= eV;
169 theShells[i] *= eV;
170 const std::size_t lsize = nNeu[i - nIso] > 0 ? nNeu[i - nIso] : 1;
171 theLegendre[i - nIso] = new G4ParticleHPLegendreTable[lsize];
172 theLegendreManager.Init(aDataFile);
173 for (ii = 0; ii < nNeu[i - nIso]; ++ii) {
174 theLegendre[i - nIso][ii].Init(aDataFile);
175 }
176 }
177 else if (tabulationType == 2) {
178 aDataFile >> theGammas[i] >> theShells[i] >> nNeu[i - nIso];
179 theGammas[i] *= eV;
180 theShells[i] *= eV;
181 const std::size_t asize = nNeu[i - nIso] > 0 ? nNeu[i - nIso] : 1;
182 theAngular[i - nIso] = new G4ParticleHPAngularP[asize];
183 for (ii = 0; ii < nNeu[i - nIso]; ++ii) {
184 theAngular[i - nIso][ii].Init(aDataFile);
185 }
186 }
187 else {
188 G4cout << "tabulation type: tabulationType" << G4endl;
190 __FILE__, __LINE__, "cannot deal with this tabulation type for angular distributions.");
191 }
192 }
193
194 if (!vct_gammas_par.empty()) {
195 // Reordering cross section data to corrsponding distribution data
196 for (i = 0; i < nDiscrete; ++i) {
197 for (G4int j = 0; j < nDiscrete; ++j) {
198 // Checking gamma and shell to identification
199 if (theGammas[i] == vct_gammas_par[j] && theShells[i] == vct_shells_par[j]) {
200 isPrimary[i] = vct_primary_par[j];
201 disType[i] = vct_distype_par[j];
202 thePartialXsec[i] = (*(vct_pXS_par[j]));
203 }
204 }
205 }
206 // Garbage collection
207 for (auto it = vct_pXS_par.cbegin(); it != vct_pXS_par.cend(); ++it) {
208 delete *it;
209 }
210 }
211 }
212}
213
214void G4ParticleHPPhotonDist::InitEnergies(std::istream& aDataFile)
215{
216 G4int i, energyDistributionsNeeded = 0;
217 for (i = 0; i < nDiscrete; ++i) {
218 if (disType[i] == 1) energyDistributionsNeeded = 1;
219 }
220 if (energyDistributionsNeeded == 0) return;
221 aDataFile >> nPartials;
222 const std::size_t dsize = nPartials > 0 ? nPartials : 1;
223 distribution = new G4int[dsize];
224 probs = new G4ParticleHPVector[dsize];
225 partials = new G4ParticleHPPartial*[dsize];
226 G4int nen;
227 G4int dummy;
228 for (i = 0; i < nPartials; ++i) {
229 aDataFile >> dummy;
230 probs[i].Init(aDataFile, eV);
231 aDataFile >> nen;
232 partials[i] = new G4ParticleHPPartial(nen);
233 partials[i]->InitInterpolation(aDataFile);
234 partials[i]->Init(aDataFile);
235 }
236}
237
238void G4ParticleHPPhotonDist::InitPartials(std::istream& aDataFile, G4ParticleHPVector* theXsec)
239{
240 if (theXsec != nullptr) theReactionXsec = theXsec;
241
242 aDataFile >> nDiscrete >> targetMass;
243 if (nDiscrete != 1) {
244 theTotalXsec.Init(aDataFile, eV);
245 }
246 const std::size_t dsize = nDiscrete > 0 ? nDiscrete : 1;
247 theGammas = new G4double[dsize];
248 theShells = new G4double[dsize];
249 isPrimary = new G4int[dsize];
250 disType = new G4int[dsize];
251 thePartialXsec = new G4ParticleHPVector[dsize];
252 for (std::size_t i = 0; i < dsize; ++i) {
253 aDataFile >> theGammas[i] >> theShells[i] >> isPrimary[i] >> disType[i];
254 theGammas[i] *= eV;
255 theShells[i] *= eV;
256 thePartialXsec[i].Init(aDataFile, eV);
257 }
258}
259
261{
262 // the partial cross-section case is not all in this yet.
263 if (actualMult.Get() == nullptr) {
264 actualMult.Get() = new std::vector<G4int>(nDiscrete);
265 }
266 G4int i, ii, iii;
267 G4int nSecondaries = 0;
268 auto thePhotons = new G4ReactionProductVector;
269
270 if (repFlag == 1) {
271 G4double current = 0;
272 for (i = 0; i < nDiscrete; ++i) {
273 current = theYield[i].GetY(anEnergy);
274 actualMult.Get()->at(i) = (G4int)G4Poisson(current); // max cut-off still missing @@@
275 if (nDiscrete == 1 && current < 1.0001) {
276 actualMult.Get()->at(i) = static_cast<G4int>(current);
277 if (current < 1) {
278 actualMult.Get()->at(i) = 0;
279 if (G4UniformRand() < current) actualMult.Get()->at(i) = 1;
280 }
281 }
282 nSecondaries += actualMult.Get()->at(i);
283 }
284 for (i = 0; i < nSecondaries; ++i) {
285 auto theOne = new G4ReactionProduct;
286 theOne->SetDefinition(G4Gamma::Gamma());
287 thePhotons->push_back(theOne);
288 }
289
290 G4int count = 0;
291
292 if (nDiscrete == 1 && nPartials == 1) {
293 if (actualMult.Get()->at(0) > 0) {
294 if (disType[0] == 1) {
295 // continuum
296 G4ParticleHPVector* temp;
297 temp = partials[0]->GetY(anEnergy); //@@@ look at, seems fishy
298 G4double maximumE = temp->GetX(temp->GetVectorLength() - 1); // This is an assumption.
299
300 std::vector<G4double> photons_e_best(actualMult.Get()->at(0), 0.0);
301 G4double best = DBL_MAX;
302 G4int maxTry = 1000;
303 for (G4int j = 0; j < maxTry; ++j) {
304 std::vector<G4double> photons_e(actualMult.Get()->at(0), 0.0);
305 for (auto it = photons_e.begin(); it < photons_e.end(); ++it) {
306 *it = temp->Sample();
307 }
308
309 if (std::accumulate(photons_e.cbegin(), photons_e.cend(), 0.0) > maximumE) {
310 if (std::accumulate(photons_e.cbegin(), photons_e.cend(), 0.0) < best)
311 photons_e_best = photons_e;
312 continue;
313 }
314 G4int iphot = 0;
315 for (auto it = photons_e.cbegin(); it < photons_e.cend(); ++it) {
316 thePhotons->operator[](iphot)->SetKineticEnergy(
317 *it); // Replace index count, which was not incremented,
318 // with iphot, which is, as per Artem Zontikov,
319 // bug report 2167
320 ++iphot;
321 }
322
323 break;
324 }
325 delete temp;
326 }
327 else {
328 // discrete
329 thePhotons->operator[](count)->SetKineticEnergy(energy[i]);
330 }
331 ++count;
332 if (count > nSecondaries)
333 throw G4HadronicException(__FILE__, __LINE__,
334 "G4ParticleHPPhotonDist::GetPhotons inconsistency");
335 }
336 }
337 else { // nDiscrete != 1 or nPartials != 1
338 for (i = 0; i < nDiscrete; ++i) {
339 for (ii = 0; ii < actualMult.Get()->at(i); ++ii) {
340 if (disType[i] == 1) {
341 // continuum
342 G4double sum = 0, run = 0;
343 for (iii = 0; iii < nPartials; ++iii)
344 sum += probs[iii].GetY(anEnergy);
345 G4double random = G4UniformRand();
346 G4int theP = 0;
347 for (iii = 0; iii < nPartials; ++iii) {
348 run += probs[iii].GetY(anEnergy);
349 theP = iii;
350 if (random < run / sum) break;
351 }
352
353 if (theP == nPartials) theP = nPartials - 1; // das sortiert J aus.
354 sum = 0;
355 G4ParticleHPVector* temp;
356 temp = partials[theP]->GetY(anEnergy); //@@@ look at, seems fishy
357 G4double eGamm = temp->Sample();
358 thePhotons->operator[](count)->SetKineticEnergy(eGamm);
359 delete temp;
360 }
361 else {
362 // discrete
363 thePhotons->operator[](count)->SetKineticEnergy(energy[i]);
364 }
365 ++count;
366 if (count > nSecondaries)
367 throw G4HadronicException(__FILE__, __LINE__,
368 "G4ParticleHPPhotonDist::GetPhotons inconsistency");
369 }
370 }
371 }
372
373 // now do the angular distributions...
374 if (isoFlag == 1) {
375 for (i = 0; i < nSecondaries; ++i) {
376 G4double costheta = 2. * G4UniformRand() - 1;
377 G4double theta = std::acos(costheta);
378 G4double phi = twopi * G4UniformRand();
379 G4double sinth = std::sin(theta);
380 G4double en = thePhotons->operator[](i)->GetTotalEnergy();
381 G4ThreeVector temp(en * sinth * std::cos(phi), en * sinth * std::sin(phi),
382 en * std::cos(theta));
383 thePhotons->operator[](i)->SetMomentum(temp);
384 }
385 }
386 else {
387 for (i = 0; i < nSecondaries; ++i) {
388 G4double currentEnergy = thePhotons->operator[](i)->GetTotalEnergy();
389 for (ii = 0; ii < nDiscrete2; ++ii) {
390 if (std::abs(currentEnergy - theGammas[ii]) < 0.1 * keV) break;
391 }
392 if (ii == nDiscrete2)
393 --ii; // fix for what seems an (file12 vs file 14) inconsistency found in the ENDF 7N14
394 // data. @@
395 if (ii < nIso) {
396 // isotropic distribution
397 //
398 // Fix Bugzilla report #1745
399 // G4double theta = pi*G4UniformRand();
400 G4double costheta = 2. * G4UniformRand() - 1;
401 G4double theta = std::acos(costheta);
402 G4double phi = twopi * G4UniformRand();
403 G4double sinth = std::sin(theta);
404 G4double en = thePhotons->operator[](i)->GetTotalEnergy();
405 // DHW G4ThreeVector tempVector(en*sinth*std::cos(phi), en*sinth*std::sin(phi),
406 // en*std::cos(theta) );
407 G4ThreeVector tempVector(en * sinth * std::cos(phi), en * sinth * std::sin(phi),
408 en * costheta);
409 thePhotons->operator[](i)->SetMomentum(tempVector);
410 }
411 else if (tabulationType == 1) {
412 // legendre polynomials
413 G4int it(0);
414 for (iii = 0; iii < nNeu[ii - nIso]; ++iii) // find the neutron energy
415 {
416 it = iii;
417 if (theLegendre[ii - nIso][iii].GetEnergy() > anEnergy) break;
418 }
420 aStore.SetCoeff(1, &(theLegendre[ii - nIso][it]));
421 if (it > 0) {
422 aStore.SetCoeff(0, &(theLegendre[ii - nIso][it - 1]));
423 }
424 else {
425 aStore.SetCoeff(0, &(theLegendre[ii - nIso][it]));
426 }
427 G4double cosTh = aStore.SampleMax(anEnergy);
428 G4double theta = std::acos(cosTh);
429 G4double phi = twopi * G4UniformRand();
430 G4double sinth = std::sin(theta);
431 G4double en = thePhotons->operator[](i)->GetTotalEnergy();
432 G4ThreeVector tempVector(en * sinth * std::cos(phi), en * sinth * std::sin(phi),
433 en * std::cos(theta));
434 thePhotons->operator[](i)->SetMomentum(tempVector);
435 }
436 else {
437 // tabulation of probabilities.
438 G4int it(0);
439 for (iii = 0; iii < nNeu[ii - nIso]; ++iii) // find the neutron energy
440 {
441 it = iii;
442 if (theAngular[ii - nIso][iii].GetEnergy() > anEnergy) break;
443 }
444 G4double costh = theAngular[ii - nIso][it].GetCosTh(); // no interpolation yet @@
445 G4double theta = std::acos(costh);
446 G4double phi = twopi * G4UniformRand();
447 G4double sinth = std::sin(theta);
448 G4double en = thePhotons->operator[](i)->GetTotalEnergy();
449 G4ThreeVector tmpVector(en * sinth * std::cos(phi), en * sinth * std::sin(phi),
450 en * costh);
451 thePhotons->operator[](i)->SetMomentum(tmpVector);
452 }
453 }
454 }
455 }
456 else if (repFlag == 2) {
457 auto running = new G4double[nGammaEnergies];
458 running[0] = theTransitionProbabilities[0];
459 for (i = 1; i < nGammaEnergies; ++i) {
460 running[i] = running[i - 1] + theTransitionProbabilities[i];
461 }
462 G4double random = G4UniformRand();
463 G4int it = 0;
464 for (i = 0; i < nGammaEnergies; ++i) {
465 it = i;
466 if (random < running[i] / running[nGammaEnergies - 1]) break;
467 }
468 delete[] running;
469 G4double totalEnergy = theBaseEnergy - theLevelEnergies[it];
470 auto theOne = new G4ReactionProduct;
471 theOne->SetDefinition(G4Gamma::Gamma());
472 random = G4UniformRand();
473 if (theInternalConversionFlag == 2 && random > thePhotonTransitionFraction[it]) {
474 theOne->SetDefinition(G4Electron::Electron());
475 // Bug reported Chao Zhang ([email protected]), Dongming Mei([email protected]) Feb. 25,
476 // 2009 But never enter at least with G4NDL3.13
477 totalEnergy +=
478 G4Electron::Electron()->GetPDGMass(); // proposed correction: add this line for electron
479 }
480 theOne->SetTotalEnergy(totalEnergy);
481 if (isoFlag == 1) {
482 G4double costheta = 2. * G4UniformRand() - 1;
483 G4double theta = std::acos(costheta);
484 G4double phi = twopi * G4UniformRand();
485 G4double sinth = std::sin(theta);
486 // Bug reported Chao Zhang ([email protected]), Dongming Mei([email protected]) Feb. 25,
487 // 2009 G4double en = theOne->GetTotalEnergy();
488 G4double en = theOne->GetTotalMomentum();
489 // But never cause real effect at least with G4NDL3.13 TK
490 G4ThreeVector temp(en * sinth * std::cos(phi), en * sinth * std::sin(phi),
491 en * std::cos(theta));
492 theOne->SetMomentum(temp);
493 }
494 else {
495 G4double currentEnergy = theOne->GetTotalEnergy();
496 for (ii = 0; ii < nDiscrete2; ++ii) {
497 if (std::abs(currentEnergy - theGammas[ii]) < 0.1 * keV) break;
498 }
499 if (ii == nDiscrete2)
500 --ii; // fix for what seems an (file12 vs file 14) inconsistency found in the ENDF 7N14
501 // data. @@
502 if (ii < nIso) {
503 // Bug reported Chao Zhang ([email protected]), Dongming Mei([email protected]) Feb. 25,
504 // 2009
505 // isotropic distribution
506 // G4double theta = pi*G4UniformRand();
507 G4double theta = std::acos(2. * G4UniformRand() - 1.);
508 // But this is alos never cause real effect at least with G4NDL3.13 TK not repFlag == 2 AND
509 // isoFlag != 1
510 G4double phi = twopi * G4UniformRand();
511 G4double sinth = std::sin(theta);
512 // Bug reported Chao Zhang ([email protected]), Dongming Mei([email protected]) Feb. 25,
513 // 2009 G4double en = theOne->GetTotalEnergy();
514 G4double en = theOne->GetTotalMomentum();
515 // But never cause real effect at least with G4NDL3.13 TK
516 G4ThreeVector tempVector(en * sinth * std::cos(phi), en * sinth * std::sin(phi),
517 en * std::cos(theta));
518 theOne->SetMomentum(tempVector);
519 }
520 else if (tabulationType == 1) {
521 // legendre polynomials
522 G4int itt(0);
523 for (iii = 0; iii < nNeu[ii - nIso]; ++iii) // find the neutron energy
524 {
525 itt = iii;
526 if (theLegendre[ii - nIso][iii].GetEnergy() > anEnergy) break;
527 }
529 aStore.SetCoeff(1, &(theLegendre[ii - nIso][itt]));
530 // aStore.SetCoeff(0, &(theLegendre[ii-nIso][it-1]));
531 // TKDB 110512
532 if (itt > 0) {
533 aStore.SetCoeff(0, &(theLegendre[ii - nIso][itt - 1]));
534 }
535 else {
536 aStore.SetCoeff(0, &(theLegendre[ii - nIso][itt]));
537 }
538 G4double cosTh = aStore.SampleMax(anEnergy);
539 G4double theta = std::acos(cosTh);
540 G4double phi = twopi * G4UniformRand();
541 G4double sinth = std::sin(theta);
542 // Bug reported Chao Zhang ([email protected]), Dongming Mei([email protected]) Feb. 25,
543 // 2009 G4double en = theOne->GetTotalEnergy();
544 G4double en = theOne->GetTotalMomentum();
545 // But never cause real effect at least with G4NDL3.13 TK
546 G4ThreeVector tempVector(en * sinth * std::cos(phi), en * sinth * std::sin(phi),
547 en * std::cos(theta));
548 theOne->SetMomentum(tempVector);
549 }
550 else {
551 // tabulation of probabilities.
552 G4int itt(0);
553 for (iii = 0; iii < nNeu[ii - nIso]; ++iii) // find the neutron energy
554 {
555 itt = iii;
556 if (theAngular[ii - nIso][iii].GetEnergy() > anEnergy) break;
557 }
558 G4double costh = theAngular[ii - nIso][itt].GetCosTh(); // no interpolation yet @@
559 G4double theta = std::acos(costh);
560 G4double phi = twopi * G4UniformRand();
561 G4double sinth = std::sin(theta);
562 // Bug reported Chao Zhang ([email protected]), Dongming Mei([email protected]) Feb. 25,
563 // 2009 G4double en = theOne->GetTotalEnergy();
564 G4double en = theOne->GetTotalMomentum();
565 // But never cause real effect at least with G4NDL3.13 TK
566 G4ThreeVector tmpVector(en * sinth * std::cos(phi), en * sinth * std::sin(phi), en * costh);
567 theOne->SetMomentum(tmpVector);
568 }
569 }
570 thePhotons->push_back(theOne);
571 }
572 else if (repFlag == 0) {
573 if (thePartialXsec == nullptr) {
574 return thePhotons;
575 }
576
577 // Partial Case
578
579 auto theOne = new G4ReactionProduct;
580 theOne->SetDefinition(G4Gamma::Gamma());
581 thePhotons->push_back(theOne);
582
583 // Energy
584
585 G4double sum = 0.0;
586 std::vector<G4double> dif(nDiscrete, 0.0);
587 for (G4int j = 0; j < nDiscrete; ++j) {
588 G4double x = thePartialXsec[j].GetXsec(anEnergy); // x in barn
589 if (x > 0) {
590 sum += x;
591 }
592 dif[j] = sum;
593 }
594
595 G4double rand = G4UniformRand();
596
597 G4int iphoton = 0;
598 for (G4int j = 0; j < nDiscrete; ++j) {
599 G4double y = rand * sum;
600 if (dif[j] > y) {
601 iphoton = j;
602 break;
603 }
604 }
605
606 // Statistically suppress the photon according to reaction cross section
607 // Fix proposed by Artem Zontikov, Bug report #1824
608 if (theReactionXsec != nullptr) {
609 if (thePartialXsec[iphoton].GetXsec(anEnergy) / theReactionXsec->GetXsec(anEnergy)
610 < G4UniformRand())
611 {
612 delete thePhotons;
613 thePhotons = nullptr;
614 return thePhotons;
615 }
616 }
617
618 // Angle
619 G4double cosTheta = 0.0; // mu
620
621 if (isoFlag == 1) {
622 // Isotropic Case
623
624 cosTheta = 2. * G4UniformRand() - 1;
625 }
626 else {
627 if (iphoton < nIso) {
628 // still Isotropic
629
630 cosTheta = 2. * G4UniformRand() - 1;
631 }
632 else {
633 if (tabulationType == 1) {
634 // Legendre polynomials
635
636 G4int iangle = 0;
637 for (G4int j = 0; j < nNeu[iphoton - nIso]; ++j) {
638 iangle = j;
639 if (theLegendre[iphoton - nIso][j].GetEnergy() > anEnergy) break;
640 }
641
643 aStore.SetCoeff(1, &(theLegendre[iphoton - nIso][iangle]));
644 aStore.SetCoeff(0, &(theLegendre[iphoton - nIso][iangle - 1]));
645
646 cosTheta = aStore.SampleMax(anEnergy);
647 }
648 else if (tabulationType == 2) {
649 // tabulation of probabilities.
650
651 G4int iangle = 0;
652 for (G4int j = 0; j < nNeu[iphoton - nIso]; ++j) {
653 iangle = j;
654 if (theAngular[iphoton - nIso][j].GetEnergy() > anEnergy) break;
655 }
656 cosTheta = theAngular[iphoton - nIso][iangle].GetCosTh();
657 // no interpolation yet @@
658 }
659 }
660 }
661
662 // Set
663 G4double phi = twopi * G4UniformRand();
664 G4double theta = std::acos(cosTheta);
665 G4double sinTheta = std::sin(theta);
666
667 G4double photonE = theGammas[iphoton];
668 G4ThreeVector direction(sinTheta * std::cos(phi), sinTheta * std::sin(phi), cosTheta);
669 G4ThreeVector photonP = photonE * direction;
670 thePhotons->operator[](0)->SetMomentum(photonP);
671 }
672 else {
673 delete thePhotons;
674 thePhotons = nullptr; // no gamma data available; some work needed @@@@@@@
675 }
676 return thePhotons;
677}
G4long G4Poisson(G4double mean)
Definition G4Poisson.hh:50
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
value_type & Get() const
Definition G4Cache.hh:315
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
void Init(G4int aScheme, G4int aRange)
void Init(std::istream &aDataFile)
void SetCoeff(G4int i, G4int l, G4double coeff)
G4double SampleMax(G4double energy)
void Init(std::istream &aDataFile)
void Init(std::istream &aDataFile)
G4double GetY(G4int i, G4int j)
void InitInterpolation(G4int i, std::istream &aDataFile)
void InitEnergies(std::istream &aDataFile)
G4ReactionProductVector * GetPhotons(G4double anEnergy)
void InitAngular(std::istream &aDataFile)
void InitPartials(std::istream &aDataFile, G4ParticleHPVector *theXsec=nullptr)
G4bool InitMean(std::istream &aDataFile)
G4double GetXsec(G4int i)
G4double GetY(G4double x)
G4double GetX(G4int i) const
void Init(std::istream &aDataFile, G4int total, G4double ux=1., G4double uy=1.)
G4int GetVectorLength() const
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
#define DBL_MAX
Definition templates.hh:62