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